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