]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliMCInfo.cxx
coverity fix
[u/mrichter/AliRoot.git] / PWGPP / TPC / 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   for (Int_t i = 0; i<32; i++) { fDig[i] = 0; }
252   Reset();
253 }
254 ////////////////////////////////////////////////////////////////////////
255 AliTPCdigitRow & AliTPCdigitRow::operator=(const AliTPCdigitRow &digOld)
256 {
257   for (Int_t i = 0; i<32; i++) fDig[i] = digOld.fDig[i];
258   return (*this);
259 }
260 ////////////////////////////////////////////////////////////////////////
261 void AliTPCdigitRow::SetRow(Int_t row) 
262 {
263   //
264   // set bit mask for given row
265   //
266   if (row >= 8*32) {
267     //    cerr<<"AliTPCdigitRow::SetRow: index "<<row<<" out of bounds."<<endl;
268     return;
269   }
270   Int_t iC = row/8;
271   Int_t iB = row%8;
272   SETBIT(fDig[iC],iB);
273 }
274
275 ////////////////////////////////////////////////////////////////////////
276 Bool_t AliTPCdigitRow::TestRow(Int_t row) const
277 {
278 //
279 // return kTRUE if row is on
280 //
281   Int_t iC = row/8;
282   Int_t iB = row%8;
283   return TESTBIT(fDig[iC],iB);
284 }
285 ////////////////////////////////////////////////////////////////////////
286 Int_t AliTPCdigitRow::RowsOn(Int_t upto) const
287 {
288 //
289 // returns number of rows with a digit  
290 // count only rows less equal row number upto
291 //
292   Int_t total = 0;
293   for (Int_t i = 0; i<32; i++) {
294     for (Int_t j = 0; j < 8; j++) {
295       if (i*8+j > upto) return total;
296       if (TESTBIT(fDig[i],j))  total++;
297     }
298   }
299   return total;
300 }
301 ////////////////////////////////////////////////////////////////////////
302 void AliTPCdigitRow::Reset()
303 {
304 //
305 // resets all rows to zero
306 //
307   for (Int_t i = 0; i<32; i++) {
308     //fDig[i] <<= 8;
309     fDig[i] = 0;
310   }
311 }
312 ////////////////////////////////////////////////////////////////////////
313 Int_t AliTPCdigitRow::Last() const
314 {
315 //
316 // returns the last row number with a digit
317 // returns -1 if now digits 
318 //
319   for (Int_t i = 32-1; i>=0; i--) {
320     for (Int_t j = 7; j >= 0; j--) {
321       if TESTBIT(fDig[i],j) return i*8+j;
322     }
323   }
324   return -1;
325 }
326 ////////////////////////////////////////////////////////////////////////
327 Int_t AliTPCdigitRow::First() const
328 {
329 //
330 // returns the first row number with a digit
331 // returns -1 if now digits 
332 //
333   for (Int_t i = 0; i<32; i++) {
334     for (Int_t j = 0; j < 8; j++) {
335       if (TESTBIT(fDig[i],j)) return i*8+j;
336     }
337   }
338   return -1;
339 }
340
341
342 void AliMCInfo::Clear(Option_t* ){
343   //
344   // clear the content
345   //
346   if (fTPCReferences) fTPCReferences->Clear();
347   if (fITSReferences) fITSReferences->Clear();
348   if (fTRDReferences) fTRDReferences->Clear();
349   if (fTOFReferences) fTOFReferences->Clear();
350   for (Int_t i=0;i<4;i++) fVDist[i]=0;
351   for (Int_t i=0;i<3;i++) fDecayCoord[i]=0;
352   fTRdecay.SetTrack(-1);
353   fCharge = 0;
354   fEventNr=-1;
355   fPrimPart=-1;               // index of primary particle in TreeH
356   fMass=0;                   // mass of the particle
357   fCharge=0;                 // charge of the particle
358   fLabel=0;                  // track label
359   fEventNr=0;                // event number
360   fMCtracks=0;               // indication of how many times the track is retuturned back
361   fPdg=0;                        //pdg code
362   fTPCdecay=0;                  //indicates decay in TPC
363   fRowsWithDigitsInn=0;          // number of rows with digits in the inner sectors
364   fRowsWithDigits=0;             // number of rows with digits in the outer sectors
365   fRowsTrackLength=0;            // last - first row with digit
366   fTPCtrackLength=0;           // distance between first and last track reference
367   fPrim=0;                     // theoretical dedx in tpc according particle momenta and mass
368   fNTPCRef=0;                    // tpc references counter
369   fNITSRef=0;                    // ITS references counter
370   fNTRDRef=0;                    // TRD references counter
371   fNTOFRef=0;                    // TOF references counter
372   //
373   fNTPCRefOut=0;                    // tpc references counter - out
374   fNITSRefOut=0;                    // ITS references counter - out
375   fNTRDRefOut=0;                    // TRD references counter - out
376   fNTOFRefOut=0;                    // TOF references counter - out
377 }
378
379
380
381 void AliMCInfo::CalcTPCrows(TClonesArray * runArrayTR){
382   //
383   // Calculates the numebr of the track references for detectors
384   // In case of changing direction - curling tracks - the counter is not increasing
385   //
386   // Rough calculation 
387   // of the first and last point in the TPC  
388   //
389   fNTPCRef = 0;
390   fNITSRef = 0;
391   fNTRDRef = 0;
392   fNTOFRef = 0;  
393   fNTPCRefOut = 0;
394   fNITSRefOut = 0;
395   fNTRDRefOut = 0;
396   fNTOFRefOut = 0;  
397   //
398   //
399   Float_t tpcminRadius=TMath::Max(90., fParticle.R());
400   Float_t tpcmaxRadius=TMath::Min(250.,fParticle.R());  
401   Int_t nover=-1;
402
403   if (runArrayTR->GetEntriesFast()>0){  
404     nover=0;
405     AliTrackReference *ref0 = (AliTrackReference*)runArrayTR->At(0);    
406     Float_t dir = ((ref0->X()*ref0->Px()+ref0->Y()*ref0->Py())>0)? 1.:-1.; //inside or outside
407     //
408     //    
409     for (Int_t iTrackRef = 0; iTrackRef < runArrayTR->GetEntriesFast(); iTrackRef++) {
410       //
411       AliTrackReference *ref = (AliTrackReference*)runArrayTR->At(iTrackRef);  
412       Float_t newdirection = ((ref->X()*ref->Px()+ref->Y()*ref->Py())>0)? 1.:-1.; //inside or outside
413       //
414       //
415       if (dir*newdirection<0.0) {
416         nover++;
417         dir = newdirection;
418       }
419       //
420       if (ref->DetectorId()== AliTrackReference::kTRD){
421         if (nover==0) {
422           tpcmaxRadius =250;
423           fNTRDRefOut++;
424         }
425         fNTRDRef++;
426       }
427       if (ref->DetectorId()== AliTrackReference::kITS){ 
428         if (nover==0) {
429           fNITSRefOut++;
430           tpcminRadius =90;
431         }
432         fNITSRef++;
433       }
434       if (ref->DetectorId()== AliTrackReference::kTOF){
435         fNTOFRef++;
436         if (nover==0) {
437           fNTOFRefOut++;
438           tpcmaxRadius =250;
439         }
440       }
441       //
442       if (ref->DetectorId()== AliTrackReference::kTPC){
443         fNTPCRef++;
444         if (nover==0) {
445           fNTPCRefOut++;
446           fTrackRefOut = (*ref);
447           if (ref->R()>tpcmaxRadius) tpcmaxRadius = ref->R();
448           if (ref->R()<tpcminRadius) tpcminRadius = ref->R();
449         }
450       }
451       if (ref->DetectorId()== AliTrackReference::kDisappeared){
452         if (nover==0 &&TMath::Abs(ref->Z())<250){
453           tpcmaxRadius = TMath::Min(Float_t(250.),ref->R());
454         }
455       }
456     }
457   }  // if we have track refs
458   fTPCtrackLength = tpcmaxRadius-tpcminRadius;
459   fMCtracks=nover+1;
460 }
461
462 //
463 //
464 //
465
466
467
468 void AliMCInfo::Update(TParticle * part, TClonesArray * runArrayTR, Double_t pvertex[4], Int_t label){
469   //
470   //
471   //
472   Clear();
473   fLabel=label;
474   fParticle=(*part);
475   AliMCInfo *info = this;
476   fVDist[0]=part->Vx()-pvertex[0];
477   fVDist[1]=part->Vy()-pvertex[1];
478   fVDist[2]=part->Vz()-pvertex[2];
479   fVDist[3]=part->R()-pvertex[3];
480   //
481   //
482   //
483   for (Int_t iTrackRef = 0; iTrackRef < runArrayTR->GetEntriesFast(); iTrackRef++) {
484     AliTrackReference *trackRef = (AliTrackReference*)runArrayTR->At(iTrackRef);         
485     if (trackRef->DetectorId()== AliTrackReference::kTPC){     
486       TClonesArray & arr = *(info->fTPCReferences);
487       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
488     }
489     if (trackRef->DetectorId()== AliTrackReference::kITS){
490       TClonesArray & arr = *(info->fITSReferences);
491       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
492     }
493     if (trackRef->DetectorId()== AliTrackReference::kTRD){
494       TClonesArray & arr = *(info->fTRDReferences);
495       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
496     }
497     if (trackRef->DetectorId()== AliTrackReference::kTOF){
498       TClonesArray & arr = *(info->fTOFReferences);
499       new (arr[arr.GetEntriesFast()]) AliTrackReference(*trackRef);     
500     }
501     //
502     // decay case
503     //
504     if (trackRef->DetectorId()== AliTrackReference::kDisappeared){
505       info->fTRdecay = *trackRef;       
506     }
507   }
508   //
509   Update();
510   CalcTPCrows(runArrayTR);
511 }
512