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