]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSReconstructor.cxx
Coding conventions
[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 ---
bb5c37a5 25#include "TGeoManager.h"
26#include "TGeoMatrix.h"
00cfce1d 27
d15a28e7 28// --- Standard library ---
364de5c6 29
d15a28e7 30// --- AliRoot header files ---
9a2cdbdf 31#include "AliLog.h"
b3006690 32#include "AliAltroMapping.h"
af885e0f 33#include "AliESDEvent.h"
33ba95c3 34#include "AliESDCaloCluster.h"
a5e00016 35#include "AliESDCaloCells.h"
f444a19f 36#include "AliPHOSReconstructor.h"
7acf6008 37#include "AliPHOSClusterizerv1.h"
7acf6008 38#include "AliPHOSTrackSegmentMakerv1.h"
39#include "AliPHOSPIDv1.h"
23904d16 40#include "AliPHOSTracker.h"
d22dd3b4 41#include "AliRawReader.h"
d3aa2291 42#include "AliPHOSCalibData.h"
7e88424f 43#include "AliCDBEntry.h"
44#include "AliCDBManager.h"
64df000d 45#include "AliPHOSTrigger.h"
46#include "AliPHOSGeometry.h"
9a2cdbdf 47#include "AliPHOSDigit.h"
48#include "AliPHOSTrackSegment.h"
49#include "AliPHOSEmcRecPoint.h"
50#include "AliPHOSRecParticle.h"
379c5c09 51#include "AliPHOSRawFitterv0.h"
52#include "AliPHOSRawFitterv1.h"
53#include "AliPHOSRawFitterv2.h"
f4a5c5fd 54#include "AliPHOSRawFitterv3.h"
de0cfd46 55#include "AliPHOSRawFitterv4.h"
9a2cdbdf 56#include "AliPHOSRawDigiProducer.h"
57#include "AliPHOSPulseGenerator.h"
e1aec4f9 58#include "AliPHOSTriggerRawDigit.h"
59#include "AliPHOSTriggerRawDigiProducer.h"
60#include "AliPHOSTriggerParameters.h"
e957fea8 61
f444a19f 62ClassImp(AliPHOSReconstructor)
d15a28e7 63
2e60107f 64Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
6483babc 65TClonesArray* AliPHOSReconstructor::fgDigitsArray = 0; // Array of PHOS digits
66TObjArray* AliPHOSReconstructor::fgEMCRecPoints = 0; // Array of EMC rec.points
d3aa2291 67AliPHOSCalibData * AliPHOSReconstructor::fgCalibData = 0 ;
e1aec4f9 68TClonesArray* AliPHOSReconstructor::fgTriggerDigits = 0; // Array of PHOS trigger digits
2e60107f 69
d15a28e7 70//____________________________________________________________________________
9a2cdbdf 71AliPHOSReconstructor::AliPHOSReconstructor() :
771123c7 72 fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL),fTmpDigLG(NULL)
d15a28e7 73{
b2a60966 74 // ctor
771123c7 75 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
76 fClusterizer = new AliPHOSClusterizerv1 (fGeom);
77 fTSM = new AliPHOSTrackSegmentMakerv1(fGeom);
78 fPID = new AliPHOSPIDv1 (fGeom);
79 fTmpDigLG = new TClonesArray("AliPHOSDigit",100);
80 fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
81 fgEMCRecPoints = new TObjArray(100) ;
d3aa2291 82 if (!fgCalibData)
83 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
e1aec4f9 84
85 fgTriggerDigits = new TClonesArray("AliPHOSTriggerRawDigit",100);
86
1d76a1f4 87 AliInfo(Form("PHOS bad channel map contains %d bad channel(s).\n",
88 fgCalibData->GetNumOfEmcBadChannels()));
d3aa2291 89
e68222ce 90}
6ad0bfa0 91
0379a13e 92//____________________________________________________________________________
771123c7 93AliPHOSReconstructor::~AliPHOSReconstructor()
0379a13e 94{
95 // dtor
9a2cdbdf 96 delete fGeom;
8d8258f6 97 delete fClusterizer;
dcab1c7e 98 delete fTSM;
99 delete fPID;
771123c7 100 delete fTmpDigLG;
6483babc 101 delete fgDigitsArray;
102 delete fgEMCRecPoints;
e1aec4f9 103 delete fgTriggerDigits;
0379a13e 104}
7acf6008 105
7acf6008 106//____________________________________________________________________________
9a2cdbdf 107void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
35293055 108{
9a2cdbdf 109 // 'single-event' local reco method called by AliReconstruction;
dfe0be07 110 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
af885e0f 111 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
dfe0be07 112 // the global tracking.
9a2cdbdf 113
7e88424f 114 fClusterizer->InitParameters();
8d8258f6 115 fClusterizer->SetInput(digitsTree);
116 fClusterizer->SetOutput(clustersTree);
a68156e6 117 if ( Debug() )
8d8258f6 118 fClusterizer->Digits2Clusters("deb all") ;
a68156e6 119 else
8d8258f6 120 fClusterizer->Digits2Clusters("") ;
a68156e6 121}
122
123//____________________________________________________________________________
9a2cdbdf 124void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
125 AliESDEvent* esd) const
a68156e6 126{
9a2cdbdf 127 // This method produces PHOS rec-particles,
128 // then it creates AliESDtracks out of them and
129 // write tracks to the ESD
130
e68222ce 131
9a2cdbdf 132 // do current event; the loop over events is done by AliReconstruction::Run()
dcab1c7e 133 fTSM->SetESD(esd) ;
134 fTSM->SetInput(clustersTree);
9a2cdbdf 135 if ( Debug() )
dcab1c7e 136 fTSM->Clusters2TrackSegments("deb all") ;
9a2cdbdf 137 else
dcab1c7e 138 fTSM->Clusters2TrackSegments("") ;
35293055 139
dcab1c7e 140 fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ;
141 fPID->SetESD(esd) ;
dfe0be07 142 if ( Debug() )
dcab1c7e 143 fPID->TrackSegments2RecParticles("deb all") ;
dfe0be07 144 else
dcab1c7e 145 fPID->TrackSegments2RecParticles("") ;
7acf6008 146
dcab1c7e 147 TClonesArray *recParticles = fPID->GetRecParticles();
148 Int_t nOfRecParticles = recParticles->GetEntriesFast();
a5fa6165 149
9a2cdbdf 150 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
a5e00016 151
152 // Read digits array
153
154 TBranch *branch = digitsTree->GetBranch("PHOS");
155 if (!branch) {
156 AliError("can't get the branch with the PHOS digits !");
157 return;
158 }
6483babc 159 branch->SetAddress(&fgDigitsArray);
a5e00016 160 branch->GetEntry(0);
161
162 // Get the clusters array
163
164 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
165 if (!emcbranch) {
166 AliError("can't get the branch with the PHOS EMC clusters !");
167 return;
168 }
169
6483babc 170 emcbranch->SetAddress(&fgEMCRecPoints);
a5e00016 171 emcbranch->GetEntry(0);
26d0974e 172
173
e1aec4f9 174 // Trigger
e1aec4f9 175
26d0974e 176 TBranch *tbranch = digitsTree->GetBranch("TPHOS");
177 if (tbranch) {
178
179 tbranch->SetAddress(&fgTriggerDigits);
180 tbranch->GetEntry(0);
181
182 AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("PHOS");
183
184 if (trgESD) {
185 trgESD->Allocate(fgTriggerDigits->GetEntriesFast());
e1aec4f9 186
26d0974e 187 for (Int_t i = 0; i < fgTriggerDigits->GetEntriesFast(); i++) {
188 AliPHOSTriggerRawDigit* tdig = (AliPHOSTriggerRawDigit*)fgTriggerDigits->At(i);
189
190 Int_t mod,modX,modZ;
191 tdig->GetModXZ(mod,modX,modZ);
192
193 const Int_t relId[4] = {5-mod,0,modX+1,modZ+1};
194 Int_t absId;
195
196 fGeom->RelToAbsNumbering(relId,absId);
197 trgESD->Add(mod,absId,tdig->GetAmp(),0.,(Int_t*)NULL,0,0,0);
198 }
199 }
e1aec4f9 200 }
e1aec4f9 201
fd3fd391 202// //#########Calculate trigger and set trigger info###########
a5e00016 203
fd3fd391 204// AliPHOSTrigger tr ;
205// // tr.SetPatchSize(1);//create 4x4 patches
206// tr.SetSimulation(kFALSE);
207// tr.Trigger(fgDigitsArray);
64df000d 208
fd3fd391 209// Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
210// Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
211// Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
212// Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
213
214// Int_t iSM2x2 = tr.Get2x2SuperModule();
215// Int_t iSMnxn = tr.GetnxnSuperModule();
216// Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
217// Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
218// Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
219// Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
220
221// AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",
222// maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
223// AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
224// maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
225
226// // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
227// Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
228// Int_t iAbsId2x2 =-1;
229// Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
230// Int_t iAbsIdnxn =-1;
231// TVector3 pos2x2(-1,-1,-1);
232// TVector3 posnxn(-1,-1,-1);
233// fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
234// fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
235// fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
236// fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
237
238// TArrayF triggerPosition(6);
239// triggerPosition[0] = pos2x2(0) ;
240// triggerPosition[1] = pos2x2(1) ;
241// triggerPosition[2] = pos2x2(2) ;
242// triggerPosition[3] = posnxn(0) ;
243// triggerPosition[4] = posnxn(1) ;
244// triggerPosition[5] = posnxn(2) ;
245
246// TArrayF triggerAmplitudes(4);
247// triggerAmplitudes[0] = maxAmp2x2 ;
248// triggerAmplitudes[1] = ampOutOfPatch2x2 ;
249// triggerAmplitudes[2] = maxAmpnxn ;
250// triggerAmplitudes[3] = ampOutOfPatchnxn ;
251
252// //esd->SetPHOSTriggerCells(triggerPosition);
253// esd->AddPHOSTriggerPosition(triggerPosition);
254// esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
64df000d 255
e68222ce 256
a5e00016 257 //########################################
258 //############# Fill CaloCells ###########
259 //########################################
260
6483babc 261 Int_t nDigits = fgDigitsArray->GetEntries();
a5e00016 262 Int_t idignew = 0 ;
263 AliDebug(1,Form("%d digits",nDigits));
264
265 const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
266 AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
267 phsCells.CreateContainer(nDigits);
268 phsCells.SetType(AliESDCaloCells::kPHOSCell);
269
270 // Add to CaloCells only EMC digits with non-zero energy
271 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
6483babc 272 const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
12dd7f10 273 if(dig->GetId() <= knEMC &&
274 Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
6f47f50d 275 phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
276 CalibrateT(dig->GetTime(),dig->GetId()));
a5e00016 277 idignew++;
278 }
e68222ce 279 }
a5e00016 280 phsCells.SetNumberOfCells(idignew);
281 phsCells.Sort();
e68222ce 282
a5e00016 283 //########################################
284 //############## Fill CaloClusters #######
285 //########################################
f4e0dd30 286
dfe0be07 287 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
09828f52 288 AliPHOSRecParticle *rp = static_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
dfe0be07 289 if (Debug())
290 rp->Print();
25ed816e 291 // Get track segment and EMC rec.point associated with this rec.particle
dcab1c7e 292 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
a5e00016 293 ->At(rp->GetPHOSTSIndex()));
9a2cdbdf 294
6483babc 295 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
25ed816e 296 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
e68222ce 297
85c60a8e 298 Float_t xyz[3];
8013c27e 299 for (Int_t ixyz=0; ixyz<3; ixyz++)
300 xyz[ixyz] = rp->GetPos()[ixyz];
dd7ee508 301
302 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
77ea1c6f 303
a5e00016 304 // Create cell lists
305
306 Int_t cellMult = emcRP->GetDigitsMultiplicity();
307 Int_t *digitsList = emcRP->GetDigitsList();
308 Float_t *rpElist = emcRP->GetEnergiesList() ;
309 UShort_t *absIdList = new UShort_t[cellMult];
310 Double_t *fracList = new Double_t[cellMult];
311
312 for (Int_t iCell=0; iCell<cellMult; iCell++) {
6483babc 313 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
a5e00016 314 absIdList[iCell] = (UShort_t)(digit->GetId());
315 if (digit->GetEnergy() > 0)
d3aa2291 316 fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
a5e00016 317 else
318 fracList[iCell] = 0;
25ed816e 319 }
77ea1c6f 320
7592dfc4 321 //Primaries
322 Int_t primMult = 0;
4dd59c4a 323 Int_t *primList = emcRP->GetPrimaries(primMult);
a5e00016 324
25312a8e 325 Float_t energy=0.;
7e88424f 326 if (GetRecoParam()->EMCEcore2ESD())
48c5db5b 327 energy = emcRP->GetCoreEnergy();
328 else
329 energy = rp->Energy();
25312a8e 330 //Apply nonlinearity correction
331 if(GetRecoParam()->GetEMCEnergyCorrectionOn())
332 energy=CorrectNonlinearity(energy) ;
48c5db5b 333
7592dfc4 334 // fills the ESDCaloCluster
c8fe2783 335 ec->SetType(AliVCluster::kPHOSNeutral);
a5e00016 336 ec->SetPosition(xyz); //rec.point position in MARS
48c5db5b 337 ec->SetE(energy); //total or core particle energy
c8fe2783 338 ec->SetDispersion(emcRP->GetDispersion()); //cluster dispersion
339 ec->SetPID(rp->GetPID()) ; //array of particle identification
25ed816e 340 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
341 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
342 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
a5e00016 343 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
95bf21ad 344 ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z"));
c8fe2783 345 ec->SetChi2(-1); //not yet implemented
6f47f50d 346 ec->SetTOF(emcRP->GetTime()); //Time of flight - already calibrated in EMCRecPoint
78902954 347
a5e00016 348 //Cells contributing to clusters
349 ec->SetNCells(cellMult);
350 ec->SetCellsAbsId(absIdList);
351 ec->SetCellsAmplitudeFraction(fracList);
7592dfc4 352
3744a6cf 353 //Distance to the nearest bad crystal
354 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
7592dfc4 355
356 //Array of MC indeces
4dd59c4a 357 TArrayI arrayPrim(primMult,primList);
7592dfc4 358 ec->AddLabels(arrayPrim);
64df000d 359
e44c41e9 360 //Matched ESD track
361 TArrayI arrayTrackMatched(1);
362 arrayTrackMatched[0]= ts->GetTrackIndex();
363 ec->AddTracksMatched(arrayTrackMatched);
364
95bf21ad 365 Int_t index = esd->AddCaloCluster(ec);
366
367 //Set pointer to this cluster in ESD track
368 Int_t nt=esd->GetNumberOfTracks();
369 for (Int_t itr=0; itr<nt; itr++) {
370 AliESDtrack *esdTrack=esd->GetTrack(itr);
371 if(!esdTrack->IsPHOS())
372 continue ;
373 if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
374 esdTrack->SetPHOScluster(index) ;
375//no garatie that only one track matched this cluster
376// break ;
377 }
378 }
379
9a2cdbdf 380 delete ec;
091d150a 381 delete [] fracList;
382 delete [] absIdList;
e68222ce 383 }
771123c7 384 fgDigitsArray ->Clear();
385 fgEMCRecPoints->Clear("C");
386 recParticles ->Clear();
bb5c37a5 387
388 //Store PHOS misalignment matrixes
389 FillMisalMatrixes(esd) ;
390
dd7ee508 391}
392
e68222ce 393//____________________________________________________________________________
d76c31f4 394AliTracker* AliPHOSReconstructor::CreateTracker() const
dd7ee508 395{
9a2cdbdf 396 // creates the PHOS tracker
397 return new AliPHOSTracker();
398}
dd7ee508 399
1267d56d 400//____________________________________________________________________________
9a2cdbdf 401void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
402{
403 // Converts raw data to
404 // PHOS digits
405 // Works on a single-event basis
9a2cdbdf 406 rawReader->Reset() ;
dd7ee508 407
379c5c09 408 AliPHOSRawFitterv0 * fitter ;
77ea1c6f 409
7e88424f 410 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
b3006690 411 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
412
47ba0d5d 413 AliAltroMapping *mapping[20];
414 for(Int_t i = 0; i < 20; i++) {
b3006690 415 mapping[i] = (AliAltroMapping*)maps->At(i);
416 }
417
adc447ab 418 if (strcmp(GetRecoParam()->EMCFitterVersion(),"v0")==0)
419 fitter=new AliPHOSRawFitterv0();
420 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
379c5c09 421 fitter=new AliPHOSRawFitterv1();
422 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
423 fitter=new AliPHOSRawFitterv2();
f4a5c5fd 424 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0)
425 fitter=new AliPHOSRawFitterv3();
7e88424f 426 else
de0cfd46 427 fitter=new AliPHOSRawFitterv4();
77ea1c6f 428
379c5c09 429 fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
430 fitter->SetAmpOffset (GetRecoParam()->GetGlobalAltroOffset());
431 fitter->SetAmpThreshold (GetRecoParam()->GetGlobalAltroThreshold());
8be3a30a 432
9a2cdbdf 433 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
434 digits->SetName("DIGITS");
435 Int_t bufsize = 32000;
436 digitsTree->Branch("PHOS", &digits, bufsize);
dd7ee508 437
379c5c09 438 AliPHOSRawDigiProducer rdp(rawReader,mapping);
8be3a30a 439
379c5c09 440 rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
441 rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
442 rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
771123c7 443 rdp.MakeDigits(digits,fTmpDigLG,fitter);
77ea1c6f 444
379c5c09 445 delete fitter ;
7592dfc4 446
e1aec4f9 447 TClonesArray *tdigits = new TClonesArray("AliPHOSTriggerRawDigit",1);
448 tdigits->SetName("TDIGITS");
449 digitsTree->Branch("TPHOS", &tdigits, bufsize);
450
451 rawReader->Reset();
452 AliPHOSTriggerRawDigiProducer tdp(rawReader);
892b88d6 453
454 AliPHOSTriggerParameters* parameters = (AliPHOSTriggerParameters*)AliPHOSRecoParam::GetTriggerParameters();
e1aec4f9 455
456 tdp.SetTriggerParameters(parameters);
457 tdp.ProcessEvent(tdigits);
458
379c5c09 459 if (AliLog::GetGlobalDebugLevel() == 1) {
460 Int_t modMax=-111;
461 Int_t colMax=-111;
462 Int_t rowMax=-111;
463 Float_t eMax=-333;
464 //!!!for debug!!!
465
466 Int_t relId[4];
467 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
468 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
469 if(digit->GetEnergy()>eMax) {
470 fGeom->AbsToRelNumbering(digit->GetId(),relId);
471 eMax=digit->GetEnergy();
472 modMax=relId[0];
473 rowMax=relId[2];
474 colMax=relId[3];
475 }
dd7ee508 476 }
379c5c09 477
478 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
479 modMax,colMax,rowMax,eMax));
dd7ee508 480 }
e1aec4f9 481
9a2cdbdf 482 digitsTree->Fill();
e1aec4f9 483
e68222ce 484 digits->Delete();
485 delete digits;
e1aec4f9 486
487 tdigits->Delete();
488 delete tdigits;
d7d3b67b 489}
d3aa2291 490//==================================================================================
491Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
492 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
493
494 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
495
496 //Determine rel.position of the cell absolute ID
497 Int_t relId[4];
498 geom->AbsToRelNumbering(absId,relId);
499 Int_t module=relId[0];
500 Int_t row =relId[2];
501 Int_t column=relId[3];
502 if(relId[1]){ //CPV
503 Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
504 return amp*calibration ;
505 }
506 else{ //EMC
507 Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
508 return amp*calibration ;
509 }
510}
bb5c37a5 511//==================================================================================
6f47f50d 512Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId)const{
513 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
514
515 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
516
517 //Determine rel.position of the cell absolute ID
518 Int_t relId[4];
519 geom->AbsToRelNumbering(absId,relId);
520 Int_t module=relId[0];
521 Int_t row =relId[2];
522 Int_t column=relId[3];
523 if(relId[1]){ //CPV
524 return 0. ;
525 }
526 else{ //EMC
527 time += fgCalibData->GetTimeShiftEmc(module,column,row);
528 return time ;
529 }
530}
531//==================================================================================
bb5c37a5 532void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
533 //Store PHOS matrixes in ESD Header
534
535 //Check, if matrixes was already stored
536 for(Int_t mod=0 ;mod<5; mod++){
537 if(esd->GetPHOSMatrix(mod)!=0)
538 return ;
539 }
540
541 //Create and store matrixes
542 if(!gGeoManager){
543 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
544 return ;
545 }
546 //Note, that owner of copied marixes will be header
547 char path[255] ;
548 TGeoHMatrix * m ;
549 for(Int_t mod=0; mod<5; mod++){
3da0f212 550 snprintf(path,255,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
bb5c37a5 551 if (gGeoManager->cd(path)){
552 m = gGeoManager->GetCurrentMatrix() ;
553 esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
554 }
555 else{
556 esd->SetPHOSMatrix(NULL,mod) ;
557 }
558 }
559
560}
25312a8e 561//==================================================================================
562Float_t AliPHOSReconstructor::CorrectNonlinearity(Float_t en){
d3aa2291 563
81046c73 564 //For backward compatibility, if no RecoParameters found
565 if(!GetRecoParam()){
566 return 0.0241+1.0504*en+0.000249*en*en ;
567 }
568
25312a8e 569 if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"NoCorrection")==0){
570 return en ;
571 }
572 if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Gustavo2005")==0){
573 const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
574 return par[0]+par[1]*en + par[2]*en*en ;
575 }
576 if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Henrik2010")==0){
577 const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
578 return en*(par[0]+par[1]*TMath::Exp(-en*par[2]))*(1.+par[3]*TMath::Exp(-en*par[4]))*(1.+par[6]/(en*en+par[5])) ;
579 }
580 //For backward compatibility
581 if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"")==0){
582 return 0.0241+1.0504*en+0.000249*en*en ;
583 }
584 return en ;
585}
d3aa2291 586
e1aec4f9 587void AliPHOSReconstructor::readTRUParameters(AliPHOSTriggerParameters* parameters) const
588{
589 //Read trigger parameters.
590
591 TString path(gSystem->Getenv("ALICE_ROOT"));
592 path += "/PHOS/macros/Trigger/OCDB/";
593
594 for (Int_t mod = 2; mod < 5; ++mod) { // module
595 for (Int_t tru = 0; tru < 4; tru++) { // tru row
596 for (Int_t branch = 0; branch < 2; branch++) { // branch
597
598 // Open the Appropriate pedestal file
599 TString fileName = path;
600 fileName += "pedestal_m";
601 fileName = fileName += mod;
602 fileName+="_r";
603 fileName+=tru;
604 fileName+="_b";
605 fileName+=branch;
606 fileName+=".dat";
607 std::ifstream instream;
608 instream.open(fileName.Data());
609
610 // Read pedestals from file
611 if( ! instream.is_open() )
612 Printf("E-TRUPedestals: could not open %s", fileName.Data());
613 else
614 {
615 Int_t ped[112];
616
04f0835e 617 char ch_s[36];
e1aec4f9 618 char *ch_s_p = ch_s;
619 //Int_t nlines = 0;
620
621 Int_t t_ped_0 =0;
622 Int_t t_ped_1 =0;
623 Int_t t_ped_2 =0;
624
625 for(Int_t n=0; n<112; n++)
626 {
627 instream.getline(ch_s_p,36);
628 if (ch_s_p[23]<='9' && ch_s_p[23]>='0')
629 {
630 t_ped_0 = ch_s_p[23]-'0';
631 }
632 else if (ch_s_p[23]>='A' && ch_s_p[23]<='Z')
633 {
634 t_ped_0 = ch_s_p[23]-'A'+10;
635
636 }
637
638 if (ch_s_p[22]<='9' && ch_s_p[22]>='0')
639 {
640 t_ped_1 = ch_s_p[22]-'0';
641 }
642 else if (ch_s_p[22]<='Z' && ch_s_p[22]>='A')
643 {
644 t_ped_1 = ch_s_p[22]-'A'+10;
645 }
646
647 if (ch_s_p[21]<='9' && ch_s_p[21]>='0')
648 {
649 t_ped_2 = ch_s_p[21]-'0';
650 }
651 else if (ch_s_p[21]<='Z' && ch_s_p[21]>='A')
652 {
653 t_ped_2 = ch_s_p[21]-'A'+10;
654 }
655
656 ped[n]=t_ped_2*256+t_ped_1*16+t_ped_0;
657
658
659 }
660 for (Int_t xrow = 0; xrow < 8; xrow++){
661 for (Int_t zcol = 0; zcol < 14; zcol++){
662 Int_t pedestal = ped[zcol*8+xrow];
663
664 if( pedestal < 612 && pedestal > 412 ) // resonable
665 parameters->SetTRUPedestal(pedestal, mod, tru, branch, xrow, zcol);
666 else // unresonable
667 continue;
668 }
669 }
670 } // Ends read of pedestals from branch from file.
671 instream.close();
672 }// end branch
673 }// end tru
674
675 }// end for mod
676}
677
678
679
680