/[PAMELA software]/calo/flight/FUTILITIES/macros/FCaloMIP.cxx
ViewVC logotype

Annotation of /calo/flight/FUTILITIES/macros/FCaloMIP.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1 - (hide annotations) (download)
Wed Mar 22 15:07:47 2006 UTC (18 years, 8 months ago) by mocchiut
Branch point for: FUTILITIES, MAIN
Initial revision

1 mocchiut 1.1 //
2     // Plot the energy [MIP], nstrip and qtot distributions - Emiliano Mocchiutti
3     //
4     // FCaloMIP.cxx version 1.01 (2006-06-03)
5     //
6     // The only input needed is the path to the directory created by YODA for the data file you want to analyze.
7     //
8     // Changelog:
9     //
10     // 1.00 - 1.01 (2006-06-03): Flight version, read unique YODA file.
11     //
12     // 0.00 - 1.00 (2006-06-03): Clone of CaloMIP.c v2r12
13     //
14     #include <iostream>
15     #include <sstream>
16     #include <fstream>
17     #include <TTree.h>
18     #include <TClassEdit.h>
19     #include <TObject.h>
20     #include <TList.h>
21     #include <TSystem.h>
22     #include <TSystemDirectory.h>
23     #include <TString.h>
24     #include <TFile.h>
25     #include <TClass.h>
26     #include <TCanvas.h>
27     #include <TH1.h>
28     #include <TH1F.h>
29     #include <TH2D.h>
30     #include <TLatex.h>
31     #include <TPad.h>
32     #include <TPaveLabel.h>
33     #include <TChain.h>
34     #include <TStyle.h>
35     #include <TF1.h>
36     //
37     extern bool existfile(TString);
38     extern void stringcopy(TString&, const TString&, Int_t, Int_t);
39     TF1 *langaufit(TH1F *, Double_t *, Double_t *, Double_t *, Double_t *, Double_t *, Double_t *, Double_t *, Int_t *, TString);
40     Double_t langaupro(Double_t *);
41     //
42     #include <caloclassesfun.h>
43     //
44     using namespace std;
45    
46     void FCaloMIP(const TString filename, Int_t view = 0, Int_t plane = 0, Int_t strip = 0, Int_t fromevent = 0, Int_t toevent = 0, TString outDir = "", TString saveas = "gif"){
47     gStyle->SetOptStat(1111);
48     const char* startingdir = gSystem->WorkingDirectory();
49     if ( !strcmp(outDir.Data(),"") ) outDir = startingdir;
50     if ( !existfile(filename) ){
51     printf(" %s :no such file or directory \n",filename.Data());
52     printf("\n ERROR: No calorimeter LEVEL1 file! \n\n Run FCaloLEVEL1 to generate level1 data. \n\n");
53     return;
54     };
55     //
56     TFile *File = new TFile(filename.Data());
57     TTree *otr = (TTree*)File->Get("CaloLevel1");
58     if ( !otr ){
59     printf(" %s are you sure? I need a LEVEL1 file! \n",filename.Data());
60     return;
61     };
62     //
63     CalorimeterLevel1 *calo = new CalorimeterLevel1();
64     otr->SetBranchAddress("Event", &calo);
65     //
66     Long64_t nevents = otr->GetEntries();
67     //
68     // input limits
69     //
70     Int_t minview;
71     Int_t maxview;
72     if ( view > 2 || view < 0 ) {
73     printf("You can choose view = 1 (X-views), view = 2 (Y-views) or view = 0 (both)\n");
74     return;
75     };
76     if ( view == 0 ) {
77     minview = 0;
78     maxview = 2;
79     } else {
80     minview = view - 1;
81     maxview = view;
82     };
83    
84     if ( plane > 22 || plane < 0 ) {
85     printf("You can choose plane between 1 and 22 or plane = 0 (all)\n");
86     return;
87     };
88     Int_t minplane;
89     Int_t maxplane;
90     if ( plane == 0 ) {
91     minplane = 0;
92     maxplane = 22;
93     } else {
94     minplane = plane - 1;
95     maxplane = plane;
96     };
97     Int_t minstrip;
98     Int_t maxstrip;
99     if ( strip > 96 || strip < 0 ) {
100     printf("You can choose strip between 0 and 96 or strip = 0 (all)\n");
101     return;
102     };
103     if ( strip == 0 ) {
104     minstrip = 0;
105     maxstrip = 96;
106     } else {
107     minstrip = strip - 1;
108     maxstrip = strip;
109     };
110    
111     if ( fromevent > toevent && toevent ){
112     printf("It must be fromevent < toevent \n");
113     return;
114     };
115    
116    
117     if ( fromevent > nevents+1 || fromevent < 0 ) {
118     printf("You can choose fromevent between 0 (all) and %i \n",(int)nevents+1);
119     return;
120     };
121    
122     if ( toevent > nevents+1 || toevent < 0 ) {
123     printf("You can choose toevent between 0 (all) and %i \n",(int)nevents+1);
124     return;
125     };
126     Int_t minevent;
127     Int_t maxevent;
128     if ( fromevent == 0 ) {
129     minevent = 0;
130     maxevent = nevents;
131     } else {
132     minevent = fromevent - 1;
133     if ( toevent > 0 ){
134     maxevent = toevent;
135     } else {
136     maxevent = fromevent - 1;
137     };
138     };
139    
140     //
141     // Book the histograms:
142     //
143     //
144     TH1F *gmip;
145     gmip = new TH1F("the MIP","",55,-0.19,1.9);
146     //gmip->SetBit(TH1F::kCanRebin);
147     //
148     TH1F *diffbase;
149     diffbase = new TH1F("diffbase","",55,-0.19,0.19);
150     diffbase->SetBit(TH1F::kCanRebin);
151     //
152     TH1F *Qtot;
153     Qtot = new TH1F("qtot","",70,0.,60.);
154     Qtot->SetBit(TH1F::kCanRebin);
155     //
156     TH1F *Nstrip;
157     Nstrip = new TH1F("nstrip","",70,0.,60.);
158     Nstrip->SetBit(TH1F::kCanRebin);
159     //
160     // run over all the events
161     //
162     //
163     Int_t tot2 = 0;
164     Int_t baseDIFFfull = 0;
165     printf("\n Processed events: \n\n");
166     //
167     for (Int_t i = minevent; i < maxevent; i++){
168     if ( i%1000 == 0 && i > 0 ) printf(" %iK \n",i/1000);
169    
170     //
171     // read from the header of the event the absolute time at which it was recorded
172     //
173     otr->GetEntry(i);
174    
175     //
176     // run over views and planes
177     //
178     for (Int_t l = minview; l < maxview; l++){
179     for (Int_t m = minplane; m < maxplane; m++){
180     Int_t pre = -1;
181     for (Int_t n = minstrip; n < maxstrip; n++){
182     //
183     if ( n%16 == 0 ) pre++;
184     if ( n == minstrip && pre == -1 ){
185     for (Int_t nm = minstrip; pre == -1 ; nm--){
186     if ( nm%16 == 0 ) pre++;
187     };
188     };
189     if ( pre == -1 ) {
190     printf("Unable to determine the preamplifier number. Quit.\n");
191     return;
192     };
193     //
194     if ( calo->nobase ) tot2++;
195     if ( calo->diffbas[l][m][pre] != 0. ) {
196     // printf("differenza %f %i %i %i %i %i \n",calo->diffbas[l][m][pre],i,l,m,pre,n);
197     baseDIFFfull = 1;
198     diffbase->Fill(calo->diffbas[l][m][pre]);
199     // return;
200     };
201     if ( calo->estrip[l][m][n] != 0. ) gmip->Fill(calo->estrip[l][m][n]);
202     };
203     };
204     };
205     if ( calo->nstrip > 0. ) Nstrip->Fill(calo->nstrip);
206     if ( calo->qtot > 0. ) Qtot->Fill(calo->qtot);
207     };
208     if ( tot2 ){
209     printf("WARNING: events with no baselines (event lost): %i\n",tot2);
210     };
211     //
212     // figures:
213     //
214     // THE MIP
215     TCanvas *figura = new TCanvas("The_calorimeter_MIP", "The_calorimeter_MIP", 1100, 900);
216     figura->SetFillColor(10);
217     figura->Range(0,0,100,100);
218     figura->cd();
219     gPad->SetLogy();
220     gmip->SetXTitle("MIP");
221     gmip->SetYTitle("Number of events");
222     gPad->SetTicks();
223     printf("\n Fitting...\n\n");
224     // Setting fit range and start values
225     Double_t fr[2];
226     Double_t sv[4], pllo[4], plhi[4], fp[4], fpe[4];
227     fr[0] = 0.7;
228     fr[1] = 2.0;
229     pllo[0]=0.01; pllo[1]=0.5; pllo[2]=1.0; pllo[3]=0.05;
230     plhi[0]=5.0; plhi[1]=5.0; plhi[2]=1000000.0; plhi[3]=1.0;
231     sv[0]=1.; sv[1]=1.0; sv[2]=500.0; sv[3]=0.2;
232     Double_t chisqr;
233     Int_t ndf;
234     TF1 *fitsnr = langaufit(gmip,fr,sv,pllo,plhi,fp,fpe,&chisqr,&ndf,"QR0");
235     //
236     Double_t SNRPeak = langaupro(fp);
237     printf("\n Conversion factor: %f \n",SNRPeak);
238     //
239     // gStyle->SetOptStat(1111);
240     gStyle->SetOptFit(111);
241     printf("\n ...done!\n");
242     gmip->Draw();
243     fitsnr->Draw("lsame");
244     // Differences between DSP and calculated baselines in full mode
245     TCanvas *base = 0;
246     if ( baseDIFFfull ){
247     base = new TCanvas("WARNING", "WARNING_full mode,_differences_between_calculated_and_DSP_baselines!!", 1100, 900);
248     base->SetFillColor(10);
249     base->Range(0,0,100,100);
250     base->cd();
251     gStyle->SetOptStat(0);
252     gPad->SetLogy();
253     diffbase->Draw();
254     gStyle->SetOptStat(1111);
255     };
256     // nstrip and qtot WARNING: only GOOD strips considered!
257     TCanvas *figura2 = new TCanvas("nstrip_and_qtot", "nstrip_and_qtot", 1100, 900);
258     figura2->SetFillColor(10);
259     figura2->Range(0,0,100,100);
260     TPad *pd1;
261     TPad *pd2;
262     pd1 = new TPad("pd1","This is pad1",0.02,0.02,0.49,0.98,10);
263     pd2 = new TPad("pd2","This is pad2",0.51,0.02,0.98,0.98,10);
264     figura2->cd();
265     gStyle->SetOptStat(1111);
266     pd1->Range(0,0,100,100);
267     pd2->Range(0,0,100,100);
268     pd1->SetGrid();
269     pd2->SetGrid();
270     pd1->SetTicks();
271     pd2->SetTicks();
272     pd1->Draw();
273     pd2->Draw();
274     pd1->cd();
275     gPad->SetLogy();
276     Nstrip->SetXTitle("Number of hit (E>0.7 MIP)");
277     Nstrip->SetYTitle("Number of events");
278     Nstrip->Draw();
279     pd2->cd();
280     gPad->SetLogy();
281     Qtot->SetXTitle("MIP");
282     Qtot->SetYTitle("Number of events");
283     Qtot->Draw();
284     printf("\n");
285     //
286     const string fil = (const char*)filename;
287     Int_t posiz = fil.find("dw_");
288     if ( posiz == -1 ) posiz = fil.find("DW_");
289     Int_t posiz2 = posiz+13;
290     TString file2;
291     stringcopy(file2,filename,posiz,posiz2);
292     //
293     const char *figrec = file2;
294     const char *outdir = outDir;
295     stringstream figsave;
296     const char *format = saveas;
297     figsave.str("");
298     figsave << outdir << "/" ;
299     figsave << figrec << "_mip1.";
300     figsave << format;
301     //
302     figura->SaveAs(figsave.str().c_str());
303     figsave.str("");
304     figsave << outdir << "/" ;
305     figsave << figrec << "_mip2.";
306     figsave << format;
307     //
308     figura2->SaveAs(figsave.str().c_str());
309     if ( baseDIFFfull ){
310     figsave.str("");
311     figsave << outdir << "/" ;
312     figsave << figrec << "_mip3.";
313     figsave << format;
314     //
315     base->SaveAs(figsave.str().c_str());
316     };
317     //
318     gSystem->ChangeDirectory(startingdir);
319     }

  ViewVC Help
Powered by ViewVC 1.1.23