Logo ROOT   6.07/09
Reference Guide
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 Geometry_classes
14 
15 Base class describing materials.
16 
17 \image html geom_material.jpg
18 */
19 
20 #include "Riostream.h"
21 #include "TMath.h"
22 #include "TObjArray.h"
23 #include "TStyle.h"
24 #include "TList.h"
25 #include "TGeoManager.h"
26 #include "TGeoExtension.h"
27 #include "TGeoMaterial.h"
28 
29 // statics and globals
30 
32 
33 ////////////////////////////////////////////////////////////////////////////////
34 /// Default constructor
35 
37  :TNamed(), TAttFill(),
38  fIndex(0),
39  fA(0.),
40  fZ(0.),
41  fDensity(0.),
42  fRadLen(0.),
43  fIntLen(0.),
44  fTemperature(0.),
45  fPressure(0.),
46  fState(kMatStateUndefined),
47  fShader(NULL),
48  fCerenkov(NULL),
49  fElement(NULL),
50  fUserExtension(0),
51  fFWExtension(0)
52 {
53  SetUsed(kFALSE);
54  fIndex = -1;
55  fTemperature = STP_temperature;
56  fPressure = STP_pressure;
57  fState = kMatStateUndefined;
58 }
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// constructor
62 
64  :TNamed(name, ""), TAttFill(),
65  fIndex(0),
66  fA(0.),
67  fZ(0.),
68  fDensity(0.),
69  fRadLen(0.),
70  fIntLen(0.),
71  fTemperature(0.),
72  fPressure(0.),
73  fState(kMatStateUndefined),
74  fShader(NULL),
75  fCerenkov(NULL),
76  fElement(NULL),
77  fUserExtension(0),
78  fFWExtension(0)
79 {
80  fName = fName.Strip();
81  SetUsed(kFALSE);
82  fIndex = -1;
86 
87  if (!gGeoManager) {
88  gGeoManager = new TGeoManager("Geometry", "default geometry");
89  }
90  gGeoManager->AddMaterial(this);
91 }
92 
93 ////////////////////////////////////////////////////////////////////////////////
94 /// constructor
95 
97  Double_t rho, Double_t radlen, Double_t intlen)
98  :TNamed(name, ""), TAttFill(),
99  fIndex(0),
100  fA(a),
101  fZ(z),
102  fDensity(rho),
103  fRadLen(0.),
104  fIntLen(0.),
105  fTemperature(0.),
106  fPressure(0.),
108  fShader(NULL),
109  fCerenkov(NULL),
110  fElement(NULL),
111  fUserExtension(0),
112  fFWExtension(0)
113 {
114  fName = fName.Strip();
115  SetUsed(kFALSE);
116  fIndex = -1;
117  fA = a;
118  fZ = z;
119  fDensity = rho;
123  SetRadLen(radlen, intlen);
124  if (!gGeoManager) {
125  gGeoManager = new TGeoManager("Geometry", "default geometry");
126  }
127  if (fZ - Int_t(fZ) > 1E-3)
128  Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
129  if (GetElement()) GetElement()->SetUsed();
130  gGeoManager->AddMaterial(this);
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// Constructor with state, temperature and pressure.
135 
137  EGeoMaterialState state, Double_t temperature, Double_t pressure)
138  :TNamed(name, ""), TAttFill(),
139  fIndex(0),
140  fA(a),
141  fZ(z),
142  fDensity(rho),
143  fRadLen(0.),
144  fIntLen(0.),
145  fTemperature(temperature),
146  fPressure(pressure),
147  fState(state),
148  fShader(NULL),
149  fCerenkov(NULL),
150  fElement(NULL),
151  fUserExtension(0),
152  fFWExtension(0)
153 {
154  fName = fName.Strip();
155  SetUsed(kFALSE);
156  fIndex = -1;
157  SetRadLen(0,0);
158  if (!gGeoManager) {
159  gGeoManager = new TGeoManager("Geometry", "default geometry");
160  }
161  if (fZ - Int_t(fZ) > 1E-3)
162  Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
163  if (GetElement()) GetElement()->SetUsed();
164  gGeoManager->AddMaterial(this);
165 }
166 
167 ////////////////////////////////////////////////////////////////////////////////
168 /// constructor
169 
171  :TNamed(name, ""), TAttFill(),
172  fIndex(0),
173  fA(0.),
174  fZ(0.),
175  fDensity(rho),
176  fRadLen(0.),
177  fIntLen(0.),
178  fTemperature(0.),
179  fPressure(0.),
181  fShader(NULL),
182  fCerenkov(NULL),
183  fElement(elem),
184  fUserExtension(0),
185  fFWExtension(0)
186 {
187  fName = fName.Strip();
188  SetUsed(kFALSE);
189  fIndex = -1;
190  fA = elem->A();
191  fZ = elem->Z();
192  SetRadLen(0,0);
196  if (!gGeoManager) {
197  gGeoManager = new TGeoManager("Geometry", "default geometry");
198  }
199  if (fZ - Int_t(fZ) > 1E-3)
200  Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
201  if (GetElement()) GetElement()->SetUsed();
202  gGeoManager->AddMaterial(this);
203 }
204 
205 ////////////////////////////////////////////////////////////////////////////////
206 
208  TNamed(gm),
209  TAttFill(gm),
210  fIndex(gm.fIndex),
211  fA(gm.fA),
212  fZ(gm.fZ),
213  fDensity(gm.fDensity),
214  fRadLen(gm.fRadLen),
215  fIntLen(gm.fIntLen),
217  fPressure(gm.fPressure),
218  fState(gm.fState),
219  fShader(gm.fShader),
220  fCerenkov(gm.fCerenkov),
221  fElement(gm.fElement),
222  fUserExtension(gm.fUserExtension->Grab()),
223  fFWExtension(gm.fFWExtension->Grab())
224 
225 {
226  //copy constructor
227 }
228 
229 ////////////////////////////////////////////////////////////////////////////////
230 ///assignment operator
231 
233 {
234  if(this!=&gm) {
235  TNamed::operator=(gm);
236  TAttFill::operator=(gm);
237  fIndex=gm.fIndex;
238  fA=gm.fA;
239  fZ=gm.fZ;
240  fDensity=gm.fDensity;
241  fRadLen=gm.fRadLen;
242  fIntLen=gm.fIntLen;
244  fPressure=gm.fPressure;
245  fState=gm.fState;
246  fShader=gm.fShader;
247  fCerenkov=gm.fCerenkov;
248  fElement=gm.fElement;
251  }
252  return *this;
253 }
254 
255 ////////////////////////////////////////////////////////////////////////////////
256 /// Destructor
257 
259 {
262 }
263 
264 ////////////////////////////////////////////////////////////////////////////////
265 /// Connect user-defined extension to the material. The material "grabs" a copy, so
266 /// the original object can be released by the producer. Release the previously
267 /// connected extension if any.
268 ///==========================================================================
269 /// NOTE: This interface is intended for user extensions and is guaranteed not
270 /// to be used by TGeo
271 ///==========================================================================
272 
274 {
276  fUserExtension = 0;
277  if (ext) fUserExtension = ext->Grab();
278 }
279 
280 ////////////////////////////////////////////////////////////////////////////////
281 /// Connect framework defined extension to the material. The material "grabs" a copy,
282 /// so the original object can be released by the producer. Release the previously
283 /// connected extension if any.
284 ///==========================================================================
285 /// NOTE: This interface is intended for the use by TGeo and the users should
286 /// NOT connect extensions using this method
287 ///==========================================================================
288 
290 {
292  fFWExtension = 0;
293  if (ext) fFWExtension = ext->Grab();
294 }
295 
296 ////////////////////////////////////////////////////////////////////////////////
297 /// Get a copy of the user extension pointer. The user must call Release() on
298 /// the copy pointer once this pointer is not needed anymore (equivalent to
299 /// delete() after calling new())
300 
302 {
303  if (fUserExtension) return fUserExtension->Grab();
304  return 0;
305 }
306 
307 ////////////////////////////////////////////////////////////////////////////////
308 /// Get a copy of the framework extension pointer. The user must call Release() on
309 /// the copy pointer once this pointer is not needed anymore (equivalent to
310 /// delete() after calling new())
311 
313 {
314  if (fFWExtension) return fFWExtension->Grab();
315  return 0;
316 }
317 
318 ////////////////////////////////////////////////////////////////////////////////
319 /// Provide a pointer name containing uid.
320 
322 {
323  static TString name;
324  name = TString::Format("pMat%d", GetUniqueID());
325  return (char*)name.Data();
326 }
327 
328 ////////////////////////////////////////////////////////////////////////////////
329 /// Set radiation/absorption lengths. If the values are negative, their absolute value
330 /// is taken, otherwise radlen is recomputed using G3 formula.
331 
333 {
334  fRadLen = TMath::Abs(radlen);
335  fIntLen = TMath::Abs(intlen);
336  // Check for vacuum
337  if (fA<0.9 || fZ<0.9) {
338  if (radlen<-1e5 || intlen<-1e-5) {
339  Error("SetRadLen","Material %s: user values taken for vacuum: radlen=%g or intlen=%g - too small", GetName(),fRadLen, fIntLen);
340  return;
341  }
342  // Ignore positive values and take big numbers
343  if (radlen>=0) fRadLen = 1.E30;
344  if (intlen>=0) fIntLen = 1.E30;
345  return;
346  }
347  // compute radlen systematically with G3 formula for a valid material
348  if (radlen>=0) {
349  //taken grom Geant3 routine GSMATE
350  const Double_t alr2av=1.39621E-03, al183=5.20948;
352  (al183-TMath::Log(fZ)/3-TGeoMaterial::Coulomb(fZ)));
353  }
354  // Compute interaction length using the same formula as in GEANT4
355  if (intlen>=0) {
356  const Double_t cm = 1.;
357  const Double_t g = 6.2415e21; // [gram = 1E-3*joule*s*s/(m*m)]
358  const Double_t amu = 1.03642688246781065e-02; // [MeV/c^2]
359  const Double_t lambda0 = 35.*g/(cm*cm); // [g/cm^2]
360  Double_t nilinv = 0.0;
361  TGeoElement *elem = GetElement();
362  if (!elem) {
363  Fatal("SetRadLen", "Element not found for material %s", GetName());
364  return;
365  }
366  Double_t nbAtomsPerVolume = TMath::Na()*fDensity/elem->A();
367  nilinv += nbAtomsPerVolume*TMath::Power(elem->Neff(), 0.6666667);
368  nilinv *= amu/lambda0;
369  fIntLen = (nilinv<=0) ? TGeoShape::Big() : (1./nilinv);
370  }
371 }
372 
373 ////////////////////////////////////////////////////////////////////////////////
374 /// static function
375 /// Compute Coulomb correction for pair production and Brem
376 /// REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
377 /// FORMULA 2.7.17
378 
380 {
381  const Double_t alpha = 7.29927E-03;
382 
383  Double_t az = alpha*z;
384  Double_t az2 = az*az;
385  Double_t az4 = az2 * az2;
386  Double_t fp = ( 0.0083*az4 + 0.20206 + 1./(1.+az2) ) * az2;
387  Double_t fm = ( 0.0020*az4 + 0.0369 ) * az4;
388  return fp - fm;
389 }
390 
391 ////////////////////////////////////////////////////////////////////////////////
392 /// return true if the other material has the same physical properties
393 
395 {
396  if (other==this) return kTRUE;
397  if (other->IsMixture()) return kFALSE;
398  if (TMath::Abs(fA-other->GetA())>1E-3) return kFALSE;
399  if (TMath::Abs(fZ-other->GetZ())>1E-3) return kFALSE;
400  if (TMath::Abs(fDensity-other->GetDensity())>1E-6) return kFALSE;
401  if (GetCerenkovProperties() != other->GetCerenkovProperties()) return kFALSE;
402 // if (fRadLen != other->GetRadLen()) return kFALSE;
403 // if (fIntLen != other->GetIntLen()) return kFALSE;
404  return kTRUE;
405 }
406 
407 ////////////////////////////////////////////////////////////////////////////////
408 /// print characteristics of this material
409 
410 void TGeoMaterial::Print(const Option_t * /*option*/) const
411 {
412  printf("Material %s %s A=%g Z=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(),
414 }
415 
416 ////////////////////////////////////////////////////////////////////////////////
417 /// Save a primitive as a C++ statement(s) on output stream "out".
418 
419 void TGeoMaterial::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
420 {
422  char *name = GetPointerName();
423  out << "// Material: " << GetName() << std::endl;
424  out << " a = " << fA << ";" << std::endl;
425  out << " z = " << fZ << ";" << std::endl;
426  out << " density = " << fDensity << ";" << std::endl;
427  out << " radl = " << fRadLen << ";" << std::endl;
428  out << " absl = " << fIntLen << ";" << std::endl;
429 
430  out << " " << name << " = new TGeoMaterial(\"" << GetName() << "\", a,z,density,radl,absl);" << std::endl;
431  out << " " << name << "->SetIndex(" << GetIndex() << ");" << std::endl;
433 }
434 
435 ////////////////////////////////////////////////////////////////////////////////
436 /// Get some default color related to this material.
437 
439 {
440  Int_t id = 1+ gGeoManager->GetListOfMaterials()->IndexOf(this);
441  return (2+id%6);
442 }
443 
444 ////////////////////////////////////////////////////////////////////////////////
445 /// Get a pointer to the element this material is made of.
446 
448 {
449  if (fElement) return fElement;
451  return table->GetElement(Int_t(fZ));
452 }
453 ////////////////////////////////////////////////////////////////////////////////
454 /// Single interface to get element properties.
455 
457 {
458  a = fA;
459  z = fZ;
460  w = 1.;
461 }
462 
463 ////////////////////////////////////////////////////////////////////////////////
464 /// Retrieve material index in the list of materials
465 
467 {
468  if (fIndex>=0) return fIndex;
469  TList *matlist = gGeoManager->GetListOfMaterials();
470  fIndex = matlist->IndexOf(this);
471  return fIndex;
472 }
473 
474 ////////////////////////////////////////////////////////////////////////////////
475 /// Create the material representing the decay product of this material at a
476 /// given time. The precision represent the minimum cumulative branching ratio for
477 /// which decay products are still taken into account.
478 
480 {
481  TObjArray *pop = new TObjArray();
482  if (!fElement || !fElement->IsRadioNuclide()) return this;
483  FillMaterialEvolution(pop, precision);
484  Int_t ncomp = pop->GetEntriesFast();
485  if (!ncomp) return this;
486  TGeoElementRN *el;
487  Double_t *weight = new Double_t[ncomp];
488  Double_t amed = 0.;
489  Int_t i;
490  for (i=0; i<ncomp; i++) {
491  el = (TGeoElementRN *)pop->At(i);
492  weight[i] = el->Ratio()->Concentration(time) * el->A();
493  amed += weight[i];
494  }
495  Double_t rho = fDensity*amed/fA;
496  TGeoMixture *mix = 0;
497  Int_t ncomp1 = ncomp;
498  for (i=0; i<ncomp; i++) {
499  if ((weight[i]/amed)<precision) {
500  amed -= weight[i];
501  ncomp1--;
502  }
503  }
504  if (ncomp1<2) {
505  el = (TGeoElementRN *)pop->At(0);
506  delete [] weight;
507  delete pop;
508  if (ncomp1==1) return new TGeoMaterial(TString::Format("%s-evol",GetName()), el, rho);
509  return NULL;
510  }
511  mix = new TGeoMixture(TString::Format("%s-evol",GetName()), ncomp, rho);
512  for (i=0; i<ncomp; i++) {
513  weight[i] /= amed;
514  if (weight[i]<precision) continue;
515  el = (TGeoElementRN *)pop->At(i);
516  mix->AddElement(el, weight[i]);
517  }
518  delete [] weight;
519  delete pop;
520  return mix;
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// Fills a user array with all the elements deriving from the possible
525 /// decay of the top element composing the mixture. Each element contained
526 /// by <population> may be a radionuclide having a Bateman solution attached.
527 /// The precision represent the minimum cumulative branching ratio for
528 /// which decay products are still taken into account.
529 /// To visualize the time evolution of each decay product one can use:
530 /// ~~~ {.cpp}
531 /// TGeoElement *elem = population->At(index);
532 /// TGeoElementRN *elemrn = 0;
533 /// if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
534 /// ~~~
535 /// One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
536 /// of one of the decay products, N1(0) is the number of atoms of the top
537 /// element at t=0.
538 /// ~~~ {.cpp}
539 /// Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
540 /// ~~~
541 /// One can also display the time evolution of the fractional weight:
542 /// ~~~ {.cpp}
543 /// elemrn->Ratio()->Draw(option);
544 /// ~~~
545 
547 {
548  if (population->GetEntriesFast()) {
549  Error("FillMaterialEvolution", "Provide an empty array !");
550  return;
551  }
553  TGeoElement *elem;
554  TGeoElementRN *elemrn;
555  TIter next(table->GetElementsRN());
556  while ((elemrn=(TGeoElementRN*)next())) elemrn->ResetRatio();
557  elem = GetElement();
558  if (!elem) {
559  Fatal("FillMaterialEvolution", "Element not found for material %s", GetName());
560  return;
561  }
562  if (!elem->IsRadioNuclide()) {
563  population->Add(elem);
564  return;
565  }
566  elemrn = (TGeoElementRN*)elem;
567  elemrn->FillPopulation(population, precision);
568 }
569 
570 
571 /*************************************************************************
572  * TGeoMixture - mixtures of elements
573  *
574  *************************************************************************/
576 
577 ////////////////////////////////////////////////////////////////////////////////
578 /// Default constructor
579 
581 {
582  fNelements = 0;
583  fZmixture = 0;
584  fAmixture = 0;
585  fWeights = 0;
586  fNatoms = 0;
587  fElements = 0;
588 }
589 
590 ////////////////////////////////////////////////////////////////////////////////
591 /// constructor
592 
593 TGeoMixture::TGeoMixture(const char *name, Int_t /*nel*/, Double_t rho)
594  :TGeoMaterial(name)
595 {
596  fZmixture = 0;
597  fAmixture = 0;
598  fWeights = 0;
599  fNelements = 0;
600  fNatoms = 0;
601  fDensity = rho;
602  fElements = 0;
603  if (fDensity < 0) fDensity = 0.001;
604 }
605 
606 ////////////////////////////////////////////////////////////////////////////////
607 ///copy constructor
608 
609 TGeoMixture::TGeoMixture(const TGeoMixture& gm) :
610  TGeoMaterial(gm),
612  fZmixture(gm.fZmixture),
613  fAmixture(gm.fAmixture),
614  fWeights(gm.fWeights),
615  fNatoms(gm.fNatoms),
616  fElements(gm.fElements)
617 {
618 }
619 
620 ////////////////////////////////////////////////////////////////////////////////
621 ///assignment operator
622 
623 TGeoMixture& TGeoMixture::operator=(const TGeoMixture& gm)
624 {
625  if(this!=&gm) {
628  fZmixture=gm.fZmixture;
629  fAmixture=gm.fAmixture;
630  fWeights=gm.fWeights;
631  fNatoms = gm.fNatoms;
632  fElements = gm.fElements;
633  }
634  return *this;
635 }
636 
637 ////////////////////////////////////////////////////////////////////////////////
638 /// Destructor
639 
641 {
642  if (fZmixture) delete[] fZmixture;
643  if (fAmixture) delete[] fAmixture;
644  if (fWeights) delete[] fWeights;
645  if (fNatoms) delete[] fNatoms;
646  if (fElements) delete fElements;
647 }
648 
649 ////////////////////////////////////////////////////////////////////////////////
650 /// Compute effective A/Z and radiation length
651 
653 {
654  const Double_t alr2av = 1.39621E-03 , al183 =5.20948;
655  const Double_t cm = 1.;
656  const Double_t g = 6.2415e21; // [gram = 1E-3*joule*s*s/(m*m)]
657  const Double_t amu = 1.03642688246781065e-02; // [MeV/c^2]
658  const Double_t lambda0 = 35.*g/(cm*cm); // [g/cm^2]
659  Double_t radinv = 0.0;
660  Double_t nilinv = 0.0;
661  Double_t nbAtomsPerVolume;
662  fA = 0;
663  fZ = 0;
664  for (Int_t j=0;j<fNelements;j++) {
665  if (fWeights[j] <= 0) continue;
666  fA += fWeights[j]*fAmixture[j];
667  fZ += fWeights[j]*fZmixture[j];
668  nbAtomsPerVolume = TMath::Na()*fDensity*fWeights[j]/GetElement(j)->A();
669  nilinv += nbAtomsPerVolume*TMath::Power(GetElement(j)->Neff(), 0.6666667);
670  Double_t zc = fZmixture[j];
671  Double_t alz = TMath::Log(zc)/3.;
672  Double_t xinv = zc*(zc+TGeoMaterial::ScreenFactor(zc))*
673  (al183-alz-TGeoMaterial::Coulomb(zc))/fAmixture[j];
674  radinv += xinv*fWeights[j];
675  }
676  radinv *= alr2av*fDensity;
677  if (radinv > 0) fRadLen = 1/radinv;
678  // Compute interaction length
679  nilinv *= amu/lambda0;
680  fIntLen = (nilinv<=0) ? TGeoShape::Big() : (1./nilinv);
681 }
682 
683 ////////////////////////////////////////////////////////////////////////////////
684 /// add an element to the mixture using fraction by weight
685 /// Check if the element is already defined
686 
688 {
690  if (z<1 || z>table->GetNelements()-1)
691  Fatal("AddElement", "Cannot add element having Z=%d to mixture %s", (Int_t)z, GetName());
692  Int_t i;
693  for (i=0; i<fNelements; i++) {
694  if (TMath::Abs(z-fZmixture[i])<1.e-6 && TMath::Abs(a-fAmixture[i])<1.e-6) {
695  fWeights[i] += weight;
697  return;
698  }
699  }
700  if (!fNelements) {
701  fZmixture = new Double_t[1];
702  fAmixture = new Double_t[1];
703  fWeights = new Double_t[1];
704  } else {
705  Int_t nelements = fNelements+1;
706  Double_t *zmixture = new Double_t[nelements];
707  Double_t *amixture = new Double_t[nelements];
708  Double_t *weights = new Double_t[nelements];
709  for (Int_t j=0; j<fNelements; j++) {
710  zmixture[j] = fZmixture[j];
711  amixture[j] = fAmixture[j];
712  weights[j] = fWeights[j];
713  }
714  delete [] fZmixture;
715  delete [] fAmixture;
716  delete [] fWeights;
717  fZmixture = zmixture;
718  fAmixture = amixture;
719  fWeights = weights;
720  }
721 
722  fNelements++;
723  i = fNelements - 1;
724  fZmixture[i] = z;
725  fAmixture[i] = a;
726  fWeights[i] = weight;
727  if (z - Int_t(z) > 1E-3)
728  Warning("DefineElement", "Mixture %s has element defined with fractional Z=%f", GetName(), z);
729  GetElement(i)->SetDefined();
730  table->GetElement((Int_t)z)->SetDefined();
731 
732  //compute equivalent radiation length (taken from Geant3/GSMIXT)
734 }
735 
736 ////////////////////////////////////////////////////////////////////////////////
737 /// Define one component of the mixture as an existing material/mixture.
738 
740 {
741  TGeoElement *elnew, *elem;
742  Double_t a,z;
743  if (!mat->IsMixture()) {
744  elem = mat->GetBaseElement();
745  if (elem) {
746  AddElement(elem, weight);
747  } else {
748  a = mat->GetA();
749  z = mat->GetZ();
750  AddElement(a, z, weight);
751  }
752  return;
753  }
754  // The material is a mixture.
755  TGeoMixture *mix = (TGeoMixture*)mat;
756  Double_t wnew;
757  Int_t nelem = mix->GetNelements();
758  Bool_t elfound;
759  Int_t i,j;
760  // loop the elements of the daughter mixture
761  for (i=0; i<nelem; i++) {
762  elfound = kFALSE;
763  elnew = mix->GetElement(i);
764  if (!elnew) continue;
765  // check if we have the element already defined in the parent mixture
766  for (j=0; j<fNelements; j++) {
767  if (fWeights[j]<=0) continue;
768  elem = GetElement(j);
769  if (elem == elnew) {
770  // element found, compute new weight
771  fWeights[j] += weight * (mix->GetWmixt())[i];
772  elfound = kTRUE;
773  break;
774  }
775  }
776  if (elfound) continue;
777  // element not found, define it
778  wnew = weight * (mix->GetWmixt())[i];
779  AddElement(elnew, wnew);
780  }
781 }
782 
783 ////////////////////////////////////////////////////////////////////////////////
784 /// add an element to the mixture using fraction by weight
785 
787 {
788  TGeoElement *elemold;
790  if (!fElements) fElements = new TObjArray(128);
791  Bool_t exist = kFALSE;
792  // If previous elements were defined by A/Z, add corresponding TGeoElements
793  for (Int_t i=0; i<fNelements; i++) {
794  elemold = (TGeoElement*)fElements->At(i);
795  if (!elemold) fElements->AddAt(elemold = table->GetElement((Int_t)fZmixture[i]), i);
796  if (elemold == elem) exist = kTRUE;
797  }
798  if (!exist) fElements->AddAtAndExpand(elem, fNelements);
799  AddElement(elem->A(), elem->Z(), weight);
800 }
801 
802 ////////////////////////////////////////////////////////////////////////////////
803 /// Add a mixture element by number of atoms in the chemical formula.
804 
806 {
807  Int_t i,j;
808  Double_t amol;
809  TGeoElement *elemold;
811  if (!fElements) fElements = new TObjArray(128);
812  // Check if the element is already defined
813  for (i=0; i<fNelements; i++) {
814  elemold = (TGeoElement*)fElements->At(i);
815  if (!elemold) fElements->AddAt(table->GetElement((Int_t)fZmixture[i]), i);
816  else if (elemold != elem) continue;
817  if ((elem==elemold) ||
818  (TMath::Abs(elem->Z()-fZmixture[i])<1.e-6 && TMath::Abs(elem->A()-fAmixture[i])<1.e-6)) {
819  fNatoms[i] += natoms;
820  amol = 0.;
821  for (j=0; j<fNelements; j++) amol += fAmixture[j]*fNatoms[j];
822  for (j=0; j<fNelements; j++) fWeights[j] = fNatoms[j]*fAmixture[j]/amol;
824  return;
825  }
826  }
827  // New element
828  if (!fNelements) {
829  fZmixture = new Double_t[1];
830  fAmixture = new Double_t[1];
831  fWeights = new Double_t[1];
832  fNatoms = new Int_t[1];
833  } else {
834  if (!fNatoms) {
835  Fatal("AddElement", "Cannot add element by natoms in mixture %s after defining elements by weight",
836  GetName());
837  return;
838  }
839  Int_t nelements = fNelements+1;
840  Double_t *zmixture = new Double_t[nelements];
841  Double_t *amixture = new Double_t[nelements];
842  Double_t *weights = new Double_t[nelements];
843  Int_t *nnatoms = new Int_t[nelements];
844  for (j=0; j<fNelements; j++) {
845  zmixture[j] = fZmixture[j];
846  amixture[j] = fAmixture[j];
847  weights[j] = fWeights[j];
848  nnatoms[j] = fNatoms[j];
849  }
850  delete [] fZmixture;
851  delete [] fAmixture;
852  delete [] fWeights;
853  delete [] fNatoms;
854  fZmixture = zmixture;
855  fAmixture = amixture;
856  fWeights = weights;
857  fNatoms = nnatoms;
858  }
859  fNelements++;
860  Int_t iel = fNelements-1;
861  fZmixture[iel] = elem->Z();
862  fAmixture[iel] = elem->A();
863  fNatoms[iel] = natoms;
864  fElements->AddAtAndExpand(elem, iel);
865  amol = 0.;
866  for (i=0; i<fNelements; i++) {
867  if (fNatoms[i]<=0) return;
868  amol += fAmixture[i]*fNatoms[i];
869  }
870  for (i=0; i<fNelements; i++) fWeights[i] = fNatoms[i]*fAmixture[i]/amol;
871  table->GetElement(elem->Z())->SetDefined();
873 }
874 
875 ////////////////////////////////////////////////////////////////////////////////
876 /// Define the mixture element at index iel by number of atoms in the chemical formula.
877 
879 {
881  TGeoElement *elem = table->GetElement(z);
882  if (!elem) {
883  Fatal("DefineElement", "In mixture %s, element with Z=%i not found",GetName(),z);
884  return;
885  }
886  AddElement(elem, natoms);
887 }
888 
889 ////////////////////////////////////////////////////////////////////////////////
890 /// Retrieve the pointer to the element corresponding to component I.
891 
893 {
894  if (i<0 || i>=fNelements) {
895  Error("GetElement", "Mixture %s has only %d elements", GetName(), fNelements);
896  return 0;
897  }
898  TGeoElement *elem = 0;
899  if (fElements) elem = (TGeoElement*)fElements->At(i);
900  if (elem) return elem;
902  return table->GetElement(Int_t(fZmixture[i]));
903 }
904 
905 ////////////////////////////////////////////////////////////////////////////////
906 /// Get specific activity (in Bq/gram) for the whole mixture (no argument) or
907 /// for a given component.
908 
910 {
911  if (i>=0 && i<fNelements) return fWeights[i]*GetElement(i)->GetSpecificActivity();
912  Double_t sa = 0;
913  for (Int_t iel=0; iel<fNelements; iel++) {
914  sa += fWeights[iel]*GetElement(iel)->GetSpecificActivity();
915  }
916  return sa;
917 }
918 
919 ////////////////////////////////////////////////////////////////////////////////
920 /// Return true if the other material has the same physical properties
921 
923 {
924  if (other->IsEqual(this)) return kTRUE;
925  if (!other->IsMixture()) return kFALSE;
926  TGeoMixture *mix = (TGeoMixture*)other;
927  if (!mix) return kFALSE;
928  if (fNelements != mix->GetNelements()) return kFALSE;
929  if (TMath::Abs(fA-other->GetA())>1E-3) return kFALSE;
930  if (TMath::Abs(fZ-other->GetZ())>1E-3) return kFALSE;
931  if (TMath::Abs(fDensity-other->GetDensity())>1E-6) return kFALSE;
932  if (GetCerenkovProperties() != other->GetCerenkovProperties()) return kFALSE;
933 // if (fRadLen != other->GetRadLen()) return kFALSE;
934 // if (fIntLen != other->GetIntLen()) return kFALSE;
935  for (Int_t i=0; i<fNelements; i++) {
936  if (TMath::Abs(fZmixture[i]-(mix->GetZmixt())[i])>1E-3) return kFALSE;
937  if (TMath::Abs(fAmixture[i]-(mix->GetAmixt())[i])>1E-3) return kFALSE;
938  if (TMath::Abs(fWeights[i]-(mix->GetWmixt())[i])>1E-3) return kFALSE;
939  }
940  return kTRUE;
941 }
942 
943 ////////////////////////////////////////////////////////////////////////////////
944 /// print characteristics of this material
945 
946 void TGeoMixture::Print(const Option_t * /*option*/) const
947 {
948  printf("Mixture %s %s Aeff=%g Zeff=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(),
950  for (Int_t i=0; i<fNelements; i++) {
951  if (fNatoms) printf(" Element #%i : %s Z=%6.2f A=%6.2f w=%6.3f natoms=%d\n", i, GetElement(i)->GetName(),fZmixture[i],
952  fAmixture[i], fWeights[i], fNatoms[i]);
953  else printf(" Element #%i : %s Z=%6.2f A=%6.2f w=%6.3f\n", i, GetElement(i)->GetName(),fZmixture[i],
954  fAmixture[i], fWeights[i]);
955  }
956 }
957 
958 ////////////////////////////////////////////////////////////////////////////////
959 /// Save a primitive as a C++ statement(s) on output stream "out".
960 
961 void TGeoMixture::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
962 {
964  char *name = GetPointerName();
965  out << "// Mixture: " << GetName() << std::endl;
966  out << " nel = " << fNelements << ";" << std::endl;
967  out << " density = " << fDensity << ";" << std::endl;
968  out << " " << name << " = new TGeoMixture(\"" << GetName() << "\", nel,density);" << std::endl;
969  for (Int_t i=0; i<fNelements; i++) {
970  TGeoElement *el = GetElement(i);
971  out << " a = " << fAmixture[i] << "; z = "<< fZmixture[i] << "; w = " << fWeights[i] << "; // " << el->GetName() << std::endl;
972  out << " " << name << "->DefineElement(" << i << ",a,z,w);" << std::endl;
973  }
974  out << " " << name << "->SetIndex(" << GetIndex() << ");" << std::endl;
976 }
977 
978 ////////////////////////////////////////////////////////////////////////////////
979 /// Create the mixture representing the decay product of this material at a
980 /// given time. The precision represent the minimum cumulative branching ratio for
981 /// which decay products are still taken into account.
982 
984 {
985  TObjArray *pop = new TObjArray();
986  FillMaterialEvolution(pop, precision);
987  Int_t ncomp = pop->GetEntriesFast();
988  if (!ncomp) return this;
989  TGeoElement *elem;
990  TGeoElementRN *el;
991  Double_t *weight = new Double_t[ncomp];
992  Double_t amed = 0.;
993  Int_t i, j;
994  for (i=0; i<ncomp; i++) {
995  elem = (TGeoElement *)pop->At(i);
996  if (!elem->IsRadioNuclide()) {
997  j = fElements->IndexOf(elem);
998  weight[i] = fWeights[j]*fAmixture[0]/fWeights[0];
999  } else {
1000  el = (TGeoElementRN*)elem;
1001  weight[i] = el->Ratio()->Concentration(time) * el->A();
1002  }
1003  amed += weight[i];
1004  }
1005  Double_t rho = fDensity * fWeights[0] * amed/fAmixture[0];
1006  TGeoMixture *mix = 0;
1007  Int_t ncomp1 = ncomp;
1008  for (i=0; i<ncomp; i++) {
1009  if ((weight[i]/amed)<precision) {
1010  amed -= weight[i];
1011  ncomp1--;
1012  }
1013  }
1014  if (ncomp1<2) {
1015  el = (TGeoElementRN *)pop->At(0);
1016  delete [] weight;
1017  delete pop;
1018  if (ncomp1==1) return new TGeoMaterial(TString::Format("%s-evol",GetName()), el, rho);
1019  return NULL;
1020  }
1021  mix = new TGeoMixture(TString::Format("%s-evol",GetName()), ncomp, rho);
1022  for (i=0; i<ncomp; i++) {
1023  weight[i] /= amed;
1024  if (weight[i]<precision) continue;
1025  el = (TGeoElementRN *)pop->At(i);
1026  mix->AddElement(el, weight[i]);
1027  }
1028  delete [] weight;
1029  delete pop;
1030  return mix;
1031 }
1032 
1033 ////////////////////////////////////////////////////////////////////////////////
1034 /// Fills a user array with all the elements deriving from the possible
1035 /// decay of the top elements composing the mixture. Each element contained
1036 /// by <population> may be a radionuclide having a Bateman solution attached.
1037 /// The precision represent the minimum cumulative branching ratio for
1038 /// which decay products are still taken into account.
1039 /// To visualize the time evolution of each decay product one can use:
1040 /// ~~~ {.cpp}
1041 /// TGeoElement *elem = population->At(index);
1042 /// TGeoElementRN *elemrn = 0;
1043 /// if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
1044 /// ~~~
1045 /// One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
1046 /// of one of the decay products, N1(0) is the number of atoms of the first top
1047 /// element at t=0.
1048 /// ~~~ {.cpp}
1049 /// Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
1050 /// ~~~
1051 /// One can also display the time evolution of the fractional weight:
1052 /// ~~~ {.cpp}
1053 /// elemrn->Ratio()->Draw(option);
1054 /// ~~~
1055 
1057 {
1058  if (population->GetEntriesFast()) {
1059  Error("FillMaterialEvolution", "Provide an empty array !");
1060  return;
1061  }
1063  TGeoElement *elem;
1064  TGeoElementRN *elemrn;
1065  TIter next(table->GetElementsRN());
1066  while ((elemrn=(TGeoElementRN*)next())) elemrn->ResetRatio();
1067  Double_t factor;
1068  for (Int_t i=0; i<fNelements; i++) {
1069  elem = GetElement(i);
1070  if (!elem->IsRadioNuclide()) {
1071  population->Add(elem);
1072  continue;
1073  }
1074  elemrn = (TGeoElementRN*)elem;
1075  factor = fWeights[i]*fAmixture[0]/(fWeights[0]*fAmixture[i]);
1076  elemrn->FillPopulation(population, precision, factor);
1077  }
1078 }
1079 
1080 ////////////////////////////////////////////////////////////////////////////////
1081 /// static function
1082 /// Compute screening factor for pair production and Bremsstrahlung
1083 /// REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
1084 /// FORMULA 2.7.22
1085 
1087 {
1088  const Double_t al183= 5.20948 , al1440 = 7.27239;
1089  Double_t alz = TMath::Log(z)/3.;
1090  Double_t factor = (al1440 - 2*alz) / (al183 - alz - TGeoMaterial::Coulomb(z));
1091  return factor;
1092 }
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
Int_t * fNatoms
Definition: TGeoMaterial.h:153
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...
An array of TObjects.
Definition: TObjArray.h:39
The manager class for any TGeo geometry.
Definition: TGeoManager.h:38
TList * GetListOfMaterials() const
Definition: TGeoManager.h:461
void DefineElement(Int_t iel, Double_t a, Double_t z, Double_t weight)
Definition: TGeoMaterial.h:199
Table of elements.
Definition: TGeoElement.h:367
Double_t * fWeights
Definition: TGeoMaterial.h:152
Double_t Log(Double_t x)
Definition: TMath.h:526
TGeoExtension * fFWExtension
Transient user-defined extension to materials.
Definition: TGeoMaterial.h:64
Double_t Concentration(Double_t time) const
Find concentration of the element at a given time.
TObject * fCerenkov
Definition: TGeoMaterial.h:61
const char Option_t
Definition: RtypesCore.h:62
Int_t Z() const
Definition: TGeoElement.h:76
Int_t fNelements
Definition: TGeoMaterial.h:149
Double_t * GetZmixt() const
Definition: TGeoMaterial.h:183
TGeoElementTable * GetElementTable()
Returns material table. Creates it if not existing.
void SetRadLen(Double_t radlen, Double_t intlen=0.)
Set radiation/absorption lengths.
Basic string class.
Definition: TString.h:137
Base class describing materials.
Definition: TGeoMaterial.h:35
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Int_t GetDefaultColor() const
Get some default color related to this material.
virtual ~TGeoMaterial()
Destructor.
Int_t GetEntriesFast() const
Definition: TObjArray.h:66
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:528
void AverageProperties()
Compute effective A/Z and radiation length.
static const Double_t STP_pressure
Definition: TGeoMaterial.h:33
virtual TGeoElement * GetElement(Int_t i=0) const
Get a pointer to the element this material is made of.
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
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 elements composi...
const char * Data() const
Definition: TString.h:349
virtual Double_t GetSpecificActivity(Int_t i=-1) const
Get specific activity (in Bq/gram) for the whole mixture (no argument) or for a given component...
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:953
Double_t A() const
Definition: TGeoElement.h:79
EGeoMaterialState fState
Definition: TGeoMaterial.h:59
Fill Area Attributes class.
Definition: TAttFill.h:24
TGeoBatemanSol * Ratio() const
Definition: TGeoElement.h:195
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:2335
Double_t * GetWmixt() const
Definition: TGeoMaterial.h:185
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
void SetUserExtension(TGeoExtension *ext)
Connect user-defined extension to the material.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
virtual Bool_t IsMixture() const
Definition: TGeoMaterial.h:114
Base class for chemical elements.
Definition: TGeoElement.h:46
Double_t Neff() const
Returns effective number of nucleons.
virtual TGeoElement * GetElement(Int_t i=0) const
Retrieve the pointer to the element corresponding to component I.
TGeoElement * GetBaseElement() const
Definition: TGeoMaterial.h:97
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
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...
virtual Double_t GetSpecificActivity() const
Definition: TGeoElement.h:85
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
Definition: TGeoMaterial.h:90
A doubly linked list.
Definition: TList.h:47
Int_t GetIndex()
Retrieve material index in the list of materials.
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 fIntLen
Definition: TGeoMaterial.h:56
virtual Int_t GetNelements() const
Definition: TGeoMaterial.h:182
virtual Bool_t IsEq(const TGeoMaterial *other) const
Return true if the other material has the same physical properties.
TGeoMaterial()
Default constructor.
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition: TNamed.cxx:42
Int_t IndexOf(const TObject *obj) const
Definition: TObjArray.cxx:552
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:222
TGeoMixture()
Default constructor.
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.
Double_t fDensity
Definition: TGeoMaterial.h:54
char * GetPointerName() const
Provide a pointer name containing uid.
TGeoMixture & operator=(const TGeoMixture &)
assignment operator
TGeoExtension * GrabUserExtension() const
Get a copy of the user extension pointer.
ABC for user objects attached to TGeoVolume or TGeoNode.
Definition: TGeoExtension.h:21
Double_t fTemperature
Definition: TGeoMaterial.h:57
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
TGeoElement * fElement
Definition: TGeoMaterial.h:62
Double_t E()
Definition: TMath.h:54
void SetFWExtension(TGeoExtension *ext)
Connect framework defined extension to the material.
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
TSubString Strip(EStripType s=kTrailing, char c= ' ') const
Return a substring of self stripped at beginning and/or end.
Definition: TString.cxx:1070
TString fName
Definition: TNamed.h:36
virtual void AddAt(TObject *obj, Int_t idx)
Add object at position ids.
Definition: TObjArray.cxx:239
static const Double_t STP_temperature
Definition: TGeoMaterial.h:32
virtual TObject * GetCerenkovProperties() const
Definition: TGeoMaterial.h:102
Double_t * fZmixture
Definition: TGeoMaterial.h:150
TObjArray * GetElementsRN() const
Definition: TGeoElement.h:411
Double_t fZ
Definition: TGeoMaterial.h:53
static Double_t ScreenFactor(Double_t z)
static function Compute screening factor for pair production and Bremsstrahlung REFERENCE : EGS MANUA...
void SetUsed(Bool_t flag=kTRUE)
Definition: TGeoElement.h:92
virtual ~TGeoMixture()
Destructor.
#define ClassImp(name)
Definition: Rtypes.h:279
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:554
virtual void Release() const =0
double Double_t
Definition: RtypesCore.h:55
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
virtual Bool_t IsRadioNuclide() const
Definition: TGeoElement.h:88
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
virtual void Print(const Option_t *option="") const
print characteristics of this material
virtual void GetElementProp(Double_t &a, Double_t &z, Double_t &w, Int_t i=0)
Single interface to get element properties.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
void SetDefined(Bool_t flag=kTRUE)
Definition: TGeoElement.h:91
Double_t * GetAmixt() const
Definition: TGeoMaterial.h:184
static Double_t Big()
Definition: TGeoShape.h:90
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition: TObject.cxx:434
Double_t Na()
Definition: TMath.h:104
virtual Double_t GetZ() const
Definition: TGeoMaterial.h:91
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Double_t fPressure
Definition: TGeoMaterial.h:58
TObjArray * fElements
Definition: TGeoMaterial.h:154
Double_t fRadLen
Definition: TGeoMaterial.h:55
virtual TGeoExtension * Grab()=0
virtual Double_t GetDensity() const
Definition: TGeoMaterial.h:93
void SetUsed(Bool_t flag=kTRUE)
Definition: TGeoMaterial.h:123
TGeoExtension * fUserExtension
Definition: TGeoMaterial.h:63
Class representing a radionuclide.
Definition: TGeoElement.h:136
Double_t fA
Definition: TGeoMaterial.h:52
#define NULL
Definition: Rtypes.h:82
Int_t AddMaterial(const TGeoMaterial *material)
Add a material to the list. Returns index of the material in list.
TGeoMaterial & operator=(const TGeoMaterial &)
assignment operator
void Add(TObject *obj)
Definition: TObjArray.h:75
Int_t GetNelements() const
Definition: TGeoElement.h:415
virtual TGeoMaterial * DecayMaterial(Double_t time, Double_t precision=0.001)
Create the mixture representing the decay product of this material at a given time.
TGeoElement * GetElement(Int_t z)
Definition: TGeoElement.h:408
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
Double_t * fAmixture
Definition: TGeoMaterial.h:151
TObject * fShader
Definition: TGeoMaterial.h:60
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void Print(const Option_t *option="") const
print characteristics of this material
virtual Bool_t IsEq(const TGeoMaterial *other) const
return true if the other material has the same physical properties
char name[80]
Definition: TGX11.cxx:109
TGeoExtension * GrabFWExtension() const
Get a copy of the framework extension pointer.
void ResetRatio()
Clears the existing ratio.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:911