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

Last change on this file since 1199 was 987, checked in by garnier, 17 years ago

fichiers manquants

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