1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
20 //-- Yves Schutz (SUBATECH)
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and
22 // derived from STEER/AliReconstructor.
24 // --- ROOT system ---
25 #include "TGeoManager.h"
26 #include "TGeoMatrix.h"
28 // --- Standard library ---
30 // --- AliRoot header files ---
32 #include "AliAltroMapping.h"
33 #include "AliESDEvent.h"
34 #include "AliESDCaloCluster.h"
35 #include "AliESDCaloCells.h"
36 #include "AliPHOSReconstructor.h"
37 #include "AliPHOSClusterizerv1.h"
38 #include "AliPHOSTrackSegmentMakerv1.h"
39 #include "AliPHOSPIDv1.h"
40 #include "AliPHOSTracker.h"
41 #include "AliRawReader.h"
42 #include "AliPHOSCalibData.h"
43 #include "AliCDBEntry.h"
44 #include "AliCDBManager.h"
45 #include "AliPHOSTrigger.h"
46 #include "AliPHOSGeometry.h"
47 #include "AliPHOSDigit.h"
48 #include "AliPHOSTrackSegment.h"
49 #include "AliPHOSEmcRecPoint.h"
50 #include "AliPHOSRecParticle.h"
51 #include "AliPHOSRawFitterv0.h"
52 #include "AliPHOSRawFitterv1.h"
53 #include "AliPHOSRawFitterv2.h"
54 #include "AliPHOSRawDigiProducer.h"
55 #include "AliPHOSPulseGenerator.h"
57 ClassImp(AliPHOSReconstructor)
59 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
60 TClonesArray* AliPHOSReconstructor::fgDigitsArray = 0; // Array of PHOS digits
61 TObjArray* AliPHOSReconstructor::fgEMCRecPoints = 0; // Array of EMC rec.points
62 AliPHOSCalibData * AliPHOSReconstructor::fgCalibData = 0 ;
65 //____________________________________________________________________________
66 AliPHOSReconstructor::AliPHOSReconstructor() :
67 fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL)
70 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
71 fClusterizer = new AliPHOSClusterizerv1 (fGeom);
72 fTSM = new AliPHOSTrackSegmentMakerv1(fGeom);
73 fPID = new AliPHOSPIDv1 (fGeom);
74 fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
75 fgEMCRecPoints= new TObjArray(100) ;
77 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
81 //____________________________________________________________________________
82 AliPHOSReconstructor::~AliPHOSReconstructor()
90 delete fgEMCRecPoints;
93 //____________________________________________________________________________
94 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
96 // 'single-event' local reco method called by AliReconstruction;
97 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
98 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
99 // the global tracking.
101 fClusterizer->InitParameters();
102 fClusterizer->SetInput(digitsTree);
103 fClusterizer->SetOutput(clustersTree);
105 fClusterizer->Digits2Clusters("deb all") ;
107 fClusterizer->Digits2Clusters("") ;
110 //____________________________________________________________________________
111 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
112 AliESDEvent* esd) const
114 // This method produces PHOS rec-particles,
115 // then it creates AliESDtracks out of them and
116 // write tracks to the ESD
119 // do current event; the loop over events is done by AliReconstruction::Run()
121 fTSM->SetInput(clustersTree);
123 fTSM->Clusters2TrackSegments("deb all") ;
125 fTSM->Clusters2TrackSegments("") ;
127 fPID->SetEnergyCorrectionOn(GetRecoParam()->GetEMCEnergyCorrectionOn());
129 fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ;
132 fPID->TrackSegments2RecParticles("deb all") ;
134 fPID->TrackSegments2RecParticles("") ;
136 TClonesArray *recParticles = fPID->GetRecParticles();
137 Int_t nOfRecParticles = recParticles->GetEntriesFast();
139 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
140 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
142 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
146 TBranch *branch = digitsTree->GetBranch("PHOS");
148 AliError("can't get the branch with the PHOS digits !");
151 branch->SetAddress(&fgDigitsArray);
154 // Get the clusters array
156 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
158 AliError("can't get the branch with the PHOS EMC clusters !");
162 emcbranch->SetAddress(&fgEMCRecPoints);
163 emcbranch->GetEntry(0);
165 //#########Calculate trigger and set trigger info###########
168 // tr.SetPatchSize(1);//create 4x4 patches
169 tr.SetSimulation(kFALSE);
170 tr.Trigger(fgDigitsArray);
172 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
173 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
174 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
175 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
177 Int_t iSM2x2 = tr.Get2x2SuperModule();
178 Int_t iSMnxn = tr.GetnxnSuperModule();
179 Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
180 Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
181 Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
182 Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
184 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",
185 maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
186 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
187 maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
189 // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
190 Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
192 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
194 TVector3 pos2x2(-1,-1,-1);
195 TVector3 posnxn(-1,-1,-1);
196 fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
197 fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
198 fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
199 fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
201 TArrayF triggerPosition(6);
202 triggerPosition[0] = pos2x2(0) ;
203 triggerPosition[1] = pos2x2(1) ;
204 triggerPosition[2] = pos2x2(2) ;
205 triggerPosition[3] = posnxn(0) ;
206 triggerPosition[4] = posnxn(1) ;
207 triggerPosition[5] = posnxn(2) ;
209 TArrayF triggerAmplitudes(4);
210 triggerAmplitudes[0] = maxAmp2x2 ;
211 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
212 triggerAmplitudes[2] = maxAmpnxn ;
213 triggerAmplitudes[3] = ampOutOfPatchnxn ;
215 //esd->SetPHOSTriggerCells(triggerPosition);
216 esd->AddPHOSTriggerPosition(triggerPosition);
217 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
220 //########################################
221 //############# Fill CaloCells ###########
222 //########################################
224 Int_t nDigits = fgDigitsArray->GetEntries();
226 AliDebug(1,Form("%d digits",nDigits));
228 const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
229 AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
230 phsCells.CreateContainer(nDigits);
231 phsCells.SetType(AliESDCaloCells::kPHOSCell);
233 // Add to CaloCells only EMC digits with non-zero energy
234 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
235 const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
236 if(dig->GetId() <= knEMC &&
237 Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
238 phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()), dig->GetTime());
242 phsCells.SetNumberOfCells(idignew);
245 //########################################
246 //############## Fill CaloClusters #######
247 //########################################
249 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
250 AliPHOSRecParticle *rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
253 // Get track segment and EMC rec.point associated with this rec.particle
254 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
255 ->At(rp->GetPHOSTSIndex()));
257 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
258 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
261 for (Int_t ixyz=0; ixyz<3; ixyz++)
262 xyz[ixyz] = rp->GetPos()[ixyz];
264 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
268 Int_t cellMult = emcRP->GetDigitsMultiplicity();
269 Int_t *digitsList = emcRP->GetDigitsList();
270 Float_t *rpElist = emcRP->GetEnergiesList() ;
271 UShort_t *absIdList = new UShort_t[cellMult];
272 Double_t *fracList = new Double_t[cellMult];
274 for (Int_t iCell=0; iCell<cellMult; iCell++) {
275 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
276 absIdList[iCell] = (UShort_t)(digit->GetId());
277 if (digit->GetEnergy() > 0)
278 fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
285 Int_t *primList = emcRP->GetPrimaries(primMult);
288 if (GetRecoParam()->EMCEcore2ESD())
289 energy = emcRP->GetCoreEnergy();
291 energy = rp->Energy();
293 // fills the ESDCaloCluster
294 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
295 ec->SetPosition(xyz); //rec.point position in MARS
296 ec->SetE(energy); //total or core particle energy
297 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
298 ec->SetPid(rp->GetPID()) ; //array of particle identification
299 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
300 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
301 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
302 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
303 ec->SetClusterChi2(-1); //not yet implemented
304 ec->SetTOF(emcRP->GetTime()); //Time of flight
306 //Cells contributing to clusters
307 ec->SetNCells(cellMult);
308 ec->SetCellsAbsId(absIdList);
309 ec->SetCellsAmplitudeFraction(fracList);
311 //Distance to the nearest bad crystal
312 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
314 //Array of MC indeces
315 TArrayI arrayPrim(primMult,primList);
316 ec->AddLabels(arrayPrim);
319 TArrayI arrayTrackMatched(1);
320 arrayTrackMatched[0]= ts->GetTrackIndex();
321 ec->AddTracksMatched(arrayTrackMatched);
323 esd->AddCaloCluster(ec);
328 fgDigitsArray ->Delete();
329 fgEMCRecPoints->Delete();
330 recParticles ->Delete();
332 //Store PHOS misalignment matrixes
333 FillMisalMatrixes(esd) ;
337 //____________________________________________________________________________
338 AliTracker* AliPHOSReconstructor::CreateTracker() const
340 // creates the PHOS tracker
341 return new AliPHOSTracker();
344 //____________________________________________________________________________
345 void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
347 // Converts raw data to
349 // Works on a single-event basis
352 AliPHOSRawFitterv0 * fitter ;
354 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
355 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
357 AliAltroMapping *mapping[4];
358 for(Int_t i = 0; i < 4; i++) {
359 mapping[i] = (AliAltroMapping*)maps->At(i);
362 if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
363 fitter=new AliPHOSRawFitterv1();
364 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
365 fitter=new AliPHOSRawFitterv2();
367 fitter=new AliPHOSRawFitterv0();
369 fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
370 fitter->SetAmpOffset (GetRecoParam()->GetGlobalAltroOffset());
371 fitter->SetAmpThreshold (GetRecoParam()->GetGlobalAltroThreshold());
373 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
374 digits->SetName("DIGITS");
375 Int_t bufsize = 32000;
376 digitsTree->Branch("PHOS", &digits, bufsize);
378 AliPHOSRawDigiProducer rdp(rawReader,mapping);
380 rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
381 rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
382 rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
383 rdp.MakeDigits(digits,fitter);
387 if (AliLog::GetGlobalDebugLevel() == 1) {
395 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
396 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
397 if(digit->GetEnergy()>eMax) {
398 fGeom->AbsToRelNumbering(digit->GetId(),relId);
399 eMax=digit->GetEnergy();
406 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
407 modMax,colMax,rowMax,eMax));
414 //==================================================================================
415 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
416 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
418 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
420 //Determine rel.position of the cell absolute ID
422 geom->AbsToRelNumbering(absId,relId);
423 Int_t module=relId[0];
425 Int_t column=relId[3];
427 Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
428 return amp*calibration ;
431 Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
432 return amp*calibration ;
435 //==================================================================================
436 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
437 //Store PHOS matrixes in ESD Header
439 //Check, if matrixes was already stored
440 for(Int_t mod=0 ;mod<5; mod++){
441 if(esd->GetPHOSMatrix(mod)!=0)
445 //Create and store matrixes
447 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
450 //Note, that owner of copied marixes will be header
453 for(Int_t mod=0; mod<5; mod++){
454 sprintf(path,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
455 if (gGeoManager->cd(path)){
456 m = gGeoManager->GetCurrentMatrix() ;
457 esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
460 esd->SetPHOSMatrix(NULL,mod) ;