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