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 //////////////////////////////////////////////////////////////
19 //////////////////////////////////////////////////////////////
21 // Utility class to evaluate the material budget from
22 // a given radius to the surface of an arbitrary cylinder
23 // along radial directions from the centre:
26 // - Interaction length
29 // Geantinos are shot in the bins in the fNtheta bins in theta
30 // and fNphi bins in phi with specified rectangular limits.
31 // The statistics are accumulated per
32 // fRadMin < r < fRadMax and <0 < z < fZMax
34 // To activate this option, you can do:
35 // Root > gAlice.RunLego();
36 // or Root > .x menu.C then select menu item "RunLego"
37 // Note that when calling gAlice->RunLego, an optional list
38 // of arguments may be specified.
42 <img src="picts/alilego.gif">
46 //////////////////////////////////////////////////////////////
48 #include "TClonesArray.h"
52 #include "TVirtualMC.h"
54 #include "AliDebugVolume.h"
56 #include "AliLegoGenerator.h"
62 //_______________________________________________________________________
81 // Default constructor
85 //_______________________________________________________________________
86 AliLego::AliLego(const AliLego &lego):
111 //_______________________________________________________________________
112 AliLego::AliLego(const char *title, Int_t ntheta, Float_t thetamin,
113 Float_t thetamax, Int_t nphi, Float_t phimin, Float_t phimax,
114 Float_t rmin, Float_t rmax, Float_t zmax):
115 TNamed("Lego Generator",title),
133 // specify the angular limits and the size of the rectangular box
135 fGener = new AliLegoGenerator(ntheta, thetamin, thetamax,
136 nphi, phimin, phimax, rmin, rmax, zmax);
137 fHistRadl = new TH2F("hradl","Radiation length map",
138 ntheta,thetamin,thetamax,nphi,phimin,phimax);
139 fHistAbso = new TH2F("habso","Interaction length map",
140 ntheta,thetamin,thetamax,nphi,phimin,phimax);
141 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
142 ntheta,thetamin,thetamax,nphi,phimin,phimax);
144 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
145 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
146 fDebug = gAlice->GetDebug();
149 //_______________________________________________________________________
150 AliLego::AliLego(const char *title, AliLegoGenerator* generator):
151 TNamed("Lego Generator",title),
169 // specify the angular limits and the size of the rectangular box
172 Float_t c1min, c1max, c2min, c2max;
173 Int_t n1 = fGener->NCoor1();
174 Int_t n2 = fGener->NCoor2();
176 fGener->Coor1Range(c1min, c1max);
177 fGener->Coor2Range(c2min, c2max);
179 fHistRadl = new TH2F("hradl","Radiation length map",
180 n2, c2min, c2max, n1, c1min, c1max);
181 fHistAbso = new TH2F("habso","Interaction length map",
182 n2, c2min, c2max, n1, c1min, c1max);
183 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
184 n2, c2min, c2max, n1, c1min, c1max);
188 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
189 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
190 fDebug = gAlice->GetDebug();
193 //_______________________________________________________________________
207 //_______________________________________________________________________
208 void AliLego::BeginEvent()
211 // --- Set to 0 radiation length, absorption length and g/cm2 ---
218 if (fErrorCondition) DumpVolumes();
219 fVolumesFwd->Delete();
220 fVolumesBwd->Delete();
224 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == 0) fStepBack = 0;
228 //_______________________________________________________________________
229 void AliLego::FinishEvent()
232 // Finish the event and update the histos
235 c1 = fGener->CurCoor1();
236 c2 = fGener->CurCoor2();
237 fHistRadl->Fill(c2,c1,fTotRadl);
238 fHistAbso->Fill(c2,c1,fTotAbso);
239 fHistGcm2->Fill(c2,c1,fTotGcm2);
242 //_______________________________________________________________________
243 void AliLego::FinishRun()
246 // Store histograms in current Root file
252 // Delete histograms from memory
253 fHistRadl->Delete(); fHistRadl=0;
254 fHistAbso->Delete(); fHistAbso=0;
255 fHistGcm2->Delete(); fHistGcm2=0;
257 if (fErrorCondition) DumpVolumes();
260 //_______________________________________________________________________
261 void AliLego::Copy(TObject&) const
264 // Copy *this onto lego -- not implemented
266 Fatal("Copy","Not implemented!\n");
269 //_______________________________________________________________________
270 void AliLego::StepManager()
273 // called from AliRun::Stepmanager from gustep.
274 // Accumulate the 3 parameters step by step
277 Float_t a,z,dens,radl,absl;
280 static Float_t vect[3], dir[3];
284 id = gMC->CurrentVolID(copy);
285 vol = gMC->VolName(id);
286 Float_t step = gMC->TrackStep();
288 TLorentzVector pos, mom;
289 gMC->TrackPosition(pos);
290 gMC->TrackMomentum(mom);
293 if (gMC->IsTrackEntering()) status = 1;
294 if (gMC->IsTrackExiting()) status = 2;
297 // printf("\n volume %s %d", vol, status);
299 // Normal Forward stepping
302 // printf("\n steps fwd %d %s %d %f", fStepsForward, vol, fErrorCondition, step);
305 // store volume if different from previous
308 TClonesArray &lvols = *fVolumesFwd;
309 if (fStepsForward > 0) {
310 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsForward-1]);
311 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
313 fVolumesFwd->RemoveAt(fStepsForward);
314 fVolumesFwd->RemoveAt(fStepsForward+1);
318 new(lvols[fStepsForward++])
319 AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
323 // Get current material properties
325 gMC->CurrentMaterial(a,z,dens,radl,absl);
329 // --- See if we have to stop now
330 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
331 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
332 if (!gMC->IsNewTrack()) {
333 // Not the first step, add past contribution
339 // generate "mirror" particle flying back
341 fStepsBackward = fStepsForward;
343 Float_t pmom[3], orig[3];
344 Float_t polar[3] = {0.,0.,0.};
353 gAlice->GetMCApp()->PushTrack(1, gAlice->GetMCApp()->GetCurrentTrackNumber(),
354 0, pmom, orig, polar, 0., kPNoProcess, ntr);
357 } // not a new track !
359 if (fDebug) fStepBack = 1;
362 } // outside scoring region ?
364 // --- See how long we have to go
370 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
374 fTotAbso += step/absl;
375 fTotRadl += step/radl;
376 fTotGcm2 += step*dens;
382 // Geometry debugging
383 // Fly back and compare volume sequence
385 TClonesArray &lvols = *fVolumesBwd;
386 if (fStepsBackward < fStepsForward) {
387 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[fStepsBackward]);
388 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
390 fVolumesBwd->RemoveAt(fStepsBackward-1);
391 fVolumesBwd->RemoveAt(fStepsBackward-2);
396 // printf("\n steps %d %s %d", fStepsBackward, vol, fErrorCondition);
397 if (fStepsBackward < 0) {
403 new(lvols[fStepsBackward]) AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
405 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsBackward]);
406 if (! (tmp->IsVEqual(vol, copy)) && (!fErrorCondition))
408 printf("\n Problem at (x,y,z): %d %f %f %f, volumes: %s %s step: %f\n",
409 fStepsBackward, pos[0], pos[1], pos[2], tmp->GetName(), vol, step);
416 //_______________________________________________________________________
417 void AliLego::DumpVolumes()
420 // Dump volume sequence in case of error
422 printf("\n Dumping Volume Sequence:");
423 printf("\n ==============================================");
425 for (Int_t i = fStepsForward-1; i>=0; i--)
427 AliDebugVolume* tmp1 = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[i]);
428 AliDebugVolume* tmp2 = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[i]);
430 printf("\n Volume Fwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
432 tmp1->GetName(), tmp1->CopyNumber(), tmp1->Step(),
433 tmp1->X(), tmp1->Y(), tmp1->Z(), tmp1->Status());
434 if (tmp2 && i>= fStepsBackward)
435 printf("\n Volume Bwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
437 tmp2->GetName(), tmp2->CopyNumber(), tmp2->Step(),
438 tmp2->X(), tmp2->Y(), tmp2->Z(), tmp2->Status());
440 printf("\n ............................................................................\n");
442 printf("\n ==============================================\n");