]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALv2.cxx
remove unnecessary histogram booking, filling, storing; QA classes handle that now
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALv2.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2004, 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 /* $Id$ */
17
18 //_________________________________________________________________________
19 //*-- Implementation version v2 of EMCAL Manager class; SHASHLYK version
20 //*-- An object of this class does not produce digits
21 //*-- It is the one to use if you do want to produce outputs in TREEH 
22 //*--                  
23 //*-- Author : Aleksei Pavlinov (WSU)
24
25 // This Class not stores information on all particles prior to EMCAL entry - in order to facilitate analysis.
26 // This is done by setting fIShunt =2, and flagging all parents of particles entering the EMCAL.
27
28 #include <cassert>
29 // --- ROOT system ---
30 #include <TBrowser.h>
31 #include <TClonesArray.h>
32 #include <TH1.h>
33 #include <TH2.h>
34 #include <TParticle.h>
35 #include <TROOT.h>
36 #include <TVirtualMC.h>
37
38 // --- Standard library ---
39
40 // --- AliRoot header files ---
41 #include "AliEMCALv2.h"
42 #include "AliEMCALHit.h"
43 #include "AliEMCALGeometry.h"
44 #include "AliRun.h"
45 #include "AliHeader.h"
46 #include "AliMC.h"
47 #include "AliStack.h"
48 #include "AliPoints.h"
49 // for TRD1 case only; May 31,2006
50
51 ClassImp(AliEMCALv2)
52
53 //______________________________________________________________________
54 AliEMCALv2::AliEMCALv2()
55   : AliEMCALv1(), 
56     fHDe(0),
57     fHNhits(0)
58 {
59   // ctor
60 }
61
62 //______________________________________________________________________
63 AliEMCALv2::AliEMCALv2(const char *name, const char *title)
64   : AliEMCALv1(name,title),
65     fHDe(0),
66     fHNhits(0)
67 {
68     // Standard Creator.
69
70     fHits= new TClonesArray("AliEMCALHit",1000);
71     gAlice->GetMCApp()->AddHitList(fHits);
72
73     fNhits    = 0;
74     fIshunt   = 2; // All hits are associated with particles entering the calorimeter
75     fTimeCut  = 30e-09;
76
77     fGeometry = GetGeometry(); 
78     fHDe = fHNhits = 0;
79     //    if (gDebug>0){
80     if (1){
81       TH1::AddDirectory(0);
82       fHDe    = new TH1F("fHDe","De in EMCAL", 1000, 0., 10.);
83       fHNhits = new TH1F("fHNhits","#hits in EMCAL", 2001, -0.5, 2000.5);
84       fHistograms->Add(fHDe);
85       fHistograms->Add(fHNhits);
86       TH1::AddDirectory(1);
87     }
88 }
89
90 //______________________________________________________________________
91 AliEMCALv2::~AliEMCALv2(){
92     // dtor
93
94     if ( fHits) {
95         fHits->Delete();
96         delete fHits;
97         fHits = 0;
98     }
99 }
100
101 //______________________________________________________________________
102 void AliEMCALv2::AddHit(Int_t shunt, Int_t primary, Int_t tracknumber, Int_t iparent, Float_t ienergy, 
103                         Int_t id, Float_t * hits,Float_t * p){
104     // Add a hit to the hit list.
105     // An EMCAL hit is the sum of all hits in a tower section
106     //   originating from the same entering particle 
107     static Int_t hitCounter;
108     static AliEMCALHit *newHit, *curHit;
109     static Bool_t deja;
110
111     deja = kFALSE;
112
113     newHit = new AliEMCALHit(shunt, primary, tracknumber, iparent, ienergy, id, hits, p);
114     for ( hitCounter = fNhits-1; hitCounter >= 0 && !deja; hitCounter-- ) {
115         curHit = (AliEMCALHit*) (*fHits)[hitCounter];
116         // We add hits with the same tracknumber, while GEANT treats
117         // primaries succesively
118         if(curHit->GetPrimary() != primary) 
119           break;
120         if( *curHit == *newHit ) {
121             *curHit = *curHit + *newHit;
122             deja = kTRUE;
123             //            break; // 30-aug-04 by PAI 
124         } // end if
125     } // end for hitCounter
126     
127     if ( !deja ) {
128         new((*fHits)[fNhits]) AliEMCALHit(*newHit);
129         fNhits++;
130     }
131     //    printf(" fNhits %i \n", fNhits); 
132     delete newHit;
133 }
134
135 //______________________________________________________________________
136 void AliEMCALv2::StepManager(void){
137   // Accumulates hits as long as the track stays in a tower
138
139   // position wrt MRS and energy deposited
140   static Float_t        xyzte[5]={0.,0.,0.,0.,0.};// position wrt MRS, time and energy deposited
141   static Float_t        pmom[4]={0.,0.,0.,0.};
142   static TLorentzVector pos;  // Lorentz vector of the track current position.
143   static TLorentzVector mom;  // Lorentz vector of the track current momentum.
144   static Float_t ienergy = 0; // part->Energy();
145   static TString curVolName;
146   static int supModuleNumber, moduleNumber, yNumber, xNumber, absid;
147   static int keyGeom=1;  //real TRD1 geometry
148   static char *vn = "SCMX"; // Apr 13, 2006 - only TRD1 case now
149   static int nSMOP[7]={1,3,5,7,9,11}; // 30-mar-05
150   static int nSMON[7]={2,4,6,8,10,12};
151   static Float_t depositedEnergy=0.0; 
152
153   if(keyGeom == 0) {
154     keyGeom = 2;
155     if(gMC->VolId("PBMO")==0 || gMC->VolId("WSUC")==1) {
156       vn      = "SCMX";   // old TRD2(TRD1) or WSUC
157       keyGeom = 1;
158     }    
159     printf("AliEMCALv2::StepManager():  keyGeom %i : Sensetive volume %s \n", 
160     keyGeom, vn); 
161     if(gMC->VolId("WSUC")==1) printf(" WSUC - cosmic ray stand geometry \n");
162   }
163   Int_t tracknumber =  gAlice->GetMCApp()->GetCurrentTrackNumber();
164
165   curVolName = gMC->CurrentVolName();
166   if(curVolName.Contains(vn) || curVolName.Contains("SCX")) { // We are in a scintillator layer; SCX for 3X3
167     
168     if( ((depositedEnergy = gMC->Edep()) > 0.)  && (gMC->TrackTime() < fTimeCut)){// Track is inside a scintillator and deposits some energy
169       //       Info("StepManager "," entry %i DE %f",++ientry, depositedEnergy); // for testing
170        if (fCurPrimary==-1) 
171         fCurPrimary=gAlice->GetMCApp()->GetPrimary(tracknumber);
172
173       if (fCurParent==-1 || tracknumber != fCurTrack) {
174         // Check parentage
175         Int_t parent=tracknumber;
176         if (fCurParent != -1) {
177           while (parent != fCurParent && parent != -1) {
178             TParticle *part=gAlice->GetMCApp()->Particle(parent);
179             parent=part->GetFirstMother();
180           }
181         }
182         if (fCurParent==-1 || parent==-1) {
183           Int_t parent=tracknumber;
184           TParticle *part=gAlice->GetMCApp()->Particle(parent);
185           while (parent != -1 && fGeometry->IsInEMCAL(part->Vx(),part->Vy(),part->Vz())) {
186             parent=part->GetFirstMother();
187             if (parent!=-1) 
188               part=gAlice->GetMCApp()->Particle(parent);
189           } 
190           fCurParent=parent;
191           if (fCurParent==-1)
192             Error("StepManager","Cannot find parent");
193           else {
194             TParticle *part=gAlice->GetMCApp()->Particle(fCurParent);
195             ienergy = part->Energy(); 
196           }
197           while (parent != -1) {
198             part=gAlice->GetMCApp()->Particle(parent);
199             part->SetBit(kKeepBit);
200             parent=part->GetFirstMother();
201           }
202         }
203         fCurTrack=tracknumber;
204       }    
205       gMC->TrackPosition(pos);
206       xyzte[0] = pos[0];
207       xyzte[1] = pos[1];
208       xyzte[2] = pos[2];
209       xyzte[3] = gMC->TrackTime() ;       
210       
211       gMC->TrackMomentum(mom);
212       pmom[0] = mom[0];
213       pmom[1] = mom[1];
214       pmom[2] = mom[2];
215       pmom[3] = mom[3];
216       
217       //      if(ientry%200 > 0) return; // testing
218       supModuleNumber = moduleNumber = yNumber = xNumber = absid = 0;
219       if(keyGeom >= 1) { // TRD1 case now
220         gMC->CurrentVolOffID(4, supModuleNumber);
221         gMC->CurrentVolOffID(3, moduleNumber);
222         gMC->CurrentVolOffID(1, yNumber);
223         gMC->CurrentVolOffID(0, xNumber); // really x number now
224         if(strcmp(gMC->CurrentVolOffName(4),"SM10")==0) supModuleNumber += 10; // 13-oct-05
225         // Nov 10,2006
226         if(strcmp(gMC->CurrentVolOffName(0),vn) != 0) { // 3X3 case
227           if     (strcmp(gMC->CurrentVolOffName(0),"SCX1")==0) xNumber=1;
228           else if(strcmp(gMC->CurrentVolOffName(0),"SCX2")==0) xNumber=2;
229           else if(strcmp(gMC->CurrentVolOffName(0),"SCX3")==0) xNumber=3;
230           else Fatal("StepManager()", "Wrong name of sensetive volume in 3X3 case : %s ", gMC->CurrentVolOffName(0));
231         }
232       } else {
233         gMC->CurrentVolOffID(5, supModuleNumber);
234         gMC->CurrentVolOffID(4, moduleNumber);
235         gMC->CurrentVolOffID(1, yNumber);
236         gMC->CurrentVolOffID(0, xNumber);
237         if     (strcmp(gMC->CurrentVolOffName(5),"SMOP")==0) supModuleNumber = nSMOP[supModuleNumber-1];
238         else if(strcmp(gMC->CurrentVolOffName(5),"SMON")==0) supModuleNumber = nSMON[supModuleNumber-1];
239         else   assert(0); // something wrong
240       }
241       absid = fGeometry->GetAbsCellId(supModuleNumber-1, moduleNumber-1, yNumber-1, xNumber-1);
242     
243       if (absid < 0) {
244         printf(" supModuleNumber %i : moduleNumber %i : yNumber %i : xNumber %i \n",
245         supModuleNumber, moduleNumber, yNumber, xNumber); 
246         Fatal("StepManager()", "Wrong id : %i ", absid) ; 
247       }
248
249       Float_t lightYield =  depositedEnergy ;
250       // Apply Birk's law (copied from G3BIRK)
251
252       if (gMC->TrackCharge()!=0) { // Check
253           Float_t birkC1Mod = 0;
254         if (fBirkC0==1){ // Apply correction for higher charge states
255           if (TMath::Abs(gMC->TrackCharge())>=2) birkC1Mod = fBirkC1*7.2/12.6;
256           else                                   birkC1Mod = fBirkC1;
257         }
258
259         Float_t dedxcm;
260         if (gMC->TrackStep()>0)  dedxcm=1000.*gMC->Edep()/gMC->TrackStep();
261         else                     dedxcm=0;
262         lightYield=lightYield/(1.+birkC1Mod*dedxcm+fBirkC2*dedxcm*dedxcm);
263       } 
264
265       // use sampling fraction to get original energy --HG
266       //      xyzte[4] = lightYield * fGeometry->GetSampling();
267       xyzte[4] = lightYield; // 15-dec-04
268         
269       if (gDebug == -2) 
270       printf("#sm %2i #m %3i #x %1i #z %1i -> absid %i : xyzte[4] = %f\n",
271       supModuleNumber,moduleNumber,yNumber,xNumber,absid, xyzte[4]);
272
273       AddHit(fIshunt, fCurPrimary,tracknumber, fCurParent, ienergy, absid,  xyzte, pmom);
274     } // there is deposited energy
275   }
276 }
277
278 //_________________________________________________________________
279 void AliEMCALv2::FinishEvent()
280
281   // Calculate deposit energy and fill control histogram; 26-may-05
282   static double de=0.;
283   fHNhits->Fill(double(fHits->GetEntries()));
284   de = GetDepositEnergy(0);
285   if(fHDe) fHDe->Fill(de);
286 }
287
288 //_________________________________________________________________
289 Double_t AliEMCALv2::GetDepositEnergy(int print)
290
291   // 23-mar-05 - for testing
292   if(fHits == 0) return 0.;
293   AliEMCALHit  *hit=0;
294   Double_t de=0.;
295   for(int ih=0; ih<fHits->GetEntries(); ih++) {
296     hit = (AliEMCALHit*)fHits->UncheckedAt(ih);
297     de += hit->GetEnergy();
298   }
299   if(print>0) {
300     cout<<"AliEMCALv2::GetDepositEnergy() : fHits "<<fHits<<endl; 
301     printf(" #hits %i de %f \n", fHits->GetEntries(), de);
302     if(print>1) {
303       printf(" #primary particles %i\n", gAlice->GetHeader()->GetNprimary()); 
304     }
305   }
306   return de;
307 }
308
309 //___________________________________________________________
310 void AliEMCALv2::Browse(TBrowser* b)
311 {
312   TObject::Browse(b);
313 }
314
315 //___________________________________________________________
316 void AliEMCALv2::DrawCalorimeterCut(const char *name, int axis, double dcut)
317
318   // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
319   TString g(fGeometry->GetName());
320   g.ToUpper();
321   gMC->Gsatt("*", "seen", 0);
322
323   int fill = 1;
324
325   if(axis!=1) SetVolumeAttributes("STPL", 1, 1, fill);
326
327   int colo=4;
328   TString sn(name);
329   if(sn.Contains("SCM")) colo=5;
330   SetVolumeAttributes(name, 1, colo, fill);
331   if(g.Contains("110DEG") && sn=="SMOD") SetVolumeAttributes("SM10", 1, colo, fill);
332
333   TString st(GetTitle());
334   st += ", zcut, ";
335   st += name;
336
337   char *optShad = "on", *optHide="on";
338   double cxy=0.02;
339   if     (axis==1) {
340     dcut = 0.;
341     //    optHide = optShad = "off";
342   } else if(axis==2){
343     if(dcut==0.) SetVolumeAttributes("SCM0", 1, 5, fill); // yellow    
344   }
345
346   gMC->Gdopt("hide", optHide);
347   gMC->Gdopt("shad", optShad);
348
349   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
350   char cmd[200];
351   sprintf(cmd,"g3->Gdrawc(\"XEN1\", %i, %5.1f, 12., 8., %3.2f, %3.2f)\n", axis, dcut, cxy,cxy);
352   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
353
354   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
355   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
356 }
357
358 //___________________________________________________________
359 void AliEMCALv2::DrawSuperModuleCut(const char *name, int axis, double dcut, int fill)
360
361  // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
362   TString sn(GetGeometry()->GetName());
363   sn.ToUpper();
364   char *tit[3]={"xcut", "ycut", "zcut"};
365   if(axis<1) axis=1; if(axis>3) axis=3;
366
367   gMC->Gsatt("*", "seen", 0);
368   //  int fill = 6;
369
370   // SetVolumeAttributes("SMOD", 1, 1, fill);  // black
371   SetVolumeAttributes(name, 1, 5, fill);    // yellow 
372
373   double cxy=0.055, x0=10., y0=10.;
374   char *optShad = "on", *optHide="on";
375   SetVolumeAttributes("STPL", 1, 3, fill);  // green 
376   if     (axis==1) {
377     gMC->Gsatt("STPL", "seen", 0);
378     dcut = 0.;
379     optHide = "off";
380     optShad = "off";
381   } else if(axis==3) cxy = 0.1;
382   
383   gMC->Gdopt("hide", optHide);
384   gMC->Gdopt("shad", optShad);
385
386   TString st("Shish-Kebab, Compact, SMOD, ");
387   st += tit[axis-1];;
388
389   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
390   char cmd[200];
391   sprintf(cmd,"g3->Gdrawc(\"SMOD\", %i, %6.3f, %5.1f, %5.1f, %5.3f, %5.3f)\n",axis,dcut,x0,y0,cxy,cxy);
392   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
393
394   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
395   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
396   // hint for testing
397   printf("Begin of super module\n");
398   printf("g3->Gdrawc(\"SMOD\", 2,  0.300, 89., 10., 0.5, 0.5)\n");
399   printf("Center of super modules\n");
400   printf("g3->Gdrawc(\"SMOD\", 2,  0.300, 0., 10., 0.5, 0.5)\n");
401   printf("end of super modules\n");
402   printf("g3->Gdrawc(\"SMOD\", 2,  0.300, -70., 10., 0.5, 0.5)\n");
403 }
404
405 //___________________________________________________________
406 void AliEMCALv2::DrawTowerCut(const char *name, int axis, double dcut, int fill, char *optShad)
407
408   // Size of tower is 5.6x5.6x24.8 (25.0); cut on Z axiz
409   if(axis<1) axis=1; if(axis>3) axis=3;
410   TString mn(name); mn.ToUpper();
411   TString sn(GetGeometry()->GetName());
412
413   gMC->Gsatt("*", "seen", 0);
414   gMC->Gsatt("*", "fill", fill);
415
416   // int mainColo=5; // yellow
417   if(mn == "EMOD") {
418     SetVolumeAttributes(mn.Data(), 1, 1, fill);
419     SetVolumeAttributes("SCM0", 1, 5, fill); // yellow
420     SetVolumeAttributes("SCPA", 1, 3, fill); // green - 13-may-05
421   } else if(mn == "SCM0") { // first division 
422     SetVolumeAttributes(mn.Data(), 1, 1, fill);
423     SetVolumeAttributes("SCMY", 1, 5, fill); // yellow
424   } else if(mn == "SCMY") { // first division 
425     SetVolumeAttributes(mn.Data(), 1, 1, fill); 
426     SetVolumeAttributes("SCMX", 1, 5, fill); // yellow
427   } else if(mn == "SCMX" || mn.Contains("SCX")) {
428     SetVolumeAttributes(mn.Data(), 1, 5, fill);// yellow
429     SetVolumeAttributes("PBTI", 1, 4, fill);
430   } else {
431     printf("<W> for volume |%s| did not defined volume attributes\n", mn.Data());
432   }
433
434   TString st("Shashlyk, 2x2 mm sampling, 77 layers, ");
435   st += name;
436
437   gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
438   double cx=0.78, cy=2.;
439   if (mn.Contains("EMOD")) {
440     cx = cy = 0.5;
441   }
442   char cmd[200];
443
444   gMC->Gdopt("hide", optShad);
445   gMC->Gdopt("shad", optShad);
446
447   sprintf(cmd,"g3->Gdrawc(\"%s\", %i, %6.2f, 10., 10., %5.3f, %5.3f)\n",name, axis,dcut,cx,cy);
448   printf("\n%s\n",cmd); gROOT->ProcessLine(cmd);
449
450   sprintf(cmd,"gMC->Gdhead(1111, \"%s\")\n", st.Data());
451   printf("%s\n",cmd); gROOT->ProcessLine(cmd);
452 }
453
454 //___________________________________________________________  
455 void AliEMCALv2::DrawAlicWithHits(int mode)
456
457  // 20-sep-04; does not work now
458   static TH2F *h2;
459   if(h2==0) h2 = new TH2F("h2","test fo hits", 60,0.5,60.5, 28,0.5,28.5);
460   else      h2->Reset();
461   if(mode==0) {
462     gROOT->ProcessLine("TGeant3 *g3 = (TGeant3*)gMC");
463     gMC->Gsatt("*","seen",0);
464     gMC->Gsatt("scm0","seen",5);
465
466     gROOT->ProcessLine("g3->Gdrawc(\"ALIC\", 1,   2.0, 12., -130, 0.3, 0.3)");
467     // g3->Gdrawc("ALIC", 1,   2.0, 10., -14, 0.05, 0.05)
468   }
469
470   TClonesArray *hits = Hits();
471   Int_t nhits = hits->GetEntries(), absId, nSupMod, nModule, nIphi, nIeta, iphi, ieta;
472   AliEMCALHit *hit = 0;
473   Double_t de, des=0.;
474   for(Int_t i=0; i<nhits; i++) {
475     hit   = (AliEMCALHit*)hits->UncheckedAt(i);
476     absId = hit->GetId();
477     de    = hit->GetEnergy();
478     des += de;
479     if(fGeometry->GetCellIndex(absId, nSupMod, nModule, nIphi, nIeta)){
480       //      printf(" de %f abs id %i smod %i tower %i | cell iphi %i : ieta %i\n",
481       // de, absId, nSupMod, nModule, nIphi, nIeta); 
482       if(nSupMod==3) {
483         fGeometry->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
484         // printf(" iphi %i : ieta %i\n", iphi,ieta);
485         h2->Fill(double(ieta),double(iphi), de);
486       }
487     }
488     //    hit->Dump();
489   }
490   printf(" #hits %i : sum de %f -> %f \n", nhits, des, h2->Integral());
491   if(mode>0 && h2->Integral()>0.) h2->Draw();
492 }
493
494 //___________________________________________________________
495 void AliEMCALv2::SetVolumeAttributes(const char *name, int seen, int color, int fill)
496 {
497  /* seen=-2:volume is visible but none of its descendants;
498          -1:volume is not visible together with all its descendants;
499           0:volume is not visible;
500           1:volume is     visible. */
501   gMC->Gsatt(name, "seen", seen);
502   gMC->Gsatt(name, "colo", color);
503   gMC->Gsatt(name, "fill", fill); 
504   printf(" %s : seen %i color %i fill %i \n", name, seen, color, fill);
505
506
507 //___________________________________________________________
508 void AliEMCALv2::TestIndexTransition(int pri, int idmax)
509
510  // Test for EMCAL SHISHKEBAB geometry
511
512   Int_t nSupMod, nModule, nIphi, nIeta, idNew, nGood=0;
513   if(idmax==0) idmax = fGeometry->GetNCells();
514   for(Int_t id=1; id<=idmax; id++) {
515     if(!fGeometry->GetCellIndex(id, nSupMod, nModule, nIphi, nIeta)){
516       printf(" Wrong abs ID %i : #cells %i\n", id, fGeometry->GetNCells());
517       break;  
518     }
519     idNew = fGeometry->GetAbsCellId(nSupMod, nModule, nIphi, nIeta);
520     if(id != idNew || pri>0) {
521       printf(" ID %i : %i <- new id\n", id, idNew);
522       printf(" nSupMod   %i ",  nSupMod);
523       printf(" nModule    %i ", nModule);
524       printf(" nIphi     %i ", nIphi);
525       printf(" nIeta     %i \n", nIeta);
526      
527     } else nGood++;
528   }
529   printf(" Good decoding %i : %i <- #cells \n", nGood, fGeometry->GetNCells());
530 }