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"
55 #include "AliDebugVolume.h"
57 #include "AliLegoGenerator.h"
63 //_______________________________________________________________________
82 // Default constructor
86 //_______________________________________________________________________
87 AliLego::AliLego(const AliLego &lego):
112 //_______________________________________________________________________
113 AliLego::AliLego(const char *title, Int_t ntheta, Float_t thetamin,
114 Float_t thetamax, Int_t nphi, Float_t phimin, Float_t phimax,
115 Float_t rmin, Float_t rmax, Float_t zmax):
116 TNamed("Lego Generator",title),
134 // specify the angular limits and the size of the rectangular box
136 fGener = new AliLegoGenerator(ntheta, thetamin, thetamax,
137 nphi, phimin, phimax, rmin, rmax, zmax);
138 fHistRadl = new TH2F("hradl","Radiation length map",
139 ntheta,thetamin,thetamax,nphi,phimin,phimax);
140 fHistAbso = new TH2F("habso","Interaction length map",
141 ntheta,thetamin,thetamax,nphi,phimin,phimax);
142 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
143 ntheta,thetamin,thetamax,nphi,phimin,phimax);
145 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
146 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
147 fDebug = gAlice->GetDebug();
150 //_______________________________________________________________________
151 AliLego::AliLego(const char *title, AliLegoGenerator* generator):
152 TNamed("Lego Generator",title),
170 // specify the angular limits and the size of the rectangular box
173 Float_t c1min, c1max, c2min, c2max;
174 Int_t n1 = fGener->NCoor1();
175 Int_t n2 = fGener->NCoor2();
177 fGener->Coor1Range(c1min, c1max);
178 fGener->Coor2Range(c2min, c2max);
180 fHistRadl = new TH2F("hradl","Radiation length map",
181 n2, c2min, c2max, n1, c1min, c1max);
182 fHistAbso = new TH2F("habso","Interaction length map",
183 n2, c2min, c2max, n1, c1min, c1max);
184 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
185 n2, c2min, c2max, n1, c1min, c1max);
189 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
190 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
191 fDebug = gAlice->GetDebug();
194 //_______________________________________________________________________
208 //_______________________________________________________________________
209 void AliLego::BeginEvent()
212 // --- Set to 0 radiation length, absorption length and g/cm2 ---
219 if (fErrorCondition) DumpVolumes();
220 fVolumesFwd->Delete();
221 fVolumesBwd->Delete();
225 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == 0) fStepBack = 0;
229 //_______________________________________________________________________
230 void AliLego::FinishEvent()
233 // Finish the event and update the histos
236 c1 = fGener->CurCoor1();
237 c2 = fGener->CurCoor2();
238 fHistRadl->Fill(c2,c1,fTotRadl);
239 fHistAbso->Fill(c2,c1,fTotAbso);
240 fHistGcm2->Fill(c2,c1,fTotGcm2);
243 //_______________________________________________________________________
244 void AliLego::FinishRun()
247 // Store histograms in current Root file
253 // Delete histograms from memory
254 fHistRadl->Delete(); fHistRadl=0;
255 fHistAbso->Delete(); fHistAbso=0;
256 fHistGcm2->Delete(); fHistGcm2=0;
258 if (fErrorCondition) DumpVolumes();
261 //_______________________________________________________________________
262 void AliLego::Copy(TObject&) const
265 // Copy *this onto lego -- not implemented
267 Fatal("Copy","Not implemented!\n");
270 //_______________________________________________________________________
271 void AliLego::StepManager()
274 // called from AliRun::Stepmanager from gustep.
275 // Accumulate the 3 parameters step by step
278 Float_t a,z,dens,radl,absl;
281 static Float_t vect[3], dir[3];
285 id = gMC->CurrentVolID(copy);
286 vol = gMC->VolName(id);
287 Float_t step = gMC->TrackStep();
289 TLorentzVector pos, mom;
290 gMC->TrackPosition(pos);
291 gMC->TrackMomentum(mom);
294 if (gMC->IsTrackEntering()) status = 1;
295 if (gMC->IsTrackExiting()) status = 2;
298 // printf("\n volume %s %d", vol, status);
300 // Normal Forward stepping
303 // printf("\n steps fwd %d %s %d %f", fStepsForward, vol, fErrorCondition, step);
306 // store volume if different from previous
309 TClonesArray &lvols = *fVolumesFwd;
310 if (fStepsForward > 0) {
311 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsForward-1]);
312 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
314 fVolumesFwd->RemoveAt(fStepsForward);
315 fVolumesFwd->RemoveAt(fStepsForward+1);
319 new(lvols[fStepsForward++])
320 AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
324 // Get current material properties
326 gMC->CurrentMaterial(a,z,dens,radl,absl);
330 // --- See if we have to stop now
331 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
332 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
333 if (!gMC->IsNewTrack()) {
334 // Not the first step, add past contribution
340 // generate "mirror" particle flying back
342 fStepsBackward = fStepsForward;
344 Float_t pmom[3], orig[3];
345 Float_t polar[3] = {0.,0.,0.};
354 gAlice->GetMCApp()->PushTrack(1, gAlice->GetMCApp()->GetCurrentTrackNumber(),
355 0, pmom, orig, polar, 0., kPNoProcess, ntr);
358 } // not a new track !
360 if (fDebug) fStepBack = 1;
363 } // outside scoring region ?
365 // --- See how long we have to go
371 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
375 fTotAbso += step/absl;
376 fTotRadl += step/radl;
377 fTotGcm2 += step*dens;
383 // Geometry debugging
384 // Fly back and compare volume sequence
386 TClonesArray &lvols = *fVolumesBwd;
387 if (fStepsBackward < fStepsForward) {
388 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[fStepsBackward]);
389 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
391 fVolumesBwd->RemoveAt(fStepsBackward-1);
392 fVolumesBwd->RemoveAt(fStepsBackward-2);
397 // printf("\n steps %d %s %d", fStepsBackward, vol, fErrorCondition);
398 if (fStepsBackward < 0) {
404 new(lvols[fStepsBackward]) AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
406 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsBackward]);
407 if (! (tmp->IsVEqual(vol, copy)) && (!fErrorCondition))
409 printf("\n Problem at (x,y,z): %d %f %f %f, volumes: %s %s step: %f\n",
410 fStepsBackward, pos[0], pos[1], pos[2], tmp->GetName(), vol, step);
417 //_______________________________________________________________________
418 void AliLego::DumpVolumes()
421 // Dump volume sequence in case of error
423 printf("\n Dumping Volume Sequence:");
424 printf("\n ==============================================");
426 for (Int_t i = fStepsForward-1; i>=0; i--)
428 AliDebugVolume* tmp1 = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[i]);
429 AliDebugVolume* tmp2 = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[i]);
431 printf("\n Volume Fwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
433 tmp1->GetName(), tmp1->CopyNumber(), tmp1->Step(),
434 tmp1->X(), tmp1->Y(), tmp1->Z(), tmp1->Status());
435 if (tmp2 && i>= fStepsBackward)
436 printf("\n Volume Bwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
438 tmp2->GetName(), tmp2->CopyNumber(), tmp2->Step(),
439 tmp2->X(), tmp2->Y(), tmp2->Z(), tmp2->Status());
441 printf("\n ............................................................................\n");
443 printf("\n ==============================================\n");