/[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.1 - (hide annotations) (download)
Mon Dec 11 18:26:55 2006 UTC (17 years, 11 months ago) by pam-fi
Branch: MAIN
examples

1 pam-fi 1.1 ////////////////////////////////////////////////////////////////
2     //
3     // Histogram set for tracker-efficiency study.
4     //
5     // It contains the functions:
6     //
7     // - CreateHistos()
8     // - FillHistos(PamLevel2* event)
9     // - Save Histos()
10     //
11     // Usage (example):
12     // [] .L My-histo-set-macro.C+
13     // [] .L My-selection.C+
14     // [] .L Loop.C+
15     // [] Loop("data-dir/","file-list.txt",100000,"+ALL fillTree fillHistos","output-file.root")
16     //
17     ////////////////////////////////////////////////////////////////
18    
19     #if !defined(__CINT__) || defined(__MAKECINT__)
20    
21     #include <TString.h>
22     #include <TH1F.h>
23     #include <TH2F.h>
24     #include <TMath.h>
25     #include <TLine.h>
26     #include <TPolyMarker.h>
27     #include <TSelector.h>
28     #include <TFile.h>
29    
30     #include <stdlib.h>
31     #include <iostream>
32     using namespace std;
33    
34     #include <PamLevel2.h>
35     #include <CaloAxis.h>
36    
37     #endif
38     //======================
39     // HITOGRAMS declaration
40     //======================
41     TH1F *hresangx;
42     TH1F *hresangy;
43     //===============================================================
44     // Create histograms
45     //===============================================================
46     void CreateHistos( TFile* outf ){
47    
48     gROOT->cd();//create histos in memory
49    
50     hresangx = new TH1F("hresangx","Angular residual on the top calo plane (x)",400,-40.,40.);
51     hresangy = new TH1F("hresangy","Angular residual on the top calo plane (y)",400,-40.,40.);
52    
53     }
54    
55     //===============================================================
56     // Fill histograms (called event-by-event)
57     //===============================================================
58    
59     void FillHistos( PamLevel2* event ){
60     // ----------------------------------
61     // retrieve histos
62     // ----------------------------------
63     // hresangx = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangx"));
64     // hresangy = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangy"));
65    
66     // ----------------------------------
67     // calo variables
68     // ----------------------------------
69     CaloAxis *x_axis = new CaloAxis();
70     CaloAxis *y_axis = new CaloAxis();
71    
72     float rcil = 1.;// tolerance (cm)
73     x_axis->FitAxis(event->GetCaloLevel1(),0,rcil);
74     y_axis->FitAxis(event->GetCaloLevel1(),1,rcil);
75    
76     float qtrack = x_axis->GetQaxis()+y_axis->GetQaxis();
77     int ntrack = x_axis->GetN()+y_axis->GetN();
78    
79     float caloangx = (float)atan((double)x_axis->par[1])*180./3.1415026;
80     float caloangy = (float)atan((double)y_axis->par[1])*180./3.1415026;
81    
82     CaloStrip top = CaloStrip();
83     top.Set(0,0,48);
84     float zcalotop = top.GetZ();
85    
86     // ----------------------------------
87     // tof variables
88     // ----------------------------------
89     // get track stored by the tof
90     ToFTrkVar *tof = event->GetToFStoredTrack(-1);
91     if(!tof){
92     cout << " no ToF-track stored "<<endl;
93     return;
94     }
95    
96     // ----------------------------------
97     // tracker variables
98     // ----------------------------------
99    
100     float resangx = -999;
101     float resangy = -999;
102     if( event->GetNTracks()==1 ){
103    
104     PamTrack *track = event->GetTrack(0);
105     if( track->chi2 >0 ){
106    
107     Float_t ztraj[13];
108    
109     Int_t i=0;
110     ztraj[i++] = event->GetZTOF(11);
111     ztraj[i++] = event->GetZTOF(12);
112     ztraj[i++] = event->GetZTOF(21);
113     ztraj[i++] = event->GetZTOF(22);
114     ztraj[i++] = event->GetZTrk(1);
115     ztraj[i++] = event->GetZTrk(2);
116     ztraj[i++] = event->GetZTrk(3);
117     ztraj[i++] = event->GetZTrk(4);
118     ztraj[i++] = event->GetZTrk(5);
119     ztraj[i++] = event->GetZTrk(6);
120     ztraj[i++] = event->GetZTOF(31);
121     ztraj[i++] = event->GetZTOF(32);
122     ztraj[i++] = zcalotop;
123     Trajectory traj = Trajectory(13,ztraj);
124     traj.DoTrack2( track->al );
125    
126     resangx = caloangx - traj.thx[12];
127     resangy = caloangy - traj.thy[12];
128    
129     }
130     if(track) delete track;
131     }
132    
133     // ----------------------------------
134     // fill histograms
135     // ----------------------------------
136     hresangx->Fill(resangx);
137     hresangy->Fill(resangy);
138    
139    
140     }
141     //===============================================================
142     // Save histograms
143     //===============================================================
144     void SaveHistos( TFile * outf ){
145    
146     gROOT->cd();
147     // ----------------------------------
148     // retrieve histos
149     // ----------------------------------
150     // hresangx = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangx"));
151     // hresangy = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangy"));
152    
153     // ----------------------------------
154     // save on file
155     // ----------------------------------
156    
157     outf->cd();
158    
159     if(hresangx) hresangx->Write();
160     if(hresangy) hresangy->Write();
161     }
162    

  ViewVC Help
Powered by ViewVC 1.1.23