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 ---
30 #include "AliESDEvent.h"
31 #include "AliESDCaloCluster.h"
32 #include "AliPHOSReconstructor.h"
33 #include "AliPHOSClusterizerv1.h"
34 #include "AliPHOSTrackSegmentMakerv1.h"
35 #include "AliPHOSPIDv1.h"
36 #include "AliPHOSTracker.h"
37 #include "AliRawReader.h"
38 #include "AliPHOSTrigger.h"
39 #include "AliPHOSGeometry.h"
40 #include "AliPHOSRecoParam.h"
41 #include "AliPHOSRecoParamEmc.h"
42 #include "AliPHOSRecoParamCpv.h"
43 #include "AliPHOSDigit.h"
44 #include "AliPHOSTrackSegment.h"
45 #include "AliPHOSEmcRecPoint.h"
46 #include "AliPHOSRecParticle.h"
47 #include "AliPHOSRawDecoder.h"
48 #include "AliPHOSRawDigiProducer.h"
49 #include "AliPHOSPulseGenerator.h"
51 ClassImp(AliPHOSReconstructor)
53 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
54 AliPHOSRecoParam* AliPHOSReconstructor::fgkRecoParamEmc =0; // EMC rec. parameters
55 AliPHOSRecoParam* AliPHOSReconstructor::fgkRecoParamCpv =0; // CPV rec. parameters
57 //____________________________________________________________________________
58 AliPHOSReconstructor::AliPHOSReconstructor() :
63 if (!fgkRecoParamEmc) {
64 AliWarning("The Reconstruction parameters for EMC nonitialized - Used default one");
65 fgkRecoParamEmc = AliPHOSRecoParamEmc::GetEmcDefaultParameters();
68 if (!fgkRecoParamCpv) {
69 AliWarning("The Reconstruction parameters for CPV nonitialized - Used default one");
70 fgkRecoParamCpv = AliPHOSRecoParamCpv::GetCpvDefaultParameters();
73 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
76 //____________________________________________________________________________
77 AliPHOSReconstructor::~AliPHOSReconstructor()
83 //____________________________________________________________________________
84 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
86 // 'single-event' local reco method called by AliReconstruction;
87 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
88 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
89 // the global tracking.
91 AliPHOSClusterizerv1 clu(fGeom);
92 clu.SetInput(digitsTree);
93 clu.SetOutput(clustersTree);
95 clu.Digits2Clusters("deb all") ;
97 clu.Digits2Clusters("") ;
100 //____________________________________________________________________________
101 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
102 AliESDEvent* esd) const
104 // This method produces PHOS rec-particles,
105 // then it creates AliESDtracks out of them and
106 // write tracks to the ESD
108 AliPHOSTrackSegmentMaker *tsm = new AliPHOSTrackSegmentMakerv1(fGeom);
109 AliPHOSPID *pid = new AliPHOSPIDv1 (fGeom);
111 // do current event; the loop over events is done by AliReconstruction::Run()
113 tsm->SetInput(clustersTree);
115 tsm->Clusters2TrackSegments("deb all") ;
117 tsm->Clusters2TrackSegments("") ;
119 pid->SetInput(clustersTree, tsm->GetTrackSegments()) ;
122 pid->TrackSegments2RecParticles("deb all") ;
124 pid->TrackSegments2RecParticles("") ;
127 // This function creates AliESDtracks from AliPHOSRecParticles
129 // writes them to the ESD
131 TClonesArray *recParticles = pid->GetRecParticles();
132 Int_t nOfRecParticles = recParticles->GetEntries();
134 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
135 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
137 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
140 //#########Calculate trigger and set trigger info###########
143 // tr.SetPatchSize(1);//create 4x4 patches
146 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
147 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
148 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
149 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
151 Int_t iSM2x2 = tr.Get2x2SuperModule();
152 Int_t iSMnxn = tr.GetnxnSuperModule();
153 Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
154 Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
155 Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
156 Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
158 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
159 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
161 Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};// PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
163 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};// PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
165 TVector3 pos2x2(-1,-1,-1);
166 TVector3 posnxn(-1,-1,-1);
167 fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
168 fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
169 fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
170 fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
172 TArrayF triggerPosition(6);
173 triggerPosition[0] = pos2x2(0) ;
174 triggerPosition[1] = pos2x2(1) ;
175 triggerPosition[2] = pos2x2(2) ;
176 triggerPosition[3] = posnxn(0) ;
177 triggerPosition[4] = posnxn(1) ;
178 triggerPosition[5] = posnxn(2) ;
180 TArrayF triggerAmplitudes(4);
181 triggerAmplitudes[0] = maxAmp2x2 ;
182 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
183 triggerAmplitudes[2] = maxAmpnxn ;
184 triggerAmplitudes[3] = ampOutOfPatchnxn ;
186 //esd->SetPHOSTriggerCells(triggerPosition);
187 esd->AddPHOSTriggerPosition(triggerPosition);
188 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
190 //######################################
193 TBranch *branch = digitsTree->GetBranch("PHOS");
195 AliError("can't get the branch with the PHOS digits !");
198 TClonesArray *fDigitsArr = new TClonesArray("AliPHOSDigit",100);
199 branch->SetAddress(&fDigitsArr);
202 // Get the clusters array
203 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
205 AliError("can't get the branch with the PHOS EMC clusters !");
209 TObjArray *fEmcRecPoints = new TObjArray(100) ;
210 emcbranch->SetAddress(&fEmcRecPoints);
211 emcbranch->GetEntry(0);
214 const Float_t kBigShort = std::numeric_limits<short int>::max() - 1;
215 const Float_t nsec100 = 1e9*100.; // units of 0.01 ns
216 const Float_t gev500 = 500.; // units of GeV/500
218 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
219 AliPHOSRecParticle * rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
222 // Get track segment and EMC rec.point associated with this rec.particle
223 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(tsm->GetTrackSegments()->At(rp->GetPHOSTSIndex()));
225 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fEmcRecPoints->At(ts->GetEmcIndex()));
226 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
229 for (Int_t ixyz=0; ixyz<3; ixyz++)
230 xyz[ixyz] = rp->GetPos()[ixyz];
232 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
234 //Create digits lists
235 Int_t digitMult = emcRP->GetDigitsMultiplicity();
236 Int_t *digitsList = emcRP->GetDigitsList();
237 Short_t *amplList = new Short_t[digitMult];
238 Short_t *timeList = new Short_t[digitMult];
239 Short_t *digiList = new Short_t[digitMult];
241 // Convert Float_t* and Int_t* to Short_t* to save memory
242 for (Int_t iDigit=0; iDigit<digitMult; iDigit++) {
243 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fDigitsArr->At(digitsList[iDigit]));
245 (Short_t)(TMath::Min(digit->GetEnergy()*gev500,kBigShort)); // Energy in units of GeV/500
247 (Short_t)(TMath::Min(digit->GetTime()*nsec100,kBigShort)); // time in units of 0.01 ns
248 digiList[iDigit] = (Short_t)(digit->GetId());
253 Int_t *primInts = emcRP->GetPrimaries(primMult);
254 Short_t *primList = new Short_t[primMult];
255 for (Int_t ipr=0; ipr<primMult; ipr++)
256 primList[ipr] = (Short_t)(primInts[ipr]);
258 // fills the ESDCaloCluster
260 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
261 ec->SetPosition(xyz); //rec.point position in MARS
262 ec->SetE(rp->Energy()); //total particle energy
263 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
264 ec->SetPid(rp->GetPID()) ; //array of particle identification
265 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
266 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
267 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
268 ec->SetEmcCpvDistance(-1); //not yet implemented
269 ec->SetClusterChi2(-1); //not yet implemented
270 ec->SetM11(-1) ; //not yet implemented
273 TArrayS arrayAmpList(digitMult,amplList);
274 TArrayS arrayTimeList(digitMult,timeList);
275 TArrayS arrayIndexList(digitMult,digiList);
276 ec->AddDigitAmplitude(arrayAmpList);
277 ec->AddDigitTime(arrayTimeList);
278 ec->AddDigitIndex(arrayIndexList);
280 //Distance to the nearest bad crystal
281 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
283 //Array of MC indeces
284 TArrayS arrayPrim(primMult,primList);
285 ec->AddLabels(arrayPrim);
287 //Array of tracks uncomment when available in future
288 //TArrayS arrayTrackMatched(1);// Only one track, temporal solution.
289 //arrayTrackMatched[0]= (Short_t)(matchedTrack[iClust]);
290 //ec->AddTracksMatched(arrayTrackMatched);
292 // add the track to the esd object
293 esd->AddCaloCluster(ec);
300 fDigitsArr ->Delete();
302 fEmcRecPoints->Delete();
303 delete fEmcRecPoints;
308 //____________________________________________________________________________
309 AliTracker* AliPHOSReconstructor::CreateTracker() const
311 // creates the PHOS tracker
312 return new AliPHOSTracker();
315 //____________________________________________________________________________
316 void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
318 // Converts raw data to
320 // Works on a single-event basis
324 AliPHOSRawDecoder dc(rawReader);
325 TString option = GetOption();
326 if (option.Contains("OldRCUFormat"))
327 dc.SetOldRCUFormat(kTRUE);
329 dc.SetOldRCUFormat(kFALSE);
331 dc.SubtractPedestals(fgkRecoParamEmc->SubtractPedestals());
333 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
334 digits->SetName("DIGITS");
335 Int_t bufsize = 32000;
336 digitsTree->Branch("PHOS", &digits, bufsize);
338 AliPHOSRawDigiProducer pr;
339 pr.MakeDigits(digits,&dc);
342 for(Int_t i=0; i<digits->GetEntries(); i++) {
343 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(i);
344 digit->SetEnergy(digit->GetEnergy()/AliPHOSPulseGenerator::GeV2ADC());
347 // Clean up digits below the noise threshold
348 // Assuming the digit noise to be 4 MeV, we suppress digits within
349 // 3-sigma of the noise.
350 // This parameter should be passed via AliPHOSRecoParamEmc later
352 const Double_t emcDigitThreshold = 0.012;
353 for(Int_t i=0; i<digits->GetEntries(); i++) {
354 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(i);
355 if(digit->GetEnergy() < emcDigitThreshold)
356 digits->RemoveAt(i) ;
368 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
369 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
370 if(digit->GetEnergy()>eMax) {
371 fGeom->AbsToRelNumbering(digit->GetId(),relId);
372 eMax=digit->GetEnergy();
379 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
380 modMax,colMax,rowMax,eMax));