/[PAMELA software]/PamelaLevel2/doc/examples/My-Histos.cpp
ViewVC logotype

Annotation of /PamelaLevel2/doc/examples/My-Histos.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (hide annotations) (download)
Wed Jan 3 13:28:49 2007 UTC (17 years, 11 months ago) by pam-fi
Branch: MAIN
CVS Tags: v0r02, v2r02
Changes since 1.1: +283 -77 lines
new example with selection cuts

1 pam-fi 1.1 ////////////////////////////////////////////////////////////////
2     //
3     // Histogram set for tracker-efficiency study.
4     //
5     // It contains the functions:
6     //
7     // - CreateHistos()
8 pam-fi 1.2 // - RetrieveHistos()
9 pam-fi 1.1 // - FillHistos(PamLevel2* event)
10     // - Save Histos()
11     //
12     // Usage (example):
13     // [] .L My-histo-set-macro.C+
14     // [] .L My-selection.C+
15     // [] .L Loop.C+
16     // [] Loop("data-dir/","file-list.txt",100000,"+ALL fillTree fillHistos","output-file.root")
17     //
18     ////////////////////////////////////////////////////////////////
19    
20     #if !defined(__CINT__) || defined(__MAKECINT__)
21    
22     #include <TString.h>
23     #include <TH1F.h>
24     #include <TH2F.h>
25     #include <TMath.h>
26     #include <TLine.h>
27     #include <TPolyMarker.h>
28     #include <TSelector.h>
29     #include <TFile.h>
30    
31     #include <stdlib.h>
32     #include <iostream>
33     using namespace std;
34    
35     #include <PamLevel2.h>
36     #include <CaloAxis.h>
37    
38     #endif
39     //======================
40     // HITOGRAMS declaration
41     //======================
42 pam-fi 1.2
43     TH2F *hchi2def;
44     TH2F *hchi2rig;
45     TH2F *hdedxtofvsdedxtrk;
46     TH2F *hdedxtrkvsrig ;
47     TH2F *hdedxtofvsrig ;
48     TH2F *hbetavsrig ;
49     TH2F *hdedxtrkvsdef ;
50     TH2F *hbetavsdef ;
51     TH2F *hbetavsdef2 ;
52    
53     TH1F *herrdef;
54     TH1F *herrx[6];
55     TH1F *herry[6];
56    
57     TH1F *hdef ;
58     TH1F *hrig1 ;
59     TH1F *hrig2 ;
60     TH1F *hen1 ;
61     TH1F *hen2 ;
62    
63    
64     Int_t tot=0;
65     Int_t sel[]={0,0,0,0,0,0,0,0};
66    
67     //===============================================================
68     // Get ToF dE/dx
69     //===============================================================
70     Float_t GetdEdx( PamTrack *trk ){
71     Float_t dedx = 0.;
72     Int_t nval=0;
73     Bool_t ISNULL = false;
74     for (Int_t i=0; i<trk->npmtadc; i++){
75     if(trk->dedx[i] == 0)ISNULL=true;
76     dedx += (trk->dedx).At(i)*(Int_t)(!ISNULL);
77     nval += (Int_t)(!ISNULL);
78     };
79     //
80     if(nval)dedx=dedx/nval;
81     // cout <<" -- dedx "<<dedx<<endl;
82     return(dedx);
83     }
84     float GetZ(float beta, float dedx){
85    
86     int nspec=2;
87     float p0[]={0.904,3.373};
88     float p1[]={0.295,2.756};
89     float zz[]={1.,4.};
90    
91     if(beta>1)beta=1;
92     float f0 = p0[0]/beta+p1[0];
93     float f1 = p0[1]/beta+p1[1];
94     float z0 = zz[0];
95     float z1 = zz[1];
96     float aa = (z0-z1)/(f0-f1);
97     float bb = -aa*f0+z0;
98     if( (bb+aa*dedx)==0 )return 0;
99     return sqrt(bb+aa*dedx);
100    
101     }
102    
103    
104 pam-fi 1.1 //===============================================================
105     // Create histograms
106     //===============================================================
107     void CreateHistos( TFile* outf ){
108    
109     gROOT->cd();//create histos in memory
110    
111 pam-fi 1.2 TString hid;
112     TString title;
113     //===========================================
114     // rigidity binning
115     Int_t nbins = 200;
116     Float_t rmin = .1;
117     Float_t rmax = 1000.;
118     Float_t dr = (log10(rmax)-log10(rmin))/nbins;
119     TArrayD* arr = new TArrayD(nbins+1);
120     arr->AddAt((Double_t)rmin,0);
121     for(Int_t i=1; i<=nbins; i++){
122     Double_t previous = arr->At(i-1);
123     Double_t value=0;
124     value = (Double_t) pow(10,log10((Float_t)previous)+dr);
125     arr->AddAt(value,i);
126     };
127     Double_t *rbinpt = arr->GetArray();
128     //===========================================
129    
130     hchi2def = new TH2F("hchi2def","(chi**2)**0.25 vs deflection",100,0.,2.5,200,0.,10.);
131     hchi2rig = new TH2F("hchi2rig","#chi^{2} vs R",nbins,rbinpt,1000,0.,500.);
132     hdedxtofvsdedxtrk = new TH2F("hdedxtofvsdedxtrk","dE/dx (ToF-average) vs dE/dx (trk-average)",600,0.,40.,600,0.,40.);
133     hdedxtrkvsrig = new TH2F("hdedxtrkvsrig","dE/dx (Trk-average) vs R",nbins,rbinpt,600,0.,40.);
134     hdedxtofvsrig = new TH2F("hdedxtofvsrig","dE/dx (ToF-average) vs R",nbins,rbinpt,600,0.,40.);
135     hbetavsrig = new TH2F("hbetavsrig","beta vs R",nbins,rbinpt,200,0.4,1.5);
136     hrig1 = new TH1F("hrig1","Rigidity distribution Z=1",nbins,rbinpt);
137     hrig2 = new TH1F("hrig2","Rigidity distribution Z=2",nbins,rbinpt);
138     hen1 = new TH1F("hen1","Kin.Energy per nucleon Z=1",nbins,rbinpt);
139     hen2 = new TH1F("hen2","Kin.Energy per nucleon Z=2",nbins,rbinpt);
140    
141     hdedxtrkvsdef = new TH2F("hdedxtrkvsdef","dE/dx (Trk-average) vs deflection",500,-10.,10.,600,0.,40.);
142     hbetavsdef = new TH2F("hbetavsdef","beta vs deflection",500,-10.,10.,200,0.4,1.5);
143     hbetavsdef2 = new TH2F("hbetavsdef2","beta vs deflection (Z=1)",500,-10.,10.,200,0.4,1.5);
144     hdef = new TH1F("hdef","Deflection distribution",5000,-2.5,2.5);
145    
146     herrdef = new TH1F("herrdef","Deflection error",1000,0.,0.01);
147    
148     for(Int_t ip=0; ip<6; ip++){
149     hid.Form("herrx%i",ip);
150     title = ""; title.Form("Spatial resolution - plane %i (x)",ip);
151     herrx[ip] = new TH1F(hid.Data(),title.Data(),500,0.,0.005);
152     hid.Form("herry%i",ip);
153     title = ""; title.Form("Spatial resolution - plane %i (y)",ip);
154     herry[ip] = new TH1F(hid.Data(),title.Data(),500,0.,0.005);
155     }
156 pam-fi 1.1
157     }
158 pam-fi 1.2 //===============================================================
159     // Retrieve histograms
160     //===============================================================
161     void RetrieveHistos(){
162    
163     TString hid;
164 pam-fi 1.1
165 pam-fi 1.2 hchi2def = dynamic_cast<TH2F*>(gDirectory->FindObject("hchi2def"));
166     hchi2rig = dynamic_cast<TH2F*>(gDirectory->FindObject("hchi2rig"));
167     hdedxtofvsdedxtrk = dynamic_cast<TH2F*>(gDirectory->FindObject("hdedxtofvsdedxtrk"));
168     hdedxtrkvsrig = dynamic_cast<TH2F*>(gDirectory->FindObject("hdedxtrkvsrig"));
169     hdedxtofvsrig = dynamic_cast<TH2F*>(gDirectory->FindObject("hdedxtofvsrig"));
170     hbetavsrig = dynamic_cast<TH2F*>(gDirectory->FindObject("hbetavsrig"));
171     hdedxtrkvsdef = dynamic_cast<TH2F*>(gDirectory->FindObject("hdedxtrkvsdef"));
172     hbetavsdef = dynamic_cast<TH2F*>(gDirectory->FindObject("hbetavsdef"));
173     hbetavsdef2 = dynamic_cast<TH2F*>(gDirectory->FindObject("hbetavsdef2"));
174     hdef = dynamic_cast<TH1F*>(gDirectory->FindObject("hdef"));
175     hrig1 = dynamic_cast<TH1F*>(gDirectory->FindObject("hrig1"));
176     hrig2 = dynamic_cast<TH1F*>(gDirectory->FindObject("hrig2"));
177     hen1 = dynamic_cast<TH1F*>(gDirectory->FindObject("hen1"));
178     hen2 = dynamic_cast<TH1F*>(gDirectory->FindObject("hen2"));
179     herrdef = dynamic_cast<TH1F*>(gDirectory->FindObject("herrdef"));
180     for(Int_t ip=0; ip<6; ip++){
181     hid.Form("herrx%i",ip);
182     herrx[ip] = dynamic_cast<TH1F*>(gDirectory->FindObject(hid.Data()));
183     hid.Form("herry%i",ip);
184     herry[ip] = dynamic_cast<TH1F*>(gDirectory->FindObject(hid.Data()));
185     }
186     }
187 pam-fi 1.1 //===============================================================
188     // Fill histograms (called event-by-event)
189     //===============================================================
190    
191     void FillHistos( PamLevel2* event ){
192     // ----------------------------------
193     // retrieve histos
194     // ----------------------------------
195 pam-fi 1.2 RetrieveHistos();
196 pam-fi 1.1
197     // ----------------------------------
198 pam-fi 1.2 // single-track events
199 pam-fi 1.1 // ----------------------------------
200 pam-fi 1.2 float beta = 0;
201     float dedxtof = 0;
202     float dedxtrk = 0;
203     float rigidity = 0;
204     float deflection = 0;
205     float chi2 = 0;
206     float errdef = 0.;
207     float errx[] = {10.,10.,10.,10.,10.,10.};
208     float erry[] = {10.,10.,10.,10.,10.,10.};
209 pam-fi 1.1
210 pam-fi 1.2 tot++;
211    
212     if( event->GetNTracks()==0 ) return;
213     PamTrack *track = event->GetTrack(0);
214     //------------------------------
215     // track selection
216     //------------------------------
217     if(
218     true
219     ){
220    
221     sel[0]++;
222     //-------------------
223     // Calo variables
224     //-------------------
225    
226    
227     //-------------------
228     // ToF variables
229     //-------------------
230     beta = track->beta[12];
231     // dedxtof = GetdEdx(track);
232    
233     //-------------------
234     // tracker variables
235     //-------------------
236     dedxtrk = track->GetDEDX()/1.3; //<<<ATTENZIONE TEMPORANEO!!!!
237     rigidity = track->GetRigidity();
238     deflection = track->GetDeflection();
239     chi2 = track->chi2;
240     errdef = sqrt(track->coval[4][4]);
241     for(Int_t ip=0; ip<6; ip++){
242     if(track->XGood(ip))errx[ip]=track->resx[ip];
243     if(track->YGood(ip))erry[ip]=track->resy[ip];
244     }
245 pam-fi 1.1
246 pam-fi 1.2 //-------------------
247     // chi2 selection
248     //-------------------
249     hchi2def->Fill( fabs(deflection), pow((double)chi2,0.25) );
250     hchi2rig->Fill( rigidity, chi2 );
251     if(pow((double)chi2,0.25) > (2.5+1.85*fabs(deflection)) )return;
252     sel[1]++;
253    
254     // hdedxtofvsdedxtrk->Fill(dedxtrk,dedxtof);
255     hdedxtrkvsrig->Fill(rigidity,dedxtrk);
256     // hdedxtofvsrig->Fill(rigidity,dedxtof);
257     hbetavsrig->Fill(rigidity,beta);
258    
259     hdedxtrkvsdef->Fill(deflection,dedxtrk);
260     hbetavsdef->Fill(deflection,beta);
261    
262     //------------------------
263     // rough charge selection
264     //------------------------
265     if(dedxtrk > (6.+16./(rigidity*rigidity)))return;
266     //--------
267     // helium
268     //--------
269     if( dedxtrk > (2.8+3.3/(rigidity*rigidity)) && deflection>0 ){
270     sel[5]++;
271     hrig2->Fill(rigidity);
272     float m=3.725;
273     float a=4.;
274     float z=2.;
275     float en=(sqrt( (z*rigidity)*(z*rigidity)+m*m )-m)/a;
276     hen2->Fill(en);
277     }
278 pam-fi 1.1
279 pam-fi 1.2 //--------
280     // Z=1
281     //--------
282     if(dedxtrk > (2.8+3.3/(rigidity*rigidity)))return;
283     if(dedxtrk < (0.4+0.35/(rigidity*rigidity)))return;
284     sel[2]++;
285    
286     if(deflection>0){
287     hrig1->Fill(rigidity);
288     float m=0.938;
289     float a=1.;
290     float z=1.;
291     float en=(sqrt( (z*rigidity)*(z*rigidity)+m*m )-m)/a;
292     hen1->Fill(en);
293 pam-fi 1.1 }
294 pam-fi 1.2 //------------------------
295     // 6-planes lever-arm
296     //------------------------
297     if( !(track->XGood(0) && track->XGood(5)) )return;
298     sel[3]++;
299    
300     for(Int_t ip=0; ip<6; ip++)herrx[ip]->Fill(errx[ip]);
301     for(Int_t ip=0; ip<6; ip++)herry[ip]->Fill(erry[ip]);
302    
303     //---------------------------
304     // cut on spatial resolution
305     //---------------------------
306     if( errx[0] > 0.0004 )return;
307     if( errx[5] > 0.0004 )return;
308     sel[4]++;
309    
310     herrdef->Fill(errdef);
311     hbetavsdef2->Fill(deflection,beta);
312     hdef->Fill(deflection);
313 pam-fi 1.1 }
314    
315 pam-fi 1.2 delete track;
316    
317    
318 pam-fi 1.1 }
319     //===============================================================
320     // Save histograms
321     //===============================================================
322     void SaveHistos( TFile * outf ){
323    
324 pam-fi 1.2 TString hid;
325    
326 pam-fi 1.1 gROOT->cd();
327     // ----------------------------------
328     // retrieve histos
329     // ----------------------------------
330 pam-fi 1.2 RetrieveHistos();
331 pam-fi 1.1
332     // ----------------------------------
333     // save on file
334     // ----------------------------------
335     outf->cd();
336    
337 pam-fi 1.2 hchi2def->Write();
338     hchi2rig->Write();
339     hdedxtofvsdedxtrk->Write();
340     hdedxtrkvsrig->Write();
341     hdedxtofvsrig->Write();
342     hbetavsrig->Write();
343     hdedxtrkvsdef->Write();
344     hbetavsdef->Write();
345     hbetavsdef2->Write();
346     hdef->Write();
347     hrig1->Write();
348     hrig2->Write();
349     hen1->Write();
350     hen2->Write();
351     herrdef->Write();
352     for(Int_t ip=0; ip<6; ip++){
353     herrx[ip]->Write();
354     herry[ip]->Write();
355     }
356    
357     cout << "------------------------------------------------------------------" <<endl;
358     cout << "Pre-selected events :"<<tot<<endl;
359     if(!tot)return;
360     cout << "chi2 selection :"<<sel[1]<< " (" << ((float)sel[1]/(float)tot) <<")"<<endl;
361     cout << "Z=2 selection :"<<sel[5]<< " (" << ((float)sel[5]/(float)tot) <<")"<<endl;
362     cout << "Z=1 selection :"<<sel[2]<< " (" << ((float)sel[2]/(float)tot) <<")"<<endl;
363     cout << "(6-plane lever-arm) :"<<sel[3]<< " (" << ((float)sel[3]/(float)tot) <<")"<<endl;
364     cout << "ext.spatial res. on planes 1-6 < 4 micron :"<<sel[4]<< " (" << ((float)sel[4]/(float)tot) <<")"<<endl;
365     cout << "------------------------------------------------------------------" <<endl;
366    
367 pam-fi 1.1 }
368    

  ViewVC Help
Powered by ViewVC 1.1.23