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 "AliAltroMapping.h"
31 #include "AliESDEvent.h"
32 #include "AliESDCaloCluster.h"
33 #include "AliESDCaloCells.h"
34 #include "AliPHOSReconstructor.h"
35 #include "AliPHOSClusterizerv1.h"
36 #include "AliPHOSTrackSegmentMakerv1.h"
37 #include "AliPHOSPIDv1.h"
38 #include "AliPHOSTracker.h"
39 #include "AliRawReader.h"
40 #include "AliPHOSTrigger.h"
41 #include "AliPHOSGeometry.h"
42 #include "AliPHOSRecoParam.h"
43 #include "AliPHOSRecoParamEmc.h"
44 #include "AliPHOSRecoParamCpv.h"
45 #include "AliPHOSDigit.h"
46 #include "AliPHOSTrackSegment.h"
47 #include "AliPHOSEmcRecPoint.h"
48 #include "AliPHOSRecParticle.h"
49 #include "AliPHOSRawDecoder.h"
50 #include "AliPHOSRawDecoderv1.h"
51 #include "AliPHOSRawDecoderv2.h"
52 #include "AliPHOSRawDigiProducer.h"
53 #include "AliPHOSPulseGenerator.h"
55 ClassImp(AliPHOSReconstructor)
57 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
58 AliPHOSRecoParam* AliPHOSReconstructor::fgkRecoParamEmc =0; // EMC rec. parameters
59 AliPHOSRecoParam* AliPHOSReconstructor::fgkRecoParamCpv =0; // CPV rec. parameters
61 //____________________________________________________________________________
62 AliPHOSReconstructor::AliPHOSReconstructor() :
63 fGeom(NULL),fClusterizer(NULL)
67 if (!fgkRecoParamEmc) {
68 AliWarning("The Reconstruction parameters for EMC nonitialized - Used default one");
69 fgkRecoParamEmc = AliPHOSRecoParamEmc::GetEmcDefaultParameters();
72 if (!fgkRecoParamCpv) {
73 AliWarning("The Reconstruction parameters for CPV nonitialized - Used default one");
74 fgkRecoParamCpv = AliPHOSRecoParamCpv::GetCpvDefaultParameters();
77 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
78 fClusterizer = new AliPHOSClusterizerv1(fGeom);
81 //____________________________________________________________________________
82 AliPHOSReconstructor::~AliPHOSReconstructor()
89 //____________________________________________________________________________
90 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
92 // 'single-event' local reco method called by AliReconstruction;
93 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
94 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
95 // the global tracking.
97 fClusterizer->SetInput(digitsTree);
98 fClusterizer->SetOutput(clustersTree);
100 fClusterizer->Digits2Clusters("deb all") ;
102 fClusterizer->Digits2Clusters("") ;
105 //____________________________________________________________________________
106 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
107 AliESDEvent* esd) const
109 // This method produces PHOS rec-particles,
110 // then it creates AliESDtracks out of them and
111 // write tracks to the ESD
113 AliPHOSTrackSegmentMaker *tsm = new AliPHOSTrackSegmentMakerv1(fGeom);
114 AliPHOSPID *pid = new AliPHOSPIDv1 (fGeom);
116 // do current event; the loop over events is done by AliReconstruction::Run()
118 tsm->SetInput(clustersTree);
120 tsm->Clusters2TrackSegments("deb all") ;
122 tsm->Clusters2TrackSegments("") ;
124 pid->SetInput(clustersTree, tsm->GetTrackSegments()) ;
127 pid->TrackSegments2RecParticles("deb all") ;
129 pid->TrackSegments2RecParticles("") ;
132 // This function creates AliESDtracks from AliPHOSRecParticles
134 // writes them to the ESD
136 TClonesArray *recParticles = pid->GetRecParticles();
137 Int_t nOfRecParticles = recParticles->GetEntries();
139 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
140 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
142 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
146 TBranch *branch = digitsTree->GetBranch("PHOS");
148 AliError("can't get the branch with the PHOS digits !");
151 TClonesArray *digitsArray = new TClonesArray("AliPHOSDigit",100);
152 branch->SetAddress(&digitsArray);
155 // Get the clusters array
157 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
159 AliError("can't get the branch with the PHOS EMC clusters !");
163 TObjArray *emcRecPoints = new TObjArray(100) ;
164 emcbranch->SetAddress(&emcRecPoints);
165 emcbranch->GetEntry(0);
167 //#########Calculate trigger and set trigger info###########
170 // tr.SetPatchSize(1);//create 4x4 patches
171 tr.SetSimulation(kFALSE);
172 tr.Trigger(digitsArray);
174 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
175 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
176 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
177 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
179 Int_t iSM2x2 = tr.Get2x2SuperModule();
180 Int_t iSMnxn = tr.GetnxnSuperModule();
181 Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
182 Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
183 Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
184 Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
186 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",
187 maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
188 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
189 maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
191 // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
192 Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
194 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
196 TVector3 pos2x2(-1,-1,-1);
197 TVector3 posnxn(-1,-1,-1);
198 fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
199 fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
200 fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
201 fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
203 TArrayF triggerPosition(6);
204 triggerPosition[0] = pos2x2(0) ;
205 triggerPosition[1] = pos2x2(1) ;
206 triggerPosition[2] = pos2x2(2) ;
207 triggerPosition[3] = posnxn(0) ;
208 triggerPosition[4] = posnxn(1) ;
209 triggerPosition[5] = posnxn(2) ;
211 TArrayF triggerAmplitudes(4);
212 triggerAmplitudes[0] = maxAmp2x2 ;
213 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
214 triggerAmplitudes[2] = maxAmpnxn ;
215 triggerAmplitudes[3] = ampOutOfPatchnxn ;
217 //esd->SetPHOSTriggerCells(triggerPosition);
218 esd->AddPHOSTriggerPosition(triggerPosition);
219 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
222 //########################################
223 //############# Fill CaloCells ###########
224 //########################################
226 Int_t nDigits = digitsArray->GetEntries();
228 AliDebug(1,Form("%d digits",nDigits));
230 const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
231 AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
232 phsCells.CreateContainer(nDigits);
233 phsCells.SetType(AliESDCaloCells::kPHOSCell);
235 // Add to CaloCells only EMC digits with non-zero energy
236 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
237 const AliPHOSDigit * dig = (const AliPHOSDigit*)digitsArray->At(idig);
238 if(dig->GetId() <= knEMC && dig->GetEnergy() > 0 ){
239 //printf("i %d; id %d; amp %f; time %e\n",
240 //idignew,dig->GetId(),dig->GetEnergy(), dig->GetTime());
241 phsCells.SetCell(idignew,dig->GetId(), dig->GetEnergy(), dig->GetTime());
245 phsCells.SetNumberOfCells(idignew);
248 //########################################
249 //############## Fill CaloClusters #######
250 //########################################
252 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
253 AliPHOSRecParticle *rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
256 // Get track segment and EMC rec.point associated with this rec.particle
257 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(tsm->GetTrackSegments()
258 ->At(rp->GetPHOSTSIndex()));
260 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(ts->GetEmcIndex()));
261 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
264 for (Int_t ixyz=0; ixyz<3; ixyz++)
265 xyz[ixyz] = rp->GetPos()[ixyz];
267 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
271 Int_t cellMult = emcRP->GetDigitsMultiplicity();
272 Int_t *digitsList = emcRP->GetDigitsList();
273 Float_t *rpElist = emcRP->GetEnergiesList() ;
274 UShort_t *absIdList = new UShort_t[cellMult];
275 Double_t *fracList = new Double_t[cellMult];
277 for (Int_t iCell=0; iCell<cellMult; iCell++) {
278 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(digitsArray->At(digitsList[iCell]));
279 absIdList[iCell] = (UShort_t)(digit->GetId());
280 if (digit->GetEnergy() > 0)
281 fracList[iCell] = rpElist[iCell]/digit->GetEnergy();
288 Int_t *primList = emcRP->GetPrimaries(primMult);
290 // fills the ESDCaloCluster
291 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
292 ec->SetPosition(xyz); //rec.point position in MARS
293 ec->SetE(rp->Energy()); //total particle energy
294 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
295 ec->SetPid(rp->GetPID()) ; //array of particle identification
296 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
297 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
298 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
299 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
300 ec->SetClusterChi2(-1); //not yet implemented
301 ec->SetM11(-1) ; //not yet implemented
303 //Cells contributing to clusters
304 ec->SetNCells(cellMult);
305 ec->SetCellsAbsId(absIdList);
306 ec->SetCellsAmplitudeFraction(fracList);
308 //Distance to the nearest bad crystal
309 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
311 //Array of MC indeces
312 TArrayI arrayPrim(primMult,primList);
313 ec->AddLabels(arrayPrim);
315 //Array of tracks uncomment when available in future
316 //TArrayS arrayTrackMatched(1);// Only one track, temporal solution.
317 //arrayTrackMatched[0]= (Short_t)(matchedTrack[iClust]);
318 //ec->AddTracksMatched(arrayTrackMatched);
320 // add the track to the esd object
322 esd->AddCaloCluster(ec);
325 digitsArray ->Delete();
327 emcRecPoints->Delete();
333 //____________________________________________________________________________
334 AliTracker* AliPHOSReconstructor::CreateTracker() const
336 // creates the PHOS tracker
337 return new AliPHOSTracker();
340 //____________________________________________________________________________
341 void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
343 // Converts raw data to
345 // Works on a single-event basis
349 AliPHOSRawDecoder * dc ;
351 const TObjArray* maps = AliPHOSRecoParamEmc::GetMappings();
352 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
354 AliAltroMapping *mapping[4];
355 for(Int_t i = 0; i < 4; i++) {
356 mapping[i] = (AliAltroMapping*)maps->At(i);
359 if(strcmp(fgkRecoParamEmc->DecoderVersion(),"v1")==0)
360 dc=new AliPHOSRawDecoderv1(rawReader,mapping);
362 if(strcmp(fgkRecoParamEmc->DecoderVersion(),"v2")==0)
363 dc=new AliPHOSRawDecoderv2(rawReader,mapping);
365 dc=new AliPHOSRawDecoder(rawReader,mapping);
367 dc->SubtractPedestals(fgkRecoParamEmc->SubtractPedestals());
369 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
370 digits->SetName("DIGITS");
371 Int_t bufsize = 32000;
372 digitsTree->Branch("PHOS", &digits, bufsize);
374 AliPHOSRawDigiProducer pr(fgkRecoParamEmc,fgkRecoParamCpv);
375 pr.MakeDigits(digits,dc);
388 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
389 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
390 if(digit->GetEnergy()>eMax) {
391 fGeom->AbsToRelNumbering(digit->GetId(),relId);
392 eMax=digit->GetEnergy();
399 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
400 modMax,colMax,rowMax,eMax));