]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/GetXsection.C
0. General code clean-up, including messages, and the like.
[u/mrichter/AliRoot.git] / FMD / scripts / GetXsection.C
1 //____________________________________________________________________
2 //
3 // $Id$
4 //
5 // Script to get the various cross sections, energy loss, and ranges
6 // of a particular particle type in a particular medium. 
7 //
8 // This script should be compiled to speed it up.  
9 // 
10 // It creates a tree on the current output file, with the relevant
11 // information. 
12 //
13 // Note, that VMC _must_ be the TGeant3TGeo VMC. 
14 //
15 #include <TArrayF.h>
16 #include <TTree.h>
17 #include <TMath.h>
18 #include <iostream>
19 #include <TVirtualMC.h>
20 #include <TDatabasePDG.h>
21 #include <TString.h>
22 #include <TGeoManager.h>
23 #include <TGeoMedium.h>
24 #include <TGeant3.h>
25
26 //____________________________________________________________________
27 struct Mech 
28 {
29   char*    name;
30   char*    title;
31   char*    unit;
32   TArrayF  values;
33   int      status;
34 };
35
36
37 //____________________________________________________________________
38 void
39 GetXsection(const char* medName, const char* pdgName,
40             Int_t n=91, Float_t emin=1e-5, Float_t emax=1e4)
41 {
42   TArrayF tkine(n);
43   Float_t dp   = 1/TMath::Log10(emax/emin);
44   Float_t pmin = TMath::Log10(emin);
45   tkine[0]     = emin;
46   for (Int_t i=1; i < tkine.fN; i++) {
47     Float_t el = pmin + i * dp;
48     tkine[i]   = TMath::Power(10, el);
49   }
50   TArrayF cuts(5);
51   cuts.Reset(1e-4);
52
53   Mech mechs[] = 
54     {{ "HADF","total hadronic x-section according to FLUKA","cm^{1}",n,0},
55      { "INEF","hadronic inelastic x-section according to FLUKA","cm^{1}",n,0},
56      { "ELAF","hadronic elastic x-section according to FLUKA","cm^{1}",n,0},
57      { "HADG","total hadronic x-section according to GHEISHA","cm^{1}",n,0},
58      { "INEG","hadronic inelastic x-section according to GHEISHA","cm^{1}",n,0},
59      { "ELAG","hadronic elastic x-section according to GHEISHA","cm^{1}",n,0},
60      { "FISG","nuclear fission x-section according to GHEISHA","cm^{1}",n,0},
61      { "CAPG","neutron capture x-section according to GHEISHA","cm^{1}",n,0},
62      { "LOSS","stopping power","cm^{1}",n,0},
63      { "PHOT","photoelectric x-section","MeV/cm",n,0},
64      { "ANNI","positron annihilation x-section","cm^{1}",n,0},
65      { "COMP","Compton effect x-section","cm^{1}",n,0},
66      { "BREM","bremsstrahlung x-section","cm^{1}",n,0},
67      { "PAIR","photon and muon direct- pair x-section","cm^{1}",n,0},
68      { "DRAY","delta-rays x-section","cm^{1}",n,0},
69      { "PFIS","photo-fission x-section","cm^{1}",n,0},
70      { "RAYL","Rayleigh scattering x-section","cm^{1}",n,0},
71      { "MUNU","muon-nuclear interaction x-section","cm^{1}",n,0},
72      { "RANG","range","cm",n,0},
73      { "STEP","maximum step","cm",n,0},
74      { 0, 0, 0, 0, 0}};
75   TGeant3* mc = (TGeant3*)gMC;
76   if (!mc) {
77     std::cerr << "Couldn't get VMC" << std::endl;
78     return;
79   }
80   TGeoMedium* medium = gGeoManager->GetMedium(medName);
81   if (!medium) {
82     std::cerr << "Couldn't find medium " << medName << std::endl;
83     return;
84   }
85   Int_t medNo = medium->GetMaterial()->GetUniqueID();
86   TDatabasePDG* pdgDb = TDatabasePDG::Instance();
87   TParticlePDG* pdgP  = pdgDb->GetParticle(pdgName);
88   if (!pdgP) {
89     std::cerr << "Couldn't find particle " << pdgName << std::endl;
90     return;
91   }
92   Int_t pdgNo = pdgP->PdgCode();
93   Int_t pidNo = mc->IdFromPDG(pdgNo);
94     
95   Mech* mech = &(mechs[0]);
96   Int_t nMech = 0;
97   Int_t nOk = 0;
98   TString vars("T/F");
99   while (mech->name) {
100     cout << mech->name << ": " << mech->title << " ... " << std::flush;
101     nMech++;
102     Int_t ixst;
103     mc->Gftmat(medNo, pidNo, mech->name, n, 
104                tkine.fArray, mech->values.fArray, cuts.fArray, ixst);
105     mech->status = ixst;
106     if (ixst) {
107       nOk++;
108       vars.Append(Form(":%s", mech->name));
109       if (!strcmp("LOSS", mech->name)) {
110         for (Int_t i = 0; i < n; i++) 
111           std::cout << i << "\t" << tkine[i] << "\t" 
112                     << mech->values[i] << std::endl;
113       }
114     }
115     std::cout << (ixst ? "ok" : "failed") << std::endl;
116     mech++;
117   }
118   // TFile* file = TFile::Open(Form("xsec-%d.root", pdgNo),
119   // "RECREATE");
120   TArrayF cache(nOk+1);
121   TTree* tree = new TTree(Form("%s_%s", medName, pdgName), 
122                           Form("%s_%s", medName, pdgName));
123   tree->Branch("xsec", cache.fArray, vars.Data());
124   for (Int_t i = 0; i < n; i++) {
125     cache[0] = tkine[i];
126     Int_t k = 0;
127     for (Int_t j = 0; j < nMech; j++) {
128       if (mechs[j].status) {
129         if (!strcmp(mechs[j].name, "LOSS")) 
130           std::cout << tkine[i] << "\t" << mechs[j].values[i] << std::endl;
131         cache[k+1] = mechs[j].values[i];
132         k++;
133       }
134     }
135     std::cout << k << "\t" << (k == nOk) << std::endl;
136     tree->Fill();
137   }
138   tree->Write();
139 }
140 //____________________________________________________________________
141 //
142 // EOF
143 //