Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoMaterial.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 25/10/01
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoMaterial
13\ingroup Materials_classes
14
15Base class describing materials.
16
17## Important note about units
18Since **v6-17-02** the geometry package adopted a system of units, upon the request to support
19an in-memory material representation consistent with the one in Geant4. The adoption was done
20gradually and starting with **v6-19-02** (back-ported to **v6-18-02**) the package supports changing
21the default units to either ROOT (CGS) or Geant4 ones. In the same version the Geant4 units were
22set to be the default ones, changing the previous behavior and making material properties such
23as radiation and interaction lengths having in memory values an order of magnitude lower. This behavior
24affected versions up to **v6-25-01**, after which the default units were restored to be the ROOT ones.
25
26For users needing to restore the CGS behavior for material properties, the following sequence needs
27to be called before creating the TGeoManager instance:
28 * From **v6-18-02** to **v6-22-06**:
29```
30 TGeoUnit::setUnitType(TGeoUnit::kTGeoUnits);
31```
32
33 * From **v6-22-08** to **v6-25-01**:
34```
35 TGeoManager::LockDefaultUnits(false);
36 TGeoManager::SetDefaultUnits(kRootUnits);
37 TGeoManager::LockDefaultUnits(true);
38```
39*/
40
41#include <iostream>
42#include <limits>
43#include "TMath.h"
44#include "TObjArray.h"
45#include "TGeoElement.h"
46#include "TGeoManager.h"
47#include "TGeoExtension.h"
48#include "TGeoMaterial.h"
51#include "TGDMLMatrix.h"
52
53// statics and globals
54
56
57////////////////////////////////////////////////////////////////////////////////
58/// Default constructor.
59
61 : TNamed(),
62 TAttFill(),
63 fIndex(0),
64 fA(0.),
65 fZ(0.),
66 fDensity(0.),
67 fRadLen(0.),
68 fIntLen(0.),
69 fTemperature(0.),
70 fPressure(0.),
71 fState(kMatStateUndefined),
72 fShader(nullptr),
73 fCerenkov(nullptr),
74 fElement(nullptr),
75 fUserExtension(nullptr),
76 fFWExtension(nullptr)
77{
78 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
80 fIndex = -1;
84}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Constructor.
88///
89/// \param name material name.
90
92 : TNamed(name, ""),
93 TAttFill(),
94 fIndex(0),
95 fA(0.),
96 fZ(0.),
97 fDensity(0.),
98 fRadLen(0.),
99 fIntLen(0.),
100 fTemperature(0.),
101 fPressure(0.),
102 fState(kMatStateUndefined),
103 fShader(nullptr),
104 fCerenkov(nullptr),
105 fElement(nullptr),
106 fUserExtension(nullptr),
107 fFWExtension(nullptr)
108{
109 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
110 fName = fName.Strip();
112 fIndex = -1;
116
117 if (!gGeoManager) {
118 gGeoManager = new TGeoManager("Geometry", "default geometry");
119 }
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Constructor.
125///
126/// \param name material name.
127/// \param a atomic mass.
128/// \param z atomic number.
129/// \param rho material density in g/cm3.
130/// \param radlen
131/// \param intlen
132
134 : TNamed(name, ""),
135 TAttFill(),
136 fIndex(0),
137 fA(a),
138 fZ(z),
139 fDensity(rho),
140 fRadLen(0.),
141 fIntLen(0.),
142 fTemperature(0.),
143 fPressure(0.),
144 fState(kMatStateUndefined),
145 fShader(nullptr),
146 fCerenkov(nullptr),
147 fElement(nullptr),
148 fUserExtension(nullptr),
149 fFWExtension(nullptr)
150{
151 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
152 fName = fName.Strip();
154 fIndex = -1;
155 fA = a;
156 fZ = z;
157 fDensity = rho;
161 SetRadLen(radlen, intlen);
162 if (!gGeoManager) {
163 gGeoManager = new TGeoManager("Geometry", "default geometry");
164 }
165 if (fZ - Int_t(fZ) > 1E-3)
166 Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
167 if (GetElement())
168 GetElement()->SetUsed();
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Constructor with state, temperature and pressure.
174///
175/// \param name material name.
176/// \param a atomic mass.
177/// \param z atomic number.
178/// \param rho material density in g/cm3.
179/// \param state
180/// \param temperature
181/// \param pressure
182
184 Double_t temperature, Double_t pressure)
185 : TNamed(name, ""),
186 TAttFill(),
187 fIndex(0),
188 fA(a),
189 fZ(z),
190 fDensity(rho),
191 fRadLen(0.),
192 fIntLen(0.),
193 fTemperature(temperature),
194 fPressure(pressure),
195 fState(state),
196 fShader(nullptr),
197 fCerenkov(nullptr),
198 fElement(nullptr),
199 fUserExtension(nullptr),
200 fFWExtension(nullptr)
201{
202 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
203 fName = fName.Strip();
205 fIndex = -1;
206 SetRadLen(0, 0);
207 if (!gGeoManager) {
208 gGeoManager = new TGeoManager("Geometry", "default geometry");
209 }
210 if (fZ - Int_t(fZ) > 1E-3)
211 Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
212 if (GetElement())
213 GetElement()->SetUsed();
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// Constructor.
219///
220/// \param name material name.
221/// \param elem
222/// \param rho material density in g/cm3.
223
225 : TNamed(name, ""),
226 TAttFill(),
227 fIndex(0),
228 fA(0.),
229 fZ(0.),
230 fDensity(rho),
231 fRadLen(0.),
232 fIntLen(0.),
233 fTemperature(0.),
234 fPressure(0.),
235 fState(kMatStateUndefined),
236 fShader(nullptr),
237 fCerenkov(nullptr),
238 fElement(elem),
239 fUserExtension(nullptr),
240 fFWExtension(nullptr)
241{
242 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
243 fName = fName.Strip();
245 fIndex = -1;
246 fA = elem->A();
247 fZ = elem->Z();
248 SetRadLen(0, 0);
252 if (!gGeoManager) {
253 gGeoManager = new TGeoManager("Geometry", "default geometry");
254 }
255 if (fZ - Int_t(fZ) > 1E-3)
256 Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
257 if (GetElement())
258 GetElement()->SetUsed();
260}
261
262////////////////////////////////////////////////////////////////////////////////
263
265 : TNamed(gm),
266 TAttFill(gm),
267 fIndex(gm.fIndex),
268 fA(gm.fA),
269 fZ(gm.fZ),
270 fDensity(gm.fDensity),
271 fRadLen(gm.fRadLen),
272 fIntLen(gm.fIntLen),
273 fTemperature(gm.fTemperature),
274 fPressure(gm.fPressure),
275 fState(gm.fState),
276 fShader(gm.fShader),
277 fCerenkov(gm.fCerenkov),
278 fElement(gm.fElement),
279 fUserExtension(gm.fUserExtension->Grab()),
280 fFWExtension(gm.fFWExtension->Grab())
281
282{
283 // copy constructor
284 TGeoManager::SetDefaultUnits(TGeoManager::GetDefaultUnits()); // Ensure nobody changes the units afterwards
286 TIter next(&fProperties);
288 while ((property = (TNamed *)next()))
290}
291
292////////////////////////////////////////////////////////////////////////////////
293/// assignment operator
294
296{
297 if (this != &gm) {
299 TAttFill::operator=(gm);
300 fIndex = gm.fIndex;
301 fA = gm.fA;
302 fZ = gm.fZ;
303 fDensity = gm.fDensity;
304 fRadLen = gm.fRadLen;
305 fIntLen = gm.fIntLen;
307 fPressure = gm.fPressure;
308 fState = gm.fState;
309 fShader = gm.fShader;
310 fCerenkov = gm.fCerenkov;
311 fElement = gm.fElement;
315 TIter next(&fProperties);
317 while ((property = (TNamed *)next()))
319 }
320 return *this;
321}
322
323////////////////////////////////////////////////////////////////////////////////
324/// Destructor
325
327{
328 if (fUserExtension) {
330 fUserExtension = nullptr;
331 }
332 if (fFWExtension) {
334 fFWExtension = nullptr;
335 }
336}
337
338////////////////////////////////////////////////////////////////////////////////
339/// Connect user-defined extension to the material. The material "grabs" a copy, so
340/// the original object can be released by the producer. Release the previously
341/// connected extension if any.
342///
343/// NOTE: This interface is intended for user extensions and is guaranteed not
344/// to be used by TGeo
345
347{
348 if (fUserExtension)
350 fUserExtension = nullptr;
351 if (ext)
352 fUserExtension = ext->Grab();
353}
354
355//_____________________________________________________________________________
356const char *TGeoMaterial::GetPropertyRef(const char *property) const
357{
358 // Find reference for a given property
360 return (prop) ? prop->GetTitle() : nullptr;
361}
362
363//_____________________________________________________________________________
365{
366 // Find reference for a given property
368 if (!prop)
369 return nullptr;
370 return gGeoManager->GetGDMLMatrix(prop->GetTitle());
371}
372
373//_____________________________________________________________________________
375{
376 // Find reference for a given property
378 if (!prop)
379 return nullptr;
380 return gGeoManager->GetGDMLMatrix(prop->GetTitle());
381}
382
383//_____________________________________________________________________________
384const char *TGeoMaterial::GetConstPropertyRef(const char *property) const
385{
386 // Find reference for a given constant property
388 return (prop) ? prop->GetTitle() : nullptr;
389}
390
391//_____________________________________________________________________________
393{
394 // Find reference for a given constant property
396 if (!prop) {
397 if (err)
398 *err = kTRUE;
399 return 0.;
400 }
401 return gGeoManager->GetProperty(prop->GetTitle(), err);
402}
403
404//_____________________________________________________________________________
406{
407 // Find reference for a given constant property
409 if (!prop) {
410 if (err)
411 *err = kTRUE;
412 return 0.;
413 }
414 return gGeoManager->GetProperty(prop->GetTitle(), err);
415}
416
417//_____________________________________________________________________________
418bool TGeoMaterial::AddProperty(const char *property, const char *ref)
419{
422 Error("AddProperty", "Property %s already added to material %s", property, GetName());
423 return false;
424 }
425 fProperties.Add(new TNamed(property, ref));
426 return true;
427}
428
429//_____________________________________________________________________________
430bool TGeoMaterial::AddConstProperty(const char *property, const char *ref)
431{
434 Error("AddConstProperty", "Constant property %s already added to material %s", property, GetName());
435 return false;
436 }
438 return true;
439}
440
441////////////////////////////////////////////////////////////////////////////////
442/// Connect framework defined extension to the material. The material "grabs" a copy,
443/// so the original object can be released by the producer. Release the previously
444/// connected extension if any.
445///
446/// NOTE: This interface is intended for the use by TGeo and the users should
447/// NOT connect extensions using this method
448
450{
451 if (fFWExtension)
453 fFWExtension = nullptr;
454 if (ext)
455 fFWExtension = ext->Grab();
456}
457
458////////////////////////////////////////////////////////////////////////////////
459/// Get a copy of the user extension pointer. The user must call Release() on
460/// the copy pointer once this pointer is not needed anymore (equivalent to
461/// delete() after calling new())
462
464{
465 if (fUserExtension)
466 return fUserExtension->Grab();
467 return nullptr;
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// Get a copy of the framework extension pointer. The user must call Release() on
472/// the copy pointer once this pointer is not needed anymore (equivalent to
473/// delete() after calling new())
474
476{
477 if (fFWExtension)
478 return fFWExtension->Grab();
479 return nullptr;
480}
481
482////////////////////////////////////////////////////////////////////////////////
483/// Provide a pointer name containing uid.
484
486{
487 static TString name;
488 name.Form("pMat%d", GetUniqueID());
489 return name.Data();
490}
491
492////////////////////////////////////////////////////////////////////////////////
493/// Set radiation/absorption lengths. If the values are negative, their absolute value
494/// is taken, otherwise radlen is recomputed using G3 formula.
495
497{
498 fRadLen = TMath::Abs(radlen);
499 fIntLen = TMath::Abs(intlen);
500 // Check for vacuum
501 if (fA < 0.9 || fZ < 0.9) {
502 if (radlen < -1e5 || intlen < -1e-5) {
503 Error("SetRadLen", "Material %s: user values taken for vacuum: radlen=%g or intlen=%g - too small", GetName(),
505 return;
506 }
507 // Ignore positive values and take big numbers
508 if (radlen >= 0)
509 fRadLen = 1.E30;
510 if (intlen >= 0)
511 fIntLen = 1.E30;
512 return;
513 }
515 // compute radlen systematically with G3 formula for a valid material
516 if (radlen >= 0) {
517 // taken grom Geant3 routine GSMATE
518 constexpr Double_t alr2av = 1.39621E-03;
519 constexpr Double_t al183 = 5.20948;
520 fRadLen = fA / (alr2av * fDensity * fZ * (fZ + TGeoMaterial::ScreenFactor(fZ)) *
521 (al183 - TMath::Log(fZ) / 3 - TGeoMaterial::Coulomb(fZ)));
522 // fRadLen is in TGeo units. Apply conversion factor in requested length-units
524 }
525 // Compute interaction length using the same formula as in GEANT4
526 if (intlen >= 0) {
527 constexpr Double_t lambda0 = 35. * TGeoUnit::g / TGeoUnit::cm2; // [g/cm^2]
528 Double_t nilinv = 0.0;
529 TGeoElement *elem = GetElement();
530 if (!elem) {
531 Fatal("SetRadLen", "Element not found for material %s", GetName());
532 return;
533 }
534 Double_t nbAtomsPerVolume = TGeoUnit::Avogadro * fDensity / elem->A();
535 nilinv += nbAtomsPerVolume * TMath::Power(elem->Neff(), 0.6666667);
536 nilinv *= TGeoUnit::amu / lambda0;
537 fIntLen = (nilinv <= 0) ? TGeoShape::Big() : (1.0 / nilinv);
538 // fIntLen is in TGeo units. Apply conversion factor in requested length-units
540 }
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// static function
545/// Compute Coulomb correction for pair production and Brem
546/// REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
547/// FORMULA 2.7.17
548
550{
553 Double_t az2 = az * az;
554 Double_t az4 = az2 * az2;
555 Double_t fp = (0.0083 * az4 + 0.20206 + 1. / (1. + az2)) * az2;
556 Double_t fm = (0.0020 * az4 + 0.0369) * az4;
557 return fp - fm;
558}
559
560////////////////////////////////////////////////////////////////////////////////
561/// return true if the other material has the same physical properties
562
564{
565 if (other == this)
566 return kTRUE;
567 if (other->IsMixture())
568 return kFALSE;
569 if (TMath::Abs(fA - other->GetA()) > 1E-3)
570 return kFALSE;
571 if (TMath::Abs(fZ - other->GetZ()) > 1E-3)
572 return kFALSE;
573 if (TMath::Abs(fDensity - other->GetDensity()) > 1E-6)
574 return kFALSE;
576 return kFALSE;
577 // if (fRadLen != other->GetRadLen()) return kFALSE;
578 // if (fIntLen != other->GetIntLen()) return kFALSE;
579 return kTRUE;
580}
581
582////////////////////////////////////////////////////////////////////////////////
583/// print characteristics of this material
584
585void TGeoMaterial::Print(const Option_t * /*option*/) const
586{
587 printf("Material %s %s A=%g Z=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(), fA, fZ, fDensity,
589}
590
591////////////////////////////////////////////////////////////////////////////////
592/// Save a primitive as a C++ statement(s) on output stream "out".
593
594void TGeoMaterial::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
595{
597 return;
598 const char *name = GetPointerName();
599 out << "// Material: " << GetName() << std::endl;
600 out << " a = " << fA << ";" << std::endl;
601 out << " z = " << fZ << ";" << std::endl;
602 out << " density = " << fDensity << ";" << std::endl;
603 out << " radl = " << fRadLen << ";" << std::endl;
604 out << " absl = " << fIntLen << ";" << std::endl;
605
606 out << " auto " << name << " = new TGeoMaterial(\"" << GetName() << "\", a, z, density, radl, absl);" << std::endl;
607 out << " " << name << "->SetIndex(" << GetIndex() << ");" << std::endl;
609}
610
611////////////////////////////////////////////////////////////////////////////////
612/// Get some default color related to this material.
613
615{
616 Int_t id = 1 + gGeoManager->GetListOfMaterials()->IndexOf(this);
617 return (2 + id % 6);
618}
619
620////////////////////////////////////////////////////////////////////////////////
621/// Get a pointer to the element this material is made of.
622/// This second call is to avoid warnings to not call a virtual
623/// method from the constructor
624
626{
627 if (fElement)
628 return fElement;
630 return table->GetElement(Int_t(fZ));
631}
632
633////////////////////////////////////////////////////////////////////////////////
634/// Get a pointer to the element this material is made of.
635
637{
638 if (fElement)
639 return fElement;
641 return table->GetElement(Int_t(fZ));
642}
643
644////////////////////////////////////////////////////////////////////////////////
645/// Single interface to get element properties.
646
648{
649 a = fA;
650 z = fZ;
651 w = 1.;
652}
653
654////////////////////////////////////////////////////////////////////////////////
655/// Retrieve material index in the list of materials
656
658{
659 if (fIndex >= 0)
660 return fIndex;
662 fIndex = matlist->IndexOf(this);
663 return fIndex;
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Create the material representing the decay product of this material at a
668/// given time. The precision represent the minimum cumulative branching ratio for
669/// which decay products are still taken into account.
670
672{
673 TObjArray *pop = new TObjArray();
674 if (!fElement || !fElement->IsRadioNuclide())
675 return this;
676 FillMaterialEvolution(pop, precision);
677 Int_t ncomp = pop->GetEntriesFast();
678 if (!ncomp)
679 return this;
680 TGeoElementRN *el;
681 Double_t *weight = new Double_t[ncomp];
682 Double_t amed = 0.;
683 Int_t i;
684 for (i = 0; i < ncomp; i++) {
685 el = (TGeoElementRN *)pop->At(i);
686 weight[i] = el->Ratio()->Concentration(time) * el->A();
687 amed += weight[i];
688 }
689 Double_t rho = fDensity * amed / fA;
690 TGeoMixture *mix = nullptr;
691 Int_t ncomp1 = ncomp;
692 for (i = 0; i < ncomp; i++) {
693 if ((weight[i] / amed) < precision) {
694 amed -= weight[i];
695 ncomp1--;
696 }
697 }
698 if (ncomp1 < 2) {
699 el = (TGeoElementRN *)pop->At(0);
700 delete[] weight;
701 delete pop;
702 if (ncomp1 == 1)
703 return new TGeoMaterial(TString::Format("%s-evol", GetName()), el, rho);
704 return nullptr;
705 }
706 mix = new TGeoMixture(TString::Format("%s-evol", GetName()), ncomp, rho);
707 for (i = 0; i < ncomp; i++) {
708 weight[i] /= amed;
709 if (weight[i] < precision)
710 continue;
711 el = (TGeoElementRN *)pop->At(i);
712 mix->AddElement(el, weight[i]);
713 }
714 delete[] weight;
715 delete pop;
716 return mix;
717}
718
719////////////////////////////////////////////////////////////////////////////////
720/// Fills a user array with all the elements deriving from the possible
721/// decay of the top element composing the mixture. Each element contained
722/// by `<population>` may be a radionuclide having a Bateman solution attached.
723/// The precision represent the minimum cumulative branching ratio for
724/// which decay products are still taken into account.
725/// To visualize the time evolution of each decay product one can use:
726/// ~~~ {.cpp}
727/// TGeoElement *elem = population->At(index);
728/// TGeoElementRN *elemrn = 0;
729/// if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
730/// ~~~
731/// One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
732/// of one of the decay products, N1(0) is the number of atoms of the top
733/// element at t=0.
734/// ~~~ {.cpp}
735/// Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
736/// ~~~
737/// One can also display the time evolution of the fractional weight:
738/// ~~~ {.cpp}
739/// elemrn->Ratio()->Draw(option);
740/// ~~~
741
743{
744 if (population->GetEntriesFast()) {
745 Error("FillMaterialEvolution", "Provide an empty array !");
746 return;
747 }
749 TGeoElement *elem;
750 TGeoElementRN *elemrn;
751 TIter next(table->GetElementsRN());
752 while ((elemrn = (TGeoElementRN *)next()))
753 elemrn->ResetRatio();
754 elem = GetElement();
755 if (!elem) {
756 Fatal("FillMaterialEvolution", "Element not found for material %s", GetName());
757 return;
758 }
759 if (!elem->IsRadioNuclide()) {
760 population->Add(elem);
761 return;
762 }
763 elemrn = (TGeoElementRN *)elem;
764 elemrn->FillPopulation(population, precision);
765}
766
767/** \class TGeoMixture
768\ingroup Materials_classes
769
770Mixtures of elements.
771
772*/
773
775
776////////////////////////////////////////////////////////////////////////////////
777/// Default constructor
778
780{
781 fNelements = 0;
782 fZmixture = nullptr;
783 fAmixture = nullptr;
784 fWeights = nullptr;
785 fNatoms = nullptr;
786 fVecNbOfAtomsPerVolume = nullptr;
787 fElements = nullptr;
788}
789
790////////////////////////////////////////////////////////////////////////////////
791/// constructor
792
794{
795 fZmixture = nullptr;
796 fAmixture = nullptr;
797 fWeights = nullptr;
798 fNelements = 0;
799 fNatoms = nullptr;
800 fVecNbOfAtomsPerVolume = nullptr;
801 fDensity = rho;
802 fElements = nullptr;
803 if (fDensity < 0)
804 fDensity = 0.001;
805}
806
807////////////////////////////////////////////////////////////////////////////////
808/// Destructor
809
811{
812 if (fZmixture)
813 delete[] fZmixture;
814 if (fAmixture)
815 delete[] fAmixture;
816 if (fWeights)
817 delete[] fWeights;
818 if (fNatoms)
819 delete[] fNatoms;
821 delete[] fVecNbOfAtomsPerVolume;
822 if (fElements)
823 delete fElements;
824}
825
826////////////////////////////////////////////////////////////////////////////////
827/// Compute effective A/Z and radiation length
828
830{
831 constexpr const Double_t na = TGeoUnit::Avogadro;
832 constexpr const Double_t alr2av = 1.39621E-03;
833 constexpr const Double_t al183 = 5.20948;
834 constexpr const Double_t lambda0 = 35. * TGeoUnit::g / TGeoUnit::cm2; // [g/cm^2]
835 Double_t radinv = 0.0;
836 Double_t nilinv = 0.0;
837 Double_t nbAtomsPerVolume;
838 fA = 0;
839 fZ = 0;
840 for (Int_t j = 0; j < fNelements; j++) {
841 if (fWeights[j] <= 0)
842 continue;
843 fA += fWeights[j] * fAmixture[j];
844 fZ += fWeights[j] * fZmixture[j];
845 nbAtomsPerVolume = na * fDensity * fWeights[j] / GetElement(j)->A();
846 nilinv += nbAtomsPerVolume * TMath::Power(GetElement(j)->Neff(), 0.6666667);
847 Double_t zc = fZmixture[j];
848 Double_t alz = TMath::Log(zc) / 3.;
849 Double_t xinv =
850 zc * (zc + TGeoMaterial::ScreenFactor(zc)) * (al183 - alz - TGeoMaterial::Coulomb(zc)) / fAmixture[j];
851 radinv += xinv * fWeights[j];
852 }
853 radinv *= alr2av * fDensity;
854 fRadLen = (radinv <= 0) ? TGeoShape::Big() : 1.0 / radinv;
855 // fRadLen is in TGeo units. Apply conversion factor in requested length-units
857
858 // Compute interaction length
859 nilinv *= TGeoUnit::amu / lambda0;
860 fIntLen = (nilinv <= 0) ? TGeoShape::Big() : 1.0 / nilinv;
861 // fIntLen is in TGeo units. Apply conversion factor in requested length-units
863}
864
865////////////////////////////////////////////////////////////////////////////////
866/// add an element to the mixture using fraction by weight
867/// Check if the element is already defined
868
870{
872
873 // Check preconditions
874 if (weight < 0e0) {
875 Fatal("AddElement", "Cannot add element with negative weight %g to mixture %s", weight, GetName());
876 } else if (weight < std::numeric_limits<Double_t>::epsilon()) {
877 return;
878 } else if (z < 1 || z > table->GetNelements() - 1) {
879 Fatal("AddElement", "Cannot add element having Z=%d to mixture %s", (Int_t)z, GetName());
880 }
881 Int_t i;
882 for (i = 0; i < fNelements; i++) {
883 if (!fElements && TMath::Abs(z - fZmixture[i]) < 1.e-6 && TMath::Abs(a - fAmixture[i]) < 1.e-6) {
884 fWeights[i] += weight;
886 return;
887 }
888 }
889 if (!fNelements) {
890 fZmixture = new Double_t[1];
891 fAmixture = new Double_t[1];
892 fWeights = new Double_t[1];
893 } else {
894 Int_t nelements = fNelements + 1;
895 Double_t *zmixture = new Double_t[nelements];
896 Double_t *amixture = new Double_t[nelements];
897 Double_t *weights = new Double_t[nelements];
898 for (Int_t j = 0; j < fNelements; j++) {
899 zmixture[j] = fZmixture[j];
900 amixture[j] = fAmixture[j];
901 weights[j] = fWeights[j];
902 }
903 delete[] fZmixture;
904 delete[] fAmixture;
905 delete[] fWeights;
906 fZmixture = zmixture;
907 fAmixture = amixture;
908 fWeights = weights;
909 }
910
911 fNelements++;
912 i = fNelements - 1;
913 fZmixture[i] = z;
914 fAmixture[i] = a;
915 fWeights[i] = weight;
916 if (z - Int_t(z) > 1E-3)
917 Warning("DefineElement", "Mixture %s has element defined with fractional Z=%f", GetName(), z);
919 table->GetElement((Int_t)z)->SetDefined();
920
921 // compute equivalent radiation length (taken from Geant3/GSMIXT)
923}
924
925////////////////////////////////////////////////////////////////////////////////
926/// Define one component of the mixture as an existing material/mixture.
927
929{
930 TGeoElement *elnew, *elem;
931 Double_t a, z;
932
933 // Check preconditions
934 if (!mat) {
935 Fatal("AddElement", "Cannot add INVALID material to mixture %s", GetName());
936 } else if (weight < 0e0) {
937 Fatal("AddElement", "Cannot add material %s with negative weight %g to mixture %s", mat->GetName(), weight,
938 GetName());
939 } else if (weight < std::numeric_limits<Double_t>::epsilon()) {
940 return;
941 }
942 if (!mat->IsMixture()) {
943 elem = mat->GetBaseElement();
944 if (elem) {
945 AddElement(elem, weight);
946 } else {
947 a = mat->GetA();
948 z = mat->GetZ();
949 AddElement(a, z, weight);
950 }
951 return;
952 }
953 // The material is a mixture.
954 TGeoMixture *mix = (TGeoMixture *)mat;
955 Double_t wnew;
956 Int_t nelem = mix->GetNelements();
957 Bool_t elfound;
958 Int_t i, j;
959 // loop the elements of the daughter mixture
960 for (i = 0; i < nelem; i++) {
961 elfound = kFALSE;
962 elnew = mix->GetElement(i);
963 if (!elnew)
964 continue;
965 // check if we have the element already defined in the parent mixture
966 for (j = 0; j < fNelements; j++) {
967 if (fWeights[j] < 0e0)
968 continue;
969 elem = GetElement(j);
970 if (elem == elnew) {
971 // element found, compute new weight
972 fWeights[j] += weight * (mix->GetWmixt())[i];
973 elfound = kTRUE;
974 break;
975 }
976 }
977 if (elfound)
978 continue;
979 // element not found, define it
980 wnew = weight * (mix->GetWmixt())[i];
981 AddElement(elnew, wnew);
982 }
983}
984
985////////////////////////////////////////////////////////////////////////////////
986/// add an element to the mixture using fraction by weight
987
989{
990 TGeoElement *elemold;
992 if (!fElements)
993 fElements = new TObjArray(128);
994 Bool_t exist = kFALSE;
995
996 // Check preconditions
997 if (!elem) {
998 Fatal("AddElement", "Cannot add INVALID element to mixture %s", GetName());
999 } else if (weight < 0e0) {
1000 Fatal("AddElement", "Cannot add element %s with negative weight %g to mixture %s", elem->GetName(), weight,
1001 GetName());
1002 } else if (weight < std::numeric_limits<Double_t>::epsilon()) {
1003 return;
1004 }
1005 // If previous elements were defined by A/Z, add corresponding TGeoElements
1006 for (Int_t i = 0; i < fNelements; i++) {
1007 elemold = (TGeoElement *)fElements->At(i);
1008 if (!elemold) {
1009 // Add element with corresponding Z in the list
1010 fElements->AddAt(elemold = table->GetElement((Int_t)fZmixture[i]), i);
1011 elemold->SetDefined();
1012 }
1013 if (elemold == elem) {
1014 fWeights[i] += weight;
1015 exist = kTRUE;
1016 }
1017 }
1018 if (!exist) {
1020 AddElement(elem->A(), elem->Z(), weight);
1021 } else {
1023 }
1024}
1025
1026////////////////////////////////////////////////////////////////////////////////
1027/// Add a mixture element by number of atoms in the chemical formula.
1028
1030{
1031 Int_t i, j;
1032 Double_t amol;
1033 TGeoElement *elemold;
1035 if (!fElements)
1036 fElements = new TObjArray(128);
1037 // Check if the element is already defined
1038 for (i = 0; i < fNelements; i++) {
1039 elemold = (TGeoElement *)fElements->At(i);
1040 if (!elemold)
1041 fElements->AddAt(table->GetElement((Int_t)fZmixture[i]), i);
1042 else if (elemold != elem)
1043 continue;
1044 if ((elem == elemold) ||
1045 (TMath::Abs(elem->Z() - fZmixture[i]) < 1.e-6 && TMath::Abs(elem->A() - fAmixture[i]) < 1.e-6)) {
1046 fNatoms[i] += natoms;
1047 amol = 0.;
1048 for (j = 0; j < fNelements; j++)
1049 amol += fAmixture[j] * fNatoms[j];
1050 for (j = 0; j < fNelements; j++)
1051 fWeights[j] = fNatoms[j] * fAmixture[j] / amol;
1053 return;
1054 }
1055 }
1056 // New element
1057 if (!fNelements) {
1058 fZmixture = new Double_t[1];
1059 fAmixture = new Double_t[1];
1060 fWeights = new Double_t[1];
1061 fNatoms = new Int_t[1];
1062 } else {
1063 if (!fNatoms) {
1064 Fatal("AddElement", "Cannot add element by natoms in mixture %s after defining elements by weight", GetName());
1065 return;
1066 }
1067 Int_t nelements = fNelements + 1;
1068 Double_t *zmixture = new Double_t[nelements];
1069 Double_t *amixture = new Double_t[nelements];
1070 Double_t *weights = new Double_t[nelements];
1071 Int_t *nnatoms = new Int_t[nelements];
1072 for (j = 0; j < fNelements; j++) {
1073 zmixture[j] = fZmixture[j];
1074 amixture[j] = fAmixture[j];
1075 weights[j] = fWeights[j];
1076 nnatoms[j] = fNatoms[j];
1077 }
1078 delete[] fZmixture;
1079 delete[] fAmixture;
1080 delete[] fWeights;
1081 delete[] fNatoms;
1082 fZmixture = zmixture;
1083 fAmixture = amixture;
1084 fWeights = weights;
1085 fNatoms = nnatoms;
1086 }
1087 fNelements++;
1088 Int_t iel = fNelements - 1;
1089 fZmixture[iel] = elem->Z();
1090 fAmixture[iel] = elem->A();
1091 fNatoms[iel] = natoms;
1092 fElements->AddAtAndExpand(elem, iel);
1093 amol = 0.;
1094 for (i = 0; i < fNelements; i++) {
1095 if (fNatoms[i] <= 0)
1096 return;
1097 amol += fAmixture[i] * fNatoms[i];
1098 }
1099 for (i = 0; i < fNelements; i++)
1100 fWeights[i] = fNatoms[i] * fAmixture[i] / amol;
1101 table->GetElement(elem->Z())->SetDefined();
1103}
1104
1105////////////////////////////////////////////////////////////////////////////////
1106/// Define the mixture element at index iel by number of atoms in the chemical formula.
1107
1109{
1111 TGeoElement *elem = table->GetElement(z);
1112 if (!elem) {
1113 Fatal("DefineElement", "In mixture %s, element with Z=%i not found", GetName(), z);
1114 return;
1115 }
1116 AddElement(elem, natoms);
1117}
1118
1119////////////////////////////////////////////////////////////////////////////////
1120/// Retrieve the pointer to the element corresponding to component I.
1121
1123{
1124 if (i < 0 || i >= fNelements) {
1125 Error("GetElement", "Mixture %s has only %d elements", GetName(), fNelements);
1126 return nullptr;
1127 }
1128 TGeoElement *elem = nullptr;
1129 if (fElements)
1130 elem = (TGeoElement *)fElements->At(i);
1131 if (elem)
1132 return elem;
1134 return table->GetElement(Int_t(fZmixture[i]));
1135}
1136
1137////////////////////////////////////////////////////////////////////////////////
1138/// Get specific activity (in Bq/gram) for the whole mixture (no argument) or
1139/// for a given component.
1140
1142{
1143 if (i >= 0 && i < fNelements)
1144 return fWeights[i] * GetElement(i)->GetSpecificActivity();
1145 Double_t sa = 0;
1146 for (Int_t iel = 0; iel < fNelements; iel++) {
1147 sa += fWeights[iel] * GetElement(iel)->GetSpecificActivity();
1148 }
1149 return sa;
1150}
1151
1152////////////////////////////////////////////////////////////////////////////////
1153/// Return true if the other material has the same physical properties
1154
1156{
1157 if (other->IsEqual(this))
1158 return kTRUE;
1159 if (!other->IsMixture())
1160 return kFALSE;
1161 TGeoMixture *mix = (TGeoMixture *)other;
1162 if (!mix)
1163 return kFALSE;
1164 if (fNelements != mix->GetNelements())
1165 return kFALSE;
1166 if (TMath::Abs(fA - other->GetA()) > 1E-3)
1167 return kFALSE;
1168 if (TMath::Abs(fZ - other->GetZ()) > 1E-3)
1169 return kFALSE;
1170 if (TMath::Abs(fDensity - other->GetDensity()) > 1E-6)
1171 return kFALSE;
1173 return kFALSE;
1174 // if (fRadLen != other->GetRadLen()) return kFALSE;
1175 // if (fIntLen != other->GetIntLen()) return kFALSE;
1176 for (Int_t i = 0; i < fNelements; i++) {
1177 if (TMath::Abs(fZmixture[i] - (mix->GetZmixt())[i]) > 1E-3)
1178 return kFALSE;
1179 if (TMath::Abs(fAmixture[i] - (mix->GetAmixt())[i]) > 1E-3)
1180 return kFALSE;
1181 if (TMath::Abs(fWeights[i] - (mix->GetWmixt())[i]) > 1E-3)
1182 return kFALSE;
1183 }
1184 return kTRUE;
1185}
1186
1187////////////////////////////////////////////////////////////////////////////////
1188/// print characteristics of this material
1189
1190void TGeoMixture::Print(const Option_t * /*option*/) const
1191{
1192 printf("Mixture %s %s Aeff=%g Zeff=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(), fA, fZ,
1194 for (Int_t i = 0; i < fNelements; i++) {
1195 if (fElements && fElements->At(i)) {
1196 printf(" Element #%i : %s Z=%6.2f A=%6.2f w=%6.3f\n", i, GetElement(i)->GetName(), fZmixture[i],
1197 fAmixture[i], fWeights[i]);
1198 continue;
1199 }
1200 if (fNatoms)
1201 printf(" Element #%i : %s Z=%6.2f A=%6.2f w=%6.3f natoms=%d\n", i, GetElement(i)->GetName(), fZmixture[i],
1202 fAmixture[i], fWeights[i], fNatoms[i]);
1203 else
1204 printf(" Element #%i : %s Z=%6.2f A=%6.2f w=%6.3f\n", i, GetElement(i)->GetName(), fZmixture[i],
1205 fAmixture[i], fWeights[i]);
1206 }
1207}
1208
1209////////////////////////////////////////////////////////////////////////////////
1210/// Save a primitive as a C++ statement(s) on output stream "out".
1211
1212void TGeoMixture::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
1213{
1215 return;
1216 const char *name = GetPointerName();
1217 out << "// Mixture: " << GetName() << std::endl;
1218 out << " nel = " << fNelements << ";" << std::endl;
1219 out << " density = " << fDensity << ";" << std::endl;
1220 out << " auto " << name << " = new TGeoMixture(\"" << GetName() << "\", nel, density);" << std::endl;
1221 for (Int_t i = 0; i < fNelements; i++) {
1222 TGeoElement *el = GetElement(i);
1223 out << " a = " << fAmixture[i] << "; z = " << fZmixture[i] << "; w = " << fWeights[i] << "; // "
1224 << el->GetName() << std::endl;
1225 out << " " << name << "->DefineElement(" << i << ",a,z,w);" << std::endl;
1226 }
1227 out << " " << name << "->SetIndex(" << GetIndex() << ");" << std::endl;
1229}
1230
1231////////////////////////////////////////////////////////////////////////////////
1232/// Create the mixture representing the decay product of this material at a
1233/// given time. The precision represent the minimum cumulative branching ratio for
1234/// which decay products are still taken into account.
1235
1237{
1238 TObjArray *pop = new TObjArray();
1239 FillMaterialEvolution(pop, precision);
1240 Int_t ncomp = pop->GetEntriesFast();
1241 if (!ncomp)
1242 return this;
1243 TGeoElement *elem;
1244 TGeoElementRN *el;
1245 Double_t *weight = new Double_t[ncomp];
1246 Double_t amed = 0.;
1247 Int_t i, j;
1248 for (i = 0; i < ncomp; i++) {
1249 elem = (TGeoElement *)pop->At(i);
1250 if (!elem->IsRadioNuclide()) {
1251 j = fElements->IndexOf(elem);
1252 weight[i] = fWeights[j] * fAmixture[0] / fWeights[0];
1253 } else {
1254 el = (TGeoElementRN *)elem;
1255 weight[i] = el->Ratio()->Concentration(time) * el->A();
1256 }
1257 amed += weight[i];
1258 }
1259 Double_t rho = fDensity * fWeights[0] * amed / fAmixture[0];
1260 TGeoMixture *mix = nullptr;
1261 Int_t ncomp1 = ncomp;
1262 for (i = 0; i < ncomp; i++) {
1263 if ((weight[i] / amed) < precision) {
1264 amed -= weight[i];
1265 ncomp1--;
1266 }
1267 }
1268 if (ncomp1 < 2) {
1269 el = (TGeoElementRN *)pop->At(0);
1270 delete[] weight;
1271 delete pop;
1272 if (ncomp1 == 1)
1273 return new TGeoMaterial(TString::Format("%s-evol", GetName()), el, rho);
1274 return nullptr;
1275 }
1276 mix = new TGeoMixture(TString::Format("%s-evol", GetName()), ncomp, rho);
1277 for (i = 0; i < ncomp; i++) {
1278 weight[i] /= amed;
1279 if (weight[i] < precision)
1280 continue;
1281 el = (TGeoElementRN *)pop->At(i);
1282 mix->AddElement(el, weight[i]);
1283 }
1284 delete[] weight;
1285 delete pop;
1286 return mix;
1287}
1288
1289////////////////////////////////////////////////////////////////////////////////
1290/// Fills a user array with all the elements deriving from the possible
1291/// decay of the top elements composing the mixture. Each element contained
1292/// by `<population>` may be a radionuclide having a Bateman solution attached.
1293/// The precision represent the minimum cumulative branching ratio for
1294/// which decay products are still taken into account.
1295/// To visualize the time evolution of each decay product one can use:
1296/// ~~~ {.cpp}
1297/// TGeoElement *elem = population->At(index);
1298/// TGeoElementRN *elemrn = 0;
1299/// if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
1300/// ~~~
1301/// One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
1302/// of one of the decay products, N1(0) is the number of atoms of the first top
1303/// element at t=0.
1304/// ~~~ {.cpp}
1305/// Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
1306/// ~~~
1307/// One can also display the time evolution of the fractional weight:
1308/// ~~~ {.cpp}
1309/// elemrn->Ratio()->Draw(option);
1310/// ~~~
1311
1313{
1314 if (population->GetEntriesFast()) {
1315 Error("FillMaterialEvolution", "Provide an empty array !");
1316 return;
1317 }
1319 TGeoElement *elem;
1320 TGeoElementRN *elemrn;
1321 TIter next(table->GetElementsRN());
1322 while ((elemrn = (TGeoElementRN *)next()))
1323 elemrn->ResetRatio();
1324 Double_t factor;
1325 for (Int_t i = 0; i < fNelements; i++) {
1326 elem = GetElement(i);
1327 if (!elem->IsRadioNuclide()) {
1328 population->Add(elem);
1329 continue;
1330 }
1331 elemrn = (TGeoElementRN *)elem;
1332 factor = fWeights[i] * fAmixture[0] / (fWeights[0] * fAmixture[i]);
1333 elemrn->FillPopulation(population, precision, factor);
1334 }
1335}
1336
1337////////////////////////////////////////////////////////////////////////////////
1338/// static function
1339/// Compute screening factor for pair production and Bremsstrahlung
1340/// REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
1341/// FORMULA 2.7.22
1342
1344{
1345 const Double_t al183 = 5.20948, al1440 = 7.27239;
1346 Double_t alz = TMath::Log(z) / 3.;
1347 Double_t factor = (al1440 - 2 * alz) / (al183 - alz - TGeoMaterial::Coulomb(z));
1348 return factor;
1349}
1350
1351////////////////////////////////////////////////////////////////////////////////
1352/// Compute Derived Quantities as in Geant4
1353
1355{
1356 const Double_t Na =
1358
1360 delete[] fVecNbOfAtomsPerVolume;
1361
1363
1364 // Formula taken from G4Material.cxx L312
1365 double sumweights = 0;
1366 for (Int_t i = 0; i < fNelements; ++i) {
1367 sumweights += fWeights[i];
1368 fVecNbOfAtomsPerVolume[i] = Na * fDensity * fWeights[i] / ((TGeoElement *)fElements->At(i))->A();
1369 }
1370 if (TMath::Abs(sumweights - 1) > 0.001)
1371 Warning("ComputeDerivedQuantities", "Mixture %s: sum of weights is: %g", GetName(), sumweights);
1374}
1375
1376////////////////////////////////////////////////////////////////////////////////
1377/// Compute Radiation Length based on Geant4 formula
1378
1380{
1381 // Formula taken from G4Material.cxx L556
1382 Double_t radinv = 0.0;
1383 // GetfRadTsai is in units of cm2 due to <unit>::alpha_rcl2. Correction must be applied to end up in TGeo cm.
1385 for (Int_t i = 0; i < fNelements; ++i) {
1386 radinv += fVecNbOfAtomsPerVolume[i] * ((TGeoElement *)fElements->At(i))->GetfRadTsai() / denom;
1387 }
1388 fRadLen = (radinv <= 0.0 ? DBL_MAX : 1.0 / radinv);
1389 // fRadLen is in TGeo units. Apply conversion factor in requested length-units
1391}
1392
1393////////////////////////////////////////////////////////////////////////////////
1394/// Compute Nuclear Interaction Length based on Geant4 formula
1396{
1397 // Formula taken from G4Material.cxx L567
1398 constexpr Double_t lambda0 = 35. * TGeoUnit::g / TGeoUnit::cm2; // [g/cm^2]
1399 const Double_t twothird = 2.0 / 3.0;
1400 Double_t NILinv = 0.0;
1401 for (Int_t i = 0; i < fNelements; ++i) {
1402 Int_t Z = static_cast<Int_t>(((TGeoElement *)fElements->At(i))->Z() + 0.5);
1403 Double_t A = ((TGeoElement *)fElements->At(i))->Neff();
1404 if (1 == Z) {
1405 NILinv += fVecNbOfAtomsPerVolume[i] * A;
1406 } else {
1407 NILinv += fVecNbOfAtomsPerVolume[i] * TMath::Exp(twothird * TMath::Log(A));
1408 }
1409 }
1410 NILinv *= TGeoUnit::amu / lambda0;
1411 fIntLen = (NILinv <= 0.0 ? DBL_MAX : 1.0 / NILinv);
1412 // fIntLen is in TGeo units. Apply conversion factor in requested length-units
1414}
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h prop
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t property
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
static const Double_t STP_temperature
static const Double_t STP_pressure
Fill Area Attributes class.
Definition TAttFill.h:19
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
This class is used in the process of reading and writing the GDML "matrix" tag.
Definition TGDMLMatrix.h:33
Double_t Concentration(Double_t time) const
Find concentration of the element at a given time.
Class representing a radionuclidevoid TGeoManager::SetDefaultRootUnits() { if ( fgDefaultUnits == kRo...
void FillPopulation(TObjArray *population, Double_t precision=0.001, Double_t factor=1.)
Fills the input array with the set of RN elements resulting from the decay of this one.
TGeoBatemanSol * Ratio() const
void ResetRatio()
Clears the existing ratio.
Table of elements.
TGeoElement * GetElement(Int_t z)
TObjArray * GetElementsRN() const
Int_t GetNelements() const
Base class for chemical elements.
Definition TGeoElement.h:36
virtual Double_t GetSpecificActivity() const
Definition TGeoElement.h:79
Double_t A() const
Definition TGeoElement.h:71
void SetDefined(Bool_t flag=kTRUE)
Definition TGeoElement.h:85
virtual Bool_t IsRadioNuclide() const
Definition TGeoElement.h:82
Double_t Neff() const
Returns effective number of nucleons.
Int_t Z() const
Definition TGeoElement.h:68
void SetUsed(Bool_t flag=kTRUE)
Definition TGeoElement.h:86
ABC for user objects attached to TGeoVolume or TGeoNode.
virtual TGeoExtension * Grab()=0
virtual void Release() const =0
The manager class for any TGeo geometry.
Definition TGeoManager.h:44
static EDefaultUnits GetDefaultUnits()
TGeoElementTable * GetElementTable()
Returns material table. Creates it if not existing.
static void SetDefaultUnits(EDefaultUnits new_value)
Int_t AddMaterial(const TGeoMaterial *material)
Add a material to the list. Returns index of the material in list.
Double_t GetProperty(const char *name, Bool_t *error=nullptr) const
Get a user-defined property.
TGDMLMatrix * GetGDMLMatrix(const char *name) const
Get GDML matrix with a given name;.
TList * GetListOfMaterials() const
Base class describing materials.
Double_t GetConstProperty(const char *property, Bool_t *error=nullptr) const
void SetUserExtension(TGeoExtension *ext)
Connect user-defined extension to the material.
virtual TObject * GetCerenkovProperties() const
EGeoMaterialState fState
const char * GetPointerName() const
Provide a pointer name containing uid.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
static Double_t ScreenFactor(Double_t z)
static function Compute screening factor for pair production and Bremsstrahlung REFERENCE : EGS MANUA...
const char * GetConstPropertyRef(const char *property) const
virtual Bool_t IsMixture() const
Double_t fPressure
void SetFWExtension(TGeoExtension *ext)
Connect framework defined extension to the material.
bool AddConstProperty(const char *property, const char *ref)
Double_t fTemperature
virtual void GetElementProp(Double_t &a, Double_t &z, Double_t &w, Int_t i=0)
Single interface to get element properties.
virtual TGeoMaterial * DecayMaterial(Double_t time, Double_t precision=0.001)
Create the material representing the decay product of this material at a given time.
TList fProperties
bool AddProperty(const char *property, const char *ref)
Double_t fZ
void SetRadLen(Double_t radlen, Double_t intlen=0.)
Set radiation/absorption lengths.
TGeoElement * GetElement() const
Get a pointer to the element this material is made of.
virtual Bool_t IsEq(const TGeoMaterial *other) const
return true if the other material has the same physical properties
TGeoElement * GetBaseElement() const
const char * GetPropertyRef(const char *property) const
Double_t fA
TGDMLMatrix * GetProperty(const char *name) const
TObject * fCerenkov
Double_t fDensity
void SetUsed(Bool_t flag=kTRUE)
virtual void FillMaterialEvolution(TObjArray *population, Double_t precision=0.001)
Fills a user array with all the elements deriving from the possible decay of the top element composin...
Int_t GetIndex()
Retrieve material index in the list of materials.
Double_t fIntLen
TGeoExtension * fUserExtension
TList fConstProperties
TGeoElement * fElement
TGeoMaterial & operator=(const TGeoMaterial &)
assignment operator
TObject * fShader
TGeoExtension * GrabUserExtension() const
Get a copy of the user extension pointer.
Double_t fRadLen
TGeoExtension * GrabFWExtension() const
Get a copy of the framework extension pointer.
virtual Int_t GetDefaultColor() const
Get some default color related to this material.
TGeoMaterial()
Default constructor.
static Double_t Coulomb(Double_t z)
static function Compute Coulomb correction for pair production and Brem REFERENCE : EGS MANUAL SLAC 2...
virtual Double_t GetA() const
TGeoExtension * fFWExtension
Transient user-defined extension to materials.
void Print(const Option_t *option="") const override
print characteristics of this material
virtual Double_t GetDensity() const
virtual Double_t GetZ() const
~TGeoMaterial() override
Destructor.
Mixtures of elements.
TObjArray * fElements
void ComputeNuclearInterLength()
Compute Nuclear Interaction Length based on Geant4 formula.
~TGeoMixture() override
Destructor.
Double_t * GetZmixt() const
Bool_t IsEq(const TGeoMaterial *other) const override
Return true if the other material has the same physical properties.
void AddElement(Double_t a, Double_t z, Double_t weight)
add an element to the mixture using fraction by weight Check if the element is already defined
Double_t * fZmixture
TGeoMixture()
Default constructor.
void ComputeDerivedQuantities()
Compute Derived Quantities as in Geant4.
Int_t GetNelements() const override
TGeoMaterial * DecayMaterial(Double_t time, Double_t precision=0.001) override
Create the mixture representing the decay product of this material at a given time.
void ComputeRadiationLength()
Compute Radiation Length based on Geant4 formula.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
Double_t * fVecNbOfAtomsPerVolume
Double_t * fAmixture
void FillMaterialEvolution(TObjArray *population, Double_t precision=0.001) override
Fills a user array with all the elements deriving from the possible decay of the top elements composi...
void AverageProperties()
Compute effective A/Z and radiation length.
void Print(const Option_t *option="") const override
print characteristics of this material
Double_t * GetWmixt() const
Double_t GetSpecificActivity(Int_t i=-1) const override
Get specific activity (in Bq/gram) for the whole mixture (no argument) or for a given component.
Int_t * fNatoms
Double_t * fWeights
void DefineElement(Int_t iel, Double_t a, Double_t z, Double_t weight)
Double_t * GetAmixt() const
TGeoElement * GetElement(Int_t i=0) const override
Retrieve the pointer to the element corresponding to component I.
static Double_t Big()
Definition TGeoShape.h:87
A doubly linked list.
Definition TList.h:38
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:576
void Add(TObject *obj) override
Definition TList.h:83
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
TNamed()
Definition TNamed.h:36
TString fName
Definition TNamed.h:32
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
Int_t IndexOf(const TObject *obj) const override
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
void Add(TObject *obj) override
Definition TObjArray.h:68
virtual Bool_t IsEqual(const TObject *obj) const
Default equal comparison (objects are equal if they have the same address in memory).
Definition TObject.cxx:565
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition TObject.cxx:457
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1015
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
Basic string class.
Definition TString.h:139
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition TString.cxx:1163
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2356
static constexpr double fine_structure_const
static constexpr double cm2
static constexpr double Avogadro
static constexpr double cm
static constexpr double amu
static constexpr double cm
static constexpr double fine_structure_const
static constexpr double cm2
static constexpr double Avogadro
static constexpr double g
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:709
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:756
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123