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