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