]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
Coverity fixes.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructor.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 //--
20 //-- Yves Schutz (SUBATECH) 
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and 
22 // derived from STEER/AliReconstructor. 
23 // 
24 // --- ROOT system ---
25 #include "TGeoManager.h"
26 #include "TGeoMatrix.h"
27
28 // --- Standard library ---
29
30 // --- AliRoot header files ---
31 #include "AliLog.h"
32 #include "AliAltroMapping.h"
33 #include "AliESDEvent.h"
34 #include "AliESDCaloCluster.h"
35 #include "AliESDCaloCells.h"
36 #include "AliPHOSReconstructor.h"
37 #include "AliPHOSClusterizerv1.h"
38 #include "AliPHOSTrackSegmentMakerv1.h"
39 #include "AliPHOSPIDv1.h"
40 #include "AliPHOSTracker.h"
41 #include "AliRawReader.h"
42 #include "AliPHOSCalibData.h"
43 #include "AliCDBEntry.h"
44 #include "AliCDBManager.h"
45 #include "AliPHOSTrigger.h"
46 #include "AliPHOSGeometry.h"
47 #include "AliPHOSDigit.h"
48 #include "AliPHOSTrackSegment.h"
49 #include "AliPHOSEmcRecPoint.h"
50 #include "AliPHOSRecParticle.h"
51 #include "AliPHOSRawFitterv0.h"
52 #include "AliPHOSRawFitterv1.h"
53 #include "AliPHOSRawFitterv2.h"
54 #include "AliPHOSRawFitterv3.h"
55 #include "AliPHOSRawDigiProducer.h"
56 #include "AliPHOSPulseGenerator.h"
57
58 ClassImp(AliPHOSReconstructor)
59
60 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ; 
61 TClonesArray*     AliPHOSReconstructor::fgDigitsArray = 0;   // Array of PHOS digits
62 TObjArray*        AliPHOSReconstructor::fgEMCRecPoints = 0;   // Array of EMC rec.points
63 AliPHOSCalibData * AliPHOSReconstructor::fgCalibData  = 0 ;
64
65
66 //____________________________________________________________________________
67 AliPHOSReconstructor::AliPHOSReconstructor() :
68   fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL),fTmpDigLG(NULL)
69 {
70   // ctor
71   fGeom          = AliPHOSGeometry::GetInstance("IHEP","");
72   fClusterizer   = new AliPHOSClusterizerv1      (fGeom);
73   fTSM           = new AliPHOSTrackSegmentMakerv1(fGeom);
74   fPID           = new AliPHOSPIDv1              (fGeom);
75   fTmpDigLG      = new TClonesArray("AliPHOSDigit",100);
76   fgDigitsArray  = new TClonesArray("AliPHOSDigit",100);
77   fgEMCRecPoints = new TObjArray(100) ;
78   if (!fgCalibData)
79     fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
80
81   AliInfo(Form("PHOS bad channel map contains %d bad channel(s).\n",
82                fgCalibData->GetNumOfEmcBadChannels()));
83  
84 }
85
86 //____________________________________________________________________________
87 AliPHOSReconstructor::~AliPHOSReconstructor()
88 {
89   // dtor
90   delete fGeom;
91   delete fClusterizer;
92   delete fTSM;
93   delete fPID;
94   delete fTmpDigLG;
95   delete fgDigitsArray;
96   delete fgEMCRecPoints;
97
98
99 //____________________________________________________________________________
100 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
101 {
102   // 'single-event' local reco method called by AliReconstruction; 
103   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
104   // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by 
105   // the global tracking.
106
107   fClusterizer->InitParameters();
108   fClusterizer->SetInput(digitsTree);
109   fClusterizer->SetOutput(clustersTree);
110   if ( Debug() ) 
111     fClusterizer->Digits2Clusters("deb all") ; 
112   else 
113     fClusterizer->Digits2Clusters("") ;
114 }
115
116 //____________________________________________________________________________
117 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
118                                    AliESDEvent* esd) const
119 {
120   // This method produces PHOS rec-particles,
121   // then it creates AliESDtracks out of them and
122   // write tracks to the ESD
123
124
125   // do current event; the loop over events is done by AliReconstruction::Run()
126   fTSM->SetESD(esd) ; 
127   fTSM->SetInput(clustersTree);
128   if ( Debug() ) 
129     fTSM->Clusters2TrackSegments("deb all") ;
130   else 
131     fTSM->Clusters2TrackSegments("") ;
132   
133   fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ; 
134   fPID->SetESD(esd) ; 
135   if ( Debug() ) 
136     fPID->TrackSegments2RecParticles("deb all") ;
137   else 
138     fPID->TrackSegments2RecParticles("") ;
139
140   TClonesArray *recParticles  = fPID->GetRecParticles();
141   Int_t nOfRecParticles = recParticles->GetEntriesFast();
142   
143   AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
144   
145   // Read digits array
146
147   TBranch *branch = digitsTree->GetBranch("PHOS");
148   if (!branch) { 
149     AliError("can't get the branch with the PHOS digits !");
150     return;
151   }
152   branch->SetAddress(&fgDigitsArray);
153   branch->GetEntry(0);
154
155   // Get the clusters array
156
157   TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
158   if (!emcbranch) { 
159     AliError("can't get the branch with the PHOS EMC clusters !");
160     return;
161   }
162
163   emcbranch->SetAddress(&fgEMCRecPoints);
164   emcbranch->GetEntry(0);
165
166 //   //#########Calculate trigger and set trigger info###########
167
168 //   AliPHOSTrigger tr ;
169 //   //   tr.SetPatchSize(1);//create 4x4 patches
170 //   tr.SetSimulation(kFALSE);
171 //   tr.Trigger(fgDigitsArray);
172   
173 //   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
174 //   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
175 //   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
176 //   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
177
178 //   Int_t iSM2x2      = tr.Get2x2SuperModule();
179 //   Int_t iSMnxn      = tr.GetnxnSuperModule();
180 //   Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
181 //   Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
182 //   Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
183 //   Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
184
185 //   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  
186 //                 maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
187 //   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
188 //                 maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
189
190 //   // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
191 //   Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
192 //   Int_t iAbsId2x2 =-1;
193 //   Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
194 //   Int_t iAbsIdnxn =-1;
195 //   TVector3    pos2x2(-1,-1,-1);
196 //   TVector3    posnxn(-1,-1,-1);
197 //   fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
198 //   fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
199 //   fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
200 //   fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
201
202 //   TArrayF triggerPosition(6);
203 //   triggerPosition[0] = pos2x2(0) ;   
204 //   triggerPosition[1] = pos2x2(1) ;   
205 //   triggerPosition[2] = pos2x2(2) ;  
206 //   triggerPosition[3] = posnxn(0) ;   
207 //   triggerPosition[4] = posnxn(1) ;   
208 //   triggerPosition[5] = posnxn(2) ;  
209
210 //   TArrayF triggerAmplitudes(4);
211 //   triggerAmplitudes[0] = maxAmp2x2 ;   
212 //   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
213 //   triggerAmplitudes[2] = maxAmpnxn ;   
214 //   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
215
216 //   //esd->SetPHOSTriggerCells(triggerPosition);
217 //   esd->AddPHOSTriggerPosition(triggerPosition);
218 //   esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
219   
220
221   //########################################
222   //############# Fill CaloCells ###########
223   //########################################
224
225   Int_t nDigits = fgDigitsArray->GetEntries();
226   Int_t idignew = 0 ;
227   AliDebug(1,Form("%d digits",nDigits));
228
229   const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
230   AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
231   phsCells.CreateContainer(nDigits);
232   phsCells.SetType(AliESDCaloCells::kPHOSCell);
233
234   // Add to CaloCells only EMC digits with non-zero energy 
235   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
236     const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
237     if(dig->GetId() <= knEMC && 
238        Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
239       phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
240                                              CalibrateT(dig->GetTime(),dig->GetId()));   
241       idignew++;
242     }
243   }
244   phsCells.SetNumberOfCells(idignew);
245   phsCells.Sort();
246
247   //########################################
248   //############## Fill CaloClusters #######
249   //########################################
250
251   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
252     AliPHOSRecParticle  *rp    = static_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
253     if (Debug()) 
254       rp->Print();
255     // Get track segment and EMC rec.point associated with this rec.particle
256     AliPHOSTrackSegment *ts    = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
257                                                                     ->At(rp->GetPHOSTSIndex()));
258
259     AliPHOSEmcRecPoint  *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
260     AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
261     
262     Float_t xyz[3];
263     for (Int_t ixyz=0; ixyz<3; ixyz++) 
264       xyz[ixyz] = rp->GetPos()[ixyz];
265     
266     AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
267    
268     // Create cell lists
269
270     Int_t     cellMult   = emcRP->GetDigitsMultiplicity();
271     Int_t    *digitsList = emcRP->GetDigitsList();
272     Float_t  *rpElist    = emcRP->GetEnergiesList() ;
273     UShort_t *absIdList  = new UShort_t[cellMult];
274     Double_t *fracList   = new Double_t[cellMult];
275
276     for (Int_t iCell=0; iCell<cellMult; iCell++) {
277       AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
278       absIdList[iCell] = (UShort_t)(digit->GetId());
279       if (digit->GetEnergy() > 0)
280         fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
281       else
282         fracList[iCell] = 0;
283     }
284
285     //Primaries
286     Int_t  primMult  = 0;
287     Int_t *primList =  emcRP->GetPrimaries(primMult);
288
289     Float_t energy=0.;
290     if (GetRecoParam()->EMCEcore2ESD())
291       energy = emcRP->GetCoreEnergy();
292     else
293       energy = rp->Energy();
294     //Apply nonlinearity correction
295     if(GetRecoParam()->GetEMCEnergyCorrectionOn())
296       energy=CorrectNonlinearity(energy) ;
297
298     // fills the ESDCaloCluster
299     ec->SetType(AliVCluster::kPHOSNeutral);
300     ec->SetPosition(xyz);                       //rec.point position in MARS
301     ec->SetE(energy);                           //total or core particle energy
302     ec->SetDispersion(emcRP->GetDispersion());  //cluster dispersion
303     ec->SetPID(rp->GetPID()) ;            //array of particle identification
304     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
305     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
306     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
307     ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
308     ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z")); 
309     ec->SetChi2(-1);                     //not yet implemented
310     ec->SetTOF(emcRP->GetTime());               //Time of flight - already calibrated in EMCRecPoint
311
312     //Cells contributing to clusters
313     ec->SetNCells(cellMult);
314     ec->SetCellsAbsId(absIdList);
315     ec->SetCellsAmplitudeFraction(fracList);
316
317     //Distance to the nearest bad crystal
318     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
319   
320     //Array of MC indeces
321     TArrayI arrayPrim(primMult,primList);
322     ec->AddLabels(arrayPrim);
323     
324     //Matched ESD track
325     TArrayI arrayTrackMatched(1);
326     arrayTrackMatched[0]= ts->GetTrackIndex();
327     ec->AddTracksMatched(arrayTrackMatched);
328     
329     Int_t index = esd->AddCaloCluster(ec);
330
331     //Set pointer to this cluster in ESD track
332     Int_t nt=esd->GetNumberOfTracks();
333     for (Int_t itr=0; itr<nt; itr++) {
334       AliESDtrack *esdTrack=esd->GetTrack(itr);
335       if(!esdTrack->IsPHOS())
336         continue ;
337       if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
338         esdTrack->SetPHOScluster(index) ;
339 //no garatie that only one track matched this cluster
340 //      break ;
341       }
342     }
343  
344     delete ec;   
345     delete [] fracList;
346     delete [] absIdList;
347   }
348   fgDigitsArray ->Clear();
349   fgEMCRecPoints->Clear("C");
350   recParticles  ->Clear();
351
352   //Store PHOS misalignment matrixes
353   FillMisalMatrixes(esd) ;
354
355 }
356
357 //____________________________________________________________________________
358 AliTracker* AliPHOSReconstructor::CreateTracker() const
359 {
360   // creates the PHOS tracker
361   return new AliPHOSTracker();
362 }
363
364 //____________________________________________________________________________
365 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
366 {
367   // Converts raw data to
368   // PHOS digits
369   // Works on a single-event basis
370   rawReader->Reset() ; 
371
372   AliPHOSRawFitterv0 * fitter ;
373
374   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
375   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
376
377   AliAltroMapping *mapping[20];
378   for(Int_t i = 0; i < 20; i++) {
379     mapping[i] = (AliAltroMapping*)maps->At(i);
380   }
381
382   if      (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0) 
383     fitter=new AliPHOSRawFitterv1();
384   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0) 
385     fitter=new AliPHOSRawFitterv2();
386   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0) 
387     fitter=new AliPHOSRawFitterv3();
388   else
389     fitter=new AliPHOSRawFitterv0();
390
391   fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
392   fitter->SetAmpOffset     (GetRecoParam()->GetGlobalAltroOffset());
393   fitter->SetAmpThreshold  (GetRecoParam()->GetGlobalAltroThreshold());
394
395   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
396   digits->SetName("DIGITS");
397   Int_t bufsize = 32000;
398   digitsTree->Branch("PHOS", &digits, bufsize);
399
400   AliPHOSRawDigiProducer rdp(rawReader,mapping);
401
402   rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
403   rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
404   rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
405   rdp.MakeDigits(digits,fTmpDigLG,fitter);
406
407   delete fitter ;
408
409   if (AliLog::GetGlobalDebugLevel() == 1) {
410     Int_t modMax=-111;
411     Int_t colMax=-111;
412     Int_t rowMax=-111;
413     Float_t eMax=-333;
414     //!!!for debug!!!
415     
416     Int_t relId[4];
417     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
418       AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
419       if(digit->GetEnergy()>eMax) {
420         fGeom->AbsToRelNumbering(digit->GetId(),relId);
421         eMax=digit->GetEnergy();
422         modMax=relId[0];
423         rowMax=relId[2];
424         colMax=relId[3];
425       }
426     }
427     
428     AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
429                     modMax,colMax,rowMax,eMax));
430   }
431
432   digitsTree->Fill();
433   digits->Delete();
434   delete digits;
435 }
436 //==================================================================================
437 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
438   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
439
440   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
441
442   //Determine rel.position of the cell absolute ID
443   Int_t relId[4];
444   geom->AbsToRelNumbering(absId,relId);
445   Int_t module=relId[0];
446   Int_t row   =relId[2];
447   Int_t column=relId[3];
448   if(relId[1]){ //CPV
449     Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
450     return amp*calibration ;
451   }
452   else{ //EMC
453     Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
454     return amp*calibration ;
455   }
456 }
457 //==================================================================================
458 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId)const{
459   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
460
461   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
462
463   //Determine rel.position of the cell absolute ID
464   Int_t relId[4];
465   geom->AbsToRelNumbering(absId,relId);
466   Int_t module=relId[0];
467   Int_t row   =relId[2];
468   Int_t column=relId[3];
469   if(relId[1]){ //CPV
470     return 0. ;
471   }
472   else{ //EMC
473     time += fgCalibData->GetTimeShiftEmc(module,column,row);
474     return time ;
475   }
476 }
477 //==================================================================================
478 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
479   //Store PHOS matrixes in ESD Header
480
481   //Check, if matrixes was already stored
482   for(Int_t mod=0 ;mod<5; mod++){
483     if(esd->GetPHOSMatrix(mod)!=0)
484       return ;
485   }
486
487   //Create and store matrixes
488   if(!gGeoManager){
489     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
490     return ;
491   }
492   //Note, that owner of copied marixes will be header
493   char path[255] ;
494   TGeoHMatrix * m ;
495   for(Int_t mod=0; mod<5; mod++){
496     sprintf(path,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
497     if (gGeoManager->cd(path)){
498       m = gGeoManager->GetCurrentMatrix() ;
499       esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
500     }
501     else{
502       esd->SetPHOSMatrix(NULL,mod) ;
503     }
504   }
505
506 }
507 //==================================================================================
508 Float_t AliPHOSReconstructor::CorrectNonlinearity(Float_t en){
509
510   //For backward compatibility, if no RecoParameters found
511   if(!GetRecoParam()){
512     return 0.0241+1.0504*en+0.000249*en*en ;
513   }
514
515   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"NoCorrection")==0){
516     return en ;
517   }
518   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Gustavo2005")==0){
519     const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
520     return par[0]+par[1]*en + par[2]*en*en ;
521   }
522   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Henrik2010")==0){
523      const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
524      return en*(par[0]+par[1]*TMath::Exp(-en*par[2]))*(1.+par[3]*TMath::Exp(-en*par[4]))*(1.+par[6]/(en*en+par[5])) ;
525   }
526   //For backward compatibility
527   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"")==0){
528     return 0.0241+1.0504*en+0.000249*en*en ;
529   }
530   return en ;
531 }
532