]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
Corrected compilation options
[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;
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       //For real data always primary==-1, we use this field to mark cells with overflow
277       //in HG channel. This is not needed in MC simulations,
278       if(primary==-1 && dig->IsLG())
279         primary=-2 ;
280       phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
281                                              CalibrateT(dig->GetTime(),dig->GetId(),dig->IsLG()),
282                                              primary) ;
283       idignew++;
284     }
285   }
286   phsCells.SetNumberOfCells(idignew);
287   phsCells.Sort();
288
289   //########################################
290   //############## Fill CaloClusters #######
291   //########################################
292
293   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
294     AliPHOSRecParticle  *rp    = static_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
295     if (Debug()) 
296       rp->Print();
297     // Get track segment and EMC rec.point associated with this rec.particle
298     AliPHOSTrackSegment *ts    = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
299                                                                     ->At(rp->GetPHOSTSIndex()));
300
301     AliPHOSEmcRecPoint  *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
302     AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
303     
304     Float_t xyz[3];
305     for (Int_t ixyz=0; ixyz<3; ixyz++) 
306       xyz[ixyz] = rp->GetPos()[ixyz];
307     
308     AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
309    
310     // Create cell lists
311
312     Int_t     cellMult   = emcRP->GetDigitsMultiplicity();
313     Int_t    *digitsList = emcRP->GetDigitsList();
314     Float_t  *rpElist    = emcRP->GetEnergiesList() ;
315     UShort_t *absIdList  = new UShort_t[cellMult];
316     Double_t *fracList   = new Double_t[cellMult];
317
318     for (Int_t iCell=0; iCell<cellMult; iCell++) {
319       AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
320       absIdList[iCell] = (UShort_t)(digit->GetId());
321       if (digit->GetEnergy() > 0)
322         fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
323       else
324         fracList[iCell] = 0;
325     }
326
327     //Primaries
328     Int_t  primMult  = 0;
329     Int_t *primList =  emcRP->GetPrimaries(primMult);
330
331     Float_t energy=0.;
332     if (GetRecoParam()->EMCEcore2ESD())
333       energy = emcRP->GetCoreEnergy();
334     else
335       energy = rp->Energy();
336     //Apply nonlinearity correction
337     if(GetRecoParam()->GetEMCEnergyCorrectionOn())
338       energy=CorrectNonlinearity(energy) ;
339
340     // fills the ESDCaloCluster
341     ec->SetType(AliVCluster::kPHOSNeutral);
342     ec->SetPosition(xyz);                       //rec.point position in MARS
343     ec->SetE(energy);                           //total or core particle energy
344     ec->SetDispersion(emcRP->GetDispersion());  //cluster dispersion
345     ec->SetPID(rp->GetPID()) ;            //array of particle identification
346     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
347     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
348     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
349     ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
350     ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z")); 
351     ec->SetChi2(-1);                     //not yet implemented
352     ec->SetTOF(emcRP->GetTime());               //Time of flight - already calibrated in EMCRecPoint
353
354     //Cells contributing to clusters
355     ec->SetNCells(cellMult);
356     ec->SetCellsAbsId(absIdList);
357     ec->SetCellsAmplitudeFraction(fracList);
358
359     //Distance to the nearest bad crystal
360     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
361   
362     //Array of MC indeces
363     TArrayI arrayPrim(primMult,primList);
364     ec->AddLabels(arrayPrim);
365     
366     //Matched ESD track
367     TArrayI arrayTrackMatched(1);
368     arrayTrackMatched[0]= ts->GetTrackIndex();
369     ec->AddTracksMatched(arrayTrackMatched);
370     
371     Int_t index = esd->AddCaloCluster(ec);
372
373     //Set pointer to this cluster in ESD track
374     Int_t nt=esd->GetNumberOfTracks();
375     for (Int_t itr=0; itr<nt; itr++) {
376       AliESDtrack *esdTrack=esd->GetTrack(itr);
377       if(!esdTrack->IsPHOS())
378         continue ;
379       if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
380         esdTrack->SetPHOScluster(index) ;
381 //no garatie that only one track matched this cluster
382 //      break ;
383       }
384     }
385  
386     delete ec;   
387     delete [] fracList;
388     delete [] absIdList;
389   }
390   fgDigitsArray ->Clear("C");
391   fgEMCRecPoints->Clear("C");
392   recParticles  ->Clear();
393
394   //Store PHOS misalignment matrixes
395   FillMisalMatrixes(esd) ;
396
397 }
398
399 //____________________________________________________________________________
400 AliTracker* AliPHOSReconstructor::CreateTracker() const
401 {
402   // creates the PHOS tracker
403   return new AliPHOSTracker();
404 }
405
406 //____________________________________________________________________________
407 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
408 {
409   // Converts raw data to
410   // PHOS digits
411   // Works on a single-event basis
412   rawReader->Reset() ; 
413
414   AliPHOSRawFitterv0 * fitter ;
415
416   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
417   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
418
419   AliAltroMapping *mapping[20];
420   for(Int_t i = 0; i < 20; i++) {
421     mapping[i] = (AliAltroMapping*)maps->At(i);
422   }
423
424   if      (strcmp(GetRecoParam()->EMCFitterVersion(),"v0")==0) 
425     fitter=new AliPHOSRawFitterv0();
426   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0) 
427     fitter=new AliPHOSRawFitterv1();
428   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0) 
429     fitter=new AliPHOSRawFitterv2();
430   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0) 
431     fitter=new AliPHOSRawFitterv3();
432   else
433     fitter=new AliPHOSRawFitterv4();
434
435   fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
436   fitter->SetAmpOffset     (GetRecoParam()->GetGlobalAltroOffset());
437   fitter->SetAmpThreshold  (GetRecoParam()->GetGlobalAltroThreshold());
438
439   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
440   digits->SetName("DIGITS");
441   Int_t bufsize = 32000;
442   digitsTree->Branch("PHOS", &digits, bufsize);
443
444   AliPHOSRawDigiProducer rdp(rawReader,mapping);
445
446   rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
447   rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
448   rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
449   rdp.MakeDigits(digits,fTmpDigLG,fitter);
450
451   delete fitter ;
452
453   TClonesArray *tdigits = new TClonesArray("AliPHOSTriggerRawDigit",1);
454   tdigits->SetName("TDIGITS");
455   digitsTree->Branch("TPHOS", &tdigits, bufsize);  
456
457   rawReader->Reset();
458   AliPHOSTriggerRawDigiProducer tdp(rawReader);
459   
460   AliPHOSTriggerParameters* parameters = (AliPHOSTriggerParameters*)AliPHOSRecoParam::GetTriggerParameters();
461   
462   tdp.SetTriggerParameters(parameters);
463   tdp.ProcessEvent(tdigits);
464   
465   if (AliLog::GetGlobalDebugLevel() == 1) {
466     Int_t modMax=-111;
467     Int_t colMax=-111;
468     Int_t rowMax=-111;
469     Float_t eMax=-333;
470     //!!!for debug!!!
471     
472     Int_t relId[4];
473     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
474       AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
475       if(digit->GetEnergy()>eMax) {
476         fGeom->AbsToRelNumbering(digit->GetId(),relId);
477         eMax=digit->GetEnergy();
478         modMax=relId[0];
479         rowMax=relId[2];
480         colMax=relId[3];
481       }
482     }
483     
484     AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
485                     modMax,colMax,rowMax,eMax));
486   }
487   
488   digitsTree->Fill();
489   
490   digits->Delete();
491   delete digits;
492   
493   tdigits->Delete();
494   delete tdigits;
495 }
496 //==================================================================================
497 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
498   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
499
500   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
501
502   //Determine rel.position of the cell absolute ID
503   Int_t relId[4];
504   geom->AbsToRelNumbering(absId,relId);
505   Int_t module=relId[0];
506   Int_t row   =relId[2];
507   Int_t column=relId[3];
508   if(relId[1]){ //CPV
509     Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
510     return amp*calibration ;
511   }
512   else{ //EMC
513     Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
514     return amp*calibration ;
515   }
516 }
517 //==================================================================================
518 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId,Bool_t isLG)const{
519   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
520
521   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
522
523   //Determine rel.position of the cell absolute ID
524   Int_t relId[4];
525   geom->AbsToRelNumbering(absId,relId);
526   Int_t module=relId[0];
527   Int_t row   =relId[2];
528   Int_t column=relId[3];
529   if(relId[1]){ //CPV
530     return 0. ;
531   }
532   else{ //EMC
533     if(isLG)
534       time += fgCalibData->GetLGTimeShiftEmc(module,column,row);
535     else
536       time += fgCalibData->GetTimeShiftEmc(module,column,row);
537     return time ;
538   }
539 }
540 //==================================================================================
541 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
542   //Store PHOS matrixes in ESD Header
543
544   //Check, if matrixes was already stored
545   for(Int_t mod=0 ;mod<5; mod++){
546     if(esd->GetPHOSMatrix(mod)!=0)
547       return ;
548   }
549
550   //Create and store matrixes
551   if(!gGeoManager){
552     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
553     return ;
554   }
555   //Note, that owner of copied marixes will be header
556   char path[255] ;
557   TGeoHMatrix * m ;
558   for(Int_t mod=0; mod<5; mod++){
559     snprintf(path,255,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
560     if (gGeoManager->cd(path)){
561       m = gGeoManager->GetCurrentMatrix() ;
562       esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
563     }
564     else{
565       esd->SetPHOSMatrix(NULL,mod) ;
566     }
567   }
568
569 }
570 //==================================================================================
571 Float_t AliPHOSReconstructor::CorrectNonlinearity(Float_t en){
572
573   //For backward compatibility, if no RecoParameters found
574   if(!GetRecoParam()){
575     return 0.0241+1.0504*en+0.000249*en*en ;
576   }
577
578   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"NoCorrection")==0){
579     return en ;
580   }
581   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Gustavo2005")==0){
582     const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
583     return par[0]+par[1]*en + par[2]*en*en ;
584   }
585   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Henrik2010")==0){
586      const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
587      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])) ;
588   }
589   //For backward compatibility
590   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"")==0){
591     return 0.0241+1.0504*en+0.000249*en*en ;
592   }
593   return en ;
594 }
595
596 void AliPHOSReconstructor::readTRUParameters(AliPHOSTriggerParameters* parameters) const
597 {
598   //Read trigger parameters.
599
600   TString path(gSystem->Getenv("ALICE_ROOT"));
601   path += "/PHOS/macros/Trigger/OCDB/";
602   
603   for (Int_t mod = 2; mod < 5; ++mod) { // module
604     for (Int_t tru = 0; tru < 4; tru++) { // tru row
605       for (Int_t branch = 0; branch < 2; branch++) { // branch
606         
607         // Open the Appropriate pedestal file
608         TString fileName = path;
609         fileName += "pedestal_m";
610         fileName = fileName += mod;
611         fileName+="_r";
612         fileName+=tru;
613         fileName+="_b";
614         fileName+=branch;
615         fileName+=".dat";
616         std::ifstream instream;
617         instream.open(fileName.Data());
618         
619         // Read pedestals from file
620         if( ! instream.is_open() )
621           Printf("E-TRUPedestals: could not open %s", fileName.Data());
622         else
623           {
624             Int_t ped[112];
625             
626             char ch_s[36];
627             char *ch_s_p = ch_s;
628             //Int_t nlines = 0;
629             
630             Int_t t_ped_0 =0;
631             Int_t t_ped_1 =0;
632             Int_t t_ped_2 =0;
633             
634             for(Int_t n=0; n<112; n++)
635               {
636                 instream.getline(ch_s_p,36);
637                 if (ch_s_p[23]<='9' && ch_s_p[23]>='0')
638                   {
639                     t_ped_0 = ch_s_p[23]-'0';
640                   }
641                 else if (ch_s_p[23]>='A' && ch_s_p[23]<='Z')
642                   {
643                     t_ped_0 = ch_s_p[23]-'A'+10;
644                     
645                   }
646                   
647                 if (ch_s_p[22]<='9' && ch_s_p[22]>='0')
648                   {
649                     t_ped_1 = ch_s_p[22]-'0';
650                   }
651                 else if (ch_s_p[22]<='Z' && ch_s_p[22]>='A')
652                   {
653                     t_ped_1 = ch_s_p[22]-'A'+10;
654                   }
655                 
656                 if (ch_s_p[21]<='9' && ch_s_p[21]>='0')
657                   {
658                     t_ped_2 = ch_s_p[21]-'0';
659                   }
660                 else if (ch_s_p[21]<='Z' && ch_s_p[21]>='A')
661                   {
662                     t_ped_2 = ch_s_p[21]-'A'+10;
663                   }
664                 
665                 ped[n]=t_ped_2*256+t_ped_1*16+t_ped_0;
666                 
667                 
668               }
669             for (Int_t xrow = 0; xrow < 8; xrow++){
670               for (Int_t zcol = 0; zcol < 14; zcol++){
671                 Int_t pedestal = ped[zcol*8+xrow];
672                 
673                 if( pedestal < 612 && pedestal > 412 ) // resonable
674                   parameters->SetTRUPedestal(pedestal, mod, tru, branch, xrow, zcol);
675                 else // unresonable
676                   continue;
677               }
678             }
679           } // Ends read of pedestals from branch from file.
680         instream.close();
681       }// end branch
682     }// end tru
683     
684   }// end for mod
685 }
686
687
688
689