]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
LHAPDF veraion 5.9.1
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructor.cxx
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
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //--
20 //-- Yves Schutz (SUBATECH) 
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and 
22 // derived from STEER/AliReconstructor. 
23 // 
24 // --- ROOT system ---
25 #include "TGeoManager.h"
26 #include "TGeoMatrix.h"
27
28 // --- Standard library ---
29
30 // --- AliRoot header files ---
31 #include "AliLog.h"
32 #include "AliAltroMapping.h"
33 #include "AliESDEvent.h"
34 #include "AliESDCaloCluster.h"
35 #include "AliESDCaloCells.h"
36 #include "AliPHOSReconstructor.h"
37 #include "AliPHOSClusterizerv1.h"
38 #include "AliPHOSTrackSegmentMakerv1.h"
39 #include "AliPHOSPIDv1.h"
40 #include "AliPHOSTracker.h"
41 #include "AliRawReader.h"
42 #include "AliPHOSCalibData.h"
43 #include "AliCDBEntry.h"
44 #include "AliCDBManager.h"
45 #include "AliPHOSTrigger.h"
46 #include "AliPHOSGeometry.h"
47 #include "AliPHOSDigit.h"
48 #include "AliPHOSTrackSegment.h"
49 #include "AliPHOSEmcRecPoint.h"
50 #include "AliPHOSRecParticle.h"
51 #include "AliPHOSRawFitterv0.h"
52 #include "AliPHOSRawFitterv1.h"
53 #include "AliPHOSRawFitterv2.h"
54 #include "AliPHOSRawFitterv3.h"
55 #include "AliPHOSRawFitterv4.h"
56 #include "AliPHOSRawDigiProducer.h"
57 #include "AliPHOSPulseGenerator.h"
58 #include "AliPHOSTriggerRawDigit.h"
59 #include "AliPHOSTriggerRawDigiProducer.h"
60 #include "AliPHOSTriggerParameters.h"
61
62 ClassImp(AliPHOSReconstructor)
63
64 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ; 
65 TClonesArray*     AliPHOSReconstructor::fgDigitsArray = 0;   // Array of PHOS digits
66 TObjArray*        AliPHOSReconstructor::fgEMCRecPoints = 0;   // Array of EMC rec.points
67 AliPHOSCalibData * AliPHOSReconstructor::fgCalibData  = 0 ;
68 TClonesArray*     AliPHOSReconstructor::fgTriggerDigits = 0;   // Array of PHOS trigger digits
69
70 //____________________________________________________________________________
71 AliPHOSReconstructor::AliPHOSReconstructor() :
72   fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL),fTmpDigLG(NULL)
73 {
74   // ctor
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) ;
82   if (!fgCalibData)
83     fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
84   
85   fgTriggerDigits  = new TClonesArray("AliPHOSTriggerRawDigit",100);
86   
87   AliInfo(Form("PHOS bad channel map contains %d bad channel(s).\n",
88                fgCalibData->GetNumOfEmcBadChannels()));
89  
90 }
91
92 //____________________________________________________________________________
93 AliPHOSReconstructor::~AliPHOSReconstructor()
94 {
95   // dtor
96   //  delete fGeom; // RS: fGeom is a static object owned by AliPHOSGeometry::fGeom, don't delete
97   delete fClusterizer;
98   delete fTSM;
99   delete fPID;
100   delete fTmpDigLG;
101   delete fgDigitsArray;
102   delete fgEMCRecPoints;
103   delete fgTriggerDigits;
104
105
106 //____________________________________________________________________________
107 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
108 {
109   // 'single-event' local reco method called by AliReconstruction; 
110   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
111   // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by 
112   // the global tracking.
113
114   fClusterizer->InitParameters();
115   fClusterizer->SetInput(digitsTree);
116   fClusterizer->SetOutput(clustersTree);
117   if ( Debug() ) 
118     fClusterizer->Digits2Clusters("deb all") ; 
119   else 
120     fClusterizer->Digits2Clusters("") ;
121 }
122
123 //____________________________________________________________________________
124 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
125                                    AliESDEvent* esd) const
126 {
127   // This method produces PHOS rec-particles,
128   // then it creates AliESDtracks out of them and
129   // write tracks to the ESD
130
131
132   // do current event; the loop over events is done by AliReconstruction::Run()
133   fTSM->SetESD(esd) ; 
134   fTSM->SetInput(clustersTree);
135   if ( Debug() ) 
136     fTSM->Clusters2TrackSegments("deb all") ;
137   else 
138     fTSM->Clusters2TrackSegments("") ;
139   
140   fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ; 
141   fPID->SetESD(esd) ; 
142   if ( Debug() ) 
143     fPID->TrackSegments2RecParticles("deb all") ;
144   else 
145     fPID->TrackSegments2RecParticles("") ;
146
147   TClonesArray *recParticles  = fPID->GetRecParticles();
148   Int_t nOfRecParticles = recParticles->GetEntriesFast();
149   
150   AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
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   }
159   branch->SetAddress(&fgDigitsArray);
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
170   emcbranch->SetAddress(&fgEMCRecPoints);
171   emcbranch->GetEntry(0);
172   
173   
174   // Trigger
175   
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());
186       
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     }  
200   }
201   
202 //   //#########Calculate trigger and set trigger info###########
203
204 //   AliPHOSTrigger tr ;
205 //   //   tr.SetPatchSize(1);//create 4x4 patches
206 //   tr.SetSimulation(kFALSE);
207 //   tr.Trigger(fgDigitsArray);
208   
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);
255   
256
257   //########################################
258   //############# Fill CaloCells ###########
259   //########################################
260
261   Int_t nDigits = fgDigitsArray->GetEntries();
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++) {
272     const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
273     if(dig->GetId() <= knEMC && 
274        Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
275       Int_t primary = dig->GetPrimary(1) ;
276       phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
277                                              CalibrateT(dig->GetTime(),dig->GetId(),dig->IsLG()),
278                                              primary,0.,!dig->IsLG()) ;
279                                              
280       idignew++;
281     }
282   }
283   phsCells.SetNumberOfCells(idignew);
284   phsCells.Sort();
285
286   //########################################
287   //############## Fill CaloClusters #######
288   //########################################
289
290   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
291     AliPHOSRecParticle  *rp    = static_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
292     if (Debug()) 
293       rp->Print();
294     // Get track segment and EMC rec.point associated with this rec.particle
295     AliPHOSTrackSegment *ts    = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
296                                                                     ->At(rp->GetPHOSTSIndex()));
297
298     AliPHOSEmcRecPoint  *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
299     AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
300     
301     Float_t xyz[3];
302     for (Int_t ixyz=0; ixyz<3; ixyz++) 
303       xyz[ixyz] = rp->GetPos()[ixyz];
304     
305     AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
306    
307     // Create cell lists
308
309     Int_t     cellMult   = emcRP->GetDigitsMultiplicity();
310     Int_t    *digitsList = emcRP->GetDigitsList();
311     Float_t  *rpElist    = emcRP->GetEnergiesList() ;
312     UShort_t *absIdList  = new UShort_t[cellMult];
313     Double_t *fracList   = new Double_t[cellMult];
314
315     for (Int_t iCell=0; iCell<cellMult; iCell++) {
316       AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
317       absIdList[iCell] = (UShort_t)(digit->GetId());
318       if (digit->GetEnergy() > 0)
319         fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
320       else
321         fracList[iCell] = 0;
322     }
323
324     //Primaries
325     Int_t  primMult  = 0;
326     Int_t *primList =  emcRP->GetPrimaries(primMult);
327
328     Float_t energy=0.;
329     if (GetRecoParam()->EMCEcore2ESD())
330       energy = emcRP->GetCoreEnergy();
331     else
332       energy = rp->Energy();
333     //Apply nonlinearity correction
334     if(GetRecoParam()->GetEMCEnergyCorrectionOn())
335       energy=CorrectNonlinearity(energy) ;
336
337     // fills the ESDCaloCluster
338     ec->SetType(AliVCluster::kPHOSNeutral);
339     ec->SetPosition(xyz);                       //rec.point position in MARS
340     ec->SetE(energy);                           //total or core particle energy
341     ec->SetDispersion(emcRP->GetDispersion());  //cluster dispersion
342     ec->SetPID(rp->GetPID()) ;            //array of particle identification
343     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
344     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
345     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
346     ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
347     ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z")); 
348     ec->SetChi2(-1);                     //not yet implemented
349     ec->SetTOF(emcRP->GetTime());               //Time of flight - already calibrated in EMCRecPoint
350
351     //Cells contributing to clusters
352     ec->SetNCells(cellMult);
353     ec->SetCellsAbsId(absIdList);
354     ec->SetCellsAmplitudeFraction(fracList);
355
356     //Distance to the nearest bad crystal
357     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
358   
359     //Array of MC indeces
360     TArrayI arrayPrim(primMult,primList);
361     ec->AddLabels(arrayPrim);
362     
363     //Matched ESD track
364     TArrayI arrayTrackMatched(1);
365     arrayTrackMatched[0]= ts->GetTrackIndex();
366     ec->AddTracksMatched(arrayTrackMatched);
367     
368     Int_t index = esd->AddCaloCluster(ec);
369
370     //Set pointer to this cluster in ESD track
371     Int_t nt=esd->GetNumberOfTracks();
372     for (Int_t itr=0; itr<nt; itr++) {
373       AliESDtrack *esdTrack=esd->GetTrack(itr);
374       if(!esdTrack->IsPHOS())
375         continue ;
376       if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
377         esdTrack->SetPHOScluster(index) ;
378 //no garatie that only one track matched this cluster
379 //      break ;
380       }
381     }
382  
383     delete ec;   
384     delete [] fracList;
385     delete [] absIdList;
386   }
387   fgDigitsArray ->Clear("C");
388   fgEMCRecPoints->Clear("C");
389   recParticles  ->Clear();
390
391   //Store PHOS misalignment matrixes
392   FillMisalMatrixes(esd) ;
393
394 }
395
396 //____________________________________________________________________________
397 AliTracker* AliPHOSReconstructor::CreateTracker() const
398 {
399   // creates the PHOS tracker
400   return new AliPHOSTracker();
401 }
402
403 //____________________________________________________________________________
404 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
405 {
406   // Converts raw data to
407   // PHOS digits
408   // Works on a single-event basis
409   rawReader->Reset() ; 
410
411   AliPHOSRawFitterv0 * fitter ;
412
413   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
414   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
415
416   AliAltroMapping *mapping[20];
417   for(Int_t i = 0; i < 20; i++) {
418     mapping[i] = (AliAltroMapping*)maps->At(i);
419   }
420
421   if      (strcmp(GetRecoParam()->EMCFitterVersion(),"v0")==0) 
422     fitter=new AliPHOSRawFitterv0();
423   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0) 
424     fitter=new AliPHOSRawFitterv1();
425   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0) 
426     fitter=new AliPHOSRawFitterv2();
427   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0) 
428     fitter=new AliPHOSRawFitterv3();
429   else
430     fitter=new AliPHOSRawFitterv4();
431
432   fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
433   fitter->SetAmpOffset     (GetRecoParam()->GetGlobalAltroOffset());
434   fitter->SetAmpThreshold  (GetRecoParam()->GetGlobalAltroThreshold());
435
436   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
437   digits->SetName("DIGITS");
438   Int_t bufsize = 32000;
439   digitsTree->Branch("PHOS", &digits, bufsize);
440
441   AliPHOSRawDigiProducer rdp(rawReader,mapping);
442
443   rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
444   rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
445   rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
446   rdp.MakeDigits(digits,fTmpDigLG,fitter);
447
448   delete fitter ;
449
450   TClonesArray *tdigits = new TClonesArray("AliPHOSTriggerRawDigit",1);
451   tdigits->SetName("TDIGITS");
452   digitsTree->Branch("TPHOS", &tdigits, bufsize);  
453
454   rawReader->Reset();
455   AliPHOSTriggerRawDigiProducer tdp(rawReader);
456   
457   AliPHOSTriggerParameters* parameters = (AliPHOSTriggerParameters*)AliPHOSRecoParam::GetTriggerParameters();
458   
459   tdp.SetTriggerParameters(parameters);
460   tdp.ProcessEvent(tdigits);
461   
462   if (AliLog::GetGlobalDebugLevel() == 1) {
463     Int_t modMax=-111;
464     Int_t colMax=-111;
465     Int_t rowMax=-111;
466     Float_t eMax=-333;
467     //!!!for debug!!!
468     
469     Int_t relId[4];
470     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
471       AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
472       if(digit->GetEnergy()>eMax) {
473         fGeom->AbsToRelNumbering(digit->GetId(),relId);
474         eMax=digit->GetEnergy();
475         modMax=relId[0];
476         rowMax=relId[2];
477         colMax=relId[3];
478       }
479     }
480     
481     AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
482                     modMax,colMax,rowMax,eMax));
483   }
484   
485   digitsTree->Fill();
486   
487   digits->Delete();
488   delete digits;
489   
490   tdigits->Delete();
491   delete tdigits;
492 }
493 //==================================================================================
494 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
495   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
496
497   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
498
499   //Determine rel.position of the cell absolute ID
500   Int_t relId[4];
501   geom->AbsToRelNumbering(absId,relId);
502   Int_t module=relId[0];
503   Int_t row   =relId[2];
504   Int_t column=relId[3];
505   if(relId[1]){ //CPV
506     Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
507     return amp*calibration ;
508   }
509   else{ //EMC
510     Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
511     return amp*calibration ;
512   }
513 }
514 //==================================================================================
515 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId,Bool_t isLG)const{
516   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
517
518   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
519
520   //Determine rel.position of the cell absolute ID
521   Int_t relId[4];
522   geom->AbsToRelNumbering(absId,relId);
523   Int_t module=relId[0];
524   Int_t row   =relId[2];
525   Int_t column=relId[3];
526   if(relId[1]){ //CPV
527     return 0. ;
528   }
529   else{ //EMC
530     if(isLG)
531       time += fgCalibData->GetLGTimeShiftEmc(module,column,row);
532     else
533       time += fgCalibData->GetTimeShiftEmc(module,column,row);
534     return time ;
535   }
536 }
537 //==================================================================================
538 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
539   //Store PHOS matrixes in ESD Header
540
541   //Check, if matrixes was already stored
542   for(Int_t mod=0 ;mod<5; mod++){
543     if(esd->GetPHOSMatrix(mod)!=0)
544       return ;
545   }
546
547   //Create and store matrixes
548   if(!gGeoManager){
549     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
550     return ;
551   }
552   //Note, that owner of copied marixes will be header
553   char path[255] ;
554   TGeoHMatrix * m ;
555   for(Int_t mod=0; mod<5; mod++){
556     snprintf(path,255,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
557     if (gGeoManager->cd(path)){
558       m = gGeoManager->GetCurrentMatrix() ;
559       esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
560     }
561     else{
562       esd->SetPHOSMatrix(NULL,mod) ;
563     }
564   }
565
566 }
567 //==================================================================================
568 Float_t AliPHOSReconstructor::CorrectNonlinearity(Float_t en){
569
570   //For backward compatibility, if no RecoParameters found
571   if(!GetRecoParam()){
572     return 0.0241+1.0504*en+0.000249*en*en ;
573   }
574
575   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"NoCorrection")==0){
576     return en ;
577   }
578   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Gustavo2005")==0){
579     const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
580     return par[0]+par[1]*en + par[2]*en*en ;
581   }
582   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Henrik2010")==0){
583      const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
584      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])) ;
585   }
586   //For backward compatibility
587   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"")==0){
588     return 0.0241+1.0504*en+0.000249*en*en ;
589   }
590   return en ;
591 }
592
593 void AliPHOSReconstructor::readTRUParameters(AliPHOSTriggerParameters* parameters) const
594 {
595   //Read trigger parameters.
596
597   TString path(gSystem->Getenv("ALICE_ROOT"));
598   path += "/PHOS/macros/Trigger/OCDB/";
599   
600   for (Int_t mod = 2; mod < 5; ++mod) { // module
601     for (Int_t tru = 0; tru < 4; tru++) { // tru row
602       for (Int_t branch = 0; branch < 2; branch++) { // branch
603         
604         // Open the Appropriate pedestal file
605         TString fileName = path;
606         fileName += "pedestal_m";
607         fileName = fileName += mod;
608         fileName+="_r";
609         fileName+=tru;
610         fileName+="_b";
611         fileName+=branch;
612         fileName+=".dat";
613         std::ifstream instream;
614         instream.open(fileName.Data());
615         
616         // Read pedestals from file
617         if( ! instream.is_open() )
618           Printf("E-TRUPedestals: could not open %s", fileName.Data());
619         else
620           {
621             Int_t ped[112];
622             
623             char ch_s[36];
624             char *ch_s_p = ch_s;
625             //Int_t nlines = 0;
626             
627             Int_t t_ped_0 =0;
628             Int_t t_ped_1 =0;
629             Int_t t_ped_2 =0;
630             
631             for(Int_t n=0; n<112; n++)
632               {
633                 instream.getline(ch_s_p,36);
634                 if (ch_s_p[23]<='9' && ch_s_p[23]>='0')
635                   {
636                     t_ped_0 = ch_s_p[23]-'0';
637                   }
638                 else if (ch_s_p[23]>='A' && ch_s_p[23]<='Z')
639                   {
640                     t_ped_0 = ch_s_p[23]-'A'+10;
641                     
642                   }
643                   
644                 if (ch_s_p[22]<='9' && ch_s_p[22]>='0')
645                   {
646                     t_ped_1 = ch_s_p[22]-'0';
647                   }
648                 else if (ch_s_p[22]<='Z' && ch_s_p[22]>='A')
649                   {
650                     t_ped_1 = ch_s_p[22]-'A'+10;
651                   }
652                 
653                 if (ch_s_p[21]<='9' && ch_s_p[21]>='0')
654                   {
655                     t_ped_2 = ch_s_p[21]-'0';
656                   }
657                 else if (ch_s_p[21]<='Z' && ch_s_p[21]>='A')
658                   {
659                     t_ped_2 = ch_s_p[21]-'A'+10;
660                   }
661                 
662                 ped[n]=t_ped_2*256+t_ped_1*16+t_ped_0;
663                 
664                 
665               }
666             for (Int_t xrow = 0; xrow < 8; xrow++){
667               for (Int_t zcol = 0; zcol < 14; zcol++){
668                 Int_t pedestal = ped[zcol*8+xrow];
669                 
670                 if( pedestal < 612 && pedestal > 412 ) // resonable
671                   parameters->SetTRUPedestal(pedestal, mod, tru, branch, xrow, zcol);
672                 else // unresonable
673                   continue;
674               }
675             }
676           } // Ends read of pedestals from branch from file.
677         instream.close();
678       }// end branch
679     }// end tru
680     
681   }// end for mod
682 }
683
684
685
686