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