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

Contents of /PamelaLevel2/doc/examples/My-Histos-0.cpp

Parent Directory Parent Directory | Revision Log Revision Log


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

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 float caloangx = 0;
83 float caloangy = 0;
84
85 CaloStrip top = CaloStrip();
86 top.Set(0,0,48);
87 float zcalotop = top.GetZ();
88
89 // ----------------------------------
90 // tof variables
91 // ----------------------------------
92 // get track stored by the tof
93 // ToFTrkVar *tof = event->GetToFStoredTrack(-1);
94 // if(!tof){
95 // cout << " no ToF-track stored "<<endl;
96 // return;
97 // }
98
99 // ----------------------------------
100 // tracker variables
101 // ----------------------------------
102
103 float resangx = -999;
104 float resangy = -999;
105 if( event->GetNTracks()==1 ){
106
107 PamTrack *track = event->GetTrack(0);
108 if( track->chi2 >0 ){
109
110 Float_t ztraj[13];
111
112 Int_t i=0;
113 ztraj[i++] = event->GetZTOF(11);
114 ztraj[i++] = event->GetZTOF(12);
115 ztraj[i++] = event->GetZTOF(21);
116 ztraj[i++] = event->GetZTOF(22);
117 ztraj[i++] = event->GetZTrk(1);
118 ztraj[i++] = event->GetZTrk(2);
119 ztraj[i++] = event->GetZTrk(3);
120 ztraj[i++] = event->GetZTrk(4);
121 ztraj[i++] = event->GetZTrk(5);
122 ztraj[i++] = event->GetZTrk(6);
123 ztraj[i++] = event->GetZTOF(31);
124 ztraj[i++] = event->GetZTOF(32);
125 ztraj[i++] = zcalotop;
126 Trajectory traj = Trajectory(13,ztraj);
127 traj.DoTrack2( track->al );
128
129 resangx = caloangx - traj.thx[12];
130 resangy = caloangy - traj.thy[12];
131
132 }
133 if(track) delete track;
134 }
135
136 // ----------------------------------
137 // fill histograms
138 // ----------------------------------
139 hresangx->Fill(resangx);
140 hresangy->Fill(resangy);
141
142
143 }
144 //===============================================================
145 // Save histograms
146 //===============================================================
147 void SaveHistos( TFile * outf ){
148
149 gROOT->cd();
150 // ----------------------------------
151 // retrieve histos
152 // ----------------------------------
153 // hresangx = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangx"));
154 // hresangy = dynamic_cast<TH1F*>(gDirectory->FindObject("hresangy"));
155
156 // ----------------------------------
157 // save on file
158 // ----------------------------------
159
160 outf->cd();
161
162 if(hresangx) hresangx->Write();
163 if(hresangy) hresangy->Write();
164 }
165

  ViewVC Help
Powered by ViewVC 1.1.23