source: trunk/source/persistency/gdml/src/G4GDMLReadMaterials.cc@ 1315

Last change on this file since 1315 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 20.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id: G4GDMLReadMaterials.cc,v 1.20 2009/04/24 15:34:20 gcosmo Exp $
28// GEANT4 tag $ Name:$
29//
30// class G4GDMLReadMaterials Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLReadMaterials.hh"
37
38#include "G4Element.hh"
39#include "G4Isotope.hh"
40#include "G4Material.hh"
41#include "G4NistManager.hh"
42
43G4GDMLReadMaterials::G4GDMLReadMaterials() : G4GDMLReadDefine()
44{
45}
46
47G4GDMLReadMaterials::~G4GDMLReadMaterials()
48{
49}
50
51G4double
52G4GDMLReadMaterials::AtomRead(const xercesc::DOMElement* const atomElement)
53{
54 G4double value = 0.0;
55 G4double unit = g/mole;
56
57 const xercesc::DOMNamedNodeMap* const attributes
58 = atomElement->getAttributes();
59 XMLSize_t attributeCount = attributes->getLength();
60
61 for (XMLSize_t attribute_index=0;
62 attribute_index<attributeCount; attribute_index++)
63 {
64 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
65
66 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
67 { continue; }
68
69 const xercesc::DOMAttr* const attribute
70 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
71 const G4String attName = Transcode(attribute->getName());
72 const G4String attValue = Transcode(attribute->getValue());
73
74 if (attName=="value") { value = eval.Evaluate(attValue); } else
75 if (attName=="unit") { unit = eval.Evaluate(attValue); }
76 }
77
78 return value*unit;
79}
80
81G4int G4GDMLReadMaterials::
82CompositeRead(const xercesc::DOMElement* const compositeElement,G4String& ref)
83{
84 G4int n = 0;
85
86 const xercesc::DOMNamedNodeMap* const attributes
87 = compositeElement->getAttributes();
88 XMLSize_t attributeCount = attributes->getLength();
89
90 for (XMLSize_t attribute_index=0;
91 attribute_index<attributeCount; attribute_index++)
92 {
93 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
94
95 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
96 { continue; }
97
98 const xercesc::DOMAttr* const attribute
99 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
100 const G4String attName = Transcode(attribute->getName());
101 const G4String attValue = Transcode(attribute->getValue());
102
103 if (attName=="n") { n = eval.EvaluateInteger(attValue); } else
104 if (attName=="ref") { ref = attValue; }
105 }
106
107 return n;
108}
109
110G4double G4GDMLReadMaterials::DRead(const xercesc::DOMElement* const DElement)
111{
112 G4double value = 0.0;
113 G4double unit = g/cm3;
114
115 const xercesc::DOMNamedNodeMap* const attributes
116 = DElement->getAttributes();
117 XMLSize_t attributeCount = attributes->getLength();
118
119 for (XMLSize_t attribute_index=0;
120 attribute_index<attributeCount; attribute_index++)
121 {
122 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
123
124 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
125 { continue; }
126
127 const xercesc::DOMAttr* const attribute
128 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
129 const G4String attName = Transcode(attribute->getName());
130 const G4String attValue = Transcode(attribute->getValue());
131
132 if (attName=="value") { value = eval.Evaluate(attValue); } else
133 if (attName=="unit") { unit = eval.Evaluate(attValue); }
134 }
135
136 return value*unit;
137}
138
139G4double G4GDMLReadMaterials::PRead(const xercesc::DOMElement* const PElement)
140{
141 G4double value = STP_Pressure;
142 G4double unit = pascal;
143
144 const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
145 XMLSize_t attributeCount = attributes->getLength();
146
147 for (XMLSize_t attribute_index=0;
148 attribute_index<attributeCount; attribute_index++)
149 {
150 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
151
152 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
153 { continue; }
154
155 const xercesc::DOMAttr* const attribute
156 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
157 const G4String attName = Transcode(attribute->getName());
158 const G4String attValue = Transcode(attribute->getValue());
159
160 if (attName=="value") { value = eval.Evaluate(attValue); } else
161 if (attName=="unit") { unit = eval.Evaluate(attValue); }
162 }
163
164 return value*unit;
165}
166
167G4double G4GDMLReadMaterials::TRead(const xercesc::DOMElement* const TElement)
168{
169 G4double value = STP_Temperature;
170 G4double unit = kelvin;
171
172 const xercesc::DOMNamedNodeMap* const attributes = TElement->getAttributes();
173 XMLSize_t attributeCount = attributes->getLength();
174
175 for (XMLSize_t attribute_index=0;
176 attribute_index<attributeCount; attribute_index++)
177 {
178 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
179
180 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
181 { continue; }
182
183 const xercesc::DOMAttr* const attribute
184 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
185 const G4String attName = Transcode(attribute->getName());
186 const G4String attValue = Transcode(attribute->getValue());
187
188 if (attName=="value") { value = eval.Evaluate(attValue); } else
189 if (attName=="unit") { unit = eval.Evaluate(attValue); }
190 }
191
192 return value*unit;
193}
194
195void G4GDMLReadMaterials::
196ElementRead(const xercesc::DOMElement* const elementElement)
197{
198 G4String name;
199 G4String formula;
200 G4double a = 0.0;
201 G4double Z = 0.0;
202
203 const xercesc::DOMNamedNodeMap* const attributes
204 = elementElement->getAttributes();
205 XMLSize_t attributeCount = attributes->getLength();
206
207 for (XMLSize_t attribute_index=0;
208 attribute_index<attributeCount; attribute_index++)
209 {
210 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
211
212 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
213 { continue; }
214
215 const xercesc::DOMAttr* const attribute
216 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
217 const G4String attName = Transcode(attribute->getName());
218 const G4String attValue = Transcode(attribute->getValue());
219
220 if (attName=="name") { name = GenerateName(attValue); } else
221 if (attName=="formula") { formula = attValue; } else
222 if (attName=="Z") { Z = eval.Evaluate(attValue); }
223 }
224
225 G4int nComponents = 0;
226
227 for (xercesc::DOMNode* iter = elementElement->getFirstChild();
228 iter != 0; iter = iter->getNextSibling())
229 {
230 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
231
232 const xercesc::DOMElement* const child
233 = dynamic_cast<xercesc::DOMElement*>(iter);
234 const G4String tag = Transcode(child->getTagName());
235
236 if (tag=="atom") { a = AtomRead(child); } else
237 if (tag=="fraction") { nComponents++; }
238 }
239
240 if (nComponents>0)
241 {
242 MixtureRead(elementElement,
243 new G4Element(Strip(name),formula,nComponents));
244 }
245 else
246 {
247 new G4Element(Strip(name),formula,Z,a);
248 }
249}
250
251G4double G4GDMLReadMaterials::
252FractionRead(const xercesc::DOMElement* const fractionElement, G4String& ref)
253{
254 G4double n = 0.0;
255
256 const xercesc::DOMNamedNodeMap* const attributes
257 = fractionElement->getAttributes();
258 XMLSize_t attributeCount = attributes->getLength();
259
260 for (XMLSize_t attribute_index=0;
261 attribute_index<attributeCount; attribute_index++)
262 {
263 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
264
265 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
266 { continue; }
267
268 const xercesc::DOMAttr* const attribute
269 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
270 const G4String attName = Transcode(attribute->getName());
271 const G4String attValue = Transcode(attribute->getValue());
272
273 if (attName=="n") { n = eval.Evaluate(attValue); } else
274 if (attName=="ref") { ref = attValue; }
275 }
276
277 return n;
278}
279
280void G4GDMLReadMaterials::
281IsotopeRead(const xercesc::DOMElement* const isotopeElement)
282{
283 G4String name;
284 G4int Z = 0;
285 G4int N = 0;
286 G4double a = 0.0;
287
288 const xercesc::DOMNamedNodeMap* const attributes
289 = isotopeElement->getAttributes();
290 XMLSize_t attributeCount = attributes->getLength();
291
292 for (XMLSize_t attribute_index=0;
293 attribute_index<attributeCount;attribute_index++)
294 {
295 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
296
297 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
298 { continue; }
299
300 const xercesc::DOMAttr* const attribute
301 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
302 const G4String attName = Transcode(attribute->getName());
303 const G4String attValue = Transcode(attribute->getValue());
304
305 if (attName=="name") { name = GenerateName(attValue); } else
306 if (attName=="Z") { Z = eval.EvaluateInteger(attValue); } else
307 if (attName=="N") { N = eval.EvaluateInteger(attValue); }
308 }
309
310 for (xercesc::DOMNode* iter = isotopeElement->getFirstChild();
311 iter != 0; iter = iter->getNextSibling())
312 {
313 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
314
315 const xercesc::DOMElement* const child
316 = dynamic_cast<xercesc::DOMElement*>(iter);
317 const G4String tag = Transcode(child->getTagName());
318
319 if (tag=="atom") { a = AtomRead(child); }
320 }
321
322 new G4Isotope(Strip(name),Z,N,a);
323}
324
325void G4GDMLReadMaterials::
326MaterialRead(const xercesc::DOMElement* const materialElement)
327{
328 G4String name;
329 G4double Z = 0.0;
330 G4double a = 0.0;
331 G4double D = 0.0;
332 G4State state = kStateUndefined;
333 G4double T = STP_Temperature;
334 G4double P = STP_Pressure;
335
336 const xercesc::DOMNamedNodeMap* const attributes
337 = materialElement->getAttributes();
338 XMLSize_t attributeCount = attributes->getLength();
339
340 for (XMLSize_t attribute_index=0;
341 attribute_index<attributeCount; attribute_index++)
342 {
343 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
344
345 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
346 { continue; }
347
348 const xercesc::DOMAttr* const attribute
349 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
350 const G4String attName = Transcode(attribute->getName());
351 const G4String attValue = Transcode(attribute->getValue());
352
353 if (attName=="name") { name = GenerateName(attValue); } else
354 if (attName=="Z") { Z = eval.Evaluate(attValue); } else
355 if (attName=="state")
356 {
357 if (attValue=="solid") { state = kStateSolid; } else
358 if (attValue=="liquid") { state = kStateLiquid; } else
359 if (attValue=="gas") { state = kStateGas; }
360 }
361 }
362
363 size_t nComponents = 0;
364
365 for (xercesc::DOMNode* iter = materialElement->getFirstChild();
366 iter != 0; iter = iter->getNextSibling())
367 {
368 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
369
370 const xercesc::DOMElement* const child
371 = dynamic_cast<xercesc::DOMElement*>(iter);
372 const G4String tag = Transcode(child->getTagName());
373
374 if (tag=="atom") { a = AtomRead(child); } else
375 if (tag=="Dref") { D = GetQuantity(GenerateName(RefRead(child))); } else
376 if (tag=="Pref") { P = GetQuantity(GenerateName(RefRead(child))); } else
377 if (tag=="Tref") { T = GetQuantity(GenerateName(RefRead(child))); } else
378 if (tag=="D") { D = DRead(child); } else
379 if (tag=="P") { P = PRead(child); } else
380 if (tag=="T") { T = TRead(child); } else
381 if (tag=="fraction" || tag=="composite") { nComponents++; }
382 }
383
384 G4Material* material = 0;
385
386 if (nComponents==0)
387 {
388 material = new G4Material(Strip(name),Z,a,D,state,T,P);
389 }
390 else
391 {
392 material = new G4Material(Strip(name),D,nComponents,state,T,P);
393 MixtureRead(materialElement, material);
394 }
395
396 for (xercesc::DOMNode* iter = materialElement->getFirstChild();
397 iter != 0; iter = iter->getNextSibling())
398 {
399 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
400
401 const xercesc::DOMElement* const child
402 = dynamic_cast<xercesc::DOMElement*>(iter);
403 const G4String tag = Transcode(child->getTagName());
404
405 if (tag=="property") { PropertyRead(child,material); }
406 }
407}
408
409void G4GDMLReadMaterials::
410MixtureRead(const xercesc::DOMElement *const mixtureElement, G4Element *element)
411{
412 for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
413 iter != 0; iter = iter->getNextSibling())
414 {
415 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
416
417 const xercesc::DOMElement* const child
418 = dynamic_cast<xercesc::DOMElement*>(iter);
419 const G4String tag = Transcode(child->getTagName());
420
421 if (tag=="fraction")
422 {
423 G4String ref;
424 G4double n = FractionRead(child,ref);
425 element->AddIsotope(GetIsotope(GenerateName(ref,true)),n);
426 }
427 }
428}
429
430void G4GDMLReadMaterials::
431MixtureRead(const xercesc::DOMElement *const mixtureElement,
432 G4Material *material)
433{
434 for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
435 iter != 0; iter = iter->getNextSibling())
436 {
437 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
438
439 const xercesc::DOMElement* const child
440 = dynamic_cast<xercesc::DOMElement*>(iter);
441 const G4String tag = Transcode(child->getTagName());
442
443 if (tag=="fraction")
444 {
445 G4String ref;
446 G4double n = FractionRead(child,ref);
447
448 G4Material *materialPtr = GetMaterial(GenerateName(ref,true), false);
449 G4Element *elementPtr = GetElement(GenerateName(ref,true), false);
450
451 if (materialPtr != 0) { material->AddMaterial(materialPtr,n); } else
452 if (elementPtr != 0) { material->AddElement(elementPtr,n); }
453
454 if ((materialPtr == 0) && (elementPtr == 0))
455 {
456 G4String error_msg = "Referenced material/element '"
457 + GenerateName(ref,true) + "' was not found!";
458 G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidSetup",
459 FatalException, error_msg);
460 }
461 }
462 else if (tag=="composite")
463 {
464 G4String ref;
465 G4int n = CompositeRead(child,ref);
466
467 G4Element *elementPtr = GetElement(GenerateName(ref,true));
468 material->AddElement(elementPtr,n);
469 }
470 }
471}
472
473void G4GDMLReadMaterials::
474PropertyRead(const xercesc::DOMElement* const propertyElement,
475 G4Material* material)
476{
477 G4String name;
478 G4String ref;
479 G4GDMLMatrix matrix;
480
481 const xercesc::DOMNamedNodeMap* const attributes
482 = propertyElement->getAttributes();
483 XMLSize_t attributeCount = attributes->getLength();
484
485 for (XMLSize_t attribute_index=0;
486 attribute_index<attributeCount; attribute_index++)
487 {
488 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
489
490 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
491 { continue; }
492
493 const xercesc::DOMAttr* const attribute
494 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
495 const G4String attName = Transcode(attribute->getName());
496 const G4String attValue = Transcode(attribute->getValue());
497
498 if (attName=="name") { name = GenerateName(attValue); } else
499 if (attName=="ref") { matrix = GetMatrix(ref=attValue); }
500 }
501
502 if (matrix.GetCols() != 2)
503 {
504 G4String error_msg = "Referenced matrix '" + ref
505 + "' should have \n two columns as a property table for material: "
506 + material->GetName();
507 G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
508 FatalException, error_msg);
509 }
510 if (matrix.GetRows() == 0) { return; }
511
512 G4MaterialPropertiesTable* matprop=material->GetMaterialPropertiesTable();
513 if (!matprop)
514 {
515 matprop = new G4MaterialPropertiesTable();
516 material->SetMaterialPropertiesTable(matprop);
517 }
518 if (matrix.GetCols() == 1) // constant property assumed
519 {
520 matprop->AddConstProperty(Strip(name), matrix.Get(0,0));
521 }
522 else // build the material properties vector
523 {
524 G4MaterialPropertyVector* propvect = new G4MaterialPropertyVector(0,0,0);
525 for (size_t i=0; i<matrix.GetRows(); i++)
526 {
527 propvect->AddElement(matrix.Get(i,0),matrix.Get(i,1));
528 }
529 matprop->AddProperty(Strip(name),propvect);
530 }
531}
532
533void G4GDMLReadMaterials::
534MaterialsRead(const xercesc::DOMElement* const materialsElement)
535{
536 G4cout << "G4GDML: Reading materials..." << G4endl;
537
538 for (xercesc::DOMNode* iter = materialsElement->getFirstChild();
539 iter != 0; iter = iter->getNextSibling())
540 {
541 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
542
543 const xercesc::DOMElement* const child
544 = dynamic_cast<xercesc::DOMElement*>(iter);
545 const G4String tag = Transcode(child->getTagName());
546
547 if (tag=="define") { DefineRead(child); } else
548 if (tag=="element") { ElementRead(child); } else
549 if (tag=="isotope") { IsotopeRead(child); } else
550 if (tag=="material") { MaterialRead(child); }
551 else
552 {
553 G4String error_msg = "Unknown tag in materials: " + tag;
554 G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidSetup",
555 FatalException, error_msg);
556 }
557 }
558}
559
560G4Element* G4GDMLReadMaterials::
561GetElement(const G4String& ref, G4bool verbose) const
562{
563 G4Element* elementPtr = G4Element::GetElement(ref,false);
564
565 if (!elementPtr)
566 {
567 elementPtr = G4NistManager::Instance()->FindOrBuildElement(ref);
568 }
569
570 if (verbose && !elementPtr)
571 {
572 G4String error_msg = "Referenced element '" + ref + "' was not found!";
573 G4Exception("G4GDMLReadMaterials::GetElement()", "InvalidRead",
574 FatalException, error_msg);
575 }
576
577 return elementPtr;
578}
579
580G4Isotope* G4GDMLReadMaterials::GetIsotope(const G4String& ref,
581 G4bool verbose) const
582{
583 G4Isotope* isotopePtr = G4Isotope::GetIsotope(ref,false);
584
585 if (verbose && !isotopePtr)
586 {
587 G4String error_msg = "Referenced isotope '" + ref + "' was not found!";
588 G4Exception("G4GDMLReadMaterials::GetIsotope()", "InvalidRead",
589 FatalException, error_msg);
590 }
591
592 return isotopePtr;
593}
594
595G4Material* G4GDMLReadMaterials::GetMaterial(const G4String& ref,
596 G4bool verbose) const
597{
598 G4Material *materialPtr = G4Material::GetMaterial(ref,false);
599
600 if (!materialPtr)
601 {
602 materialPtr = G4NistManager::Instance()->FindOrBuildMaterial(ref);
603 }
604
605 if (verbose && !materialPtr)
606 {
607 G4String error_msg = "Referenced material '" + ref + "' was not found!";
608 G4Exception("G4GDMLReadMaterials::GetMaterial()", "InvalidRead",
609 FatalException, error_msg);
610 }
611
612 return materialPtr;
613}
Note: See TracBrowser for help on using the repository browser.