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