]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSReconstructor.cxx
Making online tracklets usable in offline reconstruction
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructor.cxx
CommitLineData
d15a28e7 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
b2a60966 16/* $Id$ */
17
d15a28e7 18//_________________________________________________________________________
a5e00016 19//--
20//-- Yves Schutz (SUBATECH)
dfe0be07 21// Reconstruction class. Redesigned from the old AliReconstructionner class and
22// derived from STEER/AliReconstructor.
23//
d15a28e7 24// --- ROOT system ---
00cfce1d 25
d15a28e7 26// --- Standard library ---
364de5c6 27
d15a28e7 28// --- AliRoot header files ---
9a2cdbdf 29#include "AliLog.h"
b3006690 30#include "AliAltroMapping.h"
af885e0f 31#include "AliESDEvent.h"
33ba95c3 32#include "AliESDCaloCluster.h"
a5e00016 33#include "AliESDCaloCells.h"
f444a19f 34#include "AliPHOSReconstructor.h"
7acf6008 35#include "AliPHOSClusterizerv1.h"
7acf6008 36#include "AliPHOSTrackSegmentMakerv1.h"
37#include "AliPHOSPIDv1.h"
23904d16 38#include "AliPHOSTracker.h"
d22dd3b4 39#include "AliRawReader.h"
d3aa2291 40#include "AliPHOSCalibData.h"
7e88424f 41#include "AliCDBEntry.h"
42#include "AliCDBManager.h"
64df000d 43#include "AliPHOSTrigger.h"
44#include "AliPHOSGeometry.h"
9a2cdbdf 45#include "AliPHOSDigit.h"
46#include "AliPHOSTrackSegment.h"
47#include "AliPHOSEmcRecPoint.h"
48#include "AliPHOSRecParticle.h"
49#include "AliPHOSRawDecoder.h"
77ea1c6f 50#include "AliPHOSRawDecoderv1.h"
710634f4 51#include "AliPHOSRawDecoderv2.h"
9a2cdbdf 52#include "AliPHOSRawDigiProducer.h"
53#include "AliPHOSPulseGenerator.h"
e957fea8 54
f444a19f 55ClassImp(AliPHOSReconstructor)
d15a28e7 56
2e60107f 57Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
6483babc 58TClonesArray* AliPHOSReconstructor::fgDigitsArray = 0; // Array of PHOS digits
59TObjArray* AliPHOSReconstructor::fgEMCRecPoints = 0; // Array of EMC rec.points
d3aa2291 60AliPHOSCalibData * AliPHOSReconstructor::fgCalibData = 0 ;
61
2e60107f 62
d15a28e7 63//____________________________________________________________________________
9a2cdbdf 64AliPHOSReconstructor::AliPHOSReconstructor() :
dcab1c7e 65 fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL)
d15a28e7 66{
b2a60966 67 // ctor
8d8258f6 68 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
dcab1c7e 69 fClusterizer = new AliPHOSClusterizerv1 (fGeom);
70 fTSM = new AliPHOSTrackSegmentMakerv1(fGeom);
71 fPID = new AliPHOSPIDv1 (fGeom);
6483babc 72 fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
73 fgEMCRecPoints= new TObjArray(100) ;
d3aa2291 74 if (!fgCalibData)
75 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
76
e68222ce 77}
6ad0bfa0 78
0379a13e 79//____________________________________________________________________________
80 AliPHOSReconstructor::~AliPHOSReconstructor()
81{
82 // dtor
9a2cdbdf 83 delete fGeom;
8d8258f6 84 delete fClusterizer;
dcab1c7e 85 delete fTSM;
86 delete fPID;
6483babc 87 delete fgDigitsArray;
88 delete fgEMCRecPoints;
0379a13e 89}
7acf6008 90
7acf6008 91//____________________________________________________________________________
9a2cdbdf 92void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
35293055 93{
9a2cdbdf 94 // 'single-event' local reco method called by AliReconstruction;
dfe0be07 95 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
af885e0f 96 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
dfe0be07 97 // the global tracking.
9a2cdbdf 98
7e88424f 99 fClusterizer->InitParameters();
8d8258f6 100 fClusterizer->SetInput(digitsTree);
101 fClusterizer->SetOutput(clustersTree);
a68156e6 102 if ( Debug() )
8d8258f6 103 fClusterizer->Digits2Clusters("deb all") ;
a68156e6 104 else
8d8258f6 105 fClusterizer->Digits2Clusters("") ;
a68156e6 106}
107
108//____________________________________________________________________________
9a2cdbdf 109void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
110 AliESDEvent* esd) const
a68156e6 111{
9a2cdbdf 112 // This method produces PHOS rec-particles,
113 // then it creates AliESDtracks out of them and
114 // write tracks to the ESD
115
e68222ce 116
9a2cdbdf 117 // do current event; the loop over events is done by AliReconstruction::Run()
dcab1c7e 118 fTSM->SetESD(esd) ;
119 fTSM->SetInput(clustersTree);
9a2cdbdf 120 if ( Debug() )
dcab1c7e 121 fTSM->Clusters2TrackSegments("deb all") ;
9a2cdbdf 122 else
dcab1c7e 123 fTSM->Clusters2TrackSegments("") ;
35293055 124
71994f35 125 fPID->SetEnergyCorrectionOn(GetRecoParam()->GetEMCEnergyCorrectionOn());
126
dcab1c7e 127 fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ;
128 fPID->SetESD(esd) ;
dfe0be07 129 if ( Debug() )
dcab1c7e 130 fPID->TrackSegments2RecParticles("deb all") ;
dfe0be07 131 else
dcab1c7e 132 fPID->TrackSegments2RecParticles("") ;
7acf6008 133
dcab1c7e 134 TClonesArray *recParticles = fPID->GetRecParticles();
135 Int_t nOfRecParticles = recParticles->GetEntriesFast();
a5fa6165 136
85c60a8e 137 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
dd7ee508 138 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
a5e00016 139
9a2cdbdf 140 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
a5e00016 141
142 // Read digits array
143
144 TBranch *branch = digitsTree->GetBranch("PHOS");
145 if (!branch) {
146 AliError("can't get the branch with the PHOS digits !");
147 return;
148 }
6483babc 149 branch->SetAddress(&fgDigitsArray);
a5e00016 150 branch->GetEntry(0);
151
152 // Get the clusters array
153
154 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
155 if (!emcbranch) {
156 AliError("can't get the branch with the PHOS EMC clusters !");
157 return;
158 }
159
6483babc 160 emcbranch->SetAddress(&fgEMCRecPoints);
a5e00016 161 emcbranch->GetEntry(0);
8013c27e 162
64df000d 163 //#########Calculate trigger and set trigger info###########
a5e00016 164
64df000d 165 AliPHOSTrigger tr ;
166 // tr.SetPatchSize(1);//create 4x4 patches
bfae5a5d 167 tr.SetSimulation(kFALSE);
6483babc 168 tr.Trigger(fgDigitsArray);
64df000d 169
170 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
171 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
172 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
173 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
174
64df000d 175 Int_t iSM2x2 = tr.Get2x2SuperModule();
176 Int_t iSMnxn = tr.GetnxnSuperModule();
177 Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
178 Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
179 Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
180 Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
181
a5e00016 182 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",
183 maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
184 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
185 maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
64df000d 186
bfae5a5d 187 // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
188 Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
64df000d 189 Int_t iAbsId2x2 =-1;
bfae5a5d 190 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
64df000d 191 Int_t iAbsIdnxn =-1;
192 TVector3 pos2x2(-1,-1,-1);
193 TVector3 posnxn(-1,-1,-1);
9a2cdbdf 194 fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
195 fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
196 fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
197 fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
64df000d 198
199 TArrayF triggerPosition(6);
200 triggerPosition[0] = pos2x2(0) ;
201 triggerPosition[1] = pos2x2(1) ;
202 triggerPosition[2] = pos2x2(2) ;
203 triggerPosition[3] = posnxn(0) ;
204 triggerPosition[4] = posnxn(1) ;
205 triggerPosition[5] = posnxn(2) ;
206
207 TArrayF triggerAmplitudes(4);
208 triggerAmplitudes[0] = maxAmp2x2 ;
209 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
210 triggerAmplitudes[2] = maxAmpnxn ;
211 triggerAmplitudes[3] = ampOutOfPatchnxn ;
212
213 //esd->SetPHOSTriggerCells(triggerPosition);
214 esd->AddPHOSTriggerPosition(triggerPosition);
215 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
216
e68222ce 217
a5e00016 218 //########################################
219 //############# Fill CaloCells ###########
220 //########################################
221
6483babc 222 Int_t nDigits = fgDigitsArray->GetEntries();
a5e00016 223 Int_t idignew = 0 ;
224 AliDebug(1,Form("%d digits",nDigits));
225
226 const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
227 AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
228 phsCells.CreateContainer(nDigits);
229 phsCells.SetType(AliESDCaloCells::kPHOSCell);
230
231 // Add to CaloCells only EMC digits with non-zero energy
232 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
6483babc 233 const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
12dd7f10 234 if(dig->GetId() <= knEMC &&
235 Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
d3aa2291 236 phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()), dig->GetTime());
a5e00016 237 idignew++;
238 }
e68222ce 239 }
a5e00016 240 phsCells.SetNumberOfCells(idignew);
241 phsCells.Sort();
e68222ce 242
a5e00016 243 //########################################
244 //############## Fill CaloClusters #######
245 //########################################
f4e0dd30 246
dfe0be07 247 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
a5e00016 248 AliPHOSRecParticle *rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
dfe0be07 249 if (Debug())
250 rp->Print();
25ed816e 251 // Get track segment and EMC rec.point associated with this rec.particle
dcab1c7e 252 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
a5e00016 253 ->At(rp->GetPHOSTSIndex()));
9a2cdbdf 254
6483babc 255 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
25ed816e 256 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
e68222ce 257
85c60a8e 258 Float_t xyz[3];
8013c27e 259 for (Int_t ixyz=0; ixyz<3; ixyz++)
260 xyz[ixyz] = rp->GetPos()[ixyz];
dd7ee508 261
262 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
77ea1c6f 263
a5e00016 264 // Create cell lists
265
266 Int_t cellMult = emcRP->GetDigitsMultiplicity();
267 Int_t *digitsList = emcRP->GetDigitsList();
268 Float_t *rpElist = emcRP->GetEnergiesList() ;
269 UShort_t *absIdList = new UShort_t[cellMult];
270 Double_t *fracList = new Double_t[cellMult];
271
272 for (Int_t iCell=0; iCell<cellMult; iCell++) {
6483babc 273 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
a5e00016 274 absIdList[iCell] = (UShort_t)(digit->GetId());
275 if (digit->GetEnergy() > 0)
d3aa2291 276 fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
a5e00016 277 else
278 fracList[iCell] = 0;
25ed816e 279 }
77ea1c6f 280
7592dfc4 281 //Primaries
282 Int_t primMult = 0;
4dd59c4a 283 Int_t *primList = emcRP->GetPrimaries(primMult);
a5e00016 284
48c5db5b 285 Float_t energy;
7e88424f 286 if (GetRecoParam()->EMCEcore2ESD())
48c5db5b 287 energy = emcRP->GetCoreEnergy();
288 else
289 energy = rp->Energy();
290
7592dfc4 291 // fills the ESDCaloCluster
8ada0ffe 292 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
a5e00016 293 ec->SetPosition(xyz); //rec.point position in MARS
48c5db5b 294 ec->SetE(energy); //total or core particle energy
25ed816e 295 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
e68222ce 296 ec->SetPid(rp->GetPID()) ; //array of particle identification
25ed816e 297 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
298 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
299 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
a5e00016 300 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
25ed816e 301 ec->SetClusterChi2(-1); //not yet implemented
78902954 302 ec->SetTOF(emcRP->GetTime()); //Time of flight
303
a5e00016 304 //Cells contributing to clusters
305 ec->SetNCells(cellMult);
306 ec->SetCellsAbsId(absIdList);
307 ec->SetCellsAmplitudeFraction(fracList);
7592dfc4 308
3744a6cf 309 //Distance to the nearest bad crystal
310 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
7592dfc4 311
312 //Array of MC indeces
4dd59c4a 313 TArrayI arrayPrim(primMult,primList);
7592dfc4 314 ec->AddLabels(arrayPrim);
64df000d 315
e44c41e9 316 //Matched ESD track
317 TArrayI arrayTrackMatched(1);
318 arrayTrackMatched[0]= ts->GetTrackIndex();
319 ec->AddTracksMatched(arrayTrackMatched);
320
85c60a8e 321 esd->AddCaloCluster(ec);
9a2cdbdf 322 delete ec;
091d150a 323 delete [] fracList;
324 delete [] absIdList;
e68222ce 325 }
6483babc 326 fgDigitsArray ->Delete();
327 fgEMCRecPoints->Delete();
328 recParticles ->Delete();
dd7ee508 329}
330
e68222ce 331//____________________________________________________________________________
d76c31f4 332AliTracker* AliPHOSReconstructor::CreateTracker() const
dd7ee508 333{
9a2cdbdf 334 // creates the PHOS tracker
335 return new AliPHOSTracker();
336}
dd7ee508 337
1267d56d 338//____________________________________________________________________________
9a2cdbdf 339void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
340{
341 // Converts raw data to
342 // PHOS digits
343 // Works on a single-event basis
9a2cdbdf 344 rawReader->Reset() ;
dd7ee508 345
77ea1c6f 346 AliPHOSRawDecoder * dc ;
347
7e88424f 348 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
b3006690 349 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
350
351 AliAltroMapping *mapping[4];
352 for(Int_t i = 0; i < 4; i++) {
353 mapping[i] = (AliAltroMapping*)maps->At(i);
354 }
355
7e88424f 356 if (strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0)
b3006690 357 dc=new AliPHOSRawDecoderv1(rawReader,mapping);
7e88424f 358 else if (strcmp(GetRecoParam()->EMCDecoderVersion(),"v2")==0)
359 dc=new AliPHOSRawDecoderv2(rawReader,mapping);
360 else
361 dc=new AliPHOSRawDecoder(rawReader,mapping);
77ea1c6f 362
7e88424f 363 dc->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
8be3a30a 364 dc->SetAmpOffset(GetRecoParam()->GetGlobalAltroOffset());
365 dc->SetAmpThreshold(GetRecoParam()->GetGlobalAltroThreshold());
366
9a2cdbdf 367 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
368 digits->SetName("DIGITS");
369 Int_t bufsize = 32000;
370 digitsTree->Branch("PHOS", &digits, bufsize);
dd7ee508 371
8be3a30a 372 //AliPHOSRawDigiProducer pr(GetRecoParam());
373 AliPHOSRawDigiProducer pr;
374
375 pr.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
376 pr.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
377 pr.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
77ea1c6f 378 pr.MakeDigits(digits,dc);
379
380 delete dc ;
7592dfc4 381
9a2cdbdf 382 //!!!!for debug!!!
a28333c5 383/*
9a2cdbdf 384 Int_t modMax=-111;
385 Int_t colMax=-111;
386 Int_t rowMax=-111;
387 Float_t eMax=-333;
388 //!!!for debug!!!
389
390 Int_t relId[4];
391 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
392 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
393 if(digit->GetEnergy()>eMax) {
394 fGeom->AbsToRelNumbering(digit->GetId(),relId);
395 eMax=digit->GetEnergy();
396 modMax=relId[0];
397 rowMax=relId[2];
398 colMax=relId[3];
dd7ee508 399 }
400 }
401
9a2cdbdf 402 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
403 modMax,colMax,rowMax,eMax));
a28333c5 404*/
9a2cdbdf 405 digitsTree->Fill();
e68222ce 406 digits->Delete();
407 delete digits;
d7d3b67b 408}
d3aa2291 409//==================================================================================
410Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
411 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
412
413 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
414
415 //Determine rel.position of the cell absolute ID
416 Int_t relId[4];
417 geom->AbsToRelNumbering(absId,relId);
418 Int_t module=relId[0];
419 Int_t row =relId[2];
420 Int_t column=relId[3];
421 if(relId[1]){ //CPV
422 Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
423 return amp*calibration ;
424 }
425 else{ //EMC
426 Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
427 return amp*calibration ;
428 }
429}
430
431