]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
ConvertDigits() use information from AliPHOSRecoParam::SubtractPedestals()
[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 "AliPHOSQualAssDataMaker.h"
31 #include "AliESDEvent.h"
32 #include "AliESDCaloCluster.h"
33 #include "AliPHOSReconstructor.h"
34 #include "AliPHOSClusterizerv1.h"
35 #include "AliPHOSTrackSegmentMakerv1.h"
36 #include "AliPHOSPIDv1.h"
37 #include "AliPHOSTracker.h"
38 #include "AliRawReader.h"
39 #include "AliPHOSTrigger.h"
40 #include "AliPHOSGeometry.h"
41 #include "AliPHOSRecoParam.h"
42 #include "AliPHOSRecoParamEmc.h"
43 #include "AliPHOSRecoParamCpv.h"
44 #include "AliPHOSDigit.h"
45 #include "AliPHOSTrackSegment.h"
46 #include "AliPHOSEmcRecPoint.h"
47 #include "AliPHOSRecParticle.h"
48 #include "AliPHOSRawDecoder.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   tsm->GetQualAssDataMaker()->Init(AliQualAss::kTRACKSEGMENTS) ; 
111   AliPHOSPID *pid = new AliPHOSPIDv1(fGeom);
112   pid->GetQualAssDataMaker()->Init(AliQualAss::kRECPARTICLES) ;    
113
114   // do current event; the loop over events is done by AliReconstruction::Run()
115   tsm->SetESD(esd) ; 
116   tsm->SetInput(clustersTree);
117   if ( Debug() ) 
118     tsm->Clusters2TrackSegments("deb all") ;
119   else 
120     tsm->Clusters2TrackSegments("") ;
121   
122   tsm->GetQualAssDataMaker()->SetData(tsm->GetTrackSegments()) ; 
123   tsm->GetQualAssDataMaker()->Exec(AliQualAss::kTRACKSEGMENTS) ; 
124
125   pid->SetInput(clustersTree, tsm->GetTrackSegments()) ; 
126   pid->SetESD(esd) ; 
127   if ( Debug() ) 
128     pid->TrackSegments2RecParticles("deb all") ;
129   else 
130     pid->TrackSegments2RecParticles("") ;
131
132   pid->GetQualAssDataMaker()->SetData(pid->GetRecParticles()) ; 
133   pid->GetQualAssDataMaker()->Exec(AliQualAss::kRECPARTICLES) ; 
134   
135   // PLEASE FIX IT. SHOULD GO TO ALIRECONSTRUCTION !!
136   //  if ( eventNumber == gime->MaxEvent()-1 ) {
137   //    fTSM->GetQualAssDataMaker()->Finish(AliQualAss::kTRACKSEGMENTS) ; 
138   //    fPID->GetQualAssDataMaker()->Finish(AliQualAss::kRECPARTICLES) ; 
139   //  }
140         
141   // This function creates AliESDtracks from AliPHOSRecParticles
142   //         and
143   // writes them to the ESD
144
145   TClonesArray *recParticles  = pid->GetRecParticles();
146   Int_t nOfRecParticles = recParticles->GetEntries();
147   
148   
149   esd->SetNumberOfPHOSClusters(nOfRecParticles) ; 
150   esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
151
152   AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
153
154
155   //#########Calculate trigger and set trigger info###########
156  
157   AliPHOSTrigger tr ;
158   //   tr.SetPatchSize(1);//create 4x4 patches
159   tr.Trigger();
160   
161   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
162   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
163   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
164   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
165
166   Int_t iSM2x2      = tr.Get2x2SuperModule();
167   Int_t iSMnxn      = tr.GetnxnSuperModule();
168   Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
169   Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
170   Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
171   Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
172
173   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
174   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",  maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
175
176   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.
177   Int_t iAbsId2x2 =-1;
178   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.
179   Int_t iAbsIdnxn =-1;
180   TVector3    pos2x2(-1,-1,-1);
181   TVector3    posnxn(-1,-1,-1);
182   fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
183   fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
184   fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
185   fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
186
187   TArrayF triggerPosition(6);
188   triggerPosition[0] = pos2x2(0) ;   
189   triggerPosition[1] = pos2x2(1) ;   
190   triggerPosition[2] = pos2x2(2) ;  
191   triggerPosition[3] = posnxn(0) ;   
192   triggerPosition[4] = posnxn(1) ;   
193   triggerPosition[5] = posnxn(2) ;  
194
195   TArrayF triggerAmplitudes(4);
196   triggerAmplitudes[0] = maxAmp2x2 ;   
197   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
198   triggerAmplitudes[2] = maxAmpnxn ;   
199   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
200
201   //esd->SetPHOSTriggerCells(triggerPosition);
202   esd->AddPHOSTriggerPosition(triggerPosition);
203   esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
204   
205   //######################################
206   
207   //Fill CaloClusters 
208   const Float_t kBigShort = std::numeric_limits<short int>::max() - 1;
209   const Float_t nsec100   = 1e9*100.; // units of 0.01 ns
210   const Float_t gev500    = 500.;     // units of GeV/500
211
212   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
213     AliPHOSRecParticle * rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
214     if (Debug()) 
215       rp->Print();
216     // Get track segment and EMC rec.point associated with this rec.particle
217     AliPHOSTrackSegment *ts    = static_cast<AliPHOSTrackSegment *>(tsm->GetTrackSegments()->At(rp->GetPHOSTSIndex()));
218
219     // Get the clusters array
220     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
221     if (!emcbranch) { 
222       AliError("can't get the branch with the PHOS EMC clusters !");
223       return;
224     }
225     TObjArray *emcRecPoints = new TObjArray(100) ;
226     emcbranch->SetAddress(&emcRecPoints);
227     emcbranch->GetEntry(0);
228
229     AliPHOSEmcRecPoint  *emcRP = static_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(ts->GetEmcIndex()));
230     AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
231         
232     Float_t xyz[3];
233     for (Int_t ixyz=0; ixyz<3; ixyz++) 
234       xyz[ixyz] = rp->GetPos()[ixyz];
235     
236     AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
237     
238     //Create digits lists
239     Int_t  digitMult  = emcRP->GetDigitsMultiplicity();
240     Int_t *digitsList = emcRP->GetDigitsList();
241     Short_t *amplList  = new Short_t[digitMult];
242     Short_t *timeList  = new Short_t[digitMult];
243     Short_t *digiList  = new Short_t[digitMult];
244
245     // Read digits array
246     TBranch *branch = digitsTree->GetBranch("PHOS");
247     if (!branch) { 
248       AliError("can't get the branch with the PHOS digits !");
249       return;
250     }
251     TClonesArray *digitsArr = new TClonesArray("AliPHOSDigit",100);
252     branch->SetAddress(&digitsArr);
253     branch->GetEntry(0);
254
255     // Convert Float_t* and Int_t* to Short_t* to save memory
256     for (Int_t iDigit=0; iDigit<digitMult; iDigit++) {
257       AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(digitsArr->At(digitsList[iDigit]));
258       amplList[iDigit] =
259         (Short_t)(TMath::Min(digit->GetEnergy()*gev500,kBigShort)); // Energy in units of GeV/500
260       timeList[iDigit] =
261         (Short_t)(TMath::Min(digit->GetTime()*nsec100,kBigShort)); // time in units of 0.01 ns
262       digiList[iDigit] = (Short_t)(digit->GetId());
263     }
264     
265     //Primaries
266     Int_t  primMult  = 0;
267     Int_t *primInts =  emcRP->GetPrimaries(primMult);
268     Short_t *primList = new Short_t[primMult];
269     for (Int_t ipr=0; ipr<primMult; ipr++) 
270       primList[ipr] = (Short_t)(primInts[ipr]);  
271     
272     // fills the ESDCaloCluster
273  
274     ec->SetPHOS(kTRUE);
275     ec->SetPosition(xyz);                 //rec.point position in MARS
276     ec->SetE(rp->Energy());         //total particle energy
277     ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
278     ec->SetPid          (rp->GetPID()) ;        //array of particle identification
279     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
280     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
281     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
282     ec->SetEmcCpvDistance(-1);                  //not yet implemented
283     ec->SetClusterChi2(-1);                     //not yet implemented
284     ec->SetM11(-1) ;                            //not yet implemented
285  
286     //Digits Lists
287     TArrayS arrayAmpList(digitMult,amplList);
288     TArrayS arrayTimeList(digitMult,timeList);
289     TArrayS arrayIndexList(digitMult,digiList);
290     ec->AddDigitAmplitude(arrayAmpList);
291     ec->AddDigitTime(arrayTimeList);
292     ec->AddDigitIndex(arrayIndexList);
293
294     //Distance to the nearest bad crystal
295     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
296   
297     //Array of MC indeces
298     TArrayS arrayPrim(primMult,primList);
299     ec->AddLabels(arrayPrim);
300
301     //Array of tracks uncomment when available in future
302     //TArrayS arrayTrackMatched(1);// Only one track, temporal solution.
303     //arrayTrackMatched[0]= (Short_t)(matchedTrack[iClust]);
304     //ec->AddTracksMatched(arrayTrackMatched);
305     
306     // add the track to the esd object
307     esd->AddCaloCluster(ec);
308     delete ec;   
309     delete [] primList;
310     delete [] amplList;
311     delete [] timeList;
312     delete [] digiList;    
313   }  
314 }
315
316 AliTracker* AliPHOSReconstructor::CreateTracker(AliRunLoader* /* runLoader */) const
317 {
318   // creates the PHOS tracker
319   return new AliPHOSTracker();
320 }
321
322 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
323 {
324   // Converts raw data to
325   // PHOS digits
326   // Works on a single-event basis
327
328   rawReader->Reset() ; 
329
330   AliPHOSRawDecoder dc(rawReader);
331   TString option = GetOption();
332   if (option.Contains("OldRCUFormat"))
333     dc.SetOldRCUFormat(kTRUE);
334   else
335     dc.SetOldRCUFormat(kFALSE);
336   
337   dc.SubtractPedestals(fgkRecoParamEmc->SubtractPedestals());
338   
339   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
340   digits->SetName("DIGITS");
341   Int_t bufsize = 32000;
342   digitsTree->Branch("PHOS", &digits, bufsize);
343
344   AliPHOSRawDigiProducer pr;
345   pr.MakeDigits(digits,&dc);
346
347   //ADC counts -> GeV
348   for(Int_t i=0; i<digits->GetEntries(); i++) {
349     AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(i);
350     digit->SetEnergy(digit->GetEnergy()/AliPHOSPulseGenerator::GeV2ADC());
351   }
352   
353   //!!!!for debug!!!
354   Int_t modMax=-111;
355   Int_t colMax=-111;
356   Int_t rowMax=-111;
357   Float_t eMax=-333;
358   //!!!for debug!!!
359
360   Int_t relId[4];
361   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
362     AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
363     if(digit->GetEnergy()>eMax) {
364       fGeom->AbsToRelNumbering(digit->GetId(),relId);
365       eMax=digit->GetEnergy();
366       modMax=relId[0];
367       rowMax=relId[2];
368       colMax=relId[3];
369     }
370   }
371
372   AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
373                   modMax,colMax,rowMax,eMax));
374
375   digitsTree->Fill();
376 }