]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/scripts/GetXsection.C
Added new library libFMDutil. This library contains utility classes that
[u/mrichter/AliRoot.git] / FMD / scripts / GetXsection.C
CommitLineData
a1f80595 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
26struct Mech
27{
28 char* name;
29 char* title;
30 char* unit;
31 TArrayF values;
32 int status;
33};
34
35
36void
37GetXsection(const char* medName, const char* pdgName,
38 Int_t n=91, Float_t emin=1e-5, Float_t emax=1e4)
39{
40 TArrayF tkine(n);
41 Float_t dp = 1/TMath::Log10(emax/emin);
42 Float_t pmin = TMath::Log10(emin);
43 tkine[0] = emin;
44 for (Int_t i=1; i < tkine.fN; i++) {
45 Float_t el = pmin + i * dp;
46 tkine[i] = TMath::Power(10, el);
47 }
48 TArrayF cuts(5);
49 cuts.Reset(1e-4);
50
51 Mech mechs[] =
52 {{ "HADF","total hadronic x-section according to FLUKA","cm^{1}",n,0},
53 { "INEF","hadronic inelastic x-section according to FLUKA","cm^{1}",n,0},
54 { "ELAF","hadronic elastic x-section according to FLUKA","cm^{1}",n,0},
55 { "HADG","total hadronic x-section according to GHEISHA","cm^{1}",n,0},
56 { "INEG","hadronic inelastic x-section according to GHEISHA","cm^{1}",n,0},
57 { "ELAG","hadronic elastic x-section according to GHEISHA","cm^{1}",n,0},
58 { "FISG","nuclear fission x-section according to GHEISHA","cm^{1}",n,0},
59 { "CAPG","neutron capture x-section according to GHEISHA","cm^{1}",n,0},
60 { "LOSS","stopping power","cm^{1}",n,0},
61 { "PHOT","photoelectric x-section","MeV/cm",n,0},
62 { "ANNI","positron annihilation x-section","cm^{1}",n,0},
63 { "COMP","Compton effect x-section","cm^{1}",n,0},
64 { "BREM","bremsstrahlung x-section","cm^{1}",n,0},
65 { "PAIR","photon and muon direct- pair x-section","cm^{1}",n,0},
66 { "DRAY","delta-rays x-section","cm^{1}",n,0},
67 { "PFIS","photo-fission x-section","cm^{1}",n,0},
68 { "RAYL","Rayleigh scattering x-section","cm^{1}",n,0},
69 { "MUNU","muon-nuclear interaction x-section","cm^{1}",n,0},
70 { "RANG","range","cm",n,0},
71 { "STEP","maximum step","cm",n,0},
72 { 0, 0, 0, 0, 0}};
73 TGeant3* mc = (TGeant3*)gMC;
74 if (!mc) {
75 std::cerr << "Couldn't get VMC" << std::endl;
76 return;
77 }
78 TGeoMedium* medium = gGeoManager->GetMedium(medName);
79 if (!medium) {
80 std::cerr << "Couldn't find medium " << medName << std::endl;
81 return;
82 }
83 Int_t medNo = medium->GetMaterial()->GetUniqueID();
84 TDatabasePDG* pdgDb = TDatabasePDG::Instance();
85 TParticlePDG* pdgP = pdgDb->GetParticle(pdgName);
86 if (!pdgP) {
87 std::cerr << "Couldn't find particle " << pdgName << std::endl;
88 return;
89 }
90 Int_t pdgNo = pdgP->PdgCode();
91 Int_t pidNo = mc->IdFromPDG(pdgNo);
92
93 Mech* mech = &(mechs[0]);
94 Int_t nMech = 0;
95 Int_t nOk = 0;
96 TString vars("T/F");
97 while (mech->name) {
98 cout << mech->name << ": " << mech->title << " ... " << std::flush;
99 nMech++;
100 Int_t ixst;
101 mc->Gftmat(medNo, pidNo, mech->name, n,
102 tkine.fArray, mech->values.fArray, cuts.fArray, ixst);
103 mech->status = ixst;
104 if (ixst) {
105 nOk++;
106 vars.Append(Form(":%s", mech->name));
107 if (!strcmp("LOSS", mech->name)) {
108 for (Int_t i = 0; i < n; i++)
109 std::cout << i << "\t" << tkine[i] << "\t"
110 << mech->values[i] << std::endl;
111 }
112 }
113 std::cout << (ixst ? "ok" : "failed") << std::endl;
114 mech++;
115 }
116 // TFile* file = TFile::Open(Form("xsec-%d.root", pdgNo),
117 // "RECREATE");
118 TArrayF cache(nOk+1);
119 TTree* tree = new TTree(Form("%s_%s", medName, pdgName),
120 Form("%s_%s", medName, pdgName));
121 tree->Branch("xsec", cache.fArray, vars.Data());
122 for (Int_t i = 0; i < n; i++) {
123 cache[0] = tkine[i];
124 Int_t k = 0;
125 for (Int_t j = 0; j < nMech; j++) {
126 if (mechs[j].status) {
127 if (!strcmp(mechs[j].name, "LOSS"))
128 std::cout << tkine[i] << "\t" << mechs[j].values[i] << std::endl;
129 cache[k+1] = mechs[j].values[i];
130 k++;
131 }
132 }
133 std::cout << k << "\t" << (k == nOk) << std::endl;
134 tree->Fill();
135 }
136 tree->Write();
137}
138//____________________________________________________________________
139//
140// EOF
141//