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 //_________________________________________________________________________
19 // Algorythm class to analyze PHOSv0 events:
20 // Construct histograms and displays them.
21 // Use the macro EditorBar.C for best access to the functionnalities
23 //*-- Author: Y. Schutz (SUBATECH)
24 //////////////////////////////////////////////////////////////////////////////
26 // --- ROOT system ---
32 #include "TParticle.h"
33 #include "TClonesArray.h"
38 // --- Standard library ---
43 // --- AliRoot header files ---
46 #include "AliPHOSAnalyze.h"
47 #include "AliPHOSClusterizerv1.h"
48 #include "AliPHOSTrackSegmentMakerv1.h"
49 #include "AliPHOSPIDv1.h"
50 #include "AliPHOSReconstructioner.h"
51 #include "AliPHOSDigit.h"
52 #include "AliPHOSTrackSegment.h"
53 #include "AliPHOSRecParticle.h"
55 ClassImp(AliPHOSAnalyze)
58 //____________________________________________________________________________
59 AliPHOSAnalyze::AliPHOSAnalyze()
61 // default ctor (useless)
66 //____________________________________________________________________________
67 AliPHOSAnalyze::AliPHOSAnalyze(Text_t * name)
69 // ctor: analyze events from root file: name
71 Bool_t ok = OpenRootFile(name) ;
73 cout << " AliPHOSAnalyze > Error opening " << name << endl ;
76 gAlice = (AliRun*) fRootFile->Get("gAlice");
77 fPHOS = (AliPHOSv0 *)gAlice->GetDetector("PHOS") ;
78 fGeom = AliPHOSGeometry::GetInstance( fPHOS->GetGeometry()->GetName(), fPHOS->GetGeometry()->GetTitle() ) ;
83 //____________________________________________________________________________
84 AliPHOSAnalyze::~AliPHOSAnalyze()
109 //____________________________________________________________________________
110 void AliPHOSAnalyze::AnalyzeOneEvent(Int_t evt)
112 // analyze one single event with id=evt
116 Bool_t ok = Init(evt) ;
119 //=========== Get the number of entries in the Digits array
121 Int_t nId = fPHOS->Digits()->GetEntries();
122 printf("AnalyzeOneEvent > Number of entries in the Digit array is %d \n",nId);
124 //=========== Do the reconstruction
126 cout << "AnalyzeOneEvent > Found " << nId << " digits in PHOS" << endl ;
129 fPHOS->Reconstruction(fRec);
131 // =========== End of reconstruction
133 // =========== Write the root file
137 // =========== Finish
139 cout << "AnalyzeOneEvent > event # " << fEvt << " processed" << endl ;
142 cout << "AnalyzeOneEvent > filed to process event # " << evt << endl ;
144 ts.Stop() ; cout << "CPU time = " << ts.CpuTime() << endl ;
145 cout << "Real time = " << ts.RealTime() << endl ;
148 //____________________________________________________________________________
149 void AliPHOSAnalyze::AnalyzeManyEvents(Int_t Nevents, Int_t module)
151 // analyzes Nevents events in a single PHOS module
153 if ( fRootFile == 0 )
154 cout << "AnalyzeManyEvents > " << "Root File not openned" << endl ;
157 //========== Get AliRun object from file
158 gAlice = (AliRun*) fRootFile->Get("gAlice") ;
159 //=========== Get the PHOS object and associated geometry from the file
160 fPHOS = (AliPHOSv0 *)gAlice->GetDetector("PHOS") ;
161 fGeom = AliPHOSGeometry::GetInstance( fPHOS->GetGeometry()->GetName(), fPHOS->GetGeometry()->GetTitle() );
162 //========== Booking Histograms
163 cout << "AnalyzeManyEvents > " << "Booking Histograms" << endl ;
167 AliPHOSDigit * digit ;
168 AliPHOSEmcRecPoint * emc ;
169 AliPHOSPpsdRecPoint * ppsd ;
170 // AliPHOSTrackSegment * tracksegment ;
171 AliPHOSRecParticle * recparticle;
172 for ( ievent=0; ievent<Nevents; ievent++)
174 if (ievent==0) cout << "AnalyzeManyEvents > " << "Starting Analyzing " << endl ;
175 //========== Create the Clusterizer
176 fClu = new AliPHOSClusterizerv1() ;
177 fClu->SetEmcEnergyThreshold(0.025) ;
178 fClu->SetEmcClusteringThreshold(0.50) ;
179 fClu->SetPpsdEnergyThreshold (0.0000002) ;
180 fClu->SetPpsdClusteringThreshold(0.0000001) ;
181 fClu->SetLocalMaxCut(0.03) ;
182 fClu->SetCalibrationParameters(0., 0.00000001) ;
183 //========== Creates the track segment maker
184 fTrs = new AliPHOSTrackSegmentMakerv1() ;
185 fTrs->UnsetUnfoldFlag() ;
186 //========== Creates the particle identifier
187 fPID = new AliPHOSPIDv1() ;
188 fPID->SetShowerProfileCuts(0.3, 1.8, 0.3, 1.8 ) ;
190 //========== Creates the Reconstructioner
191 fRec = new AliPHOSReconstructioner(fClu, fTrs, fPID) ;
192 //========== Event Number
193 if ( ( log10(ievent+1) - (Int_t)(log10(ievent+1)) ) == 0. ) cout << "AnalyzeManyEvents > " << "Event is " << ievent << endl ;
194 //=========== Connects the various Tree's for evt
195 gAlice->GetEvent(ievent);
196 //=========== Gets the Digit TTree
197 gAlice->TreeD()->GetEvent(0) ;
198 //=========== Gets the number of entries in the Digits array
199 TIter nextdigit(fPHOS->Digits()) ;
200 while( ( digit = (AliPHOSDigit *)nextdigit() ) )
202 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
203 if (fClu->IsInEmc(digit)) fhEmcDigit->Fill(fClu->Calibrate(digit->GetAmp())) ;
206 if (relid[1]<17) fhVetoDigit->Fill(fClu->Calibrate(digit->GetAmp()));
207 if (relid[1]>16) fhConvertorDigit->Fill(fClu->Calibrate(digit->GetAmp()));
210 //=========== Do the reconstruction
211 fPHOS->Reconstruction(fRec);
212 //=========== Cluster in module
213 TIter nextEmc(fPHOS->EmcClusters() ) ;
214 while((emc = (AliPHOSEmcRecPoint *)nextEmc()))
216 if ( emc->GetPHOSMod() == module )
218 fhEmcCluster->Fill( emc->GetTotalEnergy() );
219 TIter nextPpsd( fPHOS->PpsdClusters()) ;
220 while((ppsd = (AliPHOSPpsdRecPoint *)nextPpsd()))
222 if ( ppsd->GetPHOSMod() == module )
224 if (!ppsd->GetUp()) fhConvertorEmc->Fill(emc->GetTotalEnergy(),ppsd->GetTotalEnergy()) ;
229 //=========== Cluster in module PPSD Down
230 TIter nextPpsd(fPHOS->PpsdClusters() ) ;
231 while((ppsd = (AliPHOSPpsdRecPoint *)nextPpsd()))
233 if ( ppsd->GetPHOSMod() == module )
235 if (!ppsd->GetUp()) fhConvertorCluster->Fill(ppsd->GetTotalEnergy()) ;
236 if (ppsd->GetUp()) fhVetoCluster ->Fill(ppsd->GetTotalEnergy()) ;
239 //========== TRackSegments in the event
240 TIter nextRecParticle(fPHOS->RecParticles() ) ;
241 while((recparticle = (AliPHOSRecParticle *)nextRecParticle()))
243 if ( recparticle->GetPHOSTrackSegment()->GetPHOSMod() == module )
245 cout << "Particle type is " << recparticle->GetType() << endl ;
246 Int_t numberofprimaries = 0 ;
247 Int_t * listofprimaries = recparticle->GetPrimaries(numberofprimaries) ;
248 cout << "Number of primaries = " << numberofprimaries << endl ;
250 for ( index = 0 ; index < numberofprimaries ; index++)
251 cout << " primary # " << index << " = " << listofprimaries[index] << endl ;
252 switch(recparticle->GetType())
255 fhPhotonEnergy->Fill(recparticle->Energy() ) ;
256 //fhPhotonPositionX->Fill(recpart. ) ;
257 //fhPhotonPositionY->Fill(recpart. ) ;
258 cout << "PHOTON" << endl;
261 fhElectronEnergy->Fill(recparticle->Energy() ) ;
262 //fhElectronPositionX->Fill(recpart. ) ;
263 //fhElectronPositionY->Fill(recpart. ) ;
264 cout << "ELECTRON" << endl;
267 fhNeutralHadronEnergy->Fill(recparticle->Energy() ) ;
268 //fhNeutralHadronPositionX->Fill(recpart. ) ;
269 //fhNeutralHadronPositionY->Fill(recpart. ) ;
270 cout << "NEUTRAl HADRON" << endl;
273 fhNeutralEMEnergy->Fill(recparticle->Energy() ) ;
274 //fhNeutralEMPositionX->Fill(recpart. ) ;
275 //fhNeutralEMPositionY->Fill(recpart. ) ;
276 //cout << "NEUTRAL EM" << endl;
279 fhChargedHadronEnergy->Fill(recparticle->Energy() ) ;
280 //fhChargedHadronPositionX->Fill(recpart. ) ;
281 //fhChargedHadronPositionY->Fill(recpart. ) ;
282 cout << "CHARGED HADRON" << endl;
285 fhPhotonHadronEnergy->Fill(recparticle->Energy() ) ;
286 //fhPhotonHadronPositionX->Fill(recpart. ) ;
287 //fhPhotonHadronPositionY->Fill(recpart. ) ;
288 cout << "PHOTON HADRON" << endl;
293 // Deleting fClu, fTrs, fPID et fRec
305 //____________________________________________________________________________
306 void AliPHOSAnalyze::BookingHistograms()
308 // Books the histograms where the results of the analysis are stored (to be changed)
314 if (fhConvertorDigit )
315 delete fhConvertorDigit ;
317 delete fhEmcCluster ;
319 delete fhVetoCluster ;
320 if (fhConvertorCluster )
321 delete fhConvertorCluster ;
323 delete fhConvertorEmc ;
325 fhEmcDigit = new TH1F("hEmcDigit", "hEmcDigit", 1000, 0. , 25.);
326 fhVetoDigit = new TH1F("hVetoDigit", "hVetoDigit", 500, 0. , 3.e-5);
327 fhConvertorDigit = new TH1F("hConvertorDigit","hConvertorDigit", 500, 0. , 3.e-5);
328 fhEmcCluster = new TH1F("hEmcCluster", "hEmcCluster", 1000, 0. , 30.);
329 fhVetoCluster = new TH1F("hVetoCluster", "hVetoCluster", 500, 0. , 3.e-5);
330 fhConvertorCluster = new TH1F("hConvertorCluster","hConvertorCluster",500, 0. , 3.e-5);
331 fhConvertorEmc = new TH2F("hConvertorEmc", "hConvertorEmc", 200, 1. , 3., 200, 0., 3.e-5);
332 fhPhotonEnergy = new TH1F("hPhotonEnergy", "hPhotonEnergy", 1000, 0. , 30.);
333 fhElectronEnergy = new TH1F("hElectronEnergy","hElectronEnergy", 1000, 0. , 30.);
334 fhNeutralHadronEnergy = new TH1F("hNeutralHadronEnergy", "hNeutralHadronEnergy", 1000, 0. , 30.);
335 fhNeutralEMEnergy = new TH1F("hNeutralEMEnergy", "hNeutralEMEnergy", 1000, 0. , 30.);
336 fhChargedHadronEnergy = new TH1F("hChargedHadronEnergy", "hChargedHadronEnergy", 1000, 0. , 30.);
337 fhPhotonHadronEnergy = new TH1F("hPhotonHadronEnergy","hPhotonHadronEnergy",500,-80. , 80.);
338 fhPhotonPositionX = new TH1F("hPhotonPositionX","hPhotonPositionX", 500,-80. , 80.);
339 fhElectronPositionX = new TH1F("hElectronPositionX","hElectronPositionX",500,-80. , 80.);
340 fhNeutralHadronPositionX = new TH1F("hNeutralHadronPositionX","hNeutralHadronPositionX",500,-80. , 80.);
341 fhNeutralEMPositionX = new TH1F("hNeutralEMPositionX","hNeutralEMPositionX",500,-80. , 80.);
342 fhChargedHadronPositionX = new TH1F("hChargedHadronPositionX","hChargedHadronPositionX",500,-80. , 80.);
343 fhPhotonHadronPositionX = new TH1F("hPhotonHadronPositionX","hPhotonHadronPositionX",500,-80. , 80.);
344 fhPhotonPositionY = new TH1F("hPhotonPositionY","hPhotonPositionY", 500,-80. , 80.);
345 fhElectronPositionY = new TH1F("hElectronPositionY","hElectronPositionY",500,-80. , 80.);
346 fhNeutralHadronPositionY = new TH1F("hNeutralHadronPositionY","hNeutralHadronPositionY",500,-80. , 80.);
347 fhNeutralEMPositionY = new TH1F("hNeutralEMPositionY","hNeutralEMPositionY",500,-80. , 80.);
348 fhChargedHadronPositionY = new TH1F("hChargedHadronPositionY","hChargedHadronPositionY",500,-80. , 80.);
349 fhPhotonHadronPositionY = new TH1F("hPhotonHadronPositionY","hPhotonHadronPositionY",500,-80. , 80.);
353 //____________________________________________________________________________
354 Bool_t AliPHOSAnalyze::Init(Int_t evt)
356 // Do a few initializations: open the root file
357 // get the AliRun object
358 // defines the clusterizer, tracksegment maker and particle identifier
359 // sets the associated parameters
363 //========== Open galice root file
365 if ( fRootFile == 0 ) {
366 Text_t * name = new Text_t[80] ;
367 cout << "AnalyzeOneEvent > Enter file root file name : " ;
369 Bool_t ok = OpenRootFile(name) ;
371 cout << " AliPHOSAnalyze > Error opening " << name << endl ;
373 //========== Get AliRun object from file
375 gAlice = (AliRun*) fRootFile->Get("gAlice") ;
377 //=========== Get the PHOS object and associated geometry from the file
379 fPHOS = (AliPHOSv0 *)gAlice->GetDetector("PHOS") ;
380 fGeom = AliPHOSGeometry::GetInstance( fPHOS->GetGeometry()->GetName(), fPHOS->GetGeometry()->GetTitle() );
386 //========== Create the Clusterizer
388 fClu = new AliPHOSClusterizerv1() ;
389 fClu->SetEmcEnergyThreshold(0.030) ;
390 fClu->SetEmcClusteringThreshold(1.0) ;
391 fClu->SetPpsdEnergyThreshold (0.0000002) ;
392 fClu->SetPpsdClusteringThreshold(0.0000001) ;
393 fClu->SetLocalMaxCut(0.03) ;
394 fClu->SetCalibrationParameters(0., 0.00000001) ;
395 cout << "AnalyzeOneEvent > using clusterizer " << fClu->GetName() << endl ;
396 fClu->PrintParameters() ;
398 //========== Creates the track segment maker
400 fTrs = new AliPHOSTrackSegmentMakerv1() ;
401 cout << "AnalyzeOneEvent > using tack segment maker " << fTrs->GetName() << endl ;
402 fTrs->UnsetUnfoldFlag() ;
404 //========== Creates the particle identifier
406 fPID = new AliPHOSPIDv1() ;
407 cout << "AnalyzeOneEvent > using particle identifier " << fPID->GetName() << endl ;
409 //========== Creates the Reconstructioner
411 fRec = new AliPHOSReconstructioner(fClu, fTrs, fPID) ;
413 //=========== Connect the various Tree's for evt
416 cout << "AnalyzeOneEvent > Enter event number : " ;
418 cout << evt << endl ;
421 gAlice->GetEvent(evt);
423 //=========== Get the Digit TTree
425 gAlice->TreeD()->GetEvent(0) ;
433 //____________________________________________________________________________
434 void AliPHOSAnalyze::DisplayKineEvent(Int_t evt)
436 // Display particles from the Kine Tree in global Alice (theta, phi) coordinates.
437 // One PHOS module at the time.
438 // The particle type can be selected.
444 cout << "DisplayKineEvent > which module (1-5, -1: all) ? " ;
445 cin >> module ; cout << module << endl ;
448 cout << " 22 : PHOTON " << endl
449 << " (-)11 : (POSITRON)ELECTRON " << endl
450 << " (-)2112 : (ANTI)NEUTRON " << endl
451 << " -999 : Everything else " << endl ;
452 cout << "DisplayKineEvent > enter PDG particle code to display " ;
453 cin >> testparticle ; cout << testparticle << endl ;
455 Text_t histoname[80] ;
456 sprintf(histoname,"Event %d: Incident particles in module %d", evt, module) ;
458 Double_t tm, tM, pm, pM ; // min and Max theta and phi covered by module
459 fGeom->EmcModuleCoverage(module, tm, tM, pm, pM, kDegre) ;
461 Double_t theta, phi ;
462 fGeom->EmcXtalCoverage(theta, phi, kDegre) ;
464 Int_t tdim = (Int_t)( (tM - tm) / theta ) ;
465 Int_t pdim = (Int_t)( (pM - pm) / phi ) ;
472 TH2F * histoparticle = new TH2F("histoparticle", histoname,
473 pdim, pm, pM, tdim, tm, tM) ;
474 histoparticle->SetStats(kFALSE) ;
476 // Get pointers to Alice Particle TClonesArray
478 TParticle * particle;
479 TClonesArray * particlearray = gAlice->Particles();
481 Text_t canvasname[80];
482 sprintf(canvasname,"Particles incident in PHOS/EMC module # %d",module) ;
483 TCanvas * kinecanvas = new TCanvas("kinecanvas", canvasname, 650, 500) ;
487 TTree * kine = gAlice->TreeK() ;
488 Stat_t nParticles = kine->GetEntries() ;
489 cout << "DisplayKineEvent > events in kine " << nParticles << endl ;
491 // loop over particles
493 Double_t kRADDEG = 180. / TMath::Pi() ;
495 Int_t nparticlein = 0 ;
496 for (index1 = 0 ; index1 < nParticles ; index1++){
497 Int_t nparticle = particlearray->GetEntriesFast() ;
499 for( index2 = 0 ; index2 < nparticle ; index2++) {
500 particle = (TParticle*)particlearray->UncheckedAt(index2) ;
501 Int_t particletype = particle->GetPdgCode() ;
502 if (testparticle == -999 || testparticle == particletype) {
503 Double_t phi = particle->Phi() ;
504 Double_t theta = particle->Theta() ;
507 fGeom->ImpactOnEmc(theta, phi, mod, z, x) ;
508 if ( mod == module ) {
510 if (particle->Energy() > fClu->GetEmcClusteringThreshold() )
511 histoparticle->Fill(phi*kRADDEG, theta*kRADDEG, particle->Energy() ) ;
517 histoparticle->Draw("color") ;
518 TPaveText * pavetext = new TPaveText(294, 100, 300, 101);
520 sprintf(text, "Particles: %d ", nparticlein) ;
521 pavetext->AddText(text) ;
523 kinecanvas->Update();
526 //____________________________________________________________________________
527 void AliPHOSAnalyze::DisplayRecParticles()
529 // Display reconstructed particles in global Alice(theta, phi) coordinates.
530 // One PHOS module at the time.
531 // Click on symbols indicate the reconstructed particle type.
534 cout << "DisplayRecParticles > Analyze an event first ... (y/n) " ;
536 cin >> answer ; cout << answer ;
543 cout << "DisplayRecParticles > which module (1-5, -1: all) ? " ;
544 cin >> module ; cout << module << endl ;
545 Text_t histoname[80] ;
546 sprintf(histoname,"Event %d: Reconstructed particles in module %d", fEvt, module) ;
547 Double_t tm, tM, pm, pM ; // min and Max theta and phi covered by module
548 fGeom->EmcModuleCoverage(module, tm, tM, pm, pM, kDegre) ;
549 Double_t theta, phi ;
550 fGeom->EmcXtalCoverage(theta, phi, kDegre) ;
551 Int_t tdim = (Int_t)( (tM - tm) / theta ) ;
552 Int_t pdim = (Int_t)( (pM - pm) / phi ) ;
556 TH2F * histoRparticle = new TH2F("histoRparticle", histoname,
557 pdim, pm, pM, tdim, tm, tM) ;
558 histoRparticle->SetStats(kFALSE) ;
559 Text_t canvasname[80] ;
560 sprintf(canvasname, "Reconstructed particles in PHOSmodule # %d", module) ;
561 TCanvas * rparticlecanvas = new TCanvas("RparticleCanvas", canvasname, 650, 500) ;
562 RecParticlesList * rpl = fPHOS->RecParticles() ;
563 Int_t nRecParticles = rpl->GetEntries() ;
564 Int_t nRecParticlesInModule = 0 ;
565 TIter nextRecPart(rpl) ;
566 AliPHOSRecParticle * rp ;
567 cout << "DisplayRecParticles > " << nRecParticles << " reconstructed particles " << endl ;
568 Double_t kRADDEG = 180. / TMath::Pi() ;
569 while ( (rp = (AliPHOSRecParticle *)nextRecPart() ) ) {
570 AliPHOSTrackSegment * ts = rp->GetPHOSTrackSegment() ;
571 if ( ts->GetPHOSMod() == module ) {
572 Int_t numberofprimaries = 0 ;
573 Int_t * listofprimaries = rp->GetPrimaries(numberofprimaries) ;
574 cout << "Number of primaries = " << numberofprimaries << endl ;
576 for ( index = 0 ; index < numberofprimaries ; index++)
577 cout << " primary # " << index << " = " << listofprimaries[index] << endl ;
579 nRecParticlesInModule++ ;
580 Double_t theta = rp->Theta() * kRADDEG ;
581 Double_t phi = rp->Phi() * kRADDEG ;
582 Double_t energy = rp->Energy() ;
583 histoRparticle->Fill(phi, theta, energy) ;
586 histoRparticle->Draw("color") ;
588 nextRecPart.Reset() ;
589 while ( (rp = (AliPHOSRecParticle *)nextRecPart() ) ) {
590 AliPHOSTrackSegment * ts = rp->GetPHOSTrackSegment() ;
591 if ( ts->GetPHOSMod() == module )
596 sprintf(text, "reconstructed particles: %d", nRecParticlesInModule) ;
597 TPaveText * pavetext = new TPaveText(292, 100, 300, 101);
598 pavetext->AddText(text) ;
600 rparticlecanvas->Update() ;
604 //____________________________________________________________________________
605 void AliPHOSAnalyze::DisplayRecPoints()
607 // Display reconstructed points in local PHOS-module (x, z) coordinates.
608 // One PHOS module at the time.
609 // Click on symbols displays the EMC cluster, or PPSD information.
612 cout << "DisplayRecPoints > Analyze an event first ... (y/n) " ;
614 cin >> answer ; cout << answer ;
621 cout << "DisplayRecPoints > which module (1-5, -1: all) ? " ;
622 cin >> module ; cout << module << endl ;
624 Text_t canvasname[80];
625 sprintf(canvasname,"Digits in PHOS/EMC module # %d",module) ;
626 TCanvas * modulecanvas = new TCanvas("module", canvasname, 650, 500) ;
627 modulecanvas->Draw() ;
629 //=========== Creating 2d-histogram of the PHOS module
630 // a little bit junkie but is used to test Geom functinalities
632 Double_t tm, tM, pm, pM ; // min and Max theta and phi covered by module
634 fGeom->EmcModuleCoverage(module, tm, tM, pm, pM);
635 // convert angles into coordinates local to the EMC module of interest
637 Int_t emcModuleNumber ;
638 Double_t emcModulexm, emcModulezm ; // minimum local coordinate in a given EMCA module
639 Double_t emcModulexM, emcModulezM ; // maximum local coordinate in a given EMCA module
640 fGeom->ImpactOnEmc(tm, pm, emcModuleNumber, emcModulezm, emcModulexm) ;
641 fGeom->ImpactOnEmc(tM, pM, emcModuleNumber, emcModulezM, emcModulexM) ;
642 Int_t xdim = (Int_t)( ( emcModulexM - emcModulexm ) / fGeom->GetCrystalSize(0) ) ;
643 Int_t zdim = (Int_t)( ( emcModulezM - emcModulezm ) / fGeom->GetCrystalSize(2) ) ;
644 Float_t xmin = emcModulexm - fGeom->GetCrystalSize(0) ;
645 Float_t xMax = emcModulexM + fGeom->GetCrystalSize(0) ;
646 Float_t zmin = emcModulezm - fGeom->GetCrystalSize(2) ;
647 Float_t zMax = emcModulezM + fGeom->GetCrystalSize(2) ;
648 Text_t histoname[80];
649 sprintf(histoname,"Event %d: Digits and RecPoints in module %d", fEvt, module) ;
650 TH2F * hModule = new TH2F("HistoReconstructed", histoname,
651 xdim, xmin, xMax, zdim, zmin, zMax) ;
652 hModule->SetMaximum(2.0);
653 hModule->SetMinimum(0.0);
654 hModule->SetStats(kFALSE);
656 TIter next(fPHOS->Digits()) ;
657 Float_t energy, y, z;
659 Int_t relid[4]; Int_t nDigits = 0 ;
660 AliPHOSDigit * digit ;
662 // Making 2D histogram of the EMC module
663 while((digit = (AliPHOSDigit *)next()))
665 fGeom->AbsToRelNumbering(digit->GetId(), relid) ;
666 if (relid[0] == module && relid[1] == 0)
668 energy = fClu->Calibrate(digit->GetAmp()) ;
669 cout << "Energy is " << energy << " and threshold is " << fClu->GetEmcEnergyThreshold() << endl;
670 if (energy > fClu->GetEmcEnergyThreshold() ){
673 fGeom->RelPosInModule(relid,y,z) ;
674 hModule->Fill(y, z, energy) ;
678 cout <<"DrawRecPoints > Found in module "
679 << module << " " << nDigits << " digits with total energy " << etot << endl ;
680 hModule->Draw("col2") ;
682 //=========== Cluster in module
684 TClonesArray * emcRP = fPHOS->EmcClusters() ;
686 Int_t totalnClusters = 0 ;
687 Int_t nClusters = 0 ;
688 TIter nextemc(emcRP) ;
689 AliPHOSEmcRecPoint * emc ;
690 while((emc = (AliPHOSEmcRecPoint *)nextemc()))
692 // Int_t numberofprimaries ;
693 // Int_t * primariesarray = new Int_t[10] ;
694 // emc->GetPrimaries(numberofprimaries, primariesarray) ;
696 if ( emc->GetPHOSMod() == module )
699 energy = emc->GetTotalEnergy() ;
704 cout << "DrawRecPoints > Found " << totalnClusters << " EMC Clusters in PHOS" << endl ;
705 cout << "DrawRecPoints > Found in module " << module << " " << nClusters << " EMC Clusters " << endl ;
706 cout << "DrawRecPoints > total energy " << etot << endl ;
708 TPaveText * pavetext = new TPaveText(22, 80, 83, 90);
710 sprintf(text, "digits: %d; clusters: %d", nDigits, nClusters) ;
711 pavetext->AddText(text) ;
713 modulecanvas->Update();
715 //=========== Cluster in module PPSD Down
717 TClonesArray * ppsdRP = fPHOS->PpsdClusters() ;
719 TIter nextPpsd(ppsdRP) ;
720 AliPHOSPpsdRecPoint * ppsd ;
721 while((ppsd = (AliPHOSPpsdRecPoint *)nextPpsd()))
724 if ( ppsd->GetPHOSMod() == module )
727 energy = ppsd->GetEnergy() ;
729 if (!ppsd->GetUp()) ppsd->Draw("P") ;
732 cout << "DrawRecPoints > Found " << totalnClusters << " Ppsd Down Clusters in PHOS" << endl ;
733 cout << "DrawRecPoints > Found in module " << module << " " << nClusters << " Ppsd Down Clusters " << endl ;
734 cout << "DrawRecPoints > total energy " << etot << endl ;
736 //=========== Cluster in module PPSD Up
738 ppsdRP = fPHOS->PpsdClusters() ;
740 TIter nextPpsdUp(ppsdRP) ;
741 while((ppsd = (AliPHOSPpsdRecPoint *)nextPpsdUp()))
744 if ( ppsd->GetPHOSMod() == module )
747 energy = ppsd->GetEnergy() ;
749 if (ppsd->GetUp()) ppsd->Draw("P") ;
752 cout << "DrawRecPoints > Found " << totalnClusters << " Ppsd Up Clusters in PHOS" << endl ;
753 cout << "DrawRecPoints > Found in module " << module << " " << nClusters << " Ppsd Up Clusters " << endl ;
754 cout << "DrawRecPoints > total energy " << etot << endl ;
759 //____________________________________________________________________________
760 void AliPHOSAnalyze::DisplayTrackSegments()
762 // Display track segments in local PHOS-module (x, z) coordinates.
763 // One PHOS module at the time.
764 // One symbol per PHOS subsystem: EMC, upper PPSD, lower PPSD.
767 cout << "DisplayTrackSegments > Analyze an event first ... (y/n) " ;
769 cin >> answer ; cout << answer ;
776 cout << "DisplayTrackSegments > which module (1-5, -1: all) ? " ;
777 cin >> module ; cout << module << endl ;
778 //=========== Creating 2d-histogram of the PHOS module
779 // a little bit junkie but is used to test Geom functinalities
781 Double_t tm, tM, pm, pM ; // min and Max theta and phi covered by module
783 fGeom->EmcModuleCoverage(module, tm, tM, pm, pM);
784 // convert angles into coordinates local to the EMC module of interest
786 Int_t emcModuleNumber ;
787 Double_t emcModulexm, emcModulezm ; // minimum local coordinate in a given EMCA module
788 Double_t emcModulexM, emcModulezM ; // maximum local coordinate in a given EMCA module
789 fGeom->ImpactOnEmc(tm, pm, emcModuleNumber, emcModulezm, emcModulexm) ;
790 fGeom->ImpactOnEmc(tM, pM, emcModuleNumber, emcModulezM, emcModulexM) ;
791 Int_t xdim = (Int_t)( ( emcModulexM - emcModulexm ) / fGeom->GetCrystalSize(0) ) ;
792 Int_t zdim = (Int_t)( ( emcModulezM - emcModulezm ) / fGeom->GetCrystalSize(2) ) ;
793 Float_t xmin = emcModulexm - fGeom->GetCrystalSize(0) ;
794 Float_t xMax = emcModulexM + fGeom->GetCrystalSize(0) ;
795 Float_t zmin = emcModulezm - fGeom->GetCrystalSize(2) ;
796 Float_t zMax = emcModulezM + fGeom->GetCrystalSize(2) ;
797 Text_t histoname[80];
798 sprintf(histoname,"Event %d: Track Segments in module %d", fEvt, module) ;
799 TH2F * histotrack = new TH2F("histotrack", histoname,
800 xdim, xmin, xMax, zdim, zmin, zMax) ;
801 histotrack->SetStats(kFALSE);
802 Text_t canvasname[80];
803 sprintf(canvasname,"Track segments in PHOS/EMC-PPSD module # %d", module) ;
804 TCanvas * trackcanvas = new TCanvas("TrackSegmentCanvas", canvasname, 650, 500) ;
807 TrackSegmentsList * trsegl = fPHOS->TrackSegments() ;
808 AliPHOSTrackSegment * trseg ;
810 Int_t nTrackSegments = trsegl->GetEntries() ;
813 Int_t nTrackSegmentsInModule = 0 ;
814 for(index = 0; index < nTrackSegments ; index++){
815 trseg = (AliPHOSTrackSegment * )trsegl->At(index) ;
816 etot+= trseg->GetEnergy() ;
817 if ( trseg->GetPHOSMod() == module ) {
818 nTrackSegmentsInModule++ ;
823 sprintf(text, "track segments: %d", nTrackSegmentsInModule) ;
824 TPaveText * pavetext = new TPaveText(22, 80, 83, 90);
825 pavetext->AddText(text) ;
827 trackcanvas->Update() ;
828 cout << "DisplayTrackSegments > Found " << trsegl->GetEntries() << " Track segments with total energy "<< etot << endl ;
832 //____________________________________________________________________________
833 Bool_t AliPHOSAnalyze::OpenRootFile(Text_t * name)
835 // Open the root file named "name"
837 fRootFile = new TFile(name, "update") ;
838 return fRootFile->IsOpen() ;
840 //____________________________________________________________________________
841 void AliPHOSAnalyze::SavingHistograms()
843 // Saves the histograms in a root file named "name.analyzed"
845 Text_t outputname[80] ;
846 sprintf(outputname,"%s.analyzed",fRootFile->GetName());
847 TFile output(outputname,"RECREATE");
850 fhEmcDigit->Write() ;
852 fhVetoDigit->Write() ;
853 if (fhConvertorDigit )
854 fhConvertorDigit->Write() ;
856 fhEmcCluster->Write() ;
858 fhVetoCluster->Write() ;
859 if (fhConvertorCluster )
860 fhConvertorCluster->Write() ;
862 fhConvertorEmc->Write() ;
864 fhPhotonEnergy->Write() ;
865 if (fhPhotonPositionX)
866 fhPhotonPositionX->Write() ;
867 if (fhPhotonPositionY)
868 fhPhotonPositionX->Write() ;
869 if (fhElectronEnergy)
870 fhElectronEnergy->Write() ;
871 if (fhElectronPositionX)
872 fhElectronPositionX->Write() ;
873 if (fhElectronPositionY)
874 fhElectronPositionX->Write() ;
875 if (fhNeutralHadronEnergy)
876 fhNeutralHadronEnergy->Write() ;
877 if (fhNeutralHadronPositionX)
878 fhNeutralHadronPositionX->Write() ;
879 if (fhNeutralHadronPositionY)
880 fhNeutralHadronPositionX->Write() ;
881 if (fhNeutralEMEnergy)
882 fhNeutralEMEnergy->Write() ;
883 if (fhNeutralEMPositionX)
884 fhNeutralEMPositionX->Write() ;
885 if (fhNeutralEMPositionY)
886 fhNeutralEMPositionX->Write() ;
887 if (fhChargedHadronEnergy)
888 fhChargedHadronEnergy->Write() ;
889 if (fhChargedHadronPositionX)
890 fhChargedHadronPositionX->Write() ;
891 if (fhChargedHadronPositionY)
892 fhChargedHadronPositionX->Write() ;
893 if (fhPhotonHadronEnergy)
894 fhPhotonHadronEnergy->Write() ;
895 if (fhPhotonHadronPositionX)
896 fhPhotonHadronPositionX->Write() ;
897 if (fhPhotonHadronPositionY)
898 fhPhotonHadronPositionX->Write() ;