source: trunk/source/materials/src/G4Material.cc@ 1250

Last change on this file since 1250 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 22.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: G4Material.cc,v 1.42 2008/08/13 16:06:42 vnivanch Exp $
28// GEANT4 tag $Name: materials-V09-02-18 $
29//
30//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31//
32// 26-06-96, Code uses operators (+=, *=, ++, -> etc.) correctly, P. Urban
33// 10-07-96, new data members added by L.Urban
34// 12-12-96, new data members added by L.Urban
35// 20-01-97, aesthetic rearrangement. RadLength calculation modified.
36// Data members Zeff and Aeff REMOVED (i.e. passed to the Elements).
37// (local definition of Zeff in DensityEffect and FluctModel...)
38// Vacuum defined as a G4State. Mixture flag removed, M.Maire.
39// 29-01-97, State=Vacuum automatically set density=0 in the contructors.
40// Subsequent protections have been put in the calculation of
41// MeanExcEnergy, ShellCorrectionVector, DensityEffect, M.Maire.
42// 11-02-97, ComputeDensityEffect() rearranged, M.Maire.
43// 20-03-97, corrected initialization of pointers, M.Maire.
44// 28-05-98, the kState=kVacuum has been removed.
45// automatic check for a minimal density, M.Maire
46// 12-06-98, new method AddMaterial() allowing mixture of materials, M.Maire
47// 09-07-98, ionisation parameters removed from the class, M.Maire
48// 05-10-98, change names: NumDensity -> NbOfAtomsPerVolume
49// 18-11-98, new interface to SandiaTable
50// 19-01-99 enlarge tolerance on test of coherence of gas conditions
51// 19-07-99, Constructors with chemicalFormula added by V.Ivanchenko
52// 16-01-01, Nuclear interaction length, M.Maire
53// 12-03-01, G4bool fImplicitElement;
54// copy constructor and assignement operator revised (mma)
55// 03-05-01, flux.precision(prec) at begin/end of operator<<
56// 17-07-01, migration to STL. M. Verderi.
57// 14-09-01, Suppression of the data member fIndexInTable
58// 26-02-02, fIndexInTable renewed
59// 16-04-02, G4Exception put in constructor with chemical formula
60// 06-05-02, remove the check of the ideal gas state equation
61// 06-08-02, remove constructors with chemical formula (mma)
62// 22-01-04, proper STL handling of theElementVector (Hisaya)
63// 30-03-05, warning in GetMaterial(materialName)
64// 09-03-06, minor change of printout (V.Ivanchenko)
65// 10-01-07, compute fAtomVector in the case of mass fraction (V.Ivanchenko)
66// 27-07-07, improve destructor (V.Ivanchenko)
67// 18-10-07, move definition of material index to InitialisePointers (V.Ivanchenko)
68// 13-08-08, do not use fixed size arrays (V.Ivanchenko)
69//
70
71//
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74#include "G4Material.hh"
75#include "G4UnitsTable.hh"
76#include <iomanip>
77
78
79G4MaterialTable G4Material::theMaterialTable;
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
83// Constructor to create a material from scratch
84
85G4Material::G4Material(const G4String& name, G4double z,
86 G4double a, G4double density,
87 G4State state, G4double temp, G4double pressure)
88 : fName(name)
89{
90 InitializePointers();
91
92 if (density < universe_mean_density)
93 { G4cerr << "--- Warning from G4Material::G4Material()"
94 << " define a material with density=0 is not allowed. \n"
95 << " The material " << name << " will be constructed with the"
96 << " default minimal density: " << universe_mean_density/(g/cm3)
97 << "g/cm3" << G4endl;
98 density = universe_mean_density;
99 }
100
101 fDensity = density;
102 fState = state;
103 fTemp = temp;
104 fPressure = pressure;
105 fChemicalFormula = " ";
106
107 // Initialize theElementVector allocating one
108 // element corresponding to this material
109 maxNbComponents = fNumberOfComponents = fNumberOfElements = 1;
110 fArrayLength = maxNbComponents;
111 fImplicitElement = true;
112 theElementVector = new G4ElementVector();
113 theElementVector->push_back( new G4Element(name, " ", z, a));
114 fMassFractionVector = new G4double[1];
115 fMassFractionVector[0] = 1. ;
116
117 (*theElementVector)[0] -> increaseCountUse();
118
119 if (fState == kStateUndefined)
120 {
121 if (fDensity > kGasThreshold) fState = kStateSolid;
122 else fState = kStateGas;
123 }
124
125 ComputeDerivedQuantities();
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129
130// Constructor to create a material from a List of constituents
131// (elements and/or materials) added with AddElement or AddMaterial
132
133G4Material::G4Material(const G4String& name, G4double density,
134 G4int nComponents,
135 G4State state, G4double temp, G4double pressure)
136 : fName(name)
137{
138 InitializePointers();
139
140 if (density < universe_mean_density)
141 {G4cerr << "--- Warning from G4Material::G4Material()"
142 << " define a material with density=0 is not allowed. \n"
143 << " The material " << name << " will be constructed with the"
144 << " default minimal density: " << universe_mean_density/(g/cm3)
145 << "g/cm3" << G4endl;
146 density = universe_mean_density;
147 }
148
149 fDensity = density;
150 fState = state;
151 fTemp = temp;
152 fPressure = pressure;
153 fChemicalFormula = " ";
154
155 maxNbComponents = nComponents;
156 fArrayLength = maxNbComponents;
157 fNumberOfComponents = fNumberOfElements = 0;
158 fImplicitElement = false;
159 theElementVector = new G4ElementVector();
160 theElementVector->reserve(maxNbComponents);
161
162 if (fState == kStateUndefined)
163 {
164 if (fDensity > kGasThreshold) fState = kStateSolid;
165 else fState = kStateGas;
166 }
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
171// Fake default constructor - sets only member data and allocates memory
172// for usage restricted to object persistency
173
174G4Material::G4Material(__void__&)
175 : fNumberOfComponents(0), fNumberOfElements(0), theElementVector(0),
176 fImplicitElement(false), fMassFractionVector(0), fAtomsVector(0),
177 fMaterialPropertiesTable(0), fIndexInTable(0),
178 VecNbOfAtomsPerVolume(0), fIonisation(0), fSandiaTable(0)
179{
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
184// AddElement -- composition by atom count
185
186void G4Material::AddElement(G4Element* element, G4int nAtoms)
187{
188 // initialization
189 if ( fNumberOfElements == 0 ) {
190 fAtomsVector = new G4int [fArrayLength];
191 fMassFractionVector = new G4double[fArrayLength];
192 }
193
194 // filling ...
195 if ( G4int(fNumberOfElements) < maxNbComponents ) {
196 theElementVector->push_back(element);
197 fAtomsVector [fNumberOfElements] = nAtoms;
198 fNumberOfComponents = ++fNumberOfElements;
199 element->increaseCountUse();
200 } else {
201 G4cerr << "G4Material::AddElement ERROR for " << fName << " nElement= "
202 << fNumberOfElements << G4endl;
203 G4Exception
204 ("ERROR!!! - Attempt to add more than the declared number of elements.");
205 }
206 // filled.
207 if ( G4int(fNumberOfElements) == maxNbComponents ) {
208 // compute proportion by mass
209 size_t i=0;
210 G4double Zmol(0.), Amol(0.);
211 for (i=0;i<fNumberOfElements;i++) {
212 Zmol += fAtomsVector[i]*(*theElementVector)[i]->GetZ();
213 Amol += fAtomsVector[i]*(*theElementVector)[i]->GetA();
214 }
215 for (i=0;i<fNumberOfElements;i++) {
216 fMassFractionVector[i] = fAtomsVector[i]
217 *(*theElementVector)[i]->GetA()/Amol;
218 }
219
220 ComputeDerivedQuantities();
221 }
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225
226// AddElement -- composition by fraction of mass
227
228void G4Material::AddElement(G4Element* element, G4double fraction)
229{
230 // initialization
231 if (fNumberOfComponents == 0) {
232 fMassFractionVector = new G4double[fArrayLength];
233 fAtomsVector = new G4int [fArrayLength];
234 }
235
236 // filling ...
237 if (G4int(fNumberOfComponents) < maxNbComponents) {
238 size_t el = 0;
239 while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) el++;
240 if (el<fNumberOfElements) fMassFractionVector[el] += fraction;
241 else {
242 theElementVector->push_back(element);
243 fMassFractionVector[el] = fraction;
244 fNumberOfElements++;
245 element->increaseCountUse();
246 }
247 fNumberOfComponents++;
248 } else {
249 G4cerr << "G4Material::AddElement ERROR for " << fName << " nElement= "
250 << fNumberOfElements << G4endl;
251 G4Exception
252 ("ERROR!!! - Attempt to add more than the declared number of components.");
253 }
254
255 // filled.
256 if (G4int(fNumberOfComponents) == maxNbComponents) {
257
258 size_t i=0;
259 G4double Zmol(0.), Amol(0.);
260 // check sum of weights -- OK?
261 G4double wtSum(0.0);
262 for (i=0;i<fNumberOfElements;i++) {
263 wtSum += fMassFractionVector[i];
264 Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
265 Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
266 }
267 if (std::abs(1.-wtSum) > perThousand) {
268 G4cerr << "WARNING !! for " << fName << " sum of fractional masses "
269 << wtSum << " is not 1 - results may be wrong"
270 << G4endl;
271 }
272 for (i=0;i<fNumberOfElements;i++) {
273 fAtomsVector[i] =
274 G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
275 }
276
277 ComputeDerivedQuantities();
278 }
279}
280
281//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
282
283// AddMaterial -- composition by fraction of mass
284
285void G4Material::AddMaterial(G4Material* material, G4double fraction)
286{
287 // initialization
288 if (fNumberOfComponents == 0) {
289 fMassFractionVector = new G4double[fArrayLength];
290 fAtomsVector = new G4int [fArrayLength];
291 }
292
293 size_t nelm = material->GetNumberOfElements();
294
295 // arrays should be extended
296 if(nelm > 1) {
297 G4int nold = fArrayLength;
298 fArrayLength += nelm - 1;
299 G4double* v1 = new G4double[fArrayLength];
300 G4int* i1 = new G4int[fArrayLength];
301 for(G4int i=0; i<nold; i++) {
302 v1[i] = fMassFractionVector[i];
303 i1[i] = fAtomsVector[i];
304 }
305 delete [] fAtomsVector;
306 delete [] fMassFractionVector;
307 fMassFractionVector = v1;
308 fAtomsVector = i1;
309 }
310
311 // filling ...
312 if (G4int(fNumberOfComponents) < maxNbComponents) {
313 for (size_t elm=0; elm<nelm; elm++)
314 {
315 G4Element* element = (*(material->GetElementVector()))[elm];
316 size_t el = 0;
317 while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) el++;
318 if (el < fNumberOfElements) fMassFractionVector[el] += fraction
319 *(material->GetFractionVector())[elm];
320 else {
321 theElementVector->push_back(element);
322 fMassFractionVector[el] = fraction
323 *(material->GetFractionVector())[elm];
324 fNumberOfElements++;
325 element->increaseCountUse();
326 }
327 }
328 fNumberOfComponents++;
329 } else {
330 G4cerr << "G4Material::AddElement ERROR for " << fName << " nElement= "
331 << fNumberOfElements << G4endl;
332 G4Exception
333 ("ERROR!!! - Attempt to add more than the declared number of components.");
334 }
335
336 // filled.
337 if (G4int(fNumberOfComponents) == maxNbComponents) {
338 size_t i=0;
339 G4double Zmol(0.), Amol(0.);
340 // check sum of weights -- OK?
341 G4double wtSum(0.0);
342 for (i=0;i<fNumberOfElements;i++) {
343 wtSum += fMassFractionVector[i];
344 Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
345 Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
346 }
347 if (std::abs(1.-wtSum) > perThousand) {
348 G4cerr << "WARNING !! for " << fName << " sum of fractional masses "
349 << wtSum << " is not 1 - results may be wrong"
350 << G4endl;
351 }
352 for (i=0;i<fNumberOfElements;i++) {
353 fAtomsVector[i] =
354 G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
355 }
356
357 ComputeDerivedQuantities();
358 }
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
362
363void G4Material::ComputeDerivedQuantities()
364{
365 // Header routine to compute various properties of material.
366 //
367
368 // Number of atoms per volume (per element), total nb of electrons per volume
369 G4double Zi, Ai;
370 TotNbOfAtomsPerVolume = 0.;
371 if (VecNbOfAtomsPerVolume) delete [] VecNbOfAtomsPerVolume;
372 VecNbOfAtomsPerVolume = new G4double[fNumberOfElements];
373 TotNbOfElectPerVolume = 0.;
374 for (size_t i=0;i<fNumberOfElements;i++) {
375 Zi = (*theElementVector)[i]->GetZ();
376 Ai = (*theElementVector)[i]->GetA();
377 VecNbOfAtomsPerVolume[i] = Avogadro*fDensity*fMassFractionVector[i]/Ai;
378 TotNbOfAtomsPerVolume += VecNbOfAtomsPerVolume[i];
379 TotNbOfElectPerVolume += VecNbOfAtomsPerVolume[i]*Zi;
380 }
381
382 ComputeRadiationLength();
383 ComputeNuclearInterLength();
384
385 if (fIonisation) delete fIonisation;
386 fIonisation = new G4IonisParamMat(this);
387 if (fSandiaTable) delete fSandiaTable;
388 fSandiaTable = new G4SandiaTable(this);
389}
390
391//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
392
393void G4Material::ComputeRadiationLength()
394{
395 G4double radinv = 0.0 ;
396 for (size_t i=0;i<fNumberOfElements;i++) {
397 radinv += VecNbOfAtomsPerVolume[i]*((*theElementVector)[i]->GetfRadTsai());
398 }
399 fRadlen = (radinv <= 0.0 ? DBL_MAX : 1./radinv);
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403
404void G4Material::ComputeNuclearInterLength()
405{
406 const G4double lambda0 = 35*g/cm2;
407 G4double NILinv = 0.0;
408 for (size_t i=0;i<fNumberOfElements;i++) {
409 NILinv +=
410 VecNbOfAtomsPerVolume[i]*std::pow(((*theElementVector)[i]->GetN()),0.6666667);
411 }
412 NILinv *= amu/lambda0;
413 fNuclInterLen = (NILinv <= 0.0 ? DBL_MAX : 1./NILinv);
414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
417
418void G4Material::InitializePointers()
419{
420 theElementVector = 0;
421 fMassFractionVector = 0;
422 fAtomsVector = 0;
423 fMaterialPropertiesTable = 0;
424
425 VecNbOfAtomsPerVolume = 0;
426 fIonisation = 0;
427 fSandiaTable = 0;
428
429 // Store in the static Table of Materials
430 theMaterialTable.push_back(this);
431 fIndexInTable = theMaterialTable.size() - 1;
432}
433
434//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
435
436const G4MaterialTable* G4Material::GetMaterialTable()
437{
438 return &theMaterialTable;
439}
440
441//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
442
443size_t G4Material::GetNumberOfMaterials()
444{
445 return theMaterialTable.size();
446}
447
448//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
449
450G4Material* G4Material::GetMaterial(G4String materialName, G4bool warning)
451{
452 // search the material by its name
453 for (size_t J=0 ; J<theMaterialTable.size() ; J++)
454 {
455 if (theMaterialTable[J]->GetName() == materialName)
456 return theMaterialTable[J];
457 }
458
459 // the material does not exist in the table
460 if (warning) {
461 G4cout << "\n---> warning from G4Material::GetMaterial(). The material: "
462 << materialName << " does not exist in the table. Return NULL pointer."
463 << G4endl;
464 }
465 return 0;
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
469
470G4Material::G4Material(const G4Material& right)
471{
472 InitializePointers();
473 *this = right;
474}
475
476//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
477
478G4Material::~G4Material()
479{
480 // G4cout << "### Destruction of material " << fName << " started" <<G4endl;
481 if (theElementVector) delete theElementVector;
482 if (fMassFractionVector) delete [] fMassFractionVector;
483 if (fAtomsVector) delete [] fAtomsVector;
484 if (VecNbOfAtomsPerVolume) delete [] VecNbOfAtomsPerVolume;
485 if (fIonisation) delete fIonisation;
486 if (fSandiaTable) delete fSandiaTable;
487
488 // Remove this material from theMaterialTable.
489 //
490 theMaterialTable[fIndexInTable] = 0;
491}
492
493//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
494
495const G4Material& G4Material::operator=(const G4Material& right)
496{
497 if (this != &right)
498 {
499 fName = right.fName;
500 fChemicalFormula = right.fChemicalFormula;
501 fDensity = right.fDensity;
502 fState = right.fState;
503 fTemp = right.fTemp;
504 fPressure = right.fPressure;
505
506 if (fImplicitElement) delete ((*theElementVector)[0]);
507 if (theElementVector) delete theElementVector;
508 if (fMassFractionVector) delete [] fMassFractionVector;
509 if (fAtomsVector) delete [] fAtomsVector;
510
511 maxNbComponents = right.maxNbComponents;
512 fNumberOfComponents = right.fNumberOfComponents;
513 fNumberOfElements = right.fNumberOfElements;
514 fImplicitElement = right.fImplicitElement;
515
516 if (fImplicitElement) {
517 G4double z = (*right.theElementVector)[0]->GetZ();
518 G4double a = (*right.theElementVector)[0]->GetA();
519 theElementVector = new G4ElementVector(1,(G4Element*)0);
520 (*theElementVector)[0] = new G4Element(fName," ",z,a);
521 fMassFractionVector = new G4double[1];
522 fMassFractionVector[0] = 1.;
523 } else {
524 theElementVector = new G4ElementVector(fNumberOfElements,0);
525 fMassFractionVector = new G4double[fNumberOfElements];
526 for (size_t i=0; i<fNumberOfElements; i++) {
527 (*theElementVector)[i]= (*right.theElementVector)[i];
528 fMassFractionVector[i]= right.fMassFractionVector[i];
529 }
530 }
531
532 if (right.fAtomsVector) {
533 fAtomsVector = new G4int[fNumberOfElements];
534 for (size_t i=0; i<fNumberOfElements; i++)
535 fAtomsVector[i] = right.fAtomsVector[i];
536 }
537
538 fMaterialPropertiesTable = right.fMaterialPropertiesTable;
539
540 ComputeDerivedQuantities();
541 }
542 return *this;
543}
544
545//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
546
547G4int G4Material::operator==(const G4Material& right) const
548{
549 return (this == (G4Material *) &right);
550}
551
552//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
553
554G4int G4Material::operator!=(const G4Material& right) const
555{
556 return (this != (G4Material *) &right);
557}
558
559//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
560
561
562std::ostream& operator<<(std::ostream& flux, G4Material* material)
563{
564 std::ios::fmtflags mode = flux.flags();
565 flux.setf(std::ios::fixed,std::ios::floatfield);
566 G4long prec = flux.precision(3);
567
568 flux
569 << " Material: " << std::setw(8) << material->fName
570 << " " << material->fChemicalFormula << " "
571 << " density: " << std::setw(6) << std::setprecision(3)
572 << G4BestUnit(material->fDensity,"Volumic Mass")
573 << " RadL: " << std::setw(7) << std::setprecision(3)
574 << G4BestUnit(material->fRadlen,"Length")
575 << " Nucl.Int.Length: " << std::setw(7) << std::setprecision(3)
576 << G4BestUnit(material->fNuclInterLen,"Length")
577 << " Imean: " << std::setw(7) << std::setprecision(3)
578 << G4BestUnit(material->GetIonisation()->GetMeanExcitationEnergy(),"Energy");
579
580 if(material->fState == kStateGas)
581 flux
582 << " temperature: " << std::setw(6) << std::setprecision(2)
583 << (material->fTemp)/kelvin << " K"
584 << " pressure: " << std::setw(6) << std::setprecision(2)
585 << (material->fPressure)/atmosphere << " atm";
586
587 for (size_t i=0; i<material->fNumberOfElements; i++)
588 flux
589 << "\n ---> " << (*(material->theElementVector))[i]
590 << " ElmMassFraction: " << std::setw(6)<< std::setprecision(2)
591 << (material->fMassFractionVector[i])/perCent << " %"
592 << " ElmAbundance " << std::setw(6)<< std::setprecision(2)
593 << 100*(material->VecNbOfAtomsPerVolume[i])/(material->TotNbOfAtomsPerVolume)
594 << " %";
595
596 flux.precision(prec);
597 flux.setf(mode,std::ios::floatfield);
598
599 return flux;
600}
601
602//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
603
604 std::ostream& operator<<(std::ostream& flux, G4Material& material)
605{
606 flux << &material;
607 return flux;
608}
609
610//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
611
612std::ostream& operator<<(std::ostream& flux, G4MaterialTable MaterialTable)
613{
614 //Dump info for all known materials
615 flux << "\n***** Table : Nb of materials = " << MaterialTable.size()
616 << " *****\n" << G4endl;
617
618 for (size_t i=0; i<MaterialTable.size(); i++) flux << MaterialTable[i]
619 << G4endl << G4endl;
620
621 return flux;
622}
623
624//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.