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