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