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