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