]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSHistos.cxx
Classes for online analysis added
[u/mrichter/AliRoot.git] / PHOS / PHOSHistos.cxx
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 //_________________________________________________________________________
16 // Cheking PHOSHistos procedure of PHOS
17 //*-- Author : Gines MARTINEZ  SUBATECH january 2000
18 //////////////////////////////////////////////////////////////////////////////
19
20 // --- ROOT system ---
21 #include "TCanvas.h"
22 #include "TFile.h"
23 #include "TH1.h"
24 #include "TPad.h"
25 #include "TTree.h"
26
27 // --- Standard library ---
28
29 #include <cstdio>
30
31
32 // --- AliRoot header files ---
33 #include "AliRun.h"
34 #include "TFile.h"
35 #include "AliPHOSGeometry.h"
36 #include "AliPHOSv0.h"
37 #include "AliPHOSDigit.h"
38 #include "AliPHOSRecPoint.h"
39 #include "AliPHOSEmcRecPoint.h"
40 #include "AliPHOSPpsdRecPoint.h"
41 #include "AliPHOSClusterizerv1.h"
42 #include "AliPHOSReconstructor.h"
43 #include "AliPHOSTrackSegment.h"
44 #include "AliPHOSTrackSegmentMakerv1.h"
45 #include "AliPHOSPIDv1.h"
46 #include "PHOSHistos.h"
47
48
49 void PHOSHistos (Text_t* infile, Int_t nevent, Int_t Module)  
50 {
51 //========== Opening galice.root file  
52   TFile * file   = new TFile(infile); 
53 //========== Get AliRun object from file 
54   gAlice = (AliRun*) file->Get("gAlice");//=========== Gets the PHOS object and associated geometry from the file 
55
56   AliPHOSv0 * PHOS  = (AliPHOSv0 *)gAlice->GetDetector("PHOS");
57   AliPHOSGeometry * Geom = AliPHOSGeometry::GetInstance(PHOS->GetGeometry()->GetName(),PHOS->GetGeometry()->GetTitle());
58 //========== Creates the track segment maker
59   AliPHOSTrackSegmentMakerv1 * tracksegmentmaker = new AliPHOSTrackSegmentMakerv1() ;
60   //========== Creates the particle identifier
61   AliPHOSPIDv1 * particleidentifier = new AliPHOSPIDv1 ;
62   Info("PHOSHistos", "AnalyzeOneEvent > using particle identifier %s\n", particleidentifier->GetName() ) ; 
63     
64   TH1F * hEmcDigit       = new TH1F("hEmcDigit","hEmcDigit",1000,0.,5.);
65   TH1F * hVetoDigit      = new TH1F("hVetoDigit","hVetoDigit",1000,0.,3.e-5);
66   TH1F * hConvertorDigit = new TH1F("hConvertorDigit","hConvertorDigit",1000,0.,3.e-5);
67   TH1F * hEmcCluster       = new TH1F("hEmcCluster","hEmcCluster",1000,-5.,5.);
68   TH1F * hVetoCluster      = new TH1F("hVetoCluster","hVetoCluster",1000,0.,3.e-5);
69   TH1F * hConvertorCluster = new TH1F("hConvertorCluster","hConvertorCluster",1000,0.,3.e-5);
70   TH2F * hConvertorEmc = new TH2F("hConvertorEmc","hConvertorEmc",100,2.,3., 100, 0., 1.e-5);
71
72   AliPHOSDigit * digit ;
73
74 //========== Loop on events
75   Int_t ievent;
76   for(ievent=0;ievent<nevent; ievent++)
77   { 
78 //========== Creates the Clusterizer
79     AliPHOSClusterizerv1 * clusterizer = new AliPHOSClusterizerv1() ; 
80     clusterizer->SetEmcEnergyThreshold(0.01) ; 
81     clusterizer->SetEmcClusteringThreshold(0.1) ; 
82     clusterizer->SetPpsdEnergyThreshold(0.00000005) ; 
83     clusterizer->SetPpsdClusteringThreshold(0.00000005) ; 
84     clusterizer->SetLocalMaxCut(0.03) ;
85     clusterizer->SetCalibrationParameters(0., 0.00000001) ;
86
87   //========== Creates the Reconstructioner  
88     AliPHOSReconstructor * Reconstructioner = new AliPHOSReconstructor(clusterizer, tracksegmentmaker, particleidentifier) ;
89      
90     Info("PHOSHistos", "Event %d\n", ievent);
91
92     Int_t RelId[4] ;
93     //=========== Connects the various Tree's for evt
94     gAlice->GetEvent(ievent);
95     //=========== Gets the Digit TTree
96     gAlice->TreeD()->GetEvent(0) ;     
97     //=========== Gets the number of entries in the Digits array
98     //    Int_t nId = PHOS->Digits()->GetEntries();      
99     TIter next(PHOS->Digits()) ;
100     Float_t Etot=0 ;
101     Int_t nVeto=0 ;
102     Int_t nConvertor=0 ;
103
104     while( ( digit = (AliPHOSDigit *)next() ) )
105     {
106        Etot+=clusterizer->Calibrate(digit->GetAmp()) ;
107        Geom->AbsToRelNumbering(digit->GetId(), RelId) ;        
108
109        if (clusterizer->IsInEmc(digit))
110        {   hEmcDigit->Fill(clusterizer->Calibrate(digit->GetAmp())) ; }
111        else    
112        {  
113          if (RelId[1]==9) {nVeto++; hVetoDigit->Fill(clusterizer->Calibrate(digit->GetAmp()));} 
114          if (RelId[1]==25) {nConvertor++; hConvertorDigit->Fill(clusterizer->Calibrate(digit->GetAmp()));}
115        }
116     }
117    
118      PHOS->Reconstruction(Reconstructioner); 
119  
120 //=========== Cluster in Module
121      TClonesArray * EmcRP = PHOS->EmcClusters() ;
122      TIter nextemc(EmcRP) ;
123      AliPHOSEmcRecPoint * emc ;
124      Float_t Energy;
125      while((emc = (AliPHOSEmcRecPoint *)nextemc())) 
126      {
127        if ( emc->GetPHOSMod() == Module )
128        {  
129          Energy = emc->GetEnergy() ;
130          hEmcCluster->Fill(Energy); 
131          printf("Energy of the EMC cluster is %f \n",Energy);
132          TClonesArray * PpsdRP = PHOS->PpsdClusters() ;
133          TIter nextPpsd(PpsdRP) ;
134          AliPHOSPpsdRecPoint * Ppsd ;
135          Float_t Energy2;
136          while((Ppsd = (AliPHOSPpsdRecPoint *)nextPpsd())) 
137            {
138              if ( Ppsd->GetPHOSMod() == Module )
139                { 
140                  Energy2 = Ppsd->GetEnergy() ;
141                  
142                  if (!Ppsd->GetUp()) hConvertorEmc->Fill(Energy,Energy2) ;
143                }
144            }
145     
146        }
147      }
148  
149     //=========== Cluster in Module PPSD Down
150      TClonesArray * PpsdRP = PHOS->PpsdClusters() ;
151      TIter nextPpsd(PpsdRP) ;
152      AliPHOSPpsdRecPoint * Ppsd ;
153      while((Ppsd = (AliPHOSPpsdRecPoint *)nextPpsd())) 
154      {
155        if ( Ppsd->GetPHOSMod() == Module )
156        { 
157          Energy = Ppsd->GetEnergy() ;
158
159          if (!Ppsd->GetUp()) hConvertorCluster->Fill(Energy) ;
160          if (Ppsd->GetUp()) hVetoCluster->Fill(Energy) ;
161    
162        }
163      }
164   } 
165
166   TCanvas * cVetoDigit = new TCanvas("VetoDigit","VetoDigit");  
167   hVetoDigit->Draw();
168   TCanvas * cConvertorDigit = new TCanvas("ConvertorDigit","ConvertorDigit");  
169   hConvertorDigit->Draw();
170   TCanvas * cEmcDigit = new TCanvas("EmcDigit","EmcDigit");  
171   hEmcDigit->Draw();
172   TCanvas * cVetoCluster = new TCanvas("VetoCluster","VetoCluster");  
173   hVetoCluster->Draw();
174   TCanvas * cConvertorCluster = new TCanvas("ConvertorCluster","ConvertorCluster");  
175   hConvertorCluster->Draw();
176   TCanvas * cEmcCluster = new TCanvas("EmcCluster","EmcCluster");  
177   hEmcCluster->Draw();
178   TCanvas * cConvertorEmc = new TCanvas("ConvertorEmc","ConvertorEmc");  
179   hConvertorEmc->Draw("col1");
180
181
182 }
183
184