]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSReconstructor.cxx
Different mapping files per module are created by WriteMapping.C
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructor.cxx
CommitLineData
d15a28e7 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
b2a60966 16/* $Id$ */
17
d15a28e7 18//_________________________________________________________________________
a5e00016 19//--
20//-- Yves Schutz (SUBATECH)
dfe0be07 21// Reconstruction class. Redesigned from the old AliReconstructionner class and
22// derived from STEER/AliReconstructor.
23//
d15a28e7 24// --- ROOT system ---
bb5c37a5 25#include "TGeoManager.h"
26#include "TGeoMatrix.h"
00cfce1d 27
d15a28e7 28// --- Standard library ---
364de5c6 29
d15a28e7 30// --- AliRoot header files ---
9a2cdbdf 31#include "AliLog.h"
b3006690 32#include "AliAltroMapping.h"
af885e0f 33#include "AliESDEvent.h"
33ba95c3 34#include "AliESDCaloCluster.h"
a5e00016 35#include "AliESDCaloCells.h"
f444a19f 36#include "AliPHOSReconstructor.h"
7acf6008 37#include "AliPHOSClusterizerv1.h"
7acf6008 38#include "AliPHOSTrackSegmentMakerv1.h"
39#include "AliPHOSPIDv1.h"
23904d16 40#include "AliPHOSTracker.h"
d22dd3b4 41#include "AliRawReader.h"
d3aa2291 42#include "AliPHOSCalibData.h"
7e88424f 43#include "AliCDBEntry.h"
44#include "AliCDBManager.h"
64df000d 45#include "AliPHOSTrigger.h"
46#include "AliPHOSGeometry.h"
9a2cdbdf 47#include "AliPHOSDigit.h"
48#include "AliPHOSTrackSegment.h"
49#include "AliPHOSEmcRecPoint.h"
50#include "AliPHOSRecParticle.h"
379c5c09 51#include "AliPHOSRawFitterv0.h"
52#include "AliPHOSRawFitterv1.h"
53#include "AliPHOSRawFitterv2.h"
9a2cdbdf 54#include "AliPHOSRawDigiProducer.h"
55#include "AliPHOSPulseGenerator.h"
e957fea8 56
f444a19f 57ClassImp(AliPHOSReconstructor)
d15a28e7 58
2e60107f 59Bool_t AliPHOSReconstructor::fgDebug = kFALSE ;
6483babc 60TClonesArray* AliPHOSReconstructor::fgDigitsArray = 0; // Array of PHOS digits
61TObjArray* AliPHOSReconstructor::fgEMCRecPoints = 0; // Array of EMC rec.points
d3aa2291 62AliPHOSCalibData * AliPHOSReconstructor::fgCalibData = 0 ;
63
2e60107f 64
d15a28e7 65//____________________________________________________________________________
9a2cdbdf 66AliPHOSReconstructor::AliPHOSReconstructor() :
dcab1c7e 67 fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL)
d15a28e7 68{
b2a60966 69 // ctor
8d8258f6 70 fGeom = AliPHOSGeometry::GetInstance("IHEP","");
dcab1c7e 71 fClusterizer = new AliPHOSClusterizerv1 (fGeom);
72 fTSM = new AliPHOSTrackSegmentMakerv1(fGeom);
73 fPID = new AliPHOSPIDv1 (fGeom);
6483babc 74 fgDigitsArray = new TClonesArray("AliPHOSDigit",100);
75 fgEMCRecPoints= new TObjArray(100) ;
d3aa2291 76 if (!fgCalibData)
77 fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
78
e68222ce 79}
6ad0bfa0 80
0379a13e 81//____________________________________________________________________________
82 AliPHOSReconstructor::~AliPHOSReconstructor()
83{
84 // dtor
9a2cdbdf 85 delete fGeom;
8d8258f6 86 delete fClusterizer;
dcab1c7e 87 delete fTSM;
88 delete fPID;
6483babc 89 delete fgDigitsArray;
90 delete fgEMCRecPoints;
0379a13e 91}
7acf6008 92
7acf6008 93//____________________________________________________________________________
9a2cdbdf 94void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
35293055 95{
9a2cdbdf 96 // 'single-event' local reco method called by AliReconstruction;
dfe0be07 97 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
af885e0f 98 // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by
dfe0be07 99 // the global tracking.
9a2cdbdf 100
7e88424f 101 fClusterizer->InitParameters();
8d8258f6 102 fClusterizer->SetInput(digitsTree);
103 fClusterizer->SetOutput(clustersTree);
a68156e6 104 if ( Debug() )
8d8258f6 105 fClusterizer->Digits2Clusters("deb all") ;
a68156e6 106 else
8d8258f6 107 fClusterizer->Digits2Clusters("") ;
a68156e6 108}
109
110//____________________________________________________________________________
9a2cdbdf 111void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
112 AliESDEvent* esd) const
a68156e6 113{
9a2cdbdf 114 // This method produces PHOS rec-particles,
115 // then it creates AliESDtracks out of them and
116 // write tracks to the ESD
117
e68222ce 118
9a2cdbdf 119 // do current event; the loop over events is done by AliReconstruction::Run()
dcab1c7e 120 fTSM->SetESD(esd) ;
121 fTSM->SetInput(clustersTree);
9a2cdbdf 122 if ( Debug() )
dcab1c7e 123 fTSM->Clusters2TrackSegments("deb all") ;
9a2cdbdf 124 else
dcab1c7e 125 fTSM->Clusters2TrackSegments("") ;
35293055 126
71994f35 127 fPID->SetEnergyCorrectionOn(GetRecoParam()->GetEMCEnergyCorrectionOn());
128
dcab1c7e 129 fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ;
130 fPID->SetESD(esd) ;
dfe0be07 131 if ( Debug() )
dcab1c7e 132 fPID->TrackSegments2RecParticles("deb all") ;
dfe0be07 133 else
dcab1c7e 134 fPID->TrackSegments2RecParticles("") ;
7acf6008 135
dcab1c7e 136 TClonesArray *recParticles = fPID->GetRecParticles();
137 Int_t nOfRecParticles = recParticles->GetEntriesFast();
a5fa6165 138
85c60a8e 139 esd->SetNumberOfPHOSClusters(nOfRecParticles) ;
dd7ee508 140 esd->SetFirstPHOSCluster(esd->GetNumberOfCaloClusters()) ;
a5e00016 141
9a2cdbdf 142 AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
a5e00016 143
144 // Read digits array
145
146 TBranch *branch = digitsTree->GetBranch("PHOS");
147 if (!branch) {
148 AliError("can't get the branch with the PHOS digits !");
149 return;
150 }
6483babc 151 branch->SetAddress(&fgDigitsArray);
a5e00016 152 branch->GetEntry(0);
153
154 // Get the clusters array
155
156 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
157 if (!emcbranch) {
158 AliError("can't get the branch with the PHOS EMC clusters !");
159 return;
160 }
161
6483babc 162 emcbranch->SetAddress(&fgEMCRecPoints);
a5e00016 163 emcbranch->GetEntry(0);
8013c27e 164
64df000d 165 //#########Calculate trigger and set trigger info###########
a5e00016 166
64df000d 167 AliPHOSTrigger tr ;
168 // tr.SetPatchSize(1);//create 4x4 patches
bfae5a5d 169 tr.SetSimulation(kFALSE);
6483babc 170 tr.Trigger(fgDigitsArray);
64df000d 171
172 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
173 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
174 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
175 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
176
64df000d 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();
183
a5e00016 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));
64df000d 188
bfae5a5d 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};
64df000d 191 Int_t iAbsId2x2 =-1;
bfae5a5d 192 Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
64df000d 193 Int_t iAbsIdnxn =-1;
194 TVector3 pos2x2(-1,-1,-1);
195 TVector3 posnxn(-1,-1,-1);
9a2cdbdf 196 fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
197 fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
198 fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
199 fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
64df000d 200
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) ;
208
209 TArrayF triggerAmplitudes(4);
210 triggerAmplitudes[0] = maxAmp2x2 ;
211 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
212 triggerAmplitudes[2] = maxAmpnxn ;
213 triggerAmplitudes[3] = ampOutOfPatchnxn ;
214
215 //esd->SetPHOSTriggerCells(triggerPosition);
216 esd->AddPHOSTriggerPosition(triggerPosition);
217 esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
218
e68222ce 219
a5e00016 220 //########################################
221 //############# Fill CaloCells ###########
222 //########################################
223
6483babc 224 Int_t nDigits = fgDigitsArray->GetEntries();
a5e00016 225 Int_t idignew = 0 ;
226 AliDebug(1,Form("%d digits",nDigits));
227
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);
232
233 // Add to CaloCells only EMC digits with non-zero energy
234 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
6483babc 235 const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
12dd7f10 236 if(dig->GetId() <= knEMC &&
237 Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
d3aa2291 238 phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()), dig->GetTime());
a5e00016 239 idignew++;
240 }
e68222ce 241 }
a5e00016 242 phsCells.SetNumberOfCells(idignew);
243 phsCells.Sort();
e68222ce 244
a5e00016 245 //########################################
246 //############## Fill CaloClusters #######
247 //########################################
f4e0dd30 248
dfe0be07 249 for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
a5e00016 250 AliPHOSRecParticle *rp = dynamic_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
dfe0be07 251 if (Debug())
252 rp->Print();
25ed816e 253 // Get track segment and EMC rec.point associated with this rec.particle
dcab1c7e 254 AliPHOSTrackSegment *ts = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
a5e00016 255 ->At(rp->GetPHOSTSIndex()));
9a2cdbdf 256
6483babc 257 AliPHOSEmcRecPoint *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
25ed816e 258 AliESDCaloCluster *ec = new AliESDCaloCluster() ;
e68222ce 259
85c60a8e 260 Float_t xyz[3];
8013c27e 261 for (Int_t ixyz=0; ixyz<3; ixyz++)
262 xyz[ixyz] = rp->GetPos()[ixyz];
dd7ee508 263
264 AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
77ea1c6f 265
a5e00016 266 // Create cell lists
267
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];
273
274 for (Int_t iCell=0; iCell<cellMult; iCell++) {
6483babc 275 AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
a5e00016 276 absIdList[iCell] = (UShort_t)(digit->GetId());
277 if (digit->GetEnergy() > 0)
d3aa2291 278 fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
a5e00016 279 else
280 fracList[iCell] = 0;
25ed816e 281 }
77ea1c6f 282
7592dfc4 283 //Primaries
284 Int_t primMult = 0;
4dd59c4a 285 Int_t *primList = emcRP->GetPrimaries(primMult);
a5e00016 286
48c5db5b 287 Float_t energy;
7e88424f 288 if (GetRecoParam()->EMCEcore2ESD())
48c5db5b 289 energy = emcRP->GetCoreEnergy();
290 else
291 energy = rp->Energy();
292
7592dfc4 293 // fills the ESDCaloCluster
8ada0ffe 294 ec->SetClusterType(AliESDCaloCluster::kPHOSCluster);
a5e00016 295 ec->SetPosition(xyz); //rec.point position in MARS
48c5db5b 296 ec->SetE(energy); //total or core particle energy
25ed816e 297 ec->SetClusterDisp(emcRP->GetDispersion()); //cluster dispersion
e68222ce 298 ec->SetPid(rp->GetPID()) ; //array of particle identification
25ed816e 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
a5e00016 302 ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
25ed816e 303 ec->SetClusterChi2(-1); //not yet implemented
78902954 304 ec->SetTOF(emcRP->GetTime()); //Time of flight
305
a5e00016 306 //Cells contributing to clusters
307 ec->SetNCells(cellMult);
308 ec->SetCellsAbsId(absIdList);
309 ec->SetCellsAmplitudeFraction(fracList);
7592dfc4 310
3744a6cf 311 //Distance to the nearest bad crystal
312 ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal());
7592dfc4 313
314 //Array of MC indeces
4dd59c4a 315 TArrayI arrayPrim(primMult,primList);
7592dfc4 316 ec->AddLabels(arrayPrim);
64df000d 317
e44c41e9 318 //Matched ESD track
319 TArrayI arrayTrackMatched(1);
320 arrayTrackMatched[0]= ts->GetTrackIndex();
321 ec->AddTracksMatched(arrayTrackMatched);
322
85c60a8e 323 esd->AddCaloCluster(ec);
9a2cdbdf 324 delete ec;
091d150a 325 delete [] fracList;
326 delete [] absIdList;
e68222ce 327 }
6483babc 328 fgDigitsArray ->Delete();
329 fgEMCRecPoints->Delete();
330 recParticles ->Delete();
bb5c37a5 331
332 //Store PHOS misalignment matrixes
333 FillMisalMatrixes(esd) ;
334
dd7ee508 335}
336
e68222ce 337//____________________________________________________________________________
d76c31f4 338AliTracker* AliPHOSReconstructor::CreateTracker() const
dd7ee508 339{
9a2cdbdf 340 // creates the PHOS tracker
341 return new AliPHOSTracker();
342}
dd7ee508 343
1267d56d 344//____________________________________________________________________________
9a2cdbdf 345void AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
346{
347 // Converts raw data to
348 // PHOS digits
349 // Works on a single-event basis
9a2cdbdf 350 rawReader->Reset() ;
dd7ee508 351
379c5c09 352 AliPHOSRawFitterv0 * fitter ;
77ea1c6f 353
7e88424f 354 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
b3006690 355 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
356
357 AliAltroMapping *mapping[4];
358 for(Int_t i = 0; i < 4; i++) {
359 mapping[i] = (AliAltroMapping*)maps->At(i);
360 }
361
379c5c09 362 if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
363 fitter=new AliPHOSRawFitterv1();
364 else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
365 fitter=new AliPHOSRawFitterv2();
7e88424f 366 else
379c5c09 367 fitter=new AliPHOSRawFitterv0();
77ea1c6f 368
379c5c09 369 fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
370 fitter->SetAmpOffset (GetRecoParam()->GetGlobalAltroOffset());
371 fitter->SetAmpThreshold (GetRecoParam()->GetGlobalAltroThreshold());
8be3a30a 372
9a2cdbdf 373 TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
374 digits->SetName("DIGITS");
375 Int_t bufsize = 32000;
376 digitsTree->Branch("PHOS", &digits, bufsize);
dd7ee508 377
379c5c09 378 AliPHOSRawDigiProducer rdp(rawReader,mapping);
8be3a30a 379
379c5c09 380 rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
381 rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
382 rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
383 rdp.MakeDigits(digits,fitter);
77ea1c6f 384
379c5c09 385 delete fitter ;
7592dfc4 386
379c5c09 387 if (AliLog::GetGlobalDebugLevel() == 1) {
388 Int_t modMax=-111;
389 Int_t colMax=-111;
390 Int_t rowMax=-111;
391 Float_t eMax=-333;
392 //!!!for debug!!!
393
394 Int_t relId[4];
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();
400 modMax=relId[0];
401 rowMax=relId[2];
402 colMax=relId[3];
403 }
dd7ee508 404 }
379c5c09 405
406 AliDebug(1,Form("Digit with max. energy: modMax %d colMax %d rowMax %d eMax %f\n\n",
407 modMax,colMax,rowMax,eMax));
dd7ee508 408 }
409
9a2cdbdf 410 digitsTree->Fill();
e68222ce 411 digits->Delete();
412 delete digits;
d7d3b67b 413}
d3aa2291 414//==================================================================================
415Float_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
417
418 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
419
420 //Determine rel.position of the cell absolute ID
421 Int_t relId[4];
422 geom->AbsToRelNumbering(absId,relId);
423 Int_t module=relId[0];
424 Int_t row =relId[2];
425 Int_t column=relId[3];
426 if(relId[1]){ //CPV
427 Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
428 return amp*calibration ;
429 }
430 else{ //EMC
431 Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
432 return amp*calibration ;
433 }
434}
bb5c37a5 435//==================================================================================
436void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
437 //Store PHOS matrixes in ESD Header
438
439 //Check, if matrixes was already stored
440 for(Int_t mod=0 ;mod<5; mod++){
441 if(esd->GetPHOSMatrix(mod)!=0)
442 return ;
443 }
444
445 //Create and store matrixes
446 if(!gGeoManager){
447 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
448 return ;
449 }
450 //Note, that owner of copied marixes will be header
451 char path[255] ;
452 TGeoHMatrix * m ;
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) ;
458 }
459 else{
460 esd->SetPHOSMatrix(NULL,mod) ;
461 }
462 }
463
464}
d3aa2291 465
466