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