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 "AliPHOSRawFitterv3.h"
55 #include "AliPHOSRawDigiProducer.h"
56 #include "AliPHOSPulseGenerator.h"
58 ClassImp(AliPHOSReconstructor)
60 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
61 TClonesArray* AliPHOSReconstructor::fgDigitsArray = 0; // Array of PHOS digits
62 TObjArray* AliPHOSReconstructor::fgEMCRecPoints = 0; // Array of EMC rec.points
63 AliPHOSCalibData * AliPHOSReconstructor::fgCalibData = 0 ;
66 //____________________________________________________________________________
67 AliPHOSReconstructor::AliPHOSReconstructor() :
68 fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL)
71 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
72 fClusterizer = new AliPHOSClusterizerv1 (fGeom);
73 fTSM = new AliPHOSTrackSegmentMakerv1(fGeom);
74 fPID = new AliPHOSPIDv1 (fGeom);
75 fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
76 fgEMCRecPoints= new TObjArray(100) ;
78 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
82 //____________________________________________________________________________
83 AliPHOSReconstructor::~AliPHOSReconstructor()
91 delete fgEMCRecPoints;
94 //____________________________________________________________________________
95 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
97 // 'single-event' local reco method called by AliReconstruction;
98 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
99 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
100 // the global tracking.
102 fClusterizer->InitParameters();
103 fClusterizer->SetInput(digitsTree);
104 fClusterizer->SetOutput(clustersTree);
106 fClusterizer->Digits2Clusters("deb all") ;
108 fClusterizer->Digits2Clusters("") ;
111 //____________________________________________________________________________
112 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
113 AliESDEvent* esd) const
115 // This method produces PHOS rec-particles,
116 // then it creates AliESDtracks out of them and
117 // write tracks to the ESD
120 // do current event; the loop over events is done by AliReconstruction::Run()
122 fTSM->SetInput(clustersTree);
124 fTSM->Clusters2TrackSegments("deb all") ;
126 fTSM->Clusters2TrackSegments("") ;
128 fPID->SetEnergyCorrectionOn(GetRecoParam()->GetEMCEnergyCorrectionOn());
130 fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ;
133 fPID->TrackSegments2RecParticles("deb all") ;
135 fPID->TrackSegments2RecParticles("") ;
137 TClonesArray *recParticles = fPID->GetRecParticles();
138 Int_t nOfRecParticles = recParticles->GetEntriesFast();
140 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
141 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
143 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
147 TBranch *branch = digitsTree->GetBranch("PHOS");
149 AliError("can't get the branch with the PHOS digits !");
152 branch->SetAddress(&fgDigitsArray);
155 // Get the clusters array
157 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
159 AliError("can't get the branch with the PHOS EMC clusters !");
163 emcbranch->SetAddress(&fgEMCRecPoints);
164 emcbranch->GetEntry(0);
166 //#########Calculate trigger and set trigger info###########
169 // tr.SetPatchSize(1);//create 4x4 patches
170 tr.SetSimulation(kFALSE);
171 tr.Trigger(fgDigitsArray);
173 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
174 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
175 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
176 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
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();
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));
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};
193 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
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);
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) ;
210 TArrayF triggerAmplitudes(4);
211 triggerAmplitudes[0] = maxAmp2x2 ;
212 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
213 triggerAmplitudes[2] = maxAmpnxn ;
214 triggerAmplitudes[3] = ampOutOfPatchnxn ;
216 //esd->SetPHOSTriggerCells(triggerPosition);
217 esd->AddPHOSTriggerPosition(triggerPosition);
218 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
221 //########################################
222 //############# Fill CaloCells ###########
223 //########################################
225 Int_t nDigits = fgDigitsArray->GetEntries();
227 AliDebug(1,Form("%d digits",nDigits));
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);
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*)fgDigitsArray->At(idig);
237 if(dig->GetId() <= knEMC &&
238 Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
239 phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
240 CalibrateT(dig->GetTime(),dig->GetId()));
244 phsCells.SetNumberOfCells(idignew);
247 //########################################
248 //############## Fill CaloClusters #######
249 //########################################
251 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
252 AliPHOSRecParticle *rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
255 // Get track segment and EMC rec.point associated with this rec.particle
256 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
257 ->At(rp->GetPHOSTSIndex()));
259 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
260 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
263 for (Int_t ixyz=0; ixyz<3; ixyz++)
264 xyz[ixyz] = rp->GetPos()[ixyz];
266 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
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];
276 for (Int_t iCell=0; iCell<cellMult; iCell++) {
277 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
278 absIdList[iCell] = (UShort_t)(digit->GetId());
279 if (digit->GetEnergy() > 0)
280 fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
287 Int_t *primList = emcRP->GetPrimaries(primMult);
290 if (GetRecoParam()->EMCEcore2ESD())
291 energy = emcRP->GetCoreEnergy();
293 energy = rp->Energy();
295 // fills the ESDCaloCluster
296 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
297 ec->SetPosition(xyz); //rec.point position in MARS
298 ec->SetE(energy); //total or core particle energy
299 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
300 ec->SetPid(rp->GetPID()) ; //array of particle identification
301 ec->SetM02(emcRP->GetM2x()) ; //second moment M2x
302 ec->SetM20(emcRP->GetM2z()) ; //second moment M2z
303 ec->SetNExMax(emcRP->GetNExMax()); //number of local maxima
304 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
305 ec->SetClusterChi2(-1); //not yet implemented
306 ec->SetTOF(emcRP->GetTime()); //Time of flight - already calibrated in EMCRecPoint
308 //Cells contributing to clusters
309 ec->SetNCells(cellMult);
310 ec->SetCellsAbsId(absIdList);
311 ec->SetCellsAmplitudeFraction(fracList);
313 //Distance to the nearest bad crystal
314 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
316 //Array of MC indeces
317 TArrayI arrayPrim(primMult,primList);
318 ec->AddLabels(arrayPrim);
321 TArrayI arrayTrackMatched(1);
322 arrayTrackMatched[0]= ts->GetTrackIndex();
323 ec->AddTracksMatched(arrayTrackMatched);
325 esd->AddCaloCluster(ec);
330 fgDigitsArray ->Delete();
331 fgEMCRecPoints->Delete();
332 recParticles ->Delete();
334 //Store PHOS misalignment matrixes
335 FillMisalMatrixes(esd) ;
339 //____________________________________________________________________________
340 AliTracker* AliPHOSReconstructor::CreateTracker() const
342 // creates the PHOS tracker
343 return new AliPHOSTracker();
346 //____________________________________________________________________________
347 void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
349 // Converts raw data to
351 // Works on a single-event basis
354 AliPHOSRawFitterv0 * fitter ;
356 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
357 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
359 AliAltroMapping *mapping[20];
360 for(Int_t i = 0; i < 20; i++) {
361 mapping[i] = (AliAltroMapping*)maps->At(i);
364 if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
365 fitter=new AliPHOSRawFitterv1();
366 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
367 fitter=new AliPHOSRawFitterv2();
368 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0)
369 fitter=new AliPHOSRawFitterv3();
371 fitter=new AliPHOSRawFitterv0();
373 fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
374 fitter->SetAmpOffset (GetRecoParam()->GetGlobalAltroOffset());
375 fitter->SetAmpThreshold (GetRecoParam()->GetGlobalAltroThreshold());
377 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
378 digits->SetName("DIGITS");
379 Int_t bufsize = 32000;
380 digitsTree->Branch("PHOS", &digits, bufsize);
382 AliPHOSRawDigiProducer rdp(rawReader,mapping);
384 rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
385 rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
386 rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
387 rdp.MakeDigits(digits,fitter);
391 if (AliLog::GetGlobalDebugLevel() == 1) {
399 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
400 AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
401 if(digit->GetEnergy()>eMax) {
402 fGeom->AbsToRelNumbering(digit->GetId(),relId);
403 eMax=digit->GetEnergy();
410 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
411 modMax,colMax,rowMax,eMax));
418 //==================================================================================
419 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
420 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
422 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
424 //Determine rel.position of the cell absolute ID
426 geom->AbsToRelNumbering(absId,relId);
427 Int_t module=relId[0];
429 Int_t column=relId[3];
431 Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
432 return amp*calibration ;
435 Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
436 return amp*calibration ;
439 //==================================================================================
440 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId)const{
441 // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
443 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
445 //Determine rel.position of the cell absolute ID
447 geom->AbsToRelNumbering(absId,relId);
448 Int_t module=relId[0];
450 Int_t column=relId[3];
455 time += fgCalibData->GetTimeShiftEmc(module,column,row);
459 //==================================================================================
460 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
461 //Store PHOS matrixes in ESD Header
463 //Check, if matrixes was already stored
464 for(Int_t mod=0 ;mod<5; mod++){
465 if(esd->GetPHOSMatrix(mod)!=0)
469 //Create and store matrixes
471 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
474 //Note, that owner of copied marixes will be header
477 for(Int_t mod=0; mod<5; mod++){
478 sprintf(path,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
479 if (gGeoManager->cd(path)){
480 m = gGeoManager->GetCurrentMatrix() ;
481 esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
484 esd->SetPHOSMatrix(NULL,mod) ;