1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.20 2000/11/30 07:12:48 alibrary
19 Introducing new Rndm and QA classes
21 Revision 1.19 2000/10/26 14:13:05 morsch
22 - Change from coordinates theta, phi to general coordinates Coor1 and Coor2.
23 - Lego generator instance can be passed in constructor.
25 Revision 1.18 2000/10/02 21:28:14 fca
26 Removal of useless dependecies via forward declarations
28 Revision 1.17 2000/07/12 08:56:25 fca
29 Coding convention correction and warning removal
31 Revision 1.16 2000/05/26 08:35:03 fca
32 Move the check on z after z has been retrieved
34 Revision 1.15 2000/05/16 13:10:40 fca
35 New method IsNewTrack and fix for a problem in Father-Daughter relations
37 Revision 1.14 2000/04/27 10:38:21 fca
38 Correct termination of Lego Run and introduce Lego getter in AliRun
40 Revision 1.13 2000/04/26 10:17:31 fca
41 Changes in Lego for G4 compatibility
43 Revision 1.12 2000/04/07 11:12:33 fca
44 G4 compatibility changes
46 Revision 1.11 2000/03/22 13:42:26 fca
47 SetGenerator does not replace an existing generator, ResetGenerator does
49 Revision 1.10 2000/02/23 16:25:22 fca
50 AliVMC and AliGeant3 classes introduced
51 ReadEuclid moved from AliRun to AliModule
53 Revision 1.9 1999/12/03 10:54:01 fca
56 Revision 1.8 1999/10/01 09:54:33 fca
57 Correct logics for Lego StepManager
59 Revision 1.7 1999/09/29 09:24:29 fca
60 Introduction of the Copyright and cvs Log
64 //////////////////////////////////////////////////////////////
65 //////////////////////////////////////////////////////////////
67 // Utility class to evaluate the material budget from
68 // a given radius to the surface of an arbitrary cylinder
69 // along radial directions from the centre:
72 // - Interaction length
75 // Geantinos are shot in the bins in the fNtheta bins in theta
76 // and fNphi bins in phi with specified rectangular limits.
77 // The statistics are accumulated per
78 // fRadMin < r < fRadMax and <0 < z < fZMax
80 // To activate this option, you can do:
81 // Root > gAlice.RunLego();
82 // or Root > .x menu.C then select menu item "RunLego"
83 // Note that when calling gAlice->RunLego, an optional list
84 // of arguments may be specified.
88 <img src="picts/alilego.gif">
92 //////////////////////////////////////////////////////////////
97 #include "AliLegoGenerator.h"
105 //___________________________________________
109 // Default constructor
117 //___________________________________________
118 AliLego::AliLego(const char *title, Int_t ntheta, Float_t thetamin,
119 Float_t thetamax, Int_t nphi, Float_t phimin, Float_t phimax,
120 Float_t rmin, Float_t rmax, Float_t zmax)
121 : TNamed("Lego Generator",title)
124 // specify the angular limits and the size of the rectangular box
126 fGener = new AliLegoGenerator(ntheta, thetamin, thetamax,
127 nphi, phimin, phimax, rmin, rmax, zmax);
131 fHistRadl = new TH2F("hradl","Radiation length map",
132 ntheta,thetamin,thetamax,nphi,phimin,phimax);
133 fHistAbso = new TH2F("habso","Interaction length map",
134 ntheta,thetamin,thetamax,nphi,phimin,phimax);
135 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
136 ntheta,thetamin,thetamax,nphi,phimin,phimax);
139 AliLego::AliLego(const char *title, AliLegoGenerator* generator)
140 : TNamed("Lego Generator",title)
143 // specify the angular limits and the size of the rectangular box
146 Float_t c1min, c1max, c2min, c2max;
147 Int_t n1 = fGener->NCoor1();
148 Int_t n2 = fGener->NCoor2();
150 fGener->Coor1Range(c1min, c1max);
151 fGener->Coor2Range(c2min, c2max);
153 fHistRadl = new TH2F("hradl","Radiation length map",
154 n2, c2min, c2max, n1, c1min, c1max);
155 fHistAbso = new TH2F("habso","Interaction length map",
156 n2, c2min, c2max, n1, c1min, c1max);
157 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
158 n2, c2min, c2max, n1, c1min, c1max);
161 //___________________________________________
173 //___________________________________________
174 void AliLego::BeginEvent()
177 // --- Set to 0 radiation length, absorption length and g/cm2 ---
184 //___________________________________________
185 void AliLego::FinishEvent()
188 // Finish the event and update the histos
191 c1 = fGener->CurCoor1();
192 c2 = fGener->CurCoor2();
193 fHistRadl->Fill(c2,c1,fTotRadl);
194 fHistAbso->Fill(c2,c1,fTotAbso);
195 fHistGcm2->Fill(c2,c1,fTotGcm2);
198 //___________________________________________
199 void AliLego::FinishRun()
202 // Store histograms in current Root file
208 // Delete histograms from memory
209 fHistRadl->Delete(); fHistRadl=0;
210 fHistAbso->Delete(); fHistAbso=0;
211 fHistGcm2->Delete(); fHistGcm2=0;
214 //___________________________________________
215 void AliLego::Copy(AliLego &lego) const
218 // Copy *this onto lego -- not implemented
220 Fatal("Copy","Not implemented!\n");
223 //___________________________________________
224 void AliLego::StepManager()
226 // called from AliRun::Stepmanager from gustep.
227 // Accumulate the 3 parameters step by step
230 Float_t a,z,dens,radl,absl;
233 Float_t step = gMC->TrackStep();
235 Float_t vect[3], dir[3];
236 TLorentzVector pos, mom;
238 gMC->CurrentMaterial(a,z,dens,radl,absl);
242 gMC->TrackPosition(pos);
243 gMC->TrackMomentum(mom);
244 // --- See if we have to stop now
245 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
246 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
247 if (!gMC->IsNewTrack()) {
248 // Not the first step, add past contribution
253 // Int_t id = gMC->CurrentVolID(copy);
254 // char* vol = gMC->VolName(id);
256 // printf("\n %f %f %f %f %s ", fTotRadl, vect[0], vect[1], vect[2], vol);
263 // --- See how long we have to go
269 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
272 fTotAbso += step/absl;
273 fTotRadl += step/radl;
274 fTotGcm2 += step*dens;
276 // Int_t id = gMC->CurrentVolID(copy);
277 // char* vol = gMC->VolName(id);
278 // printf("\n %f %f %f %f %s %f ", fTotRadl, vect[0], vect[1], vect[2], vol, t);