]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
Fixes for not filled histograms and calculation of Dijet binning
[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)
69 {
70   // ctor
71   fGeom        = AliPHOSGeometry::GetInstance("IHEP","");
72   fClusterizer = new AliPHOSClusterizerv1      (fGeom);
73   fTSM         = new AliPHOSTrackSegmentMakerv1(fGeom);
74   fPID         = new AliPHOSPIDv1              (fGeom);
75   fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
76   fgEMCRecPoints= new TObjArray(100) ;
77   if (!fgCalibData)
78     fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
79
80   AliInfo(Form("PHOS bad channel map contains %d bad channel(s).\n",
81                fgCalibData->GetNumOfEmcBadChannels()));
82  
83 }
84
85 //____________________________________________________________________________
86   AliPHOSReconstructor::~AliPHOSReconstructor()
87 {
88   // dtor
89   delete fGeom;
90   delete fClusterizer;
91   delete fTSM;
92   delete fPID;
93   delete fgDigitsArray;
94   delete fgEMCRecPoints;
95
96
97 //____________________________________________________________________________
98 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
99 {
100   // 'single-event' local reco method called by AliReconstruction; 
101   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
102   // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by 
103   // the global tracking.
104
105   fClusterizer->InitParameters();
106   fClusterizer->SetInput(digitsTree);
107   fClusterizer->SetOutput(clustersTree);
108   if ( Debug() ) 
109     fClusterizer->Digits2Clusters("deb all") ; 
110   else 
111     fClusterizer->Digits2Clusters("") ;
112 }
113
114 //____________________________________________________________________________
115 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
116                                    AliESDEvent* esd) const
117 {
118   // This method produces PHOS rec-particles,
119   // then it creates AliESDtracks out of them and
120   // write tracks to the ESD
121
122
123   // do current event; the loop over events is done by AliReconstruction::Run()
124   fTSM->SetESD(esd) ; 
125   fTSM->SetInput(clustersTree);
126   if ( Debug() ) 
127     fTSM->Clusters2TrackSegments("deb all") ;
128   else 
129     fTSM->Clusters2TrackSegments("") ;
130   
131   fPID->SetEnergyCorrectionOn(GetRecoParam()->GetEMCEnergyCorrectionOn());
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    = dynamic_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;
290     if (GetRecoParam()->EMCEcore2ESD())
291       energy = emcRP->GetCoreEnergy();
292     else
293       energy = rp->Energy();
294
295     // fills the ESDCaloCluster
296     ec->SetType(AliVCluster::kPHOSNeutral);
297     ec->SetPosition(xyz);                       //rec.point position in MARS
298     ec->SetE(energy);                           //total or core particle energy
299     ec->SetDispersion(emcRP->GetDispersion());  //cluster dispersion
300     ec->SetPID(rp->GetPID()) ;            //array of particle identification
301     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
302     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
303     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
304     ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
305     ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z")); 
306     ec->SetChi2(-1);                     //not yet implemented
307     ec->SetTOF(emcRP->GetTime());               //Time of flight - already calibrated in EMCRecPoint
308
309     //Cells contributing to clusters
310     ec->SetNCells(cellMult);
311     ec->SetCellsAbsId(absIdList);
312     ec->SetCellsAmplitudeFraction(fracList);
313
314     //Distance to the nearest bad crystal
315     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
316   
317     //Array of MC indeces
318     TArrayI arrayPrim(primMult,primList);
319     ec->AddLabels(arrayPrim);
320     
321     //Matched ESD track
322     TArrayI arrayTrackMatched(1);
323     arrayTrackMatched[0]= ts->GetTrackIndex();
324     ec->AddTracksMatched(arrayTrackMatched);
325     
326     Int_t index = esd->AddCaloCluster(ec);
327
328     //Set pointer to this cluster in ESD track
329     Int_t nt=esd->GetNumberOfTracks();
330     for (Int_t itr=0; itr<nt; itr++) {
331       AliESDtrack *esdTrack=esd->GetTrack(itr);
332       if(!esdTrack->IsPHOS())
333         continue ;
334       if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
335         esdTrack->SetPHOScluster(index) ;
336 //no garatie that only one track matched this cluster
337 //      break ;
338       }
339     }
340  
341     delete ec;   
342     delete [] fracList;
343     delete [] absIdList;
344   }
345   fgDigitsArray ->Delete();
346   fgEMCRecPoints->Delete();
347   recParticles  ->Delete();
348
349   //Store PHOS misalignment matrixes
350   FillMisalMatrixes(esd) ;
351
352 }
353
354 //____________________________________________________________________________
355 AliTracker* AliPHOSReconstructor::CreateTracker() const
356 {
357   // creates the PHOS tracker
358   return new AliPHOSTracker();
359 }
360
361 //____________________________________________________________________________
362 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
363 {
364   // Converts raw data to
365   // PHOS digits
366   // Works on a single-event basis
367   rawReader->Reset() ; 
368
369   AliPHOSRawFitterv0 * fitter ;
370
371   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
372   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
373
374   AliAltroMapping *mapping[20];
375   for(Int_t i = 0; i < 20; i++) {
376     mapping[i] = (AliAltroMapping*)maps->At(i);
377   }
378
379   if      (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0) 
380     fitter=new AliPHOSRawFitterv1();
381   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0) 
382     fitter=new AliPHOSRawFitterv2();
383   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0) 
384     fitter=new AliPHOSRawFitterv3();
385   else
386     fitter=new AliPHOSRawFitterv0();
387
388   fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
389   fitter->SetAmpOffset     (GetRecoParam()->GetGlobalAltroOffset());
390   fitter->SetAmpThreshold  (GetRecoParam()->GetGlobalAltroThreshold());
391
392   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
393   digits->SetName("DIGITS");
394   Int_t bufsize = 32000;
395   digitsTree->Branch("PHOS", &digits, bufsize);
396
397   AliPHOSRawDigiProducer rdp(rawReader,mapping);
398
399   rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
400   rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
401   rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
402   rdp.MakeDigits(digits,fitter);
403
404   delete fitter ;
405
406   if (AliLog::GetGlobalDebugLevel() == 1) {
407     Int_t modMax=-111;
408     Int_t colMax=-111;
409     Int_t rowMax=-111;
410     Float_t eMax=-333;
411     //!!!for debug!!!
412     
413     Int_t relId[4];
414     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
415       AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
416       if(digit->GetEnergy()>eMax) {
417         fGeom->AbsToRelNumbering(digit->GetId(),relId);
418         eMax=digit->GetEnergy();
419         modMax=relId[0];
420         rowMax=relId[2];
421         colMax=relId[3];
422       }
423     }
424     
425     AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
426                     modMax,colMax,rowMax,eMax));
427   }
428
429   digitsTree->Fill();
430   digits->Delete();
431   delete digits;
432 }
433 //==================================================================================
434 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
435   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
436
437   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
438
439   //Determine rel.position of the cell absolute ID
440   Int_t relId[4];
441   geom->AbsToRelNumbering(absId,relId);
442   Int_t module=relId[0];
443   Int_t row   =relId[2];
444   Int_t column=relId[3];
445   if(relId[1]){ //CPV
446     Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
447     return amp*calibration ;
448   }
449   else{ //EMC
450     Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
451     return amp*calibration ;
452   }
453 }
454 //==================================================================================
455 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId)const{
456   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
457
458   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
459
460   //Determine rel.position of the cell absolute ID
461   Int_t relId[4];
462   geom->AbsToRelNumbering(absId,relId);
463   Int_t module=relId[0];
464   Int_t row   =relId[2];
465   Int_t column=relId[3];
466   if(relId[1]){ //CPV
467     return 0. ;
468   }
469   else{ //EMC
470     time += fgCalibData->GetTimeShiftEmc(module,column,row);
471     return time ;
472   }
473 }
474 //==================================================================================
475 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
476   //Store PHOS matrixes in ESD Header
477
478   //Check, if matrixes was already stored
479   for(Int_t mod=0 ;mod<5; mod++){
480     if(esd->GetPHOSMatrix(mod)!=0)
481       return ;
482   }
483
484   //Create and store matrixes
485   if(!gGeoManager){
486     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
487     return ;
488   }
489   //Note, that owner of copied marixes will be header
490   char path[255] ;
491   TGeoHMatrix * m ;
492   for(Int_t mod=0; mod<5; mod++){
493     sprintf(path,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
494     if (gGeoManager->cd(path)){
495       m = gGeoManager->GetCurrentMatrix() ;
496       esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
497     }
498     else{
499       esd->SetPHOSMatrix(NULL,mod) ;
500     }
501   }
502
503 }
504
505