]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSPIDv1.cxx
Minor changes to the Digitizer procedure
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv1.cxx
CommitLineData
6ad0bfa0 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
6ad0bfa0 18//_________________________________________________________________________
b2a60966 19// Implementation version v1 of the PHOS particle identifier
7acf6008 20// Particle identification based on the
21// - CPV information,
a4e98857 22// - Preshower information (in MIXT or GPS2 geometries)
7acf6008 23// - shower width.
a4e98857 24//
25// CPV or Preshower clusters should be closer in PHOS plane than fCpvEmcDistance (in cm).
26// This parameter can be set by method SetCpvtoEmcDistanceCut(Float_t cut)
7acf6008 27//
28// One can set desirable ID method by the function SetIdentificationMethod(option).
a4e98857 29// Presently the following options can be used together or separately :
7acf6008 30// - "disp": use dispersion cut on shower width
31// (width can be set by method SetDispersionCut(Float_t cut)
a4e98857 32// - "ell" : use cut on the axis of the ellipse, drawn around shower
7acf6008 33// (this cut can be changed by SetShowerProfileCut(char* formula),
34// where formula - any function of two variables f(lambda[0],lambda[1]).
35// Shower is considered as EM if f() > 0 )
36// One can visualize current cuts calling method PlotDispersionCuts().
37//
a4e98857 38// use case:
39// root [0] AliPHOSPIDv1 * p1 = new AliPHOSPIDv1("galice.root")
40// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
41// root [1] p1->SetIdentificationMethod("disp ellipse")
42// root [2] p1->ExecuteTask()
43// root [3] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","ts1")
44// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
7acf6008 45// // reading headers from file galice1.root and TrackSegments
46// // with title "ts1"
a4e98857 47// root [4] p2->SetRecParticlesBranch("rp1")
7acf6008 48// // set file name for the branch RecParticles
a4e98857 49// root [5] p2->ExecuteTask("deb all time")
7acf6008 50// // available options
51// // "deb" - prints # of reconstructed particles
52// // "deb all" - prints # and list of RecParticles
53// // "time" - prints benchmarking results
54//
55//*-- Author: Yves Schutz (SUBATECH) & Gines Martinez (SUBATECH) &
56// Dmitri Peressounko (SUBATECH & Kurchatov Institute)
57// Complitely redesined by Dmitri Peressounko, March 2001
6ad0bfa0 58
59// --- ROOT system ---
7acf6008 60#include "TROOT.h"
61#include "TTree.h"
62#include "TFile.h"
63#include "TF2.h"
64#include "TFormula.h"
65#include "TCanvas.h"
66#include "TFolder.h"
67#include "TSystem.h"
68#include "TBenchmark.h"
6ad0bfa0 69// --- Standard library ---
70
de9ec31b 71#include <iostream.h>
7acf6008 72#include <iomanip.h>
6ad0bfa0 73
74// --- AliRoot header files ---
75
7acf6008 76#include "AliRun.h"
77#include "AliGenerator.h"
78#include "AliPHOS.h"
26d4b141 79#include "AliPHOSPIDv1.h"
7acf6008 80#include "AliPHOSClusterizerv1.h"
6ad0bfa0 81#include "AliPHOSTrackSegment.h"
7acf6008 82#include "AliPHOSTrackSegmentMakerv1.h"
6ad0bfa0 83#include "AliPHOSRecParticle.h"
7b7c1533 84#include "AliPHOSGeometry.h"
85#include "AliPHOSGetter.h"
6ad0bfa0 86
26d4b141 87ClassImp( AliPHOSPIDv1)
6ad0bfa0 88
1cb7c1ee 89//____________________________________________________________________________
90AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
91{
a4e98857 92 // default ctor
7b7c1533 93 fFormula = 0 ;
94 fDispersion = 0. ;
95 fCpvEmcDistance = 0 ;
96 fHeaderFileName = "" ;
97 fTrackSegmentsTitle= "" ;
98 fRecPointsTitle = "" ;
99 fRecParticlesTitle = "" ;
100 fIDOptions = "" ;
7acf6008 101}
102
103//____________________________________________________________________________
7b7c1533 104AliPHOSPIDv1::AliPHOSPIDv1(const char * headerFile,const char * name) : AliPHOSPID(headerFile, name)
7acf6008 105{
a4e98857 106 //ctor with the indication on where to look for the track segments
dd5c4038 107
7b7c1533 108 fFormula = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;
109 fDispersion = 2.0 ;
110 fCpvEmcDistance = 3.0 ;
111
112 fHeaderFileName = GetTitle() ;
113 fTrackSegmentsTitle = GetName() ;
114 fRecPointsTitle = GetName() ;
115 fRecParticlesTitle = GetName() ;
116 fIDOptions = "" ;
117
118 TString tempo(GetName()) ;
119 tempo.Append(Version()) ;
120 SetName(tempo.Data()) ;
121
2bd5457f 122 Init() ;
7acf6008 123
124}
7b7c1533 125
7acf6008 126//____________________________________________________________________________
127AliPHOSPIDv1::~AliPHOSPIDv1()
128{
a4e98857 129 //dtor
7acf6008 130}
2bd5457f 131
7acf6008 132//____________________________________________________________________________
7b7c1533 133Bool_t AliPHOSPIDv1::ReadTrackSegments(Int_t event)
7acf6008 134{
f035f6ce 135 // Reads TrackSegments an extracts the title of the RecPoints
a4e98857 136 // branch from which TS were made of.
f035f6ce 137 // Then reads both TrackSegments and RecPoints.
7acf6008 138
f035f6ce 139 //Fist read Track Segment Branch and extract RecPointsBranch from fTSMaker
7acf6008 140
7b7c1533 141 // Get TreeR header from file
142 if(gAlice->TreeR()==0){
143 cerr << "ERROR: AliPHOSPIDv1::ReadTrackSegments -> There is no Reconstruction Tree" << endl;
144 return kFALSE;
7acf6008 145 }
7b7c1533 146
147 // Find TrackSegments
148 TBranch * tsbranch = 0;
149 TBranch * tsmakerbranch = 0;
150 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
151 TIter next(lob) ;
152 TBranch * branch = 0 ;
153 Bool_t phostsfound = kFALSE, tsmakerfound = kFALSE ;
154
155 TString taskName(GetName()) ;
156 taskName.ReplaceAll(Version(), "") ;
7acf6008 157
7b7c1533 158 while ( (branch = (TBranch*)next()) && (!phostsfound || !tsmakerfound) ) {
159 if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
160 phostsfound = kTRUE ;
161 tsbranch = branch ;
7acf6008 162
7b7c1533 163 } else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
164 tsmakerfound = kTRUE ;
165 tsmakerbranch = branch ;
7acf6008 166 }
167 }
7b7c1533 168 if ( !phostsfound || !tsmakerfound ) {
169 cerr << "WARNING: AliPHOSPIDv1::ReadTrackSegments -> TrackSegments and/or TrackSegmentMaker branch with name " << taskName.Data()
170 << " not found" << endl ;
171 return kFALSE ;
172 }
7acf6008 173
7b7c1533 174 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
175
176 TClonesArray * trackSegments = gime->TrackSegments() ;
177 trackSegments->Clear() ;
178 tsbranch->SetAddress(&trackSegments) ;
179
180 AliPHOSTrackSegmentMaker * tsmaker = 0 ;
181 tsmakerbranch->SetAddress(&tsmaker) ;
182 tsmakerbranch->GetEntry(0) ;
183 TString tsmakerName( fRecParticlesTitle ) ;
184 tsmakerName.Append(tsmaker->Version()) ;
185 tsmaker = gime->TrackSegmentMaker(tsmakerName) ;
186
187 tsbranch->GetEntry(0) ;
188 tsmakerbranch->GetEntry(0) ;
189
190 fRecPointsTitle = tsmaker->GetRecPointsBranch() ;
191
192 // Find RecPoints
193 TBranch * emcbranch = 0;
194 TBranch * cpvbranch = 0;
195 TBranch * clusterizerbranch = 0;
196 lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
197 TIter next2(lob) ;
198 branch = 0 ;
199 Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ;
200
201 while ( (branch = (TBranch*)next2()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
202 if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
203 phosemcfound = kTRUE ;
204 emcbranch = branch ;
7acf6008 205 }
7b7c1533 206
207 else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
208 phoscpvfound = kTRUE ;
209 cpvbranch = branch ;
210
211 } else if ( (strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) ) {
212 clusterizerfound = kTRUE ;
213 clusterizerbranch = branch ;
7acf6008 214 }
215 }
7b7c1533 216 if ( !phosemcfound || !phoscpvfound || !clusterizerfound ) {
217 cerr << "WARNING: AliPHOSTrackPIDv1::ReadTrackSegments -> emc(cpv)RecPoints and/or Clusterizer branch with name " << taskName.Data()
218 << " not found" << endl ;
219 return kFALSE ;
220 }
221
222 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
223 emcRecPoints->Clear() ;
224 emcbranch->SetAddress(&emcRecPoints) ;
7acf6008 225
7b7c1533 226 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
227 cpvRecPoints->Clear() ;
228 cpvbranch->SetAddress(&cpvRecPoints) ;
7acf6008 229
7b7c1533 230
231 AliPHOSClusterizer * clusterizer = 0 ;
232 clusterizerbranch->SetAddress(&clusterizer) ;
233 clusterizerbranch->GetEntry(0) ;
234 TString clusterizerName( fTrackSegmentsTitle ) ;
235 clusterizerName.Append(clusterizer->Version()) ;
236 clusterizer = gime->Clusterizer(clusterizerName) ;
237
238 emcbranch->GetEntry(0) ;
239 cpvbranch->GetEntry(0) ;
240 clusterizerbranch->GetEntry(0) ;
241
7acf6008 242 return kTRUE ;
7b7c1533 243
7acf6008 244}
7b7c1533 245
69183710 246//____________________________________________________________________________
7acf6008 247Float_t AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
69183710 248{
249 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
250
7b7c1533 251 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
69183710 252 TVector3 vecEmc ;
7acf6008 253 TVector3 vecCpv ;
254
255 emc->GetLocalPosition(vecEmc) ;
256 cpv->GetLocalPosition(vecCpv) ;
257 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
258
259 // Correct to difference in CPV and EMC position due to different distance to center.
260 // we assume, that particle moves from center
7b7c1533 261 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
262 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
7acf6008 263 dEMC = dEMC / dCPV ;
264 vecCpv = dEMC * vecCpv - vecEmc ;
265 if (Axis == "X") return vecCpv.X();
266 if (Axis == "Y") return vecCpv.Y();
267 if (Axis == "Z") return vecCpv.Z();
268 if (Axis == "R") return vecCpv.Mag();
269 }
270
271 return 100000000 ;
69183710 272}
273
6ad0bfa0 274//____________________________________________________________________________
7acf6008 275void AliPHOSPIDv1::Exec(Option_t * option)
6ad0bfa0 276{
f035f6ce 277 //Steering method
278
7b7c1533 279 if( strcmp(GetName(), "")== 0 )
7acf6008 280 Init() ;
281
7b7c1533 282 if(strstr(option,"tim"))
7acf6008 283 gBenchmark->Start("PHOSPID");
7b7c1533 284
285 if(strstr(option,"print")) {
286 Print("") ;
287 return ;
288 }
289
290 //check, if the branch with name of this" already exits?
291 TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
292 TIter next(lob) ;
293 TBranch * branch = 0 ;
294 Bool_t phospidfound = kFALSE, pidfound = kFALSE ;
295
296 TString taskName(GetName()) ;
297 taskName.ReplaceAll(Version(), "") ;
7acf6008 298
7b7c1533 299 while ( (branch = (TBranch*)next()) && (!phospidfound || !pidfound) ) {
300 if ( (strcmp(branch->GetName(), "PHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
301 phospidfound = kTRUE ;
302
303 else if ( (strcmp(branch->GetName(), "AliPHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
304 pidfound = kTRUE ;
305 }
7acf6008 306
7b7c1533 307 if ( phospidfound || pidfound ) {
308 cerr << "WARNING: AliPHOSPIDv1::Exec -> RecParticles and/or PIDtMaker branch with name "
309 << taskName.Data() << " already exits" << endl ;
310 return ;
311 }
312
313 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
314 Int_t ievent ;
315
316 for(ievent = 0; ievent < nevents; ievent++){
01a599c9 317 gAlice->GetEvent(ievent) ;
318 gAlice->SetEvent(ievent) ;
319
7b7c1533 320 if(!ReadTrackSegments(ievent)) //reads TrackSegments for event ievent
01a599c9 321 continue ;
7b7c1533 322
7acf6008 323 MakeRecParticles() ;
01a599c9 324
7b7c1533 325 WriteRecParticles(ievent);
326
7acf6008 327 if(strstr(option,"deb"))
328 PrintRecParticles(option) ;
329 }
330
331 if(strstr(option,"tim")){
332 gBenchmark->Stop("PHOSPID");
333 cout << "AliPHOSPID:" << endl ;
334 cout << " took " << gBenchmark->GetCpuTime("PHOSPID") << " seconds for PID "
7b7c1533 335 << gBenchmark->GetCpuTime("PHOSPID")/nevents << " seconds per event " << endl ;
7acf6008 336 cout << endl ;
337 }
338
339}
7b7c1533 340//____________________________________________________________________________
341void AliPHOSPIDv1::Init()
342{
343 // Make all memory allocations that are not possible in default constructor
344 // Add the PID task to the list of PHOS tasks
345
346 if ( strcmp(GetTitle(), "") == 0 )
347 SetTitle("galice.root") ;
348
349 TString taskName(GetName()) ;
350 taskName.ReplaceAll(Version(), "") ;
351
352 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data()) ;
353 if ( gime == 0 ) {
354 cerr << "ERROR: AliPHOSPIDv1::Init -> Could not obtain the Getter object !" << endl ;
355 return ;
356 }
357
358 //add Task to //YSAlice/tasks/Reconstructioner/PHOS
359 TTask * aliceRe = (TTask*)gROOT->FindObjectAny("YSAlice/tasks/Reconstructioner") ;
360 TTask * phosRe = (TTask*)aliceRe->GetListOfTasks()->FindObject("PHOS") ;
361 phosRe->Add(this) ;
362 // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
363 gime->Post(GetTitle(), "P", taskName.Data() ) ;
364
365}
366
7acf6008 367//____________________________________________________________________________
368void AliPHOSPIDv1::MakeRecParticles(){
369
b2a60966 370 // Makes a RecParticle out of a TrackSegment
6ad0bfa0 371
7b7c1533 372 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
373 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
374 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
375 TClonesArray * trackSegments = gime->TrackSegments() ;
376 TClonesArray * recParticles = gime->RecParticles() ;
01a599c9 377 recParticles->Clear();
7b7c1533 378
379 TIter next(trackSegments) ;
7acf6008 380 AliPHOSTrackSegment * ts ;
6ad0bfa0 381 Int_t index = 0 ;
09fc14a0 382 AliPHOSRecParticle * rp ;
7acf6008 383
384 Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
385 Bool_t disp = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
386
387 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
fad3e5b9 388
7b7c1533 389 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
390 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
7acf6008 391 rp->SetTraskSegment(index) ;
f035f6ce 392
7acf6008 393 AliPHOSEmcRecPoint * emc = 0 ;
394 if(ts->GetEmcIndex()>=0)
7b7c1533 395 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
fad3e5b9 396
7acf6008 397 AliPHOSRecPoint * cpv = 0 ;
398 if(ts->GetCpvIndex()>=0)
7b7c1533 399 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
fad3e5b9 400
7acf6008 401 AliPHOSRecPoint * ppsd = 0 ;
402 if(ts->GetPpsdIndex()>=0)
7b7c1533 403 ppsd= (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetPpsdIndex()) ;
7acf6008 404
405 //set momentum and energy first
406 Float_t e = emc->GetEnergy() ;
407 TVector3 dir = GetMomentumDirection(emc,cpv,ppsd) ;
408 dir.SetMag(e) ;
409
410 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
411 rp->SetCalcMass(0);
412
413 //now set type (reconstructed) of the particle
414 Int_t showerprofile = 0; // 0 narrow and 1 wide
fad3e5b9 415
7acf6008 416 if(ellips){
417 Float_t lambda[2] ;
418 emc->GetElipsAxis(lambda) ;
419 if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
420 showerprofile = 1 ; // not narrow
1cb7c1ee 421 }
fad3e5b9 422
7acf6008 423 if(disp)
424 if(emc->GetDispersion() > fDispersion )
425 showerprofile = 1 ; // not narrow
426
427
2aad621e 428 // Looking at the photon conversion detector
7acf6008 429 Int_t pcdetector= 0 ; //1 hit and 0 no hit
430 if(ppsd)
431 if(GetDistance(emc, ppsd, "R") < fCpvEmcDistance)
432 pcdetector = 1 ;
433
434 // Looking at the CPV detector
435 Int_t cpvdetector= 0 ; //1 hit and 0 no hit
436 if(cpv)
437 if(GetDistance(emc, cpv, "R") < fCpvEmcDistance)
438 cpvdetector = 1 ;
439
440 Int_t type = showerprofile + 2 * pcdetector + 4 * cpvdetector ;
09fc14a0 441 rp->SetType(type) ;
6ad0bfa0 442 index++ ;
443 }
7acf6008 444
6ad0bfa0 445}
446
09fc14a0 447//____________________________________________________________________________
7acf6008 448void AliPHOSPIDv1:: Print(Option_t * option) const
09fc14a0 449{
b2a60966 450 // Print the parameters used for the particle type identification
f035f6ce 451 cout << "=============== AliPHOSPID1 ================" << endl ;
452 cout << "Making PID "<< endl ;
453 cout << " Headers file: " << fHeaderFileName.Data() << endl ;
454 cout << " RecPoints branch title: " << fRecPointsTitle.Data() << endl ;
7b7c1533 455 cout << " TrackSegments Branch title: " << fTrackSegmentsTitle.Data() << endl ;
456 cout << " RecParticles Branch title " << fRecParticlesTitle.Data() << endl;
f035f6ce 457 cout << "with parameters: " << endl ;
458 cout << " Maximal EMC - CPV (PPSD) distance (cm) " << fCpvEmcDistance << endl ;
459 if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
460 cout << " dispersion cut " << fDispersion << endl ;
461 if(fIDOptions.Contains("ell",TString::kIgnoreCase )){
462 cout << "Eliptic cuts function: " << endl ;
463 cout << fFormula->GetTitle() << endl ;
464 }
465 cout << "============================================" << endl ;
09fc14a0 466}
467
468//____________________________________________________________________________
a4e98857 469void AliPHOSPIDv1::SetShowerProfileCut(char * formula)
470{
7acf6008 471 //set shape of the cut on the axis of ellipce, drown around shouer
472 //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
473 if(fFormula)
474 delete fFormula;
f035f6ce 475 fFormula = new TFormula("Lambda Cut",formula) ;
7acf6008 476}
477//____________________________________________________________________________
7b7c1533 478void AliPHOSPIDv1::WriteRecParticles(Int_t event)
09fc14a0 479{
7b7c1533 480
481 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
482 TClonesArray * recParticles = gime->RecParticles() ;
483
7b7c1533 484 //Make branch in TreeR for RecParticles
7acf6008 485 char * filename = 0;
486 if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){ //generating file name
487 filename = new char[strlen(gAlice->GetBaseFile())+20] ;
488 sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
489 }
490
491 TDirectory *cwd = gDirectory;
492
493 //First rp
494 Int_t bufferSize = 32000 ;
7b7c1533 495 TBranch * rpBranch = gAlice->TreeR()->Branch("PHOSRP",&recParticles,bufferSize);
496 rpBranch->SetTitle(fRecParticlesTitle.Data());
7acf6008 497 if (filename) {
498 rpBranch->SetFile(filename);
499 TIter next( rpBranch->GetListOfBranches());
761e34c0 500 TBranch * sb ;
501 while ((sb=(TBranch*)next())) {
502 sb->SetFile(filename);
7acf6008 503 }
504 cwd->cd();
505 }
506
507 //second, pid
508 Int_t splitlevel = 0 ;
509 AliPHOSPIDv1 * pid = this ;
7b7c1533 510 TBranch * pidBranch = gAlice->TreeR()->Branch("AliPHOSPID","AliPHOSPIDv1",&pid,bufferSize,splitlevel);
511 pidBranch->SetTitle(fRecParticlesTitle.Data());
7acf6008 512 if (filename) {
513 pidBranch->SetFile(filename);
514 TIter next( pidBranch->GetListOfBranches());
761e34c0 515 TBranch * sb ;
516 while ((sb=(TBranch*)next())) {
517 sb->SetFile(filename);
7acf6008 518 }
519 cwd->cd();
520 }
521
761e34c0 522 rpBranch->Fill() ;
523 pidBranch->Fill() ;
eec3ac52 524
7acf6008 525 gAlice->TreeR()->Write(0,kOverwrite) ;
526
527}
7b7c1533 528
69183710 529//____________________________________________________________________________
7acf6008 530void AliPHOSPIDv1::PlotDispersionCuts()const
69183710 531{
a4e98857 532 // produces a plot of the dispersion cut
533 TCanvas* lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
dd5c4038 534
7acf6008 535 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
536 TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
537 ell->SetMinimum(0.0000001) ;
538 ell->SetMaximum(0.001) ;
539 ell->SetLineStyle(1) ;
540 ell->SetLineWidth(2) ;
541 ell->Draw() ;
542 }
543
544 if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
545 TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
546 dsp->SetParameter(0,fDispersion) ;
547 dsp->SetMinimum(0.0000001) ;
548 dsp->SetMaximum(0.001) ;
549 dsp->SetLineStyle(1) ;
550 dsp->SetLineColor(2) ;
551 dsp->SetLineWidth(2) ;
552 dsp->SetNpx(200) ;
553 dsp->SetNpy(200) ;
554 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
555 dsp->Draw("same") ;
556 else
557 dsp->Draw() ;
558 }
559 lambdas->Update();
560}
69183710 561
7acf6008 562//____________________________________________________________________________
563TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv,AliPHOSRecPoint * ppsd)const
564{
565 // Calculates the momentum direction:
566 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
567 // 2. if a EMC RecPoint and one PPSD RecPoint, direction is given by the line through the 2 recpoints
568 // 3. if a EMC RecPoint and two PPSD RecPoints, dirrection is given by the average line through
569 // the 2 pairs of recpoints
570 // However because of the poor position resolution of PPSD the direction is always taken as if we were
571 // in case 1.
572
573 TVector3 dir(0,0,0) ;
574
575 TVector3 emcglobalpos ;
576 TMatrix dummy ;
577
578 emc->GetGlobalPosition(emcglobalpos, dummy) ;
579
69183710 580
a4e98857 581 // The following commented code becomes valid once the PPSD provides
7acf6008 582 // a reasonable position resolution, at least as good as EMC !
583 // TVector3 ppsdlglobalpos ;
584 // TVector3 ppsduglobalpos ;
585 // if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
586 // fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ;
587 // dir = emcglobalpos - ppsdlglobalpos ;
588 // if( fPpsdUpRecPoint ){ // not looks like a charged
589 // fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ;
590 // dir = ( dir + emcglobalpos - ppsduglobalpos ) * 0.5 ;
591 // }
592 // }
593 // else { // looks like a neutral
594 // dir = emcglobalpos ;
595 // }
596
597 dir = emcglobalpos ;
598 dir.SetZ( -dir.Z() ) ; // why ?
599 dir.SetMag(1.) ;
600
601 //account correction to the position of IP
602 Float_t xo,yo,zo ; //Coordinates of the origin
603 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
604 TVector3 origin(xo,yo,zo);
605 dir = dir - origin ;
606
607 return dir ;
608}
609//____________________________________________________________________________
a4e98857 610void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
611{
dd5c4038 612 // Print table of reconstructed particles
613
7b7c1533 614 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
615 TClonesArray * recParticles = gime->RecParticles() ;
616
617
01a599c9 618 cout << "AliPHOSPIDv1: event "<<gAlice->GetEvNumber() << endl ;
7b7c1533 619 cout << " found " << recParticles->GetEntriesFast() << " RecParticles " << endl ;
7acf6008 620
621 if(strstr(option,"all")) { // printing found TS
622
623 cout << " PARTICLE "
624 << " Index " << endl ;
625 // << " X "
626 // << " Y "
627 // << " Z "
628 // << " # of primaries "
629 // << " Primaries list " << endl;
630
631 Int_t index ;
7b7c1533 632 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
633 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
7acf6008 634
635 Text_t particle[11];
636 switch(rp->GetType()) {
637 case AliPHOSFastRecParticle::kNEUTRALEM:
638 strcpy( particle, "NEUTRAL_EM");
639 break;
640 case AliPHOSFastRecParticle::kNEUTRALHA:
641 strcpy(particle, "NEUTRAL_HA");
642 break;
643 case AliPHOSFastRecParticle::kGAMMA:
644 strcpy(particle, "GAMMA");
645 break ;
646 case AliPHOSFastRecParticle::kGAMMAHA:
647 strcpy(particle, "GAMMA_H");
648 break ;
649 case AliPHOSFastRecParticle::kABSURDEM:
650 strcpy(particle, "ABSURD_EM") ;
651 break ;
652 case AliPHOSFastRecParticle::kABSURDHA:
653 strcpy(particle, "ABSURD_HA") ;
654 break ;
655 case AliPHOSFastRecParticle::kELECTRON:
656 strcpy(particle, "ELECTRON") ;
657 break ;
658 case AliPHOSFastRecParticle::kCHARGEDHA:
659 strcpy(particle, "CHARGED_HA") ;
660 break ;
661 }
662
663 // Int_t * primaries;
664 // Int_t nprimaries;
665 // primaries = rp->GetPrimaries(nprimaries);
666
667 cout << setw(15) << particle << " "
668 << setw(3) << rp->GetIndexInList() << " " ;
669 // << setw(4) << nprimaries << " ";
670 // for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
671 // cout << setw(4) << primaries[iprimary] << " ";
672 cout << endl;
673 }
674 cout << "-------------------------------------------" << endl ;
675 }
676
69183710 677}
678
7acf6008 679
680