]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliMCInfo.cxx
Corrected PMT gain settings
[u/mrichter/AliRoot.git] / PWG1 / AliMCInfo.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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
17 ///////////////////////////////////////////////////////////////////////////
18 /*
19
20 Origin: marian.ivanov@cern.ch
21 Container classes with MC infomation
22
23 The AliMCInfo contains the information about the particles properties 
24 during transportation throuch ALICE Detector
25
26 The base Information :
27 TParticle - fParticle             -  properties of the particle at creation point
28 AliTrackReference - fXXXRefernces -  TClonesArray of refernces in differnt detectors
29 fNXXXRef                          -  number of the track refernces in differnt detectors
30 AliTPCdigitRow    - fTPCRow       -  the map of the hitted rows - (will be repalced by TBits)
31 fRowsWith*                        - number of rows hitted by particle
32 fMCtracks         -               - number of turn over of the track inside of the TPC
33
34 ++++
35 some additional information usable for tree draw - TO SPEED UP tree queries
36 IMPORTANT FOR PROOF FAST PROTOTYPING ANALYSIS 
37                                       
38
39
40
41 */
42
43 #if !defined(__CINT__) || defined(__MAKECINT__)
44 #include <stdio.h>
45 //ROOT includes
46 #include "Rtypes.h"
47 #include "TClonesArray.h"
48 #include "TDatabasePDG.h"
49 //ALIROOT includes
50 #include "AliTrackReference.h"
51 #include "AliMCInfo.h" 
52 #include "AliMathBase.h" 
53 #endif
54
55 //
56 // 
57
58 ClassImp(AliTPCdigitRow)
59 ClassImp(AliMCInfo)
60
61
62
63 ////////////////////////////////////////////////////////////////////////
64 AliMCInfo::AliMCInfo():
65   fTrackRef(),
66   fTrackRefOut(),
67   fTRdecay(),
68   fPrimPart(0),
69   fParticle(),
70   fMass(0),
71   fCharge(0),
72   fLabel(0),
73   fEventNr(),
74   fMCtracks(),
75   fPdg(0),
76   fTPCdecay(0),
77   fRowsWithDigitsInn(0),
78   fRowsWithDigits(0),
79   fRowsTrackLength(0),
80   fTPCtrackLength(-1),
81   fPrim(0),
82   fTPCRow(), 
83   fNTPCRef(0),                    // tpc references counter
84   fNITSRef(0),                    // ITS references counter
85   fNTRDRef(0),                    // TRD references counter
86   fNTOFRef(0),                    // TOF references counter
87   fNTPCRefOut(0),                    // tpc references counter
88   fNITSRefOut(0),                    // ITS references counter
89   fNTRDRefOut(0),                    // TRD references counter
90   fNTOFRefOut(0),                    // TOF references counter
91   fTPCReferences(0),
92   fITSReferences(0),
93   fTRDReferences(0),
94   fTOFReferences(0)
95 {
96   //
97   // Default constructor
98   //
99   fTPCReferences  = new TClonesArray("AliTrackReference",10);
100   fITSReferences  = new TClonesArray("AliTrackReference",10);
101   fTRDReferences  = new TClonesArray("AliTrackReference",10);
102   fTOFReferences  = new TClonesArray("AliTrackReference",10);
103   fTRdecay.SetTrack(-1);
104   for (Int_t i=0;i<4;i++) fVDist[i]=0;
105   for (Int_t i=0;i<3;i++) fDecayCoord[i]=0;
106   fCharge = 0;
107 }
108
109 AliMCInfo::AliMCInfo(const AliMCInfo& info):
110   TObject(),
111   fTrackRef(info.fTrackRef),
112   fTrackRefOut(info.fTrackRefOut),
113   fTRdecay(info.GetTRdecay()),
114   fPrimPart(info.fPrimPart),
115   fParticle(info.fParticle),
116   fMass(info.GetMass()),
117   fCharge(info.fCharge),
118   fLabel(info.fLabel),
119   fEventNr(info.fEventNr),
120   fMCtracks(info.fMCtracks),
121   fPdg(info.fPdg),
122   fTPCdecay(info.fTPCdecay),
123   fRowsWithDigitsInn(info.fRowsWithDigitsInn),
124   fRowsWithDigits(info.fRowsWithDigits),
125   fRowsTrackLength(info.fRowsTrackLength),
126   fTPCtrackLength(info.fTPCtrackLength),
127   fPrim(info.fPrim),
128   fTPCRow(info.fTPCRow), 
129   fNTPCRef(info.fNTPCRef),                    // tpc references counter
130   fNITSRef(info.fNITSRef),                    // ITS references counter
131   fNTRDRef(info.fNTRDRef),                    // TRD references counter
132   fNTOFRef(info.fNTOFRef),                    // TOF references counter
133   fNTPCRefOut(info.fNTPCRefOut),                    // tpc references counter
134   fNITSRefOut(info.fNITSRefOut),                    // ITS references counter
135   fNTRDRefOut(info.fNTRDRefOut),                    // TRD references counter
136   fNTOFRefOut(info.fNTOFRefOut),                    // TOF references counter
137   fTPCReferences(0),
138   fITSReferences(0),
139   fTRDReferences(0),
140   fTOFReferences(0)
141 {
142   //
143   // copy constructor
144   //
145   fTPCReferences = (TClonesArray*)info.fTPCReferences->Clone();
146   fITSReferences = (TClonesArray*)info.fITSReferences->Clone();
147   fTRDReferences = (TClonesArray*)info.fTRDReferences->Clone();
148   fTOFReferences = (TClonesArray*)info.fTOFReferences->Clone();
149   for (Int_t i=0;i<4;i++) fVDist[i]=info.fVDist[i];
150   for (Int_t i=0;i<3;i++) fDecayCoord[i]=info.fDecayCoord[i];
151
152 }
153
154
155 AliMCInfo& AliMCInfo::operator=(const AliMCInfo& info) { 
156   //
157   // Assignment operator
158   //
159   this->~AliMCInfo();
160   new (this) AliMCInfo(info);
161   return *this;
162 }
163
164
165 AliMCInfo::~AliMCInfo()
166 {
167   //
168   // Destructor of the class
169   //
170   if (fTPCReferences) {
171     delete fTPCReferences;
172   }
173   if (fITSReferences){
174     delete fITSReferences;
175   }
176   if (fTRDReferences){    
177     delete fTRDReferences;  
178   }
179   if (fTOFReferences){    
180     delete fTOFReferences;  
181   }
182
183 }
184
185
186
187
188 void AliMCInfo::Update()
189 {
190   //
191   // Update MC info
192   // Calculates some derived variables
193   //
194   //
195   fMCtracks =1;
196   //
197   // decay info
198   fTPCdecay=kFALSE;
199   if (fTRdecay.GetTrack()>0){
200     fDecayCoord[0] = fTRdecay.X();
201     fDecayCoord[1] = fTRdecay.Y();
202     fDecayCoord[2] = fTRdecay.Z();
203     if ( (fTRdecay.R()<250)&&(fTRdecay.R()>85) && (TMath::Abs(fTRdecay.Z())<250) ){
204       fTPCdecay=kTRUE;     
205     }
206     else{
207       fDecayCoord[0] = 0;
208       fDecayCoord[1] = 0;
209       fDecayCoord[2] = 0;
210     }
211   }
212   //
213   //
214   //digits information update
215   fRowsWithDigits    = fTPCRow.RowsOn();    
216   fRowsWithDigitsInn = fTPCRow.RowsOn(63); // 63 = number of inner rows
217   fRowsTrackLength   = fTPCRow.Last() - fTPCRow.First();
218   //
219   //
220   // calculate primary ionization per cm  
221   if (fParticle.GetPDG()){
222     fMass = fParticle.GetMass();  
223     fCharge = fParticle.GetPDG()->Charge();
224     if (fTPCReferences->GetEntriesFast()>0){
225       fTrackRef = *((AliTrackReference*)fTPCReferences->At(0));
226     }
227     if (fMass>0){
228       Float_t p = TMath::Sqrt(fTrackRef.Px()*fTrackRef.Px()+
229                               fTrackRef.Py()*fTrackRef.Py()+
230                               fTrackRef.Pz()*fTrackRef.Pz());    
231       if (p>0.001){
232         Float_t betagama = p /fMass;
233         fPrim = AliMathBase::BetheBlochAleph(betagama);
234       }else fPrim=0;
235     }
236   }else{
237     fMass =0;
238     fPrim =0;
239   }  
240   
241
242
243 }
244
245
246
247
248 ////////////////////////////////////////////////////////////////////////
249 AliTPCdigitRow::AliTPCdigitRow()
250 {
251   Reset();
252 }
253 ////////////////////////////////////////////////////////////////////////
254 AliTPCdigitRow & AliTPCdigitRow::operator=(const AliTPCdigitRow &digOld)
255 {
256   for (Int_t i = 0; i<32; i++) fDig[i] = digOld.fDig[i];
257   return (*this);
258 }
259 ////////////////////////////////////////////////////////////////////////
260 void AliTPCdigitRow::SetRow(Int_t row) 
261 {
262   //
263   // set bit mask for given row
264   //
265   if (row >= 8*32) {
266     //    cerr<<"AliTPCdigitRow::SetRow: index "<<row<<" out of bounds."<<endl;
267     return;
268   }
269   Int_t iC = row/8;
270   Int_t iB = row%8;
271   SETBIT(fDig[iC],iB);
272 }
273
274 ////////////////////////////////////////////////////////////////////////
275 Bool_t AliTPCdigitRow::TestRow(Int_t row) const
276 {
277 //
278 // return kTRUE if row is on
279 //
280   Int_t iC = row/8;
281   Int_t iB = row%8;
282   return TESTBIT(fDig[iC],iB);
283 }
284 ////////////////////////////////////////////////////////////////////////
285 Int_t AliTPCdigitRow::RowsOn(Int_t upto) const
286 {
287 //
288 // returns number of rows with a digit  
289 // count only rows less equal row number upto
290 //
291   Int_t total = 0;
292   for (Int_t i = 0; i<32; i++) {
293     for (Int_t j = 0; j < 8; j++) {
294       if (i*8+j > upto) return total;
295       if (TESTBIT(fDig[i],j))  total++;
296     }
297   }
298   return total;
299 }
300 ////////////////////////////////////////////////////////////////////////
301 void AliTPCdigitRow::Reset()
302 {
303 //
304 // resets all rows to zero
305 //
306   for (Int_t i = 0; i<32; i++) {
307     fDig[i] <<= 8;
308   }
309 }
310 ////////////////////////////////////////////////////////////////////////
311 Int_t AliTPCdigitRow::Last() const
312 {
313 //
314 // returns the last row number with a digit
315 // returns -1 if now digits 
316 //
317   for (Int_t i = 32-1; i>=0; i--) {
318     for (Int_t j = 7; j >= 0; j--) {
319       if TESTBIT(fDig[i],j) return i*8+j;
320     }
321   }
322   return -1;
323 }
324 ////////////////////////////////////////////////////////////////////////
325 Int_t AliTPCdigitRow::First() const
326 {
327 //
328 // returns the first row number with a digit
329 // returns -1 if now digits 
330 //
331   for (Int_t i = 0; i<32; i++) {
332     for (Int_t j = 0; j < 8; j++) {
333       if (TESTBIT(fDig[i],j)) return i*8+j;
334     }
335   }
336   return -1;
337 }
338
339
340 void AliMCInfo::Clear(Option_t* ){
341   //
342   // clear the content
343   //
344   if (fTPCReferences) fTPCReferences->Clear();
345   if (fITSReferences) fITSReferences->Clear();
346   if (fTRDReferences) fTRDReferences->Clear();
347   if (fTOFReferences) fTOFReferences->Clear();
348   for (Int_t i=0;i<4;i++) fVDist[i]=0;
349   for (Int_t i=0;i<3;i++) fDecayCoord[i]=0;
350   fTRdecay.SetTrack(-1);
351   fCharge = 0;
352   fEventNr=-1;
353   fPrimPart=-1;               // index of primary particle in TreeH
354   fMass=0;                   // mass of the particle
355   fCharge=0;                 // charge of the particle
356   fLabel=0;                  // track label
357   fEventNr=0;                // event number
358   fMCtracks=0;               // indication of how many times the track is retuturned back
359   fPdg=0;                        //pdg code
360   fTPCdecay=0;                  //indicates decay in TPC
361   fRowsWithDigitsInn=0;          // number of rows with digits in the inner sectors
362   fRowsWithDigits=0;             // number of rows with digits in the outer sectors
363   fRowsTrackLength=0;            // last - first row with digit
364   fTPCtrackLength=0;           // distance between first and last track reference
365   fPrim=0;                     // theoretical dedx in tpc according particle momenta and mass
366   fNTPCRef=0;                    // tpc references counter
367   fNITSRef=0;                    // ITS references counter
368   fNTRDRef=0;                    // TRD references counter
369   fNTOFRef=0;                    // TOF references counter
370   //
371   fNTPCRefOut=0;                    // tpc references counter - out
372   fNITSRefOut=0;                    // ITS references counter - out
373   fNTRDRefOut=0;                    // TRD references counter - out
374   fNTOFRefOut=0;                    // TOF references counter - out
375 }
376
377
378
379 void AliMCInfo::CalcTPCrows(TClonesArray * runArrayTR){
380   //
381   // Calculates the numebr of the track references for detectors
382   // In case of changing direction - curling tracks - the counter is not increasing
383   //
384   // Rough calculation 
385   // of the first and last point in the TPC  
386   //
387   fNTPCRef = 0;
388   fNITSRef = 0;
389   fNTRDRef = 0;
390   fNTOFRef = 0;  
391   fNTPCRefOut = 0;
392   fNITSRefOut = 0;
393   fNTRDRefOut = 0;
394   fNTOFRefOut = 0;  
395   //
396   //
397   Float_t tpcminRadius=TMath::Max(90., fParticle.R());
398   Float_t tpcmaxRadius=TMath::Min(250.,fParticle.R());  
399   Int_t nover=-1;
400
401   if (runArrayTR->GetEntriesFast()>0){  
402     nover=0;
403     AliTrackReference *ref0 = (AliTrackReference*)runArrayTR->At(0);    
404     Float_t dir = ((ref0->X()*ref0->Px()+ref0->Y()*ref0->Py())>0)? 1.:-1.; //inside or outside
405     //
406     //    
407     for (Int_t iTrackRef = 0; iTrackRef < runArrayTR->GetEntriesFast(); iTrackRef++) {
408       //
409       AliTrackReference *ref = (AliTrackReference*)runArrayTR->At(iTrackRef);  
410       Float_t newdirection = ((ref->X()*ref->Px()+ref->Y()*ref->Py())>0)? 1.:-1.; //inside or outside
411       //
412       //
413       if (dir*newdirection<0.0) {
414         nover++;
415         dir = newdirection;
416       }
417       //
418       if (ref->DetectorId()== AliTrackReference::kTRD){
419         if (nover==0) {
420           tpcmaxRadius =250;
421           fNTRDRefOut++;
422         }
423         fNTRDRef++;
424       }
425       if (ref->DetectorId()== AliTrackReference::kITS){ 
426         if (nover==0) {
427           fNITSRefOut++;
428           tpcminRadius =90;
429         }
430         fNITSRef++;
431       }
432       if (ref->DetectorId()== AliTrackReference::kTOF){
433         fNTOFRef++;
434         if (nover==0) {
435           fNTOFRefOut++;
436           tpcmaxRadius =250;
437         }
438       }
439       //
440       if (ref->DetectorId()== AliTrackReference::kTPC){
441         fNTPCRef++;
442         if (nover==0) {
443           fNTPCRefOut++;
444           fTrackRefOut = (*ref);
445           if (ref->R()>tpcmaxRadius) tpcmaxRadius = ref->R();
446           if (ref->R()<tpcminRadius) tpcminRadius = ref->R();
447         }
448       }
449       if (ref->DetectorId()== AliTrackReference::kDisappeared){
450         if (nover==0 &&TMath::Abs(ref->Z())<250){
451           tpcmaxRadius = TMath::Min(Float_t(250.),ref->R());
452         }
453       }
454     }
455   }  // if we have track refs
456   fTPCtrackLength = tpcmaxRadius-tpcminRadius;
457   fMCtracks=nover+1;
458 }
459
460 //
461 //
462 //
463
464
465
466 void AliMCInfo::Update(TParticle * part, TClonesArray * runArrayTR, Double_t pvertex[4], Int_t label){
467   //
468   //
469   //
470   Clear();
471   fLabel=label;
472   fParticle=(*part);
473   AliMCInfo *info = this;
474   fVDist[0]=part->Vx()-pvertex[0];
475   fVDist[1]=part->Vy()-pvertex[1];
476   fVDist[2]=part->Vz()-pvertex[2];
477   fVDist[3]=part->R()-pvertex[3];
478   //
479   //
480   //
481   for (Int_t iTrackRef = 0; iTrackRef < runArrayTR->GetEntriesFast(); iTrackRef++) {
482     AliTrackReference *trackRef = (AliTrackReference*)runArrayTR->At(iTrackRef);         
483     if (trackRef->DetectorId()== AliTrackReference::kTPC){     
484       TClonesArray & arr = *(info->fTPCReferences);
485       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
486     }
487     if (trackRef->DetectorId()== AliTrackReference::kITS){
488       TClonesArray & arr = *(info->fITSReferences);
489       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
490     }
491     if (trackRef->DetectorId()== AliTrackReference::kTRD){
492       TClonesArray & arr = *(info->fTRDReferences);
493       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
494     }
495     if (trackRef->DetectorId()== AliTrackReference::kTOF){
496       TClonesArray & arr = *(info->fTOFReferences);
497       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
498     }
499     //
500     // decay case
501     //
502     if (trackRef->DetectorId()== AliTrackReference::kDisappeared){
503       info->fTRdecay = *trackRef;       
504     }
505   }
506   //
507   Update();
508   CalcTPCrows(runArrayTR);
509 }
510