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 //_______________________________________________________________________
83 // Default constructor
87 //_______________________________________________________________________
88 AliLego::AliLego(const AliLego &lego):
114 //_______________________________________________________________________
115 AliLego::AliLego(const char *title, Int_t ntheta, Float_t thetamin,
116 Float_t thetamax, Int_t nphi, Float_t phimin, Float_t phimax,
117 Float_t rmin, Float_t rmax, Float_t zmax):
118 TNamed("Lego Generator",title),
137 // specify the angular limits and the size of the rectangular box
139 fGener = new AliLegoGenerator(ntheta, thetamin, thetamax,
140 nphi, phimin, phimax, rmin, rmax, zmax);
141 fHistRadl = new TH2F("hradl","Radiation length map",
142 ntheta,thetamin,thetamax,nphi,phimin,phimax);
143 fHistAbso = new TH2F("habso","Interaction length map",
144 ntheta,thetamin,thetamax,nphi,phimin,phimax);
145 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
146 ntheta,thetamin,thetamax,nphi,phimin,phimax);
148 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
149 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
150 fDebug = AliDebugLevel();
153 //_______________________________________________________________________
154 AliLego::AliLego(const char *title, AliLegoGenerator* generator):
155 TNamed("Lego Generator",title),
174 // specify the angular limits and the size of the rectangular box
177 Float_t c1min, c1max, c2min, c2max;
178 Int_t n1 = fGener->NCoor1();
179 Int_t n2 = fGener->NCoor2();
181 fGener->Coor1Range(c1min, c1max);
182 fGener->Coor2Range(c2min, c2max);
184 fHistRadl = new TH2F("hradl","Radiation length map",
185 n2, c2min, c2max, n1, c1min, c1max);
186 fHistAbso = new TH2F("habso","Interaction length map",
187 n2, c2min, c2max, n1, c1min, c1max);
188 fHistGcm2 = new TH2F("hgcm2","g/cm2 length map",
189 n2, c2min, c2max, n1, c1min, c1max);
193 fVolumesFwd = new TClonesArray("AliDebugVolume",1000);
194 fVolumesBwd = new TClonesArray("AliDebugVolume",1000);
195 fDebug = AliDebugLevel();
198 //_______________________________________________________________________
212 //_______________________________________________________________________
213 void AliLego::BeginEvent()
216 // --- Set to 0 radiation length, absorption length and g/cm2 ---
224 if (fErrorCondition) ToAliDebug(1, DumpVolumes());
225 fVolumesFwd->Delete();
226 fVolumesBwd->Delete();
230 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == 0) fStepBack = 0;
234 //_______________________________________________________________________
235 void AliLego::FinishEvent()
238 // Finish the event and update the histos
241 c1 = fGener->CurCoor1();
242 c2 = fGener->CurCoor2();
243 fHistRadl->Fill(c2,c1,fTotRadl);
244 fHistAbso->Fill(c2,c1,fTotAbso);
245 fHistGcm2->Fill(c2,c1,fTotGcm2);
248 //_______________________________________________________________________
249 void AliLego::FinishRun()
252 // Store histograms in current Root file
258 // Delete histograms from memory
259 fHistRadl->Delete(); fHistRadl=0;
260 fHistAbso->Delete(); fHistAbso=0;
261 fHistGcm2->Delete(); fHistGcm2=0;
263 if (fErrorCondition) ToAliError(DumpVolumes());
266 //_______________________________________________________________________
267 void AliLego::Copy(TObject&) const
270 // Copy *this onto lego -- not implemented
272 AliFatal("Not implemented!");
275 //_______________________________________________________________________
276 void AliLego::StepManager()
279 // called from AliRun::Stepmanager from gustep.
280 // Accumulate the 3 parameters step by step
283 Float_t a,z,dens,radl,absl;
286 static Float_t vect[3], dir[3];
290 id = gMC->CurrentVolID(copy);
291 vol = gMC->VolName(id);
292 Float_t step = gMC->TrackStep();
294 TLorentzVector pos, mom;
295 gMC->TrackPosition(pos);
296 gMC->TrackMomentum(mom);
299 if (gMC->IsTrackEntering()) status = 1;
300 if (gMC->IsTrackExiting()) status = 2;
303 // printf("\n volume %s %d", vol, status);
305 // Normal Forward stepping
308 // printf("\n steps fwd %d %s %d %f", fStepsForward, vol, fErrorCondition, step);
311 // store volume if different from previous
314 TClonesArray &lvols = *fVolumesFwd;
315 if (fStepsForward > 0) {
316 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsForward-1]);
317 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
319 fVolumesFwd->RemoveAt(fStepsForward);
320 fVolumesFwd->RemoveAt(fStepsForward+1);
324 new(lvols[fStepsForward++])
325 AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
329 // Get current material properties
331 gMC->CurrentMaterial(a,z,dens,radl,absl);
336 // --- See if we have to stop now
337 if (TMath::Abs(pos[2]) > fGener->ZMax() ||
338 pos[0]*pos[0] +pos[1]*pos[1] > fGener->RadMax()*fGener->RadMax()) {
339 if (!gMC->IsNewTrack()) {
340 // Not the first step, add past contribution
342 if (absl) fTotAbso += t/absl;
343 if (radl) fTotRadl += t/radl;
347 // printf("We will stop now %5d %13.3f !\n", fStopped, t);
348 // printf("%13.3f %13.3f %13.3f %13.3f %13.3f %13.3f %13.3f %s %13.3f\n",
349 // pos[2], TMath::Sqrt(pos[0] * pos[0] + pos[1] * pos[1]), step, a, z, radl, absl, gMC->CurrentVolName(), fTotRadl);
352 // generate "mirror" particle flying back
354 fStepsBackward = fStepsForward;
356 Float_t pmom[3], orig[3];
357 Float_t polar[3] = {0.,0.,0.};
366 gAlice->GetMCApp()->PushTrack(1, gAlice->GetMCApp()->GetCurrentTrackNumber(),
367 0, pmom, orig, polar, 0., kPNoProcess, ntr);
370 } // not a new track !
372 if (fDebug) fStepBack = 1;
376 } // outside scoring region ?
378 // --- See how long we have to go
384 t = fGener->PropagateCylinder(vect,dir,fGener->RadMax(),fGener->ZMax());
388 if (absl) fTotAbso += step/absl;
389 if (radl) fTotRadl += step/radl;
390 fTotGcm2 += step*dens;
391 // printf("%13.3f %13.3f %13.3f %13.3f %13.3f %13.3f %13.3f %s %13.3f\n",
392 // pos[2], TMath::Sqrt(pos[0] * pos[0] + pos[1] * pos[1]), step, a, z, radl, absl, gMC->CurrentVolName(), fTotRadl);
398 // Geometry debugging
399 // Fly back and compare volume sequence
401 TClonesArray &lvols = *fVolumesBwd;
402 if (fStepsBackward < fStepsForward) {
403 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[fStepsBackward]);
404 if (tmp->IsVEqual(vol, copy) && gMC->IsTrackEntering()) {
406 fVolumesBwd->RemoveAt(fStepsBackward-1);
407 fVolumesBwd->RemoveAt(fStepsBackward-2);
412 // printf("\n steps %d %s %d", fStepsBackward, vol, fErrorCondition);
413 if (fStepsBackward < 0) {
419 new(lvols[fStepsBackward]) AliDebugVolume(vol,copy,step,pos[0], pos[1], pos[2], status);
421 AliDebugVolume* tmp = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[fStepsBackward]);
422 if (! (tmp->IsVEqual(vol, copy)) && (!fErrorCondition))
424 AliWarning(Form("Problem at (x,y,z): %d %f %f %f, volumes: %s %s step: %f\n",
425 fStepsBackward, pos[0], pos[1], pos[2], tmp->GetName(), vol, step));
432 //_______________________________________________________________________
433 void AliLego::DumpVolumes()
436 // Dump volume sequence in case of error
438 printf("\n Dumping Volume Sequence:");
439 printf("\n ==============================================");
441 for (Int_t i = fStepsForward-1; i>=0; i--)
443 AliDebugVolume* tmp1 = dynamic_cast<AliDebugVolume*>((*fVolumesFwd)[i]);
444 AliDebugVolume* tmp2 = dynamic_cast<AliDebugVolume*>((*fVolumesBwd)[i]);
446 printf("\n Volume Fwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
448 tmp1->GetName(), tmp1->CopyNumber(), tmp1->Step(),
449 tmp1->X(), tmp1->Y(), tmp1->Z(), tmp1->Status());
450 if (tmp2 && i>= fStepsBackward)
451 printf("\n Volume Bwd: %3d: %5s (%3d) step: %12.5e (x,y,z) (%12.5e %12.5e %12.5e) status: %9s \n"
453 tmp2->GetName(), tmp2->CopyNumber(), tmp2->Step(),
454 tmp2->X(), tmp2->Y(), tmp2->Z(), tmp2->Status());
456 printf("\n ............................................................................\n");
458 printf("\n ==============================================\n");