TGeoMaterial.cxx

Go to the documentation of this file.
00001 // @(#)root/geom:$Id: TGeoMaterial.cxx 36535 2010-11-08 14:41:54Z agheata $
00002 // Author: Andrei Gheata   25/10/01
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 ////////////////////////////////////////////////////////////////////////////////
00013 //
00014 //Begin_Html
00015 /*
00016 <img src="gif/t_material.jpg">
00017 */
00018 //End_Html
00019 #include "Riostream.h"
00020 #include "TMath.h"
00021 #include "TObjArray.h"
00022 #include "TStyle.h"
00023 #include "TList.h"
00024 #include "TGeoManager.h"
00025 #include "TGeoMaterial.h"
00026 
00027 // statics and globals
00028 
00029 ClassImp(TGeoMaterial)
00030 
00031 //_____________________________________________________________________________
00032 TGeoMaterial::TGeoMaterial()
00033              :TNamed(), TAttFill(),
00034               fIndex(0),
00035               fA(0.),
00036               fZ(0.),
00037               fDensity(0.),
00038               fRadLen(0.),
00039               fIntLen(0.),
00040               fTemperature(0.),
00041               fPressure(0.),
00042               fState(kMatStateUndefined),
00043               fShader(NULL),
00044               fCerenkov(NULL),
00045               fElement(NULL)
00046 {
00047 // Default constructor
00048    SetUsed(kFALSE);
00049    fIndex    = -1;
00050    fTemperature = STP_temperature;
00051    fPressure = STP_pressure;
00052    fState = kMatStateUndefined;
00053 }
00054 
00055 //_____________________________________________________________________________
00056 TGeoMaterial::TGeoMaterial(const char *name)
00057              :TNamed(name, ""), TAttFill(),
00058               fIndex(0),
00059               fA(0.),
00060               fZ(0.),
00061               fDensity(0.),
00062               fRadLen(0.),
00063               fIntLen(0.),
00064               fTemperature(0.),
00065               fPressure(0.),
00066               fState(kMatStateUndefined),
00067               fShader(NULL),
00068               fCerenkov(NULL),
00069               fElement(NULL)
00070 {
00071 // constructor
00072    fName = fName.Strip();
00073    SetUsed(kFALSE);
00074    fIndex    = -1;
00075    fTemperature = STP_temperature;
00076    fPressure = STP_pressure;
00077    fState = kMatStateUndefined;
00078    
00079    if (!gGeoManager) {
00080       gGeoManager = new TGeoManager("Geometry", "default geometry");
00081    }
00082    gGeoManager->AddMaterial(this);
00083 }
00084 
00085 //_____________________________________________________________________________
00086 TGeoMaterial::TGeoMaterial(const char *name, Double_t a, Double_t z, 
00087                 Double_t rho, Double_t radlen, Double_t intlen)
00088              :TNamed(name, ""), TAttFill(),
00089               fIndex(0),
00090               fA(a),
00091               fZ(z),
00092               fDensity(rho),
00093               fRadLen(0.),
00094               fIntLen(0.),
00095               fTemperature(0.),
00096               fPressure(0.),
00097               fState(kMatStateUndefined),
00098               fShader(NULL),
00099               fCerenkov(NULL),
00100               fElement(NULL)
00101 {
00102 // constructor
00103    fName = fName.Strip();
00104    SetUsed(kFALSE);
00105    fIndex    = -1;
00106    fA        = a;
00107    fZ        = z;
00108    fDensity  = rho;
00109    fTemperature = STP_temperature;
00110    fPressure = STP_pressure;
00111    fState = kMatStateUndefined;
00112    SetRadLen(radlen, intlen);
00113    if (!gGeoManager) {
00114       gGeoManager = new TGeoManager("Geometry", "default geometry");
00115    }
00116    if (fZ - Int_t(fZ) > 1E-3)
00117       Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
00118    GetElement()->SetUsed();
00119    gGeoManager->AddMaterial(this);
00120 }
00121 
00122 //_____________________________________________________________________________
00123 TGeoMaterial::TGeoMaterial(const char *name, Double_t a, Double_t z, Double_t rho,
00124                 EGeoMaterialState state, Double_t temperature, Double_t pressure)
00125              :TNamed(name, ""), TAttFill(),
00126               fIndex(0),
00127               fA(a),
00128               fZ(z),
00129               fDensity(rho),
00130               fRadLen(0.),
00131               fIntLen(0.),
00132               fTemperature(temperature),
00133               fPressure(pressure),
00134               fState(state),
00135               fShader(NULL),
00136               fCerenkov(NULL),
00137               fElement(NULL)
00138 {
00139 // Constructor with state, temperature and pressure.
00140    fName = fName.Strip();
00141    SetUsed(kFALSE);
00142    fIndex    = -1;
00143    SetRadLen(0,0);
00144    if (!gGeoManager) {
00145       gGeoManager = new TGeoManager("Geometry", "default geometry");
00146    }
00147    if (fZ - Int_t(fZ) > 1E-3)
00148       Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
00149    GetElement()->SetUsed();
00150    gGeoManager->AddMaterial(this);
00151 }
00152 
00153 //_____________________________________________________________________________
00154 TGeoMaterial::TGeoMaterial(const char *name, TGeoElement *elem, Double_t rho)
00155              :TNamed(name, ""), TAttFill(),
00156               fIndex(0),
00157               fA(0.),
00158               fZ(0.),
00159               fDensity(rho),
00160               fRadLen(0.),
00161               fIntLen(0.),
00162               fTemperature(0.),
00163               fPressure(0.),
00164               fState(kMatStateUndefined),
00165               fShader(NULL),
00166               fCerenkov(NULL),
00167               fElement(elem)
00168 {
00169 // constructor
00170    fName = fName.Strip();
00171    SetUsed(kFALSE);
00172    fIndex    = -1;
00173    fA        = elem->A();
00174    fZ        = elem->Z();
00175    SetRadLen(0,0);
00176    fTemperature = STP_temperature;
00177    fPressure = STP_pressure;
00178    fState = kMatStateUndefined;
00179    if (!gGeoManager) {
00180       gGeoManager = new TGeoManager("Geometry", "default geometry");
00181    }
00182    if (fZ - Int_t(fZ) > 1E-3)
00183       Warning("ctor", "Material %s defined with fractional Z=%f", GetName(), fZ);
00184    GetElement()->SetUsed();
00185    gGeoManager->AddMaterial(this);
00186 }
00187 
00188 //_____________________________________________________________________________
00189 TGeoMaterial::TGeoMaterial(const TGeoMaterial& gm) :
00190               TNamed(gm),
00191               TAttFill(gm),
00192               fIndex(gm.fIndex),
00193               fA(gm.fA),
00194               fZ(gm.fZ),
00195               fDensity(gm.fDensity),
00196               fRadLen(gm.fRadLen),
00197               fIntLen(gm.fIntLen),
00198               fTemperature(gm.fTemperature),
00199               fPressure(gm.fPressure),
00200               fState(gm.fState),
00201               fShader(gm.fShader),
00202               fCerenkov(gm.fCerenkov),
00203               fElement(gm.fElement)
00204 { 
00205    //copy constructor
00206 }
00207 
00208 //_____________________________________________________________________________
00209 TGeoMaterial& TGeoMaterial::operator=(const TGeoMaterial& gm) 
00210 {
00211    //assignment operator
00212    if(this!=&gm) {
00213       TNamed::operator=(gm);
00214       TAttFill::operator=(gm);
00215       fIndex=gm.fIndex;
00216       fA=gm.fA;
00217       fZ=gm.fZ;
00218       fDensity=gm.fDensity;
00219       fRadLen=gm.fRadLen;
00220       fIntLen=gm.fIntLen;
00221       fTemperature=gm.fTemperature;
00222       fPressure=gm.fPressure;
00223       fState=gm.fState;
00224       fShader=gm.fShader;
00225       fCerenkov=gm.fCerenkov;
00226       fElement=gm.fElement;
00227    } 
00228    return *this;
00229 }
00230 
00231 //_____________________________________________________________________________
00232 TGeoMaterial::~TGeoMaterial()
00233 {
00234 // Destructor
00235 }
00236 
00237 //_____________________________________________________________________________
00238 char *TGeoMaterial::GetPointerName() const
00239 {
00240 // Provide a pointer name containing uid.
00241    static TString name;
00242    name = TString::Format("pMat%d", GetUniqueID());
00243    return (char*)name.Data();
00244 }    
00245 
00246 //_____________________________________________________________________________
00247 void TGeoMaterial::SetRadLen(Double_t radlen, Double_t intlen)
00248 {
00249 // Set radiation/absorbtion lengths. If the values are negative, their absolute value
00250 // is taken, otherwise radlen is recomputed using G3 formula.
00251    fRadLen = TMath::Abs(radlen);
00252    fIntLen = TMath::Abs(intlen);
00253    // Check for vacuum
00254    if (fA<0.9 || fZ<0.9) {
00255       if (radlen<-1e5 || intlen<-1e-5) {
00256          Error("SetRadLen","Material %s: user values taken for vacuum: radlen=%g or intlen=%g - too small", GetName(),fRadLen, fIntLen);
00257          return;
00258       }
00259       // Ignore positive values and take big numbers
00260       if (radlen>=0) fRadLen = 1.E30;   
00261       if (intlen>=0) fIntLen = 1.E30;   
00262       return;
00263    }
00264    // compute radlen systematically with G3 formula for a valid material
00265    if (radlen>=0) {
00266       //taken grom Geant3 routine GSMATE
00267       const Double_t alr2av=1.39621E-03, al183=5.20948;
00268       fRadLen = fA/(alr2av*fDensity*fZ*(fZ +TGeoMaterial::ScreenFactor(fZ))*
00269              (al183-TMath::Log(fZ)/3-TGeoMaterial::Coulomb(fZ)));             
00270    }
00271    // Compute interaction length using the same formula as in GEANT4
00272    if (intlen>=0) {
00273       const Double_t cm = 1.;
00274       const Double_t g = 6.2415e21; // [gram = 1E-3*joule*s*s/(m*m)]
00275       const Double_t amu = 1.03642688246781065e-02; // [MeV/c^2]
00276       const Double_t lambda0 = 35.*g/(cm*cm);  // [g/cm^2]
00277       Double_t nilinv = 0.0;
00278       TGeoElement *elem = GetElement();
00279       Double_t nbAtomsPerVolume = TMath::Na()*fDensity/elem->A();
00280       nilinv += nbAtomsPerVolume*TMath::Power(elem->Neff(), 0.6666667);
00281       nilinv *= amu/lambda0;
00282       fIntLen = (nilinv<=0) ? TGeoShape::Big() : (1./nilinv);
00283    }
00284 }   
00285 
00286 //_____________________________________________________________________________
00287 Double_t TGeoMaterial::Coulomb(Double_t z)
00288 {
00289    // static function
00290    //  Compute Coulomb correction for pair production and Brem 
00291    //  REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
00292    //                        FORMULA 2.7.17
00293    
00294    const Double_t alpha = 7.29927E-03;
00295 
00296    Double_t az    = alpha*z;
00297    Double_t az2   = az*az;
00298    Double_t az4   =   az2 * az2;
00299    Double_t fp    = ( 0.0083*az4 + 0.20206 + 1./(1.+az2) ) * az2;
00300    Double_t fm    = ( 0.0020*az4 + 0.0369  ) * az4;
00301    return fp - fm;
00302 }
00303 
00304 //_____________________________________________________________________________
00305 Bool_t TGeoMaterial::IsEq(const TGeoMaterial *other) const
00306 {
00307 // return true if the other material has the same physical properties
00308    if (other==this) return kTRUE;
00309    if (other->IsMixture()) return kFALSE;
00310    if (TMath::Abs(fA-other->GetA())>1E-3) return kFALSE;
00311    if (TMath::Abs(fZ-other->GetZ())>1E-3) return kFALSE;
00312    if (TMath::Abs(fDensity-other->GetDensity())>1E-6) return kFALSE;
00313    if (GetCerenkovProperties() != other->GetCerenkovProperties()) return kFALSE;
00314 //   if (fRadLen != other->GetRadLen()) return kFALSE;
00315 //   if (fIntLen != other->GetIntLen()) return kFALSE;
00316    return kTRUE;
00317 }
00318 
00319 //_____________________________________________________________________________
00320 void TGeoMaterial::Print(const Option_t * /*option*/) const
00321 {
00322 // print characteristics of this material
00323    printf("Material %s %s   A=%g Z=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(),
00324           fA,fZ,fDensity, fRadLen, fIntLen, fIndex);
00325 }
00326 
00327 //_____________________________________________________________________________
00328 void TGeoMaterial::SavePrimitive(ostream &out, Option_t * /*option*/ /*= ""*/)
00329 {
00330 // Save a primitive as a C++ statement(s) on output stream "out".
00331    if (TestBit(TGeoMaterial::kMatSavePrimitive)) return;
00332    char *name = GetPointerName();
00333    out << "// Material: " << GetName() << endl;
00334    out << "   a       = " << fA << ";" << endl;
00335    out << "   z       = " << fZ << ";" << endl;
00336    out << "   density = " << fDensity << ";" << endl;
00337    out << "   radl    = " << fRadLen << ";" << endl;
00338    out << "   absl    = " << fIntLen << ";" << endl;
00339    
00340    out << "   " << name << " = new TGeoMaterial(\"" << GetName() << "\", a,z,density,radl,absl);" << endl;
00341    out << "   " << name << "->SetIndex(" << GetIndex() << ");" << endl;
00342    SetBit(TGeoMaterial::kMatSavePrimitive);
00343 }
00344 
00345 //_____________________________________________________________________________
00346 Int_t TGeoMaterial::GetDefaultColor() const
00347 {
00348 // Get some default color related to this material.
00349    Int_t id = 1+ gGeoManager->GetListOfMaterials()->IndexOf(this);
00350    return (2+id%6);
00351 }
00352 
00353 //_____________________________________________________________________________
00354 TGeoElement *TGeoMaterial::GetElement(Int_t) const
00355 {
00356 // Get a pointer to the element this material is made of.
00357    if (fElement) return fElement;
00358    TGeoElementTable *table = gGeoManager->GetElementTable();
00359    return table->GetElement(Int_t(fZ));
00360 }
00361 
00362 //_____________________________________________________________________________
00363 Int_t TGeoMaterial::GetIndex()
00364 {
00365 // Retreive material index in the list of materials
00366    if (fIndex>=0) return fIndex;
00367    TList *matlist = gGeoManager->GetListOfMaterials();
00368    fIndex = matlist->IndexOf(this);
00369    return fIndex;
00370 }      
00371 
00372 //_____________________________________________________________________________
00373 TGeoMaterial *TGeoMaterial::DecayMaterial(Double_t time, Double_t precision)
00374 {
00375 // Create the material representing the decay product of this material at a
00376 // given time. The precision represent the minimum cumulative branching ratio for 
00377 // which decay products are still taken into account.
00378    TObjArray *pop = new TObjArray();
00379    if (!fElement || !fElement->IsRadioNuclide()) return this;
00380    FillMaterialEvolution(pop, precision);
00381    Int_t ncomp = pop->GetEntriesFast();
00382    if (!ncomp) return this;
00383    TGeoElementRN *el;
00384    Double_t *weight = new Double_t[ncomp];
00385    Double_t amed = 0.;
00386    Int_t i;
00387    for (i=0; i<ncomp; i++) {
00388       el = (TGeoElementRN *)pop->At(i);
00389       weight[i] = el->Ratio()->Concentration(time) * el->A();
00390       amed += weight[i];
00391    }  
00392    Double_t rho = fDensity*amed/fA;
00393    TGeoMixture *mix = 0;
00394    Int_t ncomp1 = ncomp;
00395    for (i=0; i<ncomp; i++) {
00396       if ((weight[i]/amed)<precision) {
00397          amed -= weight[i];
00398          ncomp1--;
00399       }
00400    }
00401    if (ncomp1<2) {
00402       el = (TGeoElementRN *)pop->At(0);
00403       delete [] weight;
00404       delete pop;
00405       if (ncomp1==1) return new TGeoMaterial(TString::Format("%s-evol",GetName()), el, rho);
00406       return NULL;
00407    }   
00408    mix = new TGeoMixture(TString::Format("%s-evol",GetName()), ncomp, rho);
00409    for (i=0; i<ncomp; i++) {
00410       weight[i] /= amed;
00411       if (weight[i]<precision) continue;
00412       el = (TGeoElementRN *)pop->At(i);
00413       mix->AddElement(el, weight[i]);
00414    }
00415    delete [] weight;
00416    delete pop;
00417    return mix;
00418 }      
00419 
00420 //_____________________________________________________________________________
00421 void TGeoMaterial::FillMaterialEvolution(TObjArray *population, Double_t precision)
00422 {
00423 // Fills a user array with all the elements deriving from the possible
00424 // decay of the top element composing the mixture. Each element contained
00425 // by <population> may be a radionuclide having a Bateman solution attached.
00426 // The precision represent the minimum cumulative branching ratio for 
00427 // which decay products are still taken into account.
00428 // To visualize the time evolution of each decay product one can use:
00429 //    TGeoElement *elem = population->At(index);
00430 //    TGeoElementRN *elemrn = 0;
00431 //    if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
00432 // One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
00433 // of one of the decay products, N1(0) is the number of atoms of the top
00434 // element at t=0.
00435 //    Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
00436 // One can also display the time evolution of the fractional weigth:
00437 //    elemrn->Ratio()->Draw(option);
00438    if (population->GetEntriesFast()) {
00439       Error("FillMaterialEvolution", "Provide an empty array !");
00440       return;
00441    }
00442    TGeoElementTable *table = gGeoManager->GetElementTable();
00443    TGeoElement *elem;
00444    TGeoElementRN *elemrn;
00445    TIter next(table->GetElementsRN());
00446    while ((elemrn=(TGeoElementRN*)next())) elemrn->ResetRatio();
00447    elem = GetElement();
00448    if (!elem->IsRadioNuclide()) {
00449       population->Add(elem);
00450       return;
00451    }
00452    elemrn = (TGeoElementRN*)elem;
00453    elemrn->FillPopulation(population, precision);
00454 }      
00455 
00456 
00457 /*************************************************************************
00458  * TGeoMixture - mixtures of elements 
00459  *
00460  *************************************************************************/
00461 ClassImp(TGeoMixture)
00462 
00463 //_____________________________________________________________________________
00464 TGeoMixture::TGeoMixture()
00465 {
00466 // Default constructor
00467    fNelements = 0;
00468    fZmixture  = 0;
00469    fAmixture  = 0;
00470    fWeights   = 0;
00471    fNatoms    = 0;
00472    fElements  = 0;
00473 }
00474 
00475 //_____________________________________________________________________________
00476 TGeoMixture::TGeoMixture(const char *name, Int_t /*nel*/, Double_t rho)
00477             :TGeoMaterial(name)
00478 {
00479 // constructor
00480    fZmixture   = 0;
00481    fAmixture   = 0;
00482    fWeights    = 0;
00483    fNelements  = 0;
00484    fNatoms     = 0;
00485    fDensity = rho;
00486    fElements   = 0;
00487    if (fDensity < 0) fDensity = 0.001;
00488 }
00489 
00490 //_____________________________________________________________________________
00491 TGeoMixture::TGeoMixture(const TGeoMixture& gm) :
00492   TGeoMaterial(gm),
00493   fNelements(gm.fNelements),
00494   fZmixture(gm.fZmixture),
00495   fAmixture(gm.fAmixture),
00496   fWeights(gm.fWeights),
00497   fNatoms(gm.fNatoms),
00498   fElements(gm.fElements)
00499 { 
00500    //copy constructor
00501 }
00502 
00503 //_____________________________________________________________________________
00504 TGeoMixture& TGeoMixture::operator=(const TGeoMixture& gm) 
00505 {
00506    //assignment operator
00507    if(this!=&gm) {
00508       TGeoMaterial::operator=(gm);
00509       fNelements=gm.fNelements;
00510       fZmixture=gm.fZmixture;
00511       fAmixture=gm.fAmixture;
00512       fWeights=gm.fWeights;
00513       fNatoms = gm.fNatoms;
00514       fElements = gm.fElements;
00515    } 
00516    return *this;
00517 }
00518 
00519 //_____________________________________________________________________________
00520 TGeoMixture::~TGeoMixture()
00521 {
00522 // Destructor
00523    if (fZmixture) delete[] fZmixture;
00524    if (fAmixture) delete[] fAmixture;
00525    if (fWeights)  delete[] fWeights;
00526    if (fNatoms)   delete[] fNatoms;
00527    if (fElements) delete fElements;
00528 }
00529 
00530 //_____________________________________________________________________________
00531 void TGeoMixture::AverageProperties()
00532 {
00533 // Compute effective A/Z and radiation length
00534    const Double_t alr2av = 1.39621E-03 , al183 =5.20948;
00535    const Double_t cm = 1.;
00536    const Double_t g = 6.2415e21; // [gram = 1E-3*joule*s*s/(m*m)]
00537    const Double_t amu = 1.03642688246781065e-02; // [MeV/c^2]
00538    const Double_t lambda0 = 35.*g/(cm*cm);  // [g/cm^2]
00539    Double_t radinv = 0.0;
00540    Double_t nilinv = 0.0;
00541    Double_t nbAtomsPerVolume;
00542    fA = 0;
00543    fZ = 0;
00544    for (Int_t j=0;j<fNelements;j++) {
00545       if (fWeights[j] <= 0) continue;
00546       fA += fWeights[j]*fAmixture[j];
00547       fZ += fWeights[j]*fZmixture[j];
00548       nbAtomsPerVolume = TMath::Na()*fDensity*fWeights[j]/GetElement(j)->A();
00549       nilinv += nbAtomsPerVolume*TMath::Power(GetElement(j)->Neff(), 0.6666667);
00550       Double_t zc = fZmixture[j];
00551       Double_t alz = TMath::Log(zc)/3.;
00552       Double_t xinv = zc*(zc+TGeoMaterial::ScreenFactor(zc))*
00553          (al183-alz-TGeoMaterial::Coulomb(zc))/fAmixture[j];
00554       radinv += xinv*fWeights[j];
00555    }
00556    radinv *= alr2av*fDensity;
00557    if (radinv > 0) fRadLen = 1/radinv;
00558    // Compute interaction length
00559    nilinv *= amu/lambda0;
00560    fIntLen = (nilinv<=0) ? TGeoShape::Big() : (1./nilinv);
00561 }
00562 
00563 //_____________________________________________________________________________
00564 void TGeoMixture::AddElement(Double_t a, Double_t z, Double_t weight)
00565 {
00566 // add an element to the mixture using fraction by weight
00567    // Check if the element is already defined
00568    TGeoElementTable *table = gGeoManager->GetElementTable();
00569    if (z<1 || z>table->GetNelements()-1)
00570       Fatal("AddElement", "Cannot add element having Z=%d to mixture %s", (Int_t)z, GetName());
00571    Int_t i;
00572    for (i=0; i<fNelements; i++) {
00573       if (TMath::Abs(z-fZmixture[i])<1.e-6  && TMath::Abs(a-fAmixture[i])<1.e-6) {
00574          fWeights[i] += weight;
00575          AverageProperties();
00576          return;
00577       }
00578    }      
00579    if (!fNelements) {
00580       fZmixture = new Double_t[1];
00581       fAmixture = new Double_t[1];
00582       fWeights  = new Double_t[1];
00583    } else {   
00584       Int_t nelements = fNelements+1;
00585       Double_t *zmixture = new Double_t[nelements];
00586       Double_t *amixture = new Double_t[nelements];
00587       Double_t *weights  = new Double_t[nelements];
00588       for (Int_t j=0; j<fNelements; j++) {
00589          zmixture[j] = fZmixture[j];
00590          amixture[j] = fAmixture[j];
00591          weights[j]  = fWeights[j];
00592       }
00593       delete [] fZmixture;
00594       delete [] fAmixture;
00595       delete [] fWeights;
00596       fZmixture = zmixture;
00597       fAmixture = amixture;
00598       fWeights  = weights;
00599    }       
00600    
00601    fNelements++;
00602    i = fNelements - 1;   
00603    fZmixture[i] = z;
00604    fAmixture[i] = a;
00605    fWeights[i]  = weight;
00606    if (z - Int_t(z) > 1E-3)
00607       Warning("DefineElement", "Mixture %s has element defined with fractional Z=%f", GetName(), z);
00608    GetElement(i)->SetDefined();
00609    table->GetElement((Int_t)z)->SetDefined();
00610    
00611    //compute equivalent radiation length (taken from Geant3/GSMIXT)
00612    AverageProperties();
00613 }
00614 
00615 //_____________________________________________________________________________
00616 void TGeoMixture::AddElement(TGeoMaterial *mat, Double_t weight)
00617 {
00618 // Define one component of the mixture as an existing material/mixture.
00619    TGeoElement *elnew, *elem;   
00620    Double_t a,z;
00621    if (!mat->IsMixture()) {
00622       elem = mat->GetBaseElement();
00623       if (elem) {
00624          AddElement(elem, weight);
00625       } else {   
00626          a = mat->GetA();
00627          z = mat->GetZ();
00628          AddElement(a, z, weight);
00629       }   
00630       return;
00631    }
00632    // The material is a mixture.
00633    TGeoMixture *mix = (TGeoMixture*)mat;
00634    Double_t wnew;
00635    Int_t nelem = mix->GetNelements();
00636    Bool_t elfound;
00637    Int_t i,j;
00638    // loop the elements of the daughter mixture
00639    for (i=0; i<nelem; i++) {
00640       elfound = kFALSE;
00641       elnew = mix->GetElement(i);
00642       if (!elnew) continue;
00643       // check if we have the element already defined in the parent mixture
00644       for (j=0; j<fNelements; j++) {
00645          if (fWeights[j]<=0) continue;
00646          elem = GetElement(j);
00647          if (elem == elnew) {
00648             // element found, compute new weight
00649             fWeights[j] += weight * (mix->GetWmixt())[i];
00650             elfound = kTRUE;
00651             break;
00652          }
00653       }
00654       if (elfound) continue;
00655       // element not found, define it
00656       wnew = weight * (mix->GetWmixt())[i];
00657       AddElement(elnew, wnew);
00658    }   
00659 }         
00660 
00661 //_____________________________________________________________________________
00662 void TGeoMixture::AddElement(TGeoElement *elem, Double_t weight)
00663 {
00664 // add an element to the mixture using fraction by weight
00665    TGeoElement *elemold;
00666    TGeoElementTable *table = gGeoManager->GetElementTable();
00667    if (!fElements) fElements = new TObjArray(128);
00668    Bool_t exist = kFALSE;
00669    // If previous elements were defined by A/Z, add corresponding TGeoElements
00670    for (Int_t i=0; i<fNelements; i++) {
00671       elemold = (TGeoElement*)fElements->At(i);
00672       if (!elemold) fElements->AddAt(elemold = table->GetElement((Int_t)fZmixture[i]), i);   
00673       if (elemold == elem) exist = kTRUE;
00674    }
00675    if (!exist) fElements->AddAtAndExpand(elem, fNelements);   
00676    AddElement(elem->A(), elem->Z(), weight);
00677 }   
00678 
00679 //_____________________________________________________________________________
00680 void TGeoMixture::AddElement(TGeoElement *elem, Int_t natoms)
00681 {
00682 // Add a mixture element by number of atoms in the chemical formula.
00683    Int_t i,j;
00684    Double_t amol;
00685    TGeoElement *elemold;
00686    TGeoElementTable *table = gGeoManager->GetElementTable();
00687    if (!fElements) fElements = new TObjArray(128);
00688    // Check if the element is already defined
00689    for (i=0; i<fNelements; i++) {
00690       elemold = (TGeoElement*)fElements->At(i);
00691       if (!elemold) fElements->AddAt(table->GetElement((Int_t)fZmixture[i]), i);
00692       else if (elemold != elem) continue;
00693       if ((elem==elemold) || 
00694           (TMath::Abs(elem->Z()-fZmixture[i])<1.e-6 && TMath::Abs(elem->A()-fAmixture[i])<1.e-6)) {
00695          fNatoms[i] += natoms;
00696          amol = 0.;
00697          for (j=0; j<fNelements; j++) amol += fAmixture[j]*fNatoms[j];
00698          for (j=0; j<fNelements; j++) fWeights[j] = fNatoms[j]*fAmixture[j]/amol;
00699          AverageProperties();
00700          return;
00701       }
00702    }
00703    // New element      
00704    if (!fNelements) {
00705       fZmixture = new Double_t[1];
00706       fAmixture = new Double_t[1];
00707       fWeights  = new Double_t[1];
00708       fNatoms   = new Int_t[1];
00709    } else {   
00710       if (!fNatoms) {
00711          Fatal("AddElement", "Cannot add element by natoms in mixture %s after defining elements by weight",
00712                GetName());
00713          return;
00714       }         
00715       Int_t nelements = fNelements+1;
00716       Double_t *zmixture = new Double_t[nelements];
00717       Double_t *amixture = new Double_t[nelements];
00718       Double_t *weights  = new Double_t[nelements];
00719       Int_t *nnatoms  = new Int_t[nelements];
00720       for (j=0; j<fNelements; j++) {
00721          zmixture[j] = fZmixture[j];
00722          amixture[j] = fAmixture[j];
00723          weights[j]  = fWeights[j];
00724          nnatoms[j]  = fNatoms[j];
00725       }
00726       delete [] fZmixture;
00727       delete [] fAmixture;
00728       delete [] fWeights;
00729       delete [] fNatoms;
00730       fZmixture = zmixture;
00731       fAmixture = amixture;
00732       fWeights  = weights;
00733       fNatoms   = nnatoms;
00734    }
00735    fNelements++;       
00736    Int_t iel = fNelements-1;
00737    fZmixture[iel] = elem->Z();
00738    fAmixture[iel] = elem->A();
00739    fNatoms[iel]  = natoms;
00740    fElements->AddAtAndExpand(elem, iel);
00741    amol = 0.;
00742    for (i=0; i<fNelements; i++) {
00743       if (fNatoms[i]<=0) return;
00744       amol += fAmixture[i]*fNatoms[i];
00745    }   
00746    for (i=0; i<fNelements; i++) fWeights[i] = fNatoms[i]*fAmixture[i]/amol;
00747    table->GetElement(elem->Z())->SetDefined();
00748    AverageProperties();
00749 }          
00750 
00751 //_____________________________________________________________________________
00752 void TGeoMixture::DefineElement(Int_t /*iel*/, Int_t z, Int_t natoms)
00753 {
00754 // Define the mixture element at index iel by number of atoms in the chemical formula.
00755    TGeoElementTable *table = gGeoManager->GetElementTable();
00756    TGeoElement *elem = table->GetElement(z);
00757    if (!elem) {
00758       Fatal("DefineElement", "In mixture %s, element with Z=%i not found",GetName(),z);
00759       return;
00760    }   
00761    AddElement(elem, natoms);
00762 }
00763    
00764 //_____________________________________________________________________________
00765 TGeoElement *TGeoMixture::GetElement(Int_t i) const
00766 {
00767 // Retreive the pointer to the element corresponding to component I.
00768    if (i<0 || i>=fNelements) {
00769       Error("GetElement", "Mixture %s has only %d elements", GetName(), fNelements);
00770       return 0;
00771    }   
00772    TGeoElement *elem = 0;
00773    if (fElements) elem = (TGeoElement*)fElements->At(i);
00774    if (elem) return elem;
00775    TGeoElementTable *table = gGeoManager->GetElementTable();
00776    return table->GetElement(Int_t(fZmixture[i]));
00777 }
00778 
00779 //_____________________________________________________________________________
00780 Bool_t TGeoMixture::IsEq(const TGeoMaterial *other) const
00781 {
00782 // Return true if the other material has the same physical properties
00783    if (other->IsEqual(this)) return kTRUE;
00784    if (!other->IsMixture()) return kFALSE;
00785    TGeoMixture *mix = (TGeoMixture*)other;
00786    if (!mix) return kFALSE;
00787    if (fNelements != mix->GetNelements()) return kFALSE;
00788    if (TMath::Abs(fA-other->GetA())>1E-3) return kFALSE;
00789    if (TMath::Abs(fZ-other->GetZ())>1E-3) return kFALSE;
00790    if (TMath::Abs(fDensity-other->GetDensity())>1E-6) return kFALSE;
00791    if (GetCerenkovProperties() != other->GetCerenkovProperties()) return kFALSE;
00792 //   if (fRadLen != other->GetRadLen()) return kFALSE;
00793 //   if (fIntLen != other->GetIntLen()) return kFALSE;
00794    for (Int_t i=0; i<fNelements; i++) {
00795       if (TMath::Abs(fZmixture[i]-(mix->GetZmixt())[i])>1E-3) return kFALSE;
00796       if (TMath::Abs(fAmixture[i]-(mix->GetAmixt())[i])>1E-3) return kFALSE;
00797       if (TMath::Abs(fWeights[i]-(mix->GetWmixt())[i])>1E-3) return kFALSE;
00798    }
00799    return kTRUE;
00800 }
00801 
00802 //_____________________________________________________________________________
00803 void TGeoMixture::Print(const Option_t * /*option*/) const
00804 {
00805 // print characteristics of this material
00806    printf("Mixture %s %s   Aeff=%g Zeff=%g rho=%g radlen=%g intlen=%g index=%i\n", GetName(), GetTitle(),
00807           fA,fZ,fDensity, fRadLen, fIntLen, fIndex);
00808    for (Int_t i=0; i<fNelements; i++) {
00809       if (fNatoms) printf("   Element #%i : %s  Z=%6.2f A=%6.2f w=%6.3f natoms=%d\n", i, GetElement(i)->GetName(),fZmixture[i],
00810              fAmixture[i], fWeights[i], fNatoms[i]);
00811       else printf("   Element #%i : %s  Z=%6.2f A=%6.2f w=%6.3f\n", i, GetElement(i)->GetName(),fZmixture[i],
00812              fAmixture[i], fWeights[i]);
00813    }
00814 }
00815 
00816 //_____________________________________________________________________________
00817 void TGeoMixture::SavePrimitive(ostream &out, Option_t * /*option*/ /*= ""*/)
00818 {
00819 // Save a primitive as a C++ statement(s) on output stream "out".
00820    if (TestBit(TGeoMaterial::kMatSavePrimitive)) return;
00821    char *name = GetPointerName();
00822    out << "// Mixture: " << GetName() << endl;
00823    out << "   nel     = " << fNelements << ";" << endl;
00824    out << "   density = " << fDensity << ";" << endl;
00825    out << "   " << name << " = new TGeoMixture(\"" << GetName() << "\", nel,density);" << endl;
00826    for (Int_t i=0; i<fNelements; i++) {
00827       TGeoElement *el = GetElement(i);
00828       out << "      a = " << fAmixture[i] << ";   z = "<< fZmixture[i] << ";   w = " << fWeights[i] << ";  // " << el->GetName() << endl;
00829       out << "   " << name << "->DefineElement(" << i << ",a,z,w);" << endl;
00830    }         
00831    out << "   " << name << "->SetIndex(" << GetIndex() << ");" << endl;
00832    SetBit(TGeoMaterial::kMatSavePrimitive);
00833 }
00834 
00835 //_____________________________________________________________________________
00836 TGeoMaterial *TGeoMixture::DecayMaterial(Double_t time, Double_t precision)
00837 {
00838 // Create the mixture representing the decay product of this material at a
00839 // given time. The precision represent the minimum cumulative branching ratio for 
00840 // which decay products are still taken into account.
00841    TObjArray *pop = new TObjArray();
00842    FillMaterialEvolution(pop, precision);
00843    Int_t ncomp = pop->GetEntriesFast();
00844    if (!ncomp) return this;
00845    TGeoElement *elem;
00846    TGeoElementRN *el;
00847    Double_t *weight = new Double_t[ncomp];
00848    Double_t amed = 0.;
00849    Int_t i, j;
00850    for (i=0; i<ncomp; i++) {
00851       elem = (TGeoElement *)pop->At(i);
00852       if (!elem->IsRadioNuclide()) {
00853          j = fElements->IndexOf(elem);
00854          weight[i] = fWeights[j]*fAmixture[0]/fWeights[0];
00855       } else {   
00856          el = (TGeoElementRN*)elem;
00857          weight[i] = el->Ratio()->Concentration(time) * el->A();
00858       }   
00859       amed += weight[i];
00860    }  
00861    Double_t rho = fDensity * fWeights[0] * amed/fAmixture[0];
00862    TGeoMixture *mix = 0;
00863    Int_t ncomp1 = ncomp;
00864    for (i=0; i<ncomp; i++) {
00865       if ((weight[i]/amed)<precision) {
00866          amed -= weight[i];
00867          ncomp1--;
00868       }
00869    }
00870    if (ncomp1<2) {
00871       el = (TGeoElementRN *)pop->At(0);
00872       delete [] weight;
00873       delete pop;
00874       if (ncomp1==1) return new TGeoMaterial(TString::Format("%s-evol",GetName()), el, rho);
00875       return NULL;
00876    }
00877    mix = new TGeoMixture(TString::Format("%s-evol",GetName()), ncomp, rho); 
00878    for (i=0; i<ncomp; i++) {
00879       weight[i] /= amed;
00880       if (weight[i]<precision) continue;
00881       el = (TGeoElementRN *)pop->At(i);
00882       mix->AddElement(el, weight[i]);
00883    }
00884    delete [] weight;
00885    delete pop;
00886    return mix;
00887 }      
00888 
00889 //_____________________________________________________________________________
00890 void TGeoMixture::FillMaterialEvolution(TObjArray *population, Double_t precision)
00891 {
00892 // Fills a user array with all the elements deriving from the possible
00893 // decay of the top elements composing the mixture. Each element contained
00894 // by <population> may be a radionuclide having a Bateman solution attached.
00895 // The precision represent the minimum cumulative branching ratio for 
00896 // which decay products are still taken into account.
00897 // To visualize the time evolution of each decay product one can use:
00898 //    TGeoElement *elem = population->At(index);
00899 //    TGeoElementRN *elemrn = 0;
00900 //    if (elem->IsRadioNuclide()) elemrn = (TGeoElementRN*)elem;
00901 // One can get Ni/N1(t=0) at any moment of time. Ni is the number of atoms
00902 // of one of the decay products, N1(0) is the number of atoms of the first top
00903 // element at t=0.
00904 //    Double_t fraction_weight = elemrn->Ratio()->Concentration(time);
00905 // One can also display the time evolution of the fractional weigth:
00906 //    elemrn->Ratio()->Draw(option);
00907    if (population->GetEntriesFast()) {
00908       Error("FillMaterialEvolution", "Provide an empty array !");
00909       return;
00910    }
00911    TGeoElementTable *table = gGeoManager->GetElementTable();
00912    TGeoElement *elem;
00913    TGeoElementRN *elemrn;
00914    TIter next(table->GetElementsRN());
00915    while ((elemrn=(TGeoElementRN*)next())) elemrn->ResetRatio();
00916    Double_t factor;
00917    for (Int_t i=0; i<fNelements; i++) {
00918       elem = GetElement(i);
00919       if (!elem->IsRadioNuclide()) {
00920          population->Add(elem);
00921          continue;
00922       }
00923       elemrn = (TGeoElementRN*)elem;
00924       factor = fWeights[i]*fAmixture[0]/(fWeights[0]*fAmixture[i]);
00925       elemrn->FillPopulation(population, precision, factor);
00926    }   
00927 }      
00928 
00929 //_____________________________________________________________________________
00930 Double_t TGeoMaterial::ScreenFactor(Double_t z)
00931 {
00932    // static function
00933    //  Compute screening factor for pair production and Bremstrahlung
00934    //  REFERENCE : EGS MANUAL SLAC 210 - UC32 - JUNE 78
00935    //                        FORMULA 2.7.22
00936    
00937    const Double_t al183= 5.20948 , al1440 = 7.27239;
00938    Double_t alz  = TMath::Log(z)/3.;
00939    Double_t factor = (al1440 - 2*alz) / (al183 - alz - TGeoMaterial::Coulomb(z));
00940    return factor;
00941 }

Generated on Tue Jul 5 14:12:56 2011 for ROOT_528-00b_version by  doxygen 1.5.1