]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
chi2, mean and rms checks for input samples for fitting
[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 "AliESDVZERO.h"
57 #include "AliCDBManager.h"
58 #include "AliRunLoader.h"
59 #include "AliRun.h"
60 #include "AliEMCALTriggerData.h"
61 #include "AliEMCALTriggerElectronics.h"
62 #include "AliEMCALTriggerDCSConfigDB.h"
63 #include "AliEMCALTriggerDCSConfig.h"
64
65 ClassImp(AliEMCALReconstructor) 
66
67 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0;  // EMCAL rec. parameters
68 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0;   // EMCAL raw utilities class
69 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0;   // EMCAL clusterizer class
70 TClonesArray*     AliEMCALReconstructor::fgDigitsArr = 0;  // shoud read just once at event
71 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
72 //____________________________________________________________________________
73 AliEMCALReconstructor::AliEMCALReconstructor() 
74   : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0) 
75 {
76   // ctor
77
78   fgRawUtils = new AliEMCALRawUtils;
79
80   //To make sure we match with the geometry in a simulation file,
81   //let's try to get it first.  If not, take the default geometry
82   AliRunLoader *rl = AliRunLoader::Instance();
83   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
84     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
85   } else {
86     AliInfo(Form("Using default geometry in reconstruction"));
87     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
88   }
89
90   //Get calibration parameters  
91   if(!fCalibData)
92     {
93                 AliCDBEntry *entry = (AliCDBEntry*) 
94                 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
95                 if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
96     }
97         
98   if(!fCalibData)
99                 AliFatal("Calibration parameters not found in CDB!");
100         
101   //Get calibration parameters  
102   if(!fPedestalData)
103     {
104                 AliCDBEntry *entry = (AliCDBEntry*) 
105                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
106                 if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
107     }
108         
109         if(!fPedestalData)
110                 AliFatal("Dead map not found in CDB!");
111         
112         
113   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
114   fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData); 
115         
116   if(!fGeom) AliFatal(Form("Could not get geometry!"));
117
118   AliEMCALTriggerDCSConfigDB* dcsConfigDB = AliEMCALTriggerDCSConfigDB::Instance();
119
120   const AliEMCALTriggerDCSConfig* dcsConfig = dcsConfigDB->GetTriggerDCSConfig();
121
122   if (!dcsConfig) AliFatal("No Trigger DCS Configuration from OCDB!");
123   fgTriggerProcessor = new AliEMCALTriggerElectronics( dcsConfig );
124
125
126 //____________________________________________________________________________
127 AliEMCALReconstructor::~AliEMCALReconstructor()
128 {
129   // dtor
130   delete fGeom;
131   delete fgRawUtils;
132   delete fgClusterizer;
133   delete fgTriggerProcessor;
134
135   AliCodeTimer::Instance()->Print();
136
137
138 // //____________________________________________________________________________
139 // void AliEMCALReconstructor::Init()
140 // {
141 //   // Trigger hists - Oct 24, 2007
142 //    fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
143 // }
144
145 //____________________________________________________________________________
146 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
147 {
148   // method called by AliReconstruction; 
149   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
150   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
151   // the global tracking.
152   // Works on the current event.
153
154   AliCodeTimerAuto("",0)
155
156   ReadDigitsArrayFromTree(digitsTree);
157   fgClusterizer->InitParameters();
158   fgClusterizer->SetOutput(clustersTree);
159
160   AliEMCALTriggerData* trgData = new AliEMCALTriggerData();
161         
162   Int_t bufferSize = 32000;
163
164   if (TBranch* triggerBranch = clustersTree->GetBranch("EMTRG"))
165           triggerBranch->SetAddress(&trgData);
166   else
167           clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
168
169   TClonesArray *trgDigits = new TClonesArray("AliEMCALRawDigit",1000);
170   TBranch *branchdig = digitsTree->GetBranch("EMTRG");
171   if (!branchdig) 
172   { 
173           AliError("Can't get the branch with the EMCAL trigger digits !");
174           return;
175   }
176
177   branchdig->SetAddress(&trgDigits);
178   branchdig->GetEntry(0);
179
180   Int_t v0M[2] = {0,0};
181   fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, trgData);
182         
183   trgDigits->Delete();
184   delete trgDigits; trgDigits = 0x0;
185
186   if(fgDigitsArr && fgDigitsArr->GetEntries()) {
187
188     fgClusterizer->SetInput(digitsTree);
189     
190     if(Debug())
191       fgClusterizer->Digits2Clusters("deb all") ;
192     else
193       fgClusterizer->Digits2Clusters("");
194     
195     fgClusterizer->Clear();
196
197   }
198
199   clustersTree->Fill(); 
200
201   delete trgData; trgData = 0x0;
202 }
203
204 //____________________________________________________________________________
205 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
206
207 {
208   // Conversion from raw data to
209   // EMCAL digits.
210   // Works on a single-event basis
211
212   rawReader->Reset() ; 
213
214   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
215   TClonesArray *digitsTrg = new TClonesArray("AliEMCALRawDigit", 200);
216
217   Int_t bufsize = 32000;
218   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
219   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
220
221   //must be done here because, in constructor, option is not yet known
222   fgRawUtils->SetOption(GetOption());
223
224   fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
225   fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
226   fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
227   fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
228   fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
229   fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
230   fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
231   fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
232   fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
233   fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
234         
235   fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
236
237   digitsTree->Fill();
238   digitsArr->Delete();
239   digitsTrg->Delete();
240   delete digitsArr;
241   delete digitsTrg;
242
243 }
244
245
246 //____________________________________________________________________________
247 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
248                                     AliESDEvent* esd) const
249 {
250   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
251   // and V0 
252   // Works on the current event
253   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
254   //return;
255
256   //########################################
257   //##############Fill CaloCells###############
258   //########################################
259
260   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
261   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
262   if (!branchdig) { 
263     AliError("can't get the branch with the EMCAL digits !");
264     return;
265   }
266   branchdig->SetAddress(&digits);
267   digitsTree->GetEvent(0);
268   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
269   AliDebug(1,Form("%d digits",nDigits));
270
271   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
272   emcCells.CreateContainer(nDigits);
273   emcCells.SetType(AliESDCaloCells::kEMCALCell);
274   Float_t energy = 0;
275   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
276     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
277     if(dig->GetAmplitude() > 0 ){
278           energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
279           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
280                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
281                   idignew++;
282           }
283     }
284   }
285   emcCells.SetNumberOfCells(idignew);
286   emcCells.Sort();
287
288   //------------------------------------------------------------
289   //-----------------CLUSTERS-----------------------------
290   //------------------------------------------------------------
291   clustersTree->SetBranchStatus("*",0); //disable all branches
292   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
293
294   TObjArray *clusters = new TObjArray(100);
295   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
296   branch->SetAddress(&clusters);
297   branch->GetEntry(0);
298   //clustersTree->GetEvent(0);
299
300   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
301   AliDebug(1,Form("%d clusters",nClusters));
302
303   //######################################################
304   //#######################TRACK MATCHING###############
305   //######################################################
306   //Fill list of integers, each one is index of track to which the cluster belongs.
307
308   // step 1 - initialize array of matched track indexes
309   Int_t *matchedTrack = new Int_t[nClusters];
310   for (Int_t iclus = 0; iclus < nClusters; iclus++)
311     matchedTrack[iclus] = -1;  // neg. index --> no matched track
312   
313   // step 2, change the flag for all matched clusters found in tracks
314   Int_t iemcalMatch = -1;
315   Int_t endtpc = esd->GetNumberOfTracks();
316   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
317     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
318     iemcalMatch = track->GetEMCALcluster();
319     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
320   } 
321   
322   //########################################
323   //##############Fill CaloClusters#############
324   //########################################
325   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
326     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
327     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
328     if (Debug()) clust->Print();
329     // Get information from EMCAL reconstruction points
330     Float_t xyz[3];
331     TVector3 gpos;
332     clust->GetGlobalPosition(gpos);
333     for (Int_t ixyz=0; ixyz<3; ixyz++)
334       xyz[ixyz] = gpos[ixyz];
335     Float_t elipAxis[2];
336     clust->GetElipsAxis(elipAxis);
337        //Create digits lists
338     Int_t cellMult = clust->GetMultiplicity();
339     //TArrayS digiList(digitMult);
340     Float_t *amplFloat = clust->GetEnergiesList();
341     Int_t   *digitInts = clust->GetAbsId();
342     TArrayS absIdList(cellMult);
343     TArrayD fracList(cellMult);
344
345     Int_t newCellMult = 0;
346     for (Int_t iCell=0; iCell<cellMult; iCell++) {
347       if (amplFloat[iCell] > 0) {
348       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
349       //Uncomment when unfolding is done
350       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
351       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
352       //else
353       fracList[newCellMult] = 0; 
354       newCellMult++;
355       }
356     }
357
358     absIdList.Set(newCellMult);
359     fracList.Set(newCellMult);
360     
361     if(newCellMult > 0) { // accept cluster if it has some digit
362       nClustersNew++;
363       //Primaries
364       Int_t  parentMult  = 0;
365       Int_t *parentList =  clust->GetParents(parentMult);
366       // fills the ESDCaloCluster
367       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
368       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
369       ec->SetPosition(xyz);
370       ec->SetE(clust->GetEnergy());
371                 
372           //Distance to the nearest bad crystal
373           ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
374
375       ec->SetNCells(newCellMult);
376       //Change type of list from short to ushort
377       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
378       Double_t *newFracList  = new Double_t[newCellMult];
379       for(Int_t i = 0; i < newCellMult ; i++) {
380         newAbsIdList[i]=absIdList[i];
381         newFracList[i]=fracList[i];
382       }
383       ec->SetCellsAbsId(newAbsIdList);
384       ec->SetCellsAmplitudeFraction(newFracList);
385       ec->SetClusterDisp(clust->GetDispersion());
386       ec->SetClusterChi2(-1); //not yet implemented
387       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
388       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
389       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
390       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
391       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
392       arrayTrackMatched[0]= matchedTrack[iClust];
393       ec->AddTracksMatched(arrayTrackMatched);
394
395       TArrayI arrayParents(parentMult,parentList);
396       ec->AddLabels(arrayParents);
397
398       // add the cluster to the esd object
399       esd->AddCaloCluster(ec);
400       delete ec;
401       delete [] newAbsIdList ;
402       delete [] newFracList ;
403    }
404  } // cycle on clusters
405
406  delete [] matchedTrack;
407
408  //Fill ESDCaloCluster with PID weights
409  AliEMCALPID *pid = new AliEMCALPID;
410  //pid->SetPrintInfo(kTRUE);
411  pid->SetReconstructor(kTRUE);
412  pid->RunPID(esd);
413  delete pid;
414   
415  delete digits;
416  delete clusters;
417   
418  //Store EMCAL misalignment matrixes
419  FillMisalMatrixes(esd) ;
420
421 }
422
423 //==================================================================================
424 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
425         //Store EMCAL matrixes in ESD Header
426         
427         //Check, if matrixes was already stored
428         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
429                 if(esd->GetEMCALMatrix(sm)!=0)
430                         return ;
431         }
432         
433         //Create and store matrixes
434         if(!gGeoManager){
435                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
436                 return ;
437         }
438         //Note, that owner of copied marixes will be header
439         char path[255] ;
440         TGeoHMatrix * m = 0x0;
441         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
442                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
443                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
444                 
445                 if (gGeoManager->CheckPath(path)){
446                         gGeoManager->cd(path);
447                         m = gGeoManager->GetCurrentMatrix() ;
448 //                      printf("================================================= \n");
449 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
450 //                      m->Print("");
451                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
452 //                      printf("================================================= \n");
453                 }
454                 else{
455                         esd->SetEMCALMatrix(NULL,sm) ;
456                 }
457         }
458 }
459
460
461
462 //__________________________________________________________________________
463 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
464 {
465   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
466   if(fgDigitsArr) {
467     // Clear previous digits 
468     fgDigitsArr->Delete();
469     delete fgDigitsArr;
470   }
471   // Read the digits from the input tree
472   TBranch *branch = digitsTree->GetBranch("EMCAL");
473   if (!branch) { 
474     AliError("can't get the branch with the EMCAL digits !");
475     return;
476   }
477   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
478   branch->SetAddress(&fgDigitsArr);
479   branch->GetEntry(0);
480 }
481
482