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