]>
Commit | Line | Data |
---|---|---|
4c039060 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* | |
17 | $Log$ | |
8918e700 | 18 | Revision 1.16 2000/05/26 08:35:03 fca |
19 | Move the check on z after z has been retrieved | |
20 | ||
119a6af6 | 21 | Revision 1.15 2000/05/16 13:10:40 fca |
22 | New method IsNewTrack and fix for a problem in Father-Daughter relations | |
23 | ||
a01a8b12 | 24 | Revision 1.14 2000/04/27 10:38:21 fca |
25 | Correct termination of Lego Run and introduce Lego getter in AliRun | |
26 | ||
838edcaf | 27 | Revision 1.13 2000/04/26 10:17:31 fca |
28 | Changes in Lego for G4 compatibility | |
29 | ||
dffd31ef | 30 | Revision 1.12 2000/04/07 11:12:33 fca |
31 | G4 compatibility changes | |
32 | ||
875c717b | 33 | Revision 1.11 2000/03/22 13:42:26 fca |
34 | SetGenerator does not replace an existing generator, ResetGenerator does | |
35 | ||
ee1dd322 | 36 | Revision 1.10 2000/02/23 16:25:22 fca |
37 | AliVMC and AliGeant3 classes introduced | |
38 | ReadEuclid moved from AliRun to AliModule | |
39 | ||
b13db077 | 40 | Revision 1.9 1999/12/03 10:54:01 fca |
41 | Fix lego summary | |
42 | ||
00719c1b | 43 | Revision 1.8 1999/10/01 09:54:33 fca |
44 | Correct logics for Lego StepManager | |
45 | ||
f059c84a | 46 | Revision 1.7 1999/09/29 09:24:29 fca |
47 | Introduction of the Copyright and cvs Log | |
48 | ||
4c039060 | 49 | */ |
50 | ||
fe4da5cc | 51 | ////////////////////////////////////////////////////////////// |
52 | ////////////////////////////////////////////////////////////// | |
53 | // | |
54 | // Utility class to evaluate the material budget from | |
55 | // a given radius to the surface of an arbitrary cylinder | |
56 | // along radial directions from the centre: | |
57 | // | |
58 | // - radiation length | |
59 | // - Interaction length | |
60 | // - g/cm2 | |
61 | // | |
62 | // Geantinos are shot in the bins in the fNtheta bins in theta | |
63 | // and fNphi bins in phi with specified rectangular limits. | |
64 | // The statistics are accumulated per | |
65 | // fRadMin < r < fRadMax and <0 < z < fZMax | |
66 | // | |
67 | // To activate this option, you can do: | |
68 | // Root > gAlice.RunLego(); | |
69 | // or Root > .x menu.C then select menu item "RunLego" | |
70 | // Note that when calling gAlice->RunLego, an optional list | |
71 | // of arguments may be specified. | |
72 | // | |
73 | //Begin_Html | |
74 | /* | |
1439f98e | 75 | <img src="picts/alilego.gif"> |
fe4da5cc | 76 | */ |
77 | //End_Html | |
78 | // | |
79 | ////////////////////////////////////////////////////////////// | |
80 | ||
81 | #include "TMath.h" | |
1578254f | 82 | #include "AliLego.h" |
8918e700 | 83 | #include "AliLegoGenerator.h" |
fe4da5cc | 84 | #include "AliRun.h" |
85 | #include "AliConst.h" | |
875c717b | 86 | #include "AliMC.h" |
fe4da5cc | 87 | |
88 | ClassImp(AliLego) | |
89 | ||
90 | ||
91 | //___________________________________________ | |
92 | AliLego::AliLego() | |
93 | { | |
8918e700 | 94 | // |
95 | // Default constructor | |
96 | // | |
97 | fHistRadl = 0; | |
98 | fHistAbso = 0; | |
99 | fHistGcm2 = 0; | |
100 | fHistReta = 0; | |
fe4da5cc | 101 | } |
102 | ||
103 | //___________________________________________ | |
8918e700 | 104 | AliLego::AliLego(const char *title, Int_t ntheta, Float_t themin, |
105 | Float_t themax, Int_t nphi, Float_t phimin, Float_t phimax, | |
b13db077 | 106 | Float_t rmin, Float_t rmax, Float_t zmax) |
107 | : TNamed("Lego Generator",title) | |
fe4da5cc | 108 | { |
8918e700 | 109 | // |
110 | // specify the angular limits and the size of the rectangular box | |
111 | // | |
b13db077 | 112 | fGener = new AliLegoGenerator(ntheta, themin, themax, |
113 | nphi, phimin, phimax, rmin, rmax, zmax); | |
114 | ||
ee1dd322 | 115 | gAlice->ResetGenerator(fGener); |
b13db077 | 116 | |
117 | Float_t etamin = -TMath::Log(TMath::Tan(TMath::Min((Double_t)themax*kDegrad/2,TMath::Pi()/2-1.e-10))); | |
118 | Float_t etamax = -TMath::Log(TMath::Tan(TMath::Max((Double_t)themin*kDegrad/2, 1.e-10))); | |
119 | ||
120 | fHistRadl = new TH2F("hradl","Radiation length map", | |
121 | nphi,phimin,phimax,ntheta,themin,themax); | |
122 | fHistAbso = new TH2F("habso","Interaction length map", | |
123 | nphi,phimin,phimax,ntheta,themin,themax); | |
124 | fHistGcm2 = new TH2F("hgcm2","g/cm2 length map", | |
125 | nphi,phimin,phimax,ntheta,themin,themax); | |
126 | fHistReta = new TH2F("hetar","Radiation length vs. eta", | |
127 | nphi,phimin,phimax,ntheta,etamin,etamax); | |
128 | ||
fe4da5cc | 129 | } |
130 | ||
131 | //___________________________________________ | |
132 | AliLego::~AliLego() | |
133 | { | |
8918e700 | 134 | // |
135 | // Destructor | |
136 | // | |
137 | delete fHistRadl; | |
138 | delete fHistAbso; | |
139 | delete fHistGcm2; | |
140 | delete fHistReta; | |
141 | gAlice->ResetGenerator(0); | |
142 | delete fGener; | |
fe4da5cc | 143 | } |
144 | ||
b13db077 | 145 | //___________________________________________ |
dffd31ef | 146 | void AliLego::BeginEvent() |
b13db077 | 147 | { |
8918e700 | 148 | // |
149 | // --- Set to 0 radiation length, absorption length and g/cm2 --- | |
150 | // | |
dffd31ef | 151 | fTotRadl = 0; |
152 | fTotAbso = 0; | |
153 | fTotGcm2 = 0; | |
154 | } | |
155 | ||
156 | //___________________________________________ | |
157 | void AliLego::FinishEvent() | |
158 | { | |
8918e700 | 159 | // |
160 | // Finish the event and update the histos | |
161 | // | |
dffd31ef | 162 | Double_t thed, phid, eta; |
163 | thed = fGener->CurTheta()*kRaddeg; | |
164 | phid = fGener->CurPhi()*kRaddeg; | |
165 | eta = -TMath::Log(TMath::Tan(TMath::Max( | |
b13db077 | 166 | TMath::Min((Double_t)(fGener->CurTheta())/2, |
167 | TMath::Pi()/2-1.e-10),1.e-10))); | |
168 | ||
dffd31ef | 169 | fHistRadl->Fill(phid,thed,fTotRadl); |
170 | fHistAbso->Fill(phid,thed,fTotAbso); | |
171 | fHistGcm2->Fill(phid,thed,fTotGcm2); | |
172 | fHistReta->Fill(phid,eta,fTotRadl); | |
b13db077 | 173 | } |
174 | ||
dffd31ef | 175 | //___________________________________________ |
176 | void AliLego::FinishRun() | |
177 | { | |
8918e700 | 178 | // |
179 | // Store histograms in current Root file | |
180 | // | |
dffd31ef | 181 | fHistRadl->Write(); |
182 | fHistAbso->Write(); | |
183 | fHistGcm2->Write(); | |
184 | fHistReta->Write(); | |
185 | ||
186 | // Delete histograms from memory | |
187 | fHistRadl->Delete(); fHistRadl=0; | |
188 | fHistAbso->Delete(); fHistAbso=0; | |
189 | fHistGcm2->Delete(); fHistGcm2=0; | |
190 | fHistReta->Delete(); fHistReta=0; | |
191 | ||
192 | } | |
193 | ||
8918e700 | 194 | //___________________________________________ |
195 | void AliLego::Copy(AliLego &lego) const | |
196 | { | |
197 | // | |
198 | // Copy *this onto lego -- not implemented | |
199 | // | |
200 | Fatal("Copy","Not implemented!\n"); | |
201 | } | |
dffd31ef | 202 | |
b13db077 | 203 | //___________________________________________ |
204 | void AliLego::StepManager() | |
205 | { | |
206 | // called from AliRun::Stepmanager from gustep. | |
207 | // Accumulate the 3 parameters step by step | |
208 | ||
209 | static Float_t t; | |
210 | Float_t a,z,dens,radl,absl; | |
211 | Int_t i; | |
212 | ||
213 | Float_t step = gMC->TrackStep(); | |
214 | ||
215 | Float_t vect[3], dir[3]; | |
216 | TLorentzVector pos, mom; | |
217 | ||
119a6af6 | 218 | gMC->CurrentMaterial(a,z,dens,radl,absl); |
219 | ||
a01a8b12 | 220 | if (z < 1) return; |
221 | ||
b13db077 | 222 | gMC->TrackPosition(pos); |
223 | gMC->TrackMomentum(mom); | |
b13db077 | 224 | // --- See if we have to stop now |
225 | if (TMath::Abs(pos[2]) > fGener->ZMax() || | |
226 | pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) { | |
a01a8b12 | 227 | if (!gMC->IsNewTrack()) { |
b13db077 | 228 | // Not the first step, add past contribution |
229 | fTotAbso += t/absl; | |
230 | fTotRadl += t/radl; | |
231 | fTotGcm2 += t*dens; | |
232 | } | |
233 | gMC->StopTrack(); | |
234 | return; | |
235 | } | |
236 | ||
237 | // --- See how long we have to go | |
238 | for(i=0;i<3;++i) { | |
239 | vect[i]=pos[i]; | |
240 | dir[i]=mom[i]; | |
241 | } | |
242 | ||
243 | t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax()); | |
244 | ||
245 | if(step) { | |
246 | fTotAbso += step/absl; | |
247 | fTotRadl += step/radl; | |
248 | fTotGcm2 += step*dens; | |
249 | } | |
250 | } | |
251 |