1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
20 //-- Yves Schutz (SUBATECH)
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and
22 // derived from STEER/AliReconstructor.
24 // --- ROOT system ---
26 // --- Standard library ---
28 // --- AliRoot header files ---
29 #include "AliEMCALReconstructor.h"
31 #include "AliCodeTimer.h"
32 #include "AliESDEvent.h"
33 #include "AliESDCaloCluster.h"
34 #include "AliESDCaloCells.h"
35 #include "AliESDtrack.h"
36 #include "AliEMCALLoader.h"
37 #include "AliEMCALRawUtils.h"
38 #include "AliEMCALDigit.h"
39 #include "AliEMCALClusterizerv1.h"
40 #include "AliEMCALRecPoint.h"
41 #include "AliEMCALPID.h"
42 #include "AliEMCALTrigger.h"
43 #include "AliRawReader.h"
44 #include "AliCDBEntry.h"
45 #include "AliCDBManager.h"
46 #include "AliEMCALRecParam.h"
47 #include "AliEMCALGeometry.h"
50 ClassImp(AliEMCALReconstructor)
52 AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
53 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
54 //____________________________________________________________________________
55 AliEMCALReconstructor::AliEMCALReconstructor()
56 : fDebug(kFALSE),fGeom(0)
61 fgRawUtils = new AliEMCALRawUtils;
62 fGeom = AliEMCALGeometry::GetInstance();
64 fGeom = AliEMCALGeometry::GetInstance("","");
65 if(!fGeom) AliFatal(Form("Could not get geometry!"));
70 //____________________________________________________________________________
71 AliEMCALReconstructor::AliEMCALReconstructor(const AliEMCALReconstructor & rec)
72 : AliReconstructor(rec),
79 //____________________________________________________________________________
80 AliEMCALReconstructor::~AliEMCALReconstructor()
84 AliCodeTimer::Instance()->Print();
87 //____________________________________________________________________________
88 void AliEMCALReconstructor::InitRecParam() const
90 // Check if the instance of AliEMCALRecParam exists,
91 // if not, get it from OCDB if available, otherwise create a default one
93 if (!fgkRecParam && (AliCDBManager::Instance()->IsDefaultStorageSet())) {
94 AliCDBEntry *entry = (AliCDBEntry*)
95 AliCDBManager::Instance()->Get("EMCAL/Config/RecParam");
96 if (entry) fgkRecParam = (AliEMCALRecParam*) entry->GetObject();
100 AliWarning("The Reconstruction parameters for EMCAL nonitialized - Used default one");
101 fgkRecParam = new AliEMCALRecParam;
105 //____________________________________________________________________________
106 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
108 // method called by AliReconstruction;
109 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
110 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
111 // the global tracking.
112 // Works on the current event.
116 AliEMCALClusterizerv1 clu;
117 clu.SetInput(digitsTree);
118 clu.SetOutput(clustersTree);
120 clu.Digits2Clusters("deb all") ;
122 clu.Digits2Clusters("") ;
126 //____________________________________________________________________________
127 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
130 // Conversion from raw data to
132 // Works on a single-event basis
136 TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",100);
137 Int_t bufsize = 32000;
138 digitsTree->Branch("EMCAL", &digitsArr, bufsize);
140 //must be done here because, in constructor, option is not yet known
141 fgRawUtils->SetOption(GetOption());
142 fgRawUtils->Raw2Digits(rawReader,digitsArr);
150 //____________________________________________________________________________
151 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
152 AliESDEvent* esd) const
154 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
155 // Works on the current event
156 // Creates AliESDCaloCluster from AliEMCALRecPoints
157 // and AliESDCaloCells from AliEMCALDigits
158 // Also, fills ESD with calorimeter trigger information
160 //######################################################
161 //#########Calculate trigger and set trigger info###########
162 //######################################################
165 // tr.SetPatchSize(1);//create 4x4 patches
168 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
169 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
170 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
171 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
173 Int_t iSM2x2 = tr.Get2x2SuperModule();
174 Int_t iSMnxn = tr.GetnxnSuperModule();
175 Int_t iCellPhi2x2 = tr.Get2x2CellPhi();
176 Int_t iCellPhinxn = tr.GetnxnCellPhi();
177 Int_t iCellEta2x2 = tr.Get2x2CellEta();
178 Int_t iCellEtanxn = tr.GetnxnCellEta();
180 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCellPhi2x2, iCellEta2x2));
181 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCellPhinxn, iCellEtanxn));
183 TVector3 pos2x2(-1,-1,-1);
184 TVector3 posnxn(-1,-1,-1);
186 Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iCellPhi2x2, iCellEta2x2) ;
187 Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iCellPhinxn, iCellEtanxn) ;
188 fGeom->GetGlobal(iAbsId2x2, pos2x2);
189 fGeom->GetGlobal(iAbsIdnxn, posnxn);
191 TArrayF triggerPosition(6);
192 triggerPosition[0] = pos2x2(0) ;
193 triggerPosition[1] = pos2x2(1) ;
194 triggerPosition[2] = pos2x2(2) ;
195 triggerPosition[3] = posnxn(0) ;
196 triggerPosition[4] = posnxn(1) ;
197 triggerPosition[5] = posnxn(2) ;
199 TArrayF triggerAmplitudes(4);
200 triggerAmplitudes[0] = maxAmp2x2 ;
201 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
202 triggerAmplitudes[2] = maxAmpnxn ;
203 triggerAmplitudes[3] = ampOutOfPatchnxn ;
205 esd->AddEMCALTriggerPosition(triggerPosition);
206 esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
208 //########################################
209 //##############Fill CaloCells###############
210 //########################################
212 TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
213 TBranch *branchdig = digitsTree->GetBranch("EMCAL");
215 AliError("can't get the branch with the PHOS digits !");
218 branchdig->SetAddress(&digits);
219 digitsTree->GetEvent(0);
220 Int_t nDigits = digits->GetEntries(), idignew = 0 ;
221 AliDebug(1,Form("%d digits",nDigits));
223 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
224 emcCells.CreateContainer(nDigits);
225 emcCells.SetType(AliESDCaloCells::kEMCALCell);
226 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
227 const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
228 if(dig->GetAmp() > 0 ){
229 emcCells.SetCell(idignew,dig->GetId(),dig->GetAmp(), dig->GetTime());
233 emcCells.SetNumberOfCells(idignew);
236 //------------------------------------------------------------
237 //-----------------CLUSTERS-----------------------------
238 //------------------------------------------------------------
239 TObjArray *clusters = new TObjArray(100);
240 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
241 branch->SetAddress(&clusters);
242 clustersTree->GetEvent(0);
244 Int_t nClusters = clusters->GetEntries(), nClustersNew=0;
245 AliDebug(1,Form("%d clusters",nClusters));
246 esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters
248 //######################################################
249 //#######################TRACK MATCHING###############
250 //######################################################
251 //Fill list of integers, each one is index of track to which the cluster belongs.
253 // step 1 - initialize array of matched track indexes
254 Int_t *matchedTrack = new Int_t[nClusters];
255 for (Int_t iclus = 0; iclus < nClusters; iclus++)
256 matchedTrack[iclus] = -1; // neg. index --> no matched track
258 // step 2, change the flag for all matched clusters found in tracks
259 Int_t iemcalMatch = -1;
260 Int_t endtpc = esd->GetNumberOfTracks();
261 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
262 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
263 iemcalMatch = track->GetEMCALcluster();
264 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
267 //########################################
268 //##############Fill CaloClusters############
269 //########################################
271 esd->SetNumberOfEMCALClusters(nClusters);
272 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
273 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
274 //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
275 if (Debug()) clust->Print();
276 // Get information from EMCAL reconstruction points
279 clust->GetGlobalPosition(gpos);
280 for (Int_t ixyz=0; ixyz<3; ixyz++)
281 xyz[ixyz] = gpos[ixyz];
283 clust->GetElipsAxis(elipAxis);
285 //Create digits lists
286 Int_t cellMult = clust->GetMultiplicity();
287 //TArrayS digiList(digitMult);
288 Float_t *amplFloat = clust->GetEnergiesList();
289 Int_t *digitInts = clust->GetAbsId();
290 TArrayS absIdList(cellMult);
291 //Uncomment when unfolding is done
292 //TArrayD fracList(cellMult);
294 Int_t newCellMult = 0;
295 for (Int_t iCell=0; iCell<cellMult; iCell++) {
296 if (amplFloat[iCell] > 0) {
297 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
298 //Uncomment when unfolding is done
299 //fracList[newCellMult] = amplFloat[iCell]/emcCells.GetCellAmplitude(digitInts[iCell]);
303 absIdList.Set(newCellMult);
304 //Uncomment when unfolding is done
305 //fracList.Set(newCellMult);
307 if(newCellMult > 0) { // accept cluster if it has some digit
310 Int_t parentMult = 0;
311 Int_t *parentList = clust->GetParents(parentMult);
313 // fills the ESDCaloCluster
314 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
315 ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
316 ec->SetPosition(xyz);
317 ec->SetE(clust->GetEnergy());
318 ec->SetNCells(newCellMult);
319 //Change type of list from short to ushort
320 UShort_t *newAbsIdList = new UShort_t[newCellMult];
321 //Uncomment when unfolding is done
322 //Double_t *newFracList = new Double_t[newCellMult];
323 for(Int_t i = 0; i < newCellMult ; i++) {
324 newAbsIdList[i]=absIdList[i];
325 //Uncomment when unfolding is done
326 //newFracList[i]=fracList[i];
328 ec->SetCellsAbsId(newAbsIdList);
329 //Uncomment when unfolding is done
330 //ec->SetCellsAmplitudeFraction(newFracList);
332 ec->SetClusterDisp(clust->GetDispersion());
333 ec->SetClusterChi2(-1); //not yet implemented
334 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
335 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
336 ec->SetM11(-1) ; //not yet implemented
338 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
339 arrayTrackMatched[0]= matchedTrack[iClust];
340 ec->AddTracksMatched(arrayTrackMatched);
342 TArrayI arrayParents(parentMult,parentList);
343 ec->AddLabels(arrayParents);
346 // add the cluster to the esd object
347 esd->AddCaloCluster(ec);
349 //delete [] newAbsIdList ;
350 //delete [] newFracList ;
352 } // cycle on clusters
354 delete [] matchedTrack;
356 esd->SetNumberOfEMCALClusters(nClustersNew);
357 //if(nClustersNew != nClusters)
358 //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
360 //Fill ESDCaloCluster with PID weights
361 AliEMCALPID *pid = new AliEMCALPID;
362 //pid->SetPrintInfo(kTRUE);
363 pid->SetReconstructor(kTRUE);