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