]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliLego.cxx
Containers definition
[u/mrichter/AliRoot.git] / STEER / AliLego.cxx
CommitLineData
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$
0b27b8c7 18Revision 1.18 2000/10/02 21:28:14 fca
19Removal of useless dependecies via forward declarations
20
94de3818 21Revision 1.17 2000/07/12 08:56:25 fca
22Coding convention correction and warning removal
23
8918e700 24Revision 1.16 2000/05/26 08:35:03 fca
25Move the check on z after z has been retrieved
26
119a6af6 27Revision 1.15 2000/05/16 13:10:40 fca
28New method IsNewTrack and fix for a problem in Father-Daughter relations
29
a01a8b12 30Revision 1.14 2000/04/27 10:38:21 fca
31Correct termination of Lego Run and introduce Lego getter in AliRun
32
838edcaf 33Revision 1.13 2000/04/26 10:17:31 fca
34Changes in Lego for G4 compatibility
35
dffd31ef 36Revision 1.12 2000/04/07 11:12:33 fca
37G4 compatibility changes
38
875c717b 39Revision 1.11 2000/03/22 13:42:26 fca
40SetGenerator does not replace an existing generator, ResetGenerator does
41
ee1dd322 42Revision 1.10 2000/02/23 16:25:22 fca
43AliVMC and AliGeant3 classes introduced
44ReadEuclid moved from AliRun to AliModule
45
b13db077 46Revision 1.9 1999/12/03 10:54:01 fca
47Fix lego summary
48
00719c1b 49Revision 1.8 1999/10/01 09:54:33 fca
50Correct logics for Lego StepManager
51
f059c84a 52Revision 1.7 1999/09/29 09:24:29 fca
53Introduction of the Copyright and cvs Log
54
4c039060 55*/
56
fe4da5cc 57//////////////////////////////////////////////////////////////
58//////////////////////////////////////////////////////////////
59//
60// Utility class to evaluate the material budget from
61// a given radius to the surface of an arbitrary cylinder
62// along radial directions from the centre:
63//
64// - radiation length
65// - Interaction length
66// - g/cm2
67//
68// Geantinos are shot in the bins in the fNtheta bins in theta
69// and fNphi bins in phi with specified rectangular limits.
70// The statistics are accumulated per
71// fRadMin < r < fRadMax and <0 < z < fZMax
72//
73// To activate this option, you can do:
74// Root > gAlice.RunLego();
75// or Root > .x menu.C then select menu item "RunLego"
76// Note that when calling gAlice->RunLego, an optional list
77// of arguments may be specified.
78//
79//Begin_Html
80/*
1439f98e 81<img src="picts/alilego.gif">
fe4da5cc 82*/
83//End_Html
84//
85//////////////////////////////////////////////////////////////
86
87#include "TMath.h"
94de3818 88
1578254f 89#include "AliLego.h"
8918e700 90#include "AliLegoGenerator.h"
fe4da5cc 91#include "AliRun.h"
92#include "AliConst.h"
875c717b 93#include "AliMC.h"
94de3818 94#include "TH2.h"
fe4da5cc 95
96ClassImp(AliLego)
97
98
99//___________________________________________
100AliLego::AliLego()
101{
8918e700 102 //
103 // Default constructor
104 //
105 fHistRadl = 0;
106 fHistAbso = 0;
107 fHistGcm2 = 0;
108 fHistReta = 0;
fe4da5cc 109}
110
111//___________________________________________
0b27b8c7 112AliLego::AliLego(const char *title, Int_t ntheta, Float_t thetamin,
113 Float_t thetamax, Int_t nphi, Float_t phimin, Float_t phimax,
b13db077 114 Float_t rmin, Float_t rmax, Float_t zmax)
115 : TNamed("Lego Generator",title)
fe4da5cc 116{
8918e700 117 //
118 // specify the angular limits and the size of the rectangular box
119 //
0b27b8c7 120 fGener = new AliLegoGenerator(ntheta, thetamin, thetamax,
b13db077 121 nphi, phimin, phimax, rmin, rmax, zmax);
122
b13db077 123
0b27b8c7 124
b13db077 125 fHistRadl = new TH2F("hradl","Radiation length map",
0b27b8c7 126 ntheta,thetamin,thetamax,nphi,phimin,phimax);
b13db077 127 fHistAbso = new TH2F("habso","Interaction length map",
0b27b8c7 128 ntheta,thetamin,thetamax,nphi,phimin,phimax);
b13db077 129 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
0b27b8c7 130 ntheta,thetamin,thetamax,nphi,phimin,phimax);
131}
132
133AliLego::AliLego(const char *title, AliLegoGenerator* generator)
134 : TNamed("Lego Generator",title)
135{
136 //
137 // specify the angular limits and the size of the rectangular box
138 //
139 fGener = generator;
140 Float_t c1min, c1max, c2min, c2max;
141 Int_t n1 = fGener->NCoor1();
142 Int_t n2 = fGener->NCoor2();
143
144 fGener->Coor1Range(c1min, c1max);
145 fGener->Coor2Range(c2min, c2max);
b13db077 146
0b27b8c7 147 fHistRadl = new TH2F("hradl","Radiation length map",
148 n1, c1min, c1max, n2, c2min, c2max);
149 fHistAbso = new TH2F("habso","Interaction length map",
150 n1, c1min, c1max, n2, c2min, c2max);
151 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
152 n1, c1min, c1max, n2, c2min, c2max);
fe4da5cc 153}
154
155//___________________________________________
156AliLego::~AliLego()
157{
8918e700 158 //
159 // Destructor
160 //
161 delete fHistRadl;
162 delete fHistAbso;
163 delete fHistGcm2;
8918e700 164 delete fGener;
fe4da5cc 165}
166
b13db077 167//___________________________________________
dffd31ef 168void AliLego::BeginEvent()
b13db077 169{
8918e700 170 //
171 // --- Set to 0 radiation length, absorption length and g/cm2 ---
172 //
dffd31ef 173 fTotRadl = 0;
174 fTotAbso = 0;
175 fTotGcm2 = 0;
176}
177
178//___________________________________________
179void AliLego::FinishEvent()
180{
8918e700 181 //
182 // Finish the event and update the histos
183 //
0b27b8c7 184 Double_t c1, c2;
185 c1 = fGener->CurCoor1();
186 c2 = fGener->CurCoor2();
187 fHistRadl->Fill(c1,c2,fTotRadl);
188 fHistAbso->Fill(c1,c2,fTotAbso);
189 fHistGcm2->Fill(c1,c2,fTotGcm2);
b13db077 190}
191
dffd31ef 192//___________________________________________
193void AliLego::FinishRun()
194{
8918e700 195 //
196 // Store histograms in current Root file
197 //
dffd31ef 198 fHistRadl->Write();
199 fHistAbso->Write();
200 fHistGcm2->Write();
dffd31ef 201
202 // Delete histograms from memory
203 fHistRadl->Delete(); fHistRadl=0;
204 fHistAbso->Delete(); fHistAbso=0;
205 fHistGcm2->Delete(); fHistGcm2=0;
dffd31ef 206}
207
8918e700 208//___________________________________________
209void AliLego::Copy(AliLego &lego) const
210{
211 //
212 // Copy *this onto lego -- not implemented
213 //
214 Fatal("Copy","Not implemented!\n");
215}
dffd31ef 216
b13db077 217//___________________________________________
218void AliLego::StepManager()
219{
220// called from AliRun::Stepmanager from gustep.
221// Accumulate the 3 parameters step by step
222
223 static Float_t t;
224 Float_t a,z,dens,radl,absl;
225 Int_t i;
226
227 Float_t step = gMC->TrackStep();
228
229 Float_t vect[3], dir[3];
230 TLorentzVector pos, mom;
231
119a6af6 232 gMC->CurrentMaterial(a,z,dens,radl,absl);
233
a01a8b12 234 if (z < 1) return;
235
b13db077 236 gMC->TrackPosition(pos);
237 gMC->TrackMomentum(mom);
b13db077 238// --- See if we have to stop now
239 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
240 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
a01a8b12 241 if (!gMC->IsNewTrack()) {
b13db077 242 // Not the first step, add past contribution
243 fTotAbso += t/absl;
244 fTotRadl += t/radl;
245 fTotGcm2 += t*dens;
0b27b8c7 246// Int_t copy;
247// Int_t id = gMC->CurrentVolID(copy);
248// char* vol = gMC->VolName(id);
249
250// printf("\n %f %f %f %f %s ", fTotRadl, vect[0], vect[1], vect[2], vol);
251
b13db077 252 }
253 gMC->StopTrack();
254 return;
255 }
256
257// --- See how long we have to go
258 for(i=0;i<3;++i) {
259 vect[i]=pos[i];
260 dir[i]=mom[i];
261 }
262
263 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
264
265 if(step) {
266 fTotAbso += step/absl;
267 fTotRadl += step/radl;
268 fTotGcm2 += step*dens;
0b27b8c7 269// Int_t copy;
270// Int_t id = gMC->CurrentVolID(copy);
271// char* vol = gMC->VolName(id);
272// printf("\n %f %f %f %f %s %f ", fTotRadl, vect[0], vect[1], vect[2], vol, t);
b13db077 273 }
274}
275
0b27b8c7 276
277
278
279