]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Fix from Yves, if fitting is not good, recalculate from parabola; reject bad channels...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALReconstructor.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 //-- Aleksei Pavlinov : added staf for EMCAL jet trigger 9Apr 25, 2008)
25 //                    : fgDigitsArr should read just once at event
26
27 // --- ROOT system ---
28 #include <TList.h>
29 #include <TClonesArray.h>
30 #include <TH2.h>
31 #include "TGeoManager.h"
32 #include "TGeoMatrix.h"
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliEMCALReconstructor.h"
38
39 #include "AliCodeTimer.h"
40 #include "AliESDEvent.h"
41 #include "AliESDCaloCluster.h"
42 #include "AliESDCaloCells.h"
43 #include "AliESDtrack.h"
44 #include "AliEMCALLoader.h"
45 #include "AliEMCALRawUtils.h"
46 #include "AliEMCALDigit.h"
47 #include "AliEMCALClusterizerv1.h"
48 #include "AliEMCALRecPoint.h"
49 #include "AliEMCALPID.h"
50 #include "AliEMCALTrigger.h"
51 #include "AliRawReader.h"
52 #include "AliCDBEntry.h"
53 #include "AliCDBManager.h"
54 #include "AliEMCALGeometry.h"
55 #include "AliEMCAL.h"
56 #include "AliEMCALHistoUtilities.h"
57 #include "AliESDVZERO.h"
58 #include "AliCDBManager.h"
59 #include "AliRunLoader.h"
60 #include "AliRun.h"
61
62 ClassImp(AliEMCALReconstructor) 
63
64 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0;  // EMCAL rec. parameters
65 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0;   // EMCAL raw utilities class
66 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0;   // EMCAL clusterizer class
67 TClonesArray*     AliEMCALReconstructor::fgDigitsArr = 0;  // shoud read just once at event
68 //____________________________________________________________________________
69 AliEMCALReconstructor::AliEMCALReconstructor() 
70   : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0) 
71 {
72   // ctor
73
74   fgRawUtils = new AliEMCALRawUtils;
75
76   //To make sure we match with the geometry in a simulation file,
77   //let's try to get it first.  If not, take the default geometry
78   AliRunLoader *rl = AliRunLoader::Instance();
79   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
80     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
81   } else {
82     AliInfo(Form("Using default geometry in reconstruction"));
83     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
84   }
85
86   //Get calibration parameters  
87   if(!fCalibData)
88     {
89                 AliCDBEntry *entry = (AliCDBEntry*) 
90                 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
91                 if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
92     }
93         
94   if(!fCalibData)
95                 AliFatal("Calibration parameters not found in CDB!");
96         
97   //Get calibration parameters  
98   if(!fPedestalData)
99     {
100                 AliCDBEntry *entry = (AliCDBEntry*) 
101                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
102                 if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
103     }
104         
105         if(!fPedestalData)
106                 AliFatal("Dead map not found in CDB!");
107         
108         
109   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
110   fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData); 
111         
112   if(!fGeom) AliFatal(Form("Could not get geometry!"));
113
114
115
116 //____________________________________________________________________________
117 AliEMCALReconstructor::~AliEMCALReconstructor()
118 {
119   // dtor
120   delete fGeom;
121   delete fgRawUtils;
122   delete fgClusterizer;
123         
124   AliCodeTimer::Instance()->Print();
125
126
127 //____________________________________________________________________________
128 void AliEMCALReconstructor::Init()
129 {
130   // Trigger hists - Oct 24, 2007
131   fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
132 }
133
134 //____________________________________________________________________________
135 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
136 {
137   // method called by AliReconstruction; 
138   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
139   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
140   // the global tracking.
141   // Works on the current event.
142
143   AliCodeTimerAuto("",0)
144
145   ReadDigitsArrayFromTree(digitsTree);
146   fgClusterizer->InitParameters();
147   fgClusterizer->SetOutput(clustersTree);
148         
149   if(fgDigitsArr && fgDigitsArr->GetEntries()) {
150
151     fgClusterizer->SetInput(digitsTree);
152     
153     if(Debug())
154       fgClusterizer->Digits2Clusters("deb all") ;
155     else
156       fgClusterizer->Digits2Clusters("");
157     
158     fgClusterizer->Clear();
159
160   }
161
162 }
163
164 //____________________________________________________________________________
165 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
166
167 {
168   // Conversion from raw data to
169   // EMCAL digits.
170   // Works on a single-event basis
171
172   rawReader->Reset() ; 
173
174   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
175   Int_t bufsize = 32000;
176   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
177
178   //must be done here because, in constructor, option is not yet known
179   fgRawUtils->SetOption(GetOption());
180
181   fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
182   fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
183   fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
184   fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
185   fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
186
187   fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData);
188
189   digitsTree->Fill();
190   //digitsArr->Delete(); //Do not delete digits array are not created here.
191   //delete digitsArr;
192
193 }
194
195
196 //____________________________________________________________________________
197 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
198                                     AliESDEvent* esd) const
199 {
200   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
201   // and V0 
202   // Works on the current event
203   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
204   //return;
205
206   //######################################################
207   //#########Calculate trigger and set trigger info###########
208   //######################################################
209  
210   AliEMCALTrigger tr;
211   //   tr.SetPatchSize(1);  // create 4x4 patches
212   tr.SetSimulation(kFALSE); // Reconstruction mode
213   tr.SetDigitsList(fgDigitsArr);
214   // Get VZERO total multiplicity for jet trigger simulation 
215   // The simulation of jey trigger will be incorrect if no VZERO data 
216   // at ESD
217   AliESDVZERO* vZero = esd->GetVZEROData();
218   if(vZero) {
219     tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
220   }
221   //
222   tr.Trigger();
223
224   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
225   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
226   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
227   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
228
229   Int_t iSM2x2      = tr.Get2x2SuperModule();
230   Int_t iSMnxn      = tr.GetnxnSuperModule();
231   Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
232   Int_t iModulePhinxn = tr.GetnxnModulePhi();
233   Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
234   Int_t iModuleEtanxn = tr.GetnxnModuleEta();
235
236   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
237   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",  maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
238
239   TVector3    pos2x2(-1,-1,-1);
240   TVector3    posnxn(-1,-1,-1);
241
242   Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
243   Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
244   fGeom->GetGlobal(iAbsId2x2, pos2x2);
245   fGeom->GetGlobal(iAbsIdnxn, posnxn);
246   //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
247   
248   TArrayF triggerPosition(6);
249   triggerPosition[0] = pos2x2(0) ;   
250   triggerPosition[1] = pos2x2(1) ;   
251   triggerPosition[2] = pos2x2(2) ;  
252   triggerPosition[3] = posnxn(0) ;   
253   triggerPosition[4] = posnxn(1) ;   
254   triggerPosition[5] = posnxn(2) ;
255   //printf(" triggerPosition ");
256   //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
257
258   TArrayF triggerAmplitudes(4);
259   triggerAmplitudes[0] = maxAmp2x2 ;   
260   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
261   triggerAmplitudes[2] = maxAmpnxn ;   
262   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
263   //printf("\n triggerAmplitudes ");
264   //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
265   //printf("\n");
266   //tr.Print("");
267   //
268   // Trigger jet staff
269   //
270   if(tr.GetNJetThreshold()>0) {
271     // Jet phi/eta
272     Int_t n0 = triggerPosition.GetSize();
273     const TH2F *hpatch = tr.GetJetMatrixE();
274     triggerPosition.Set(n0 + 2);
275     for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);   
276     // Add jet ampitudes
277     n0 = triggerAmplitudes.GetSize();
278     triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
279     Double_t *ampJet = tr.GetL1JetThresholds();
280     for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
281       triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
282     }
283   }
284   esd->AddEMCALTriggerPosition(triggerPosition);
285   esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
286   // Fill trigger hists
287   AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
288
289   //########################################
290   //##############Fill CaloCells###############
291   //########################################
292
293   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
294   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
295   if (!branchdig) { 
296     AliError("can't get the branch with the EMCAL digits !");
297     return;
298   }
299   branchdig->SetAddress(&digits);
300   digitsTree->GetEvent(0);
301   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
302   AliDebug(1,Form("%d digits",nDigits));
303
304   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
305   emcCells.CreateContainer(nDigits);
306   emcCells.SetType(AliESDCaloCells::kEMCALCell);
307   Float_t energy = 0;
308   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
309     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
310     if(dig->GetAmp() > 0 ){
311           energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmp(),dig->GetId());
312           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
313                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
314                   idignew++;
315           }
316     }
317   }
318   emcCells.SetNumberOfCells(idignew);
319   emcCells.Sort();
320
321   //------------------------------------------------------------
322   //-----------------CLUSTERS-----------------------------
323   //------------------------------------------------------------
324   TObjArray *clusters = new TObjArray(100);
325   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
326   branch->SetAddress(&clusters);
327   clustersTree->GetEvent(0);
328
329   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
330   AliDebug(1,Form("%d clusters",nClusters));
331   esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters 
332
333
334   //######################################################
335   //#######################TRACK MATCHING###############
336   //######################################################
337   //Fill list of integers, each one is index of track to which the cluster belongs.
338
339   // step 1 - initialize array of matched track indexes
340   Int_t *matchedTrack = new Int_t[nClusters];
341   for (Int_t iclus = 0; iclus < nClusters; iclus++)
342     matchedTrack[iclus] = -1;  // neg. index --> no matched track
343   
344   // step 2, change the flag for all matched clusters found in tracks
345   Int_t iemcalMatch = -1;
346   Int_t endtpc = esd->GetNumberOfTracks();
347   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
348     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
349     iemcalMatch = track->GetEMCALcluster();
350     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
351   } 
352   
353   //########################################
354   //##############Fill CaloClusters#############
355   //########################################
356   esd->SetNumberOfEMCALClusters(nClusters);
357   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
358     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
359     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
360     if (Debug()) clust->Print();
361     // Get information from EMCAL reconstruction points
362     Float_t xyz[3];
363     TVector3 gpos;
364     clust->GetGlobalPosition(gpos);
365     for (Int_t ixyz=0; ixyz<3; ixyz++)
366       xyz[ixyz] = gpos[ixyz];
367     Float_t elipAxis[2];
368     clust->GetElipsAxis(elipAxis);
369        //Create digits lists
370     Int_t cellMult = clust->GetMultiplicity();
371     //TArrayS digiList(digitMult);
372     Float_t *amplFloat = clust->GetEnergiesList();
373     Int_t   *digitInts = clust->GetAbsId();
374     TArrayS absIdList(cellMult);
375     TArrayD fracList(cellMult);
376
377     Int_t newCellMult = 0;
378     for (Int_t iCell=0; iCell<cellMult; iCell++) {
379       if (amplFloat[iCell] > 0) {
380       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
381       //Uncomment when unfolding is done
382       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
383       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
384       //else
385       fracList[newCellMult] = 0; 
386       newCellMult++;
387       }
388     }
389
390     absIdList.Set(newCellMult);
391     fracList.Set(newCellMult);
392     
393     if(newCellMult > 0) { // accept cluster if it has some digit
394       nClustersNew++;
395       //Primaries
396       Int_t  parentMult  = 0;
397       Int_t *parentList =  clust->GetParents(parentMult);
398       // fills the ESDCaloCluster
399       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
400       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
401       ec->SetPosition(xyz);
402       ec->SetE(clust->GetEnergy());
403                 
404           //Distance to the nearest bad crystal
405           ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
406
407       ec->SetNCells(newCellMult);
408       //Change type of list from short to ushort
409       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
410       Double_t *newFracList  = new Double_t[newCellMult];
411       for(Int_t i = 0; i < newCellMult ; i++) {
412         newAbsIdList[i]=absIdList[i];
413         newFracList[i]=fracList[i];
414       }
415       ec->SetCellsAbsId(newAbsIdList);
416       ec->SetCellsAmplitudeFraction(newFracList);
417       ec->SetClusterDisp(clust->GetDispersion());
418       ec->SetClusterChi2(-1); //not yet implemented
419       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
420       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
421       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
422       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
423       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
424       arrayTrackMatched[0]= matchedTrack[iClust];
425       ec->AddTracksMatched(arrayTrackMatched);
426
427       TArrayI arrayParents(parentMult,parentList);
428       ec->AddLabels(arrayParents);
429
430       // add the cluster to the esd object
431       esd->AddCaloCluster(ec);
432       delete ec;
433       delete [] newAbsIdList ;
434       delete [] newFracList ;
435    }
436  } // cycle on clusters
437
438  delete [] matchedTrack;
439
440  esd->SetNumberOfEMCALClusters(nClustersNew);
441  //if(nClustersNew != nClusters)
442  //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
443
444  //Fill ESDCaloCluster with PID weights
445   AliEMCALPID *pid = new AliEMCALPID;
446   //pid->SetPrintInfo(kTRUE);
447   pid->SetReconstructor(kTRUE);
448   pid->RunPID(esd);
449   delete pid;
450   
451   //delete digits;
452   delete clusters;
453   
454   // printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew); 
455
456   //Store EMCAL misalignment matrixes
457   FillMisalMatrixes(esd) ;
458
459 }
460
461 //==================================================================================
462 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
463         //Store EMCAL matrixes in ESD Header
464         
465         //Check, if matrixes was already stored
466         for(Int_t sm = 0 ; sm < 12; sm++){
467                 if(esd->GetEMCALMatrix(sm)!=0)
468                         return ;
469         }
470         
471         //Create and store matrixes
472         if(!gGeoManager){
473                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
474                 return ;
475         }
476         //Note, that owner of copied marixes will be header
477         char path[255] ;
478         TGeoHMatrix * m = 0x0;
479         for(Int_t sm = 0; sm < 12; sm++){
480                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
481                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
482                 
483                 if (gGeoManager->cd(path)){
484                         m = gGeoManager->GetCurrentMatrix() ;
485                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
486                 }
487                 else{
488                         esd->SetEMCALMatrix(NULL,sm) ;
489                 }
490         }
491 }
492
493
494
495 //__________________________________________________________________________
496 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
497 {
498   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
499   if(fgDigitsArr) {
500     // Clear previous digits 
501     fgDigitsArr->Delete();
502     delete fgDigitsArr;
503   }
504   // Read the digits from the input tree
505   TBranch *branch = digitsTree->GetBranch("EMCAL");
506   if (!branch) { 
507     AliError("can't get the branch with the EMCAL digits !");
508     return;
509   }
510   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
511   branch->SetAddress(&fgDigitsArr);
512   branch->GetEntry(0);
513 }
514