/[PAMELA software]/quicklook/tof/src/TofScan.cpp
ViewVC logotype

Annotation of /quicklook/tof/src/TofScan.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.4 - (hide annotations) (download)
Mon Oct 8 12:28:54 2007 UTC (17 years, 2 months ago) by pam-de
Branch: MAIN
Changes since 1.3: +440 -7 lines
HV-status check added to TOFScan13

1 campana 1.1 /**
2     * TOFScan
3     * Author Nagni
4     * Version 1.2
5     * Modified by G.De Rosa
6     * Date 27 Apr 2006
7 campana 1.3 * Modified by G.De Rosa
8     * Date 03 Jul 2006
9 pam-de 1.4 * Modified by W. Menn to select helium particles for PMT gain check
10     * Date 09 Aug 2007
11     * Last version 08 Oct 2007
12 campana 1.1 *
13     * Description:
14     * Describe the performance of the TOF.
15     *
16     * Parameters:
17     * TString base - the path to the root directory for the specific Pamela unpack session
18     * TString outDirectory - the path where to save the output image (Default = base)
19     * TString format - the format which will be used fo rsave the produced images (Default = "gif")
20     */
21    
22 campana 1.3 #include <TROOT.h>
23     #include <TH1.h>
24     #include <TFile.h>
25     #include <TObjArray.h>
26 campana 1.1 #include <TString.h>
27     #include <TObjString.h>
28     #include <TTree.h>
29     #include <TBranch.h>
30     #include <TGraph.h>
31     #include <TStyle.h>
32     #include <TH2S.h>
33     #include <TPaveText.h>
34     #include <TCanvas.h>
35     #include <physics/tof/TofEvent.h>
36    
37 pam-de 1.4 #include <iostream>
38     #include <fstream>
39 campana 1.1
40     using namespace std;
41    
42     void TofScan(TString base, TString outDirectory = "", TString format = ""){
43    
44     std::stringstream sst;
45     if (outDirectory == "") outDirectory = base.Data();
46     TString filename = ((TObjString*)base.Tokenize('/')->Last())->GetString();
47    
48     TFile *file =new TFile(base.Data()) ;
49     if (!file){
50     printf("file not Found \n");
51     return;
52     }
53    
54     TTree *PhysicsTr = (TTree*)file->Get("Physics");
55     TBranch *TofBr = PhysicsTr->GetBranch("Tof");
56     pamela::tof::TofEvent *tofEvent = 0;
57     PhysicsTr->SetBranchAddress("Tof", &tofEvent);
58    
59     Long64_t nevents = TofBr->GetEntries();
60     if (nevents <= 0) {
61     printf("nevents = %llu \n", nevents);
62     file->Close();
63     return;
64     }
65    
66     /*
67     * Array to convert hdc/adc to the real Photomultiplier
68     * The array rows definitions are:
69     * tof[0][] = chxxA (strip or channel xxA)
70     * tof[1][] = hbxxA (halfboard xxA)
71     * tof[2][] = chxxB (strip or channel xxB)
72     * tof[3][] = hbxxB (halfboard xxB)
73     *
74     * Each single row is a sequence of photomultipliers in this shape
75     * - The elements from 0 to 7 correspond to S11_1->S11_8
76     * - The elements from 8 to 13 correspond to S12_1->S12_6
77     * - The elements from 14 to 15 correspond to S21_1->S21_2
78     * - The elements from 16 to 17 correspond to S22_1->S22_2
79     * - The elements from 18 to 20 correspond to S31_1->S31_3
80     * - The elements from 21 to 23 correspond to S32_1->S32_3
81     *
82     * Example:
83     * -------> the tdc of the S12_3B photomultiplier correspond to tdc[(tof[2][10])][(tof[3][10])]
84     * -------> the tdc of the S31_3A photomultiplier correspond to tdc[(tof[0][20])][(tof[1][20])]
85     */
86     short tof[4][24] = {
87     {4, 4, 4, 4, 1, 1, 2, 2, 3, 3, 3, 3, 3, 3, 1, 1, 1, 1, 2, 3, 3, 3, 3, 4},
88     {1, 3, 5, 7, 10, 12, 2, 4, 2, 4, 6, 8, 10, 12, 1, 5, 3, 9, 7, 9, 11, 1, 5, 9},
89     {2, 2, 2, 2, 1, 1, 1, 1, 4, 4, 4, 4, 4, 4, 2, 1, 2, 1, 2, 2, 2, 3, 3, 4},
90     {6, 8, 12, 10, 8, 6, 4, 2, 12, 10, 8, 6, 4, 2, 9, 7, 11, 11, 5, 3, 1, 3, 7, 11}
91     };
92    
93     TString photoS[48] = {
94     "S11_1A", "S11_1B", "S11_2A", "S11_2B", "S11_3A", "S11_3B", "S11_4A", "S11_4B",
95     "S11_5A", "S11_5B", "S11_6A", "S11_6B", "S11_7A", "S11_7B", "S11_8A", "S11_8B",
96     "S12_1A", "S12_1B", "S12_2A", "S12_2B", "S12_3A", "S12_3B", "S12_4A", "S12_4B", "S12_5A", "S12_5B", "S12_6A", "S12_6B",
97     "S21_1A", "S21_1B", "S21_2A", "S21_2B",
98     "S22_1A", "S22_1B", "S22_2A", "S22_2B",
99     "S31_1A", "S31_1B", "S31_2A", "S31_2B", "S31_3A", "S31_3B",
100     "S32_1A", "S32_1B", "S32_2A", "S32_2B", "S32_3A", "S32_3B"
101     };
102    
103 campana 1.3 const Int_t nh = 48;
104     TH1F *htdc[nh];
105     TH1F *hadc[nh];
106    
107     TObjArray *hhtdc = new TObjArray(nh);
108     TObjArray *hhadc = new TObjArray(nh);
109     char tdcname[48]="";
110     char adcname[48]="";
111 pam-de 1.4
112     char htitle[50];
113     TH1F *adche[48];
114     for(int i=0;i<48;i++) {
115     sprintf(htitle, "adche_%d",(i+1));
116     adche[i] = new TH1F(htitle,htitle,100,0.,1500.);
117     }
118    
119    
120     Float_t adca[48]; // vector with adc values according to "ind"=pmt_id
121     Float_t tdca[48]; // the same for tdc
122 campana 1.3
123 campana 1.1 int j = 0;
124     int k = 0;
125     int z = 0;
126     int ch = 0;
127     int hb = 0;
128 campana 1.3 int ind =0;
129    
130 pam-de 1.4 int heevent =0;
131    
132     // upper and lower limits for the helium selection
133     Float_t A_l[24]={200,190,300,210,220,200,210,60, 60, 120,220,120,160,50, 300,200,120,250,350,300,350,250,280,300};
134     Float_t A_h[24]={550,490,800,600,650,600,600,260,200,380,620,380,550,200,850,560,400,750,900,800,880,800,750,800};
135    
136     // The k1 constants for the beta calculation, only for S1-S3
137     // k2 constant is taken to be the standard 2D/c
138     Float_t k1[72] = {50,59.3296,28.4328,-26.0818,5.91253,-19.588,-9.26316,24.7544,2.32465,
139     -50.5058,-15.3195,-39.1443,-91.2546,-58.6243,-84.5641,-63.1516,-32.2091,
140     -58.3358,13.8084,45.5322,33.2416,-11.5313,51.3271,75,-14.1141,
141     42.8466,15.1794,-63.6672,-6.07739,-32.164,-41.771,10.5274,-9.46096,
142     -81.7404,-28.783,-52.7167,-127.394,-69.6166,-93.4655,-98.9543,-42.863,
143     -67.8244,-19.3238,31.1221,8.7319,-43.1627,5.55573,-14.4078,-83.4466,
144     -47.4647,-77.8379,-108.222,-75.986,-101.297,-96.0205,-63.1881,-90.1372,
145     -22.7347,8.31409,-19.6912,-7.49008,23.6979,-1.66677,1.81556,34.4668,
146     6.23693,-100,-59.5861,-90.9159,-141.639,-89.2521,-112.881} ;
147    
148     //-------------------------------------------------------------------
149    
150 campana 1.3
151     for (int i=0; i < nevents; i++){
152    
153 campana 1.1 TofBr->GetEntry(i);
154 campana 1.3
155 campana 1.1 k = 0;
156     while (k < 24){
157     j = 0;
158     while (j < 2){
159     ch = tof[2*j][k] - 1;
160     hb = tof[2*j + 1][k] - 1;
161 campana 1.3 ind = 2*k + j;
162    
163     if(i==0){
164     sprintf(tdcname,"TDChist%4.4d",ind);
165     sprintf(adcname,"ADChist%4.4d",ind);
166    
167     htdc[ind] = new TH1F(tdcname,tdcname,409,0,4096);
168     hadc[ind] = new TH1F(adcname,adcname,409,0,4096);
169    
170     hhtdc->Add(htdc[ind]);
171     hhadc->Add(hadc[ind]);
172     }
173    
174     htdc[ind]->Fill(tofEvent->tdc[ch][hb]);
175     hadc[ind]->Fill(tofEvent->adc[ch][hb]);
176 pam-de 1.4 tdca[ind]=tofEvent->tdc[ch][hb];
177     adca[ind]=tofEvent->adc[ch][hb];
178 campana 1.1 j++;
179     }
180     k++;
181     }
182 pam-de 1.4
183     //============ calculate beta and select helium ====================
184    
185     // find hitted paddle by looking for ADC values on both sides
186     // since we looking for helium this gives decent results
187    
188     Int_t tof11_i,tof12_i,tof21_i,tof22_i,tof31_i,tof32_i;
189     Float_t a1,a2;
190     Int_t jj;
191    
192     // reset values
193     tof11_i = -1;
194     tof12_i = -1;
195     tof21_i = -1;
196     tof22_i = -1;
197     tof31_i = -1;
198     tof32_i = -1;
199    
200     for(jj=0; jj<8; jj++){
201     a1 = adca[2*jj];
202     a2 = adca[2*jj+1];
203     if ((a1 < 3000) && (a2 < 3000)) tof11_i = jj;
204     }
205     for(jj=0; jj<6; jj++){
206     a1 = adca[16+2*jj];
207     a2 = adca[16+2*jj+1];
208     if ((a1 < 3000) && (a2 < 3000)) tof12_i = jj;
209     }
210     for(jj=0; jj<2; jj++){
211     a1 = adca[28+2*jj];
212     a2 = adca[28+2*jj+1];
213     if ((a1 < 3000) && (a2 < 3000)) tof21_i = jj;
214     }
215     for(jj=0; jj<2; jj++){
216     a1 = adca[32+2*jj];
217     a2 = adca[32+2*jj+1];
218     if ((a1 < 3000) && (a2 < 3000)) tof22_i = jj;
219     }
220     for(jj=0; jj<3; jj++){
221     a1 = adca[36+2*jj];
222     a2 = adca[36+2*jj+1];
223     if ((a1 < 3000) && (a2 < 3000)) tof31_i = jj;
224     }
225     for(jj=0; jj<3; jj++){
226     a1 = adca[42+2*jj];
227     a2 = adca[42+2*jj+1];
228     if ((a1 < 3000) && (a2 < 3000)) tof32_i = jj;
229     }
230    
231    
232     //----------------------------------------------------------------
233    
234     Float_t zin[6] = {53.74, 53.04, 23.94, 23.44, -23.49, -24.34};
235     Float_t c1,c2,xhelp,xhelp1,xhelp2,ds,dist,F;
236     Float_t sw,sxw,beta_mean_tof,w_i;
237     Float_t theta,x1,x2,y1,y2,dx,dy,dr;
238     Int_t ihelp;
239     Int_t ipmt[4];
240     Float_t time[4];
241     Float_t beta1[4];
242    
243     // Only use events with: S11 and S12 and S31 and S32
244    
245     if ( (tof11_i>-1) && (tof12_i>-1) && (tof31_i>-1) && (tof32_i>-1) ) {
246    
247     // calculate zenith angle theta using the locations of the hitted paddles
248    
249    
250     Float_t tof11_x[8] = {-17.85,-12.75,-7.65,-2.55,2.55,7.65,12.75,17.85};
251     Float_t tof12_y[6] = { -13.75,-8.25,-2.75,2.75,8.25,13.75};
252     // Float_t tof21_y[2] = { 3.75,-3.75};
253     // Float_t tof22_x[2] = { -4.5,4.5};
254     Float_t tof31_x[3] = { -6.0,0.,6.0};
255     Float_t tof32_y[3] = { -5.0,0.0,5.0};
256    
257     // S11 8 paddles 33.0 x 5.1 cm
258     // S12 6 paddles 40.8 x 5.5 cm
259     // S21 2 paddles 18.0 x 7.5 cm
260     // S22 2 paddles 15.0 x 9.0 cm
261     // S31 3 paddles 15.0 x 6.0 cm
262     // S32 3 paddles 18.0 x 5.0 cm
263    
264     x1 = 0.;
265     x2 = 0.;
266     y1 = 0.;
267     y2 = 0.;
268    
269     x1 = tof11_x[tof11_i] ;
270     y1 = tof12_y[tof12_i] ;
271     x2 = tof31_x[tof31_i] ;
272     y2 = tof32_y[tof32_i] ;
273    
274     theta=0.;
275     dx=0.;
276     dy=0.;
277     dr=0.;
278    
279     dx = x1-x2;
280     dy = y1-y2;
281     dr = sqrt(dx*dx+dy*dy);
282     theta = atan(dr/77.5);
283    
284    
285     beta_mean_tof=100.;
286    
287     for (Int_t jj=0; jj< 4; jj++) beta1[jj] = 100. ;
288    
289    
290     //----------------------------------------------------------------
291     //--------- S1 - S3 ---------------------------------------------
292     //----------------------------------------------------------------
293    
294     //--------- S11 - S31 -------------------------------------------
295    
296     if ((tof11_i>-1)&&(tof31_i>-1)) {
297    
298     dist = zin[0] - zin[4];
299     c2 = (2.*0.01*dist)/(3.E08*50.E-12);
300     F = 1./cos(theta);
301    
302     ipmt[0] = (tof11_i)*2;
303     ipmt[1] = (tof11_i)*2+1;
304     ipmt[2] = 36+(tof31_i)*2;
305     ipmt[3] = 36+(tof31_i)*2+1;
306    
307     for (Int_t jj=0; jj< 4; jj++) time[jj] = tdca[(ipmt[jj])] ;
308    
309     if ((time[0]<4095)&&(time[1]<4095)&&(time[2]<4095)&&(time[3]<4095)) {
310     xhelp1 = time[0] + time[1] ;
311     xhelp2 = time[2] + time[3] ;
312     ds = xhelp1-xhelp2;
313     ihelp=0+(tof11_i)*3+tof31_i ;
314     c1 = k1[ihelp] ;
315     beta1[0] = c2*F/(ds-c1);
316     }
317     }
318    
319     //--------- S11 - S32 -------------------------------------------
320    
321     if ((tof11_i>-1)&&(tof32_i>-1)) {
322    
323     dist = zin[0] - zin[5];
324     F = 1./cos(theta);
325     c2 = (2.*0.01*dist)/(3.E08*50.E-12);
326    
327     ipmt[0] = (tof11_i)*2;
328     ipmt[1] = (tof11_i)*2+1;
329     ipmt[2] = 42+(tof32_i)*2;
330     ipmt[3] = 42+(tof32_i)*2+1;
331    
332     for (Int_t jj=0; jj< 4; jj++) time[jj] = tdca[(ipmt[jj])] ;
333    
334     if ((time[0]<4095)&&(time[1]<4095)&&(time[2]<4095)&&(time[3]<4095)) {
335     xhelp1 = time[0] + time[1] ;
336     xhelp2 = time[2] + time[3] ;
337     ds = xhelp1-xhelp2;
338     ihelp=24+(tof11_i)*3+tof32_i ;
339     c1 = k1[ihelp] ;
340     beta1[1] = c2*F/(ds-c1);
341     }
342     }
343    
344     //--------- S12 - S31 -------------------------------------------
345    
346     if ((tof12_i>-1)&&(tof31_i>-1)) {
347    
348     dist = zin[1] - zin[4];
349     F = 1./cos(theta);
350     c2 = (2.*0.01*dist)/(3.E08*50.E-12);
351    
352     ipmt[0] = 16+(tof12_i)*2;
353     ipmt[1] = 16+(tof12_i)*2+1;
354     ipmt[2] = 36+(tof31_i)*2;
355     ipmt[3] = 36+(tof31_i)*2+1;
356    
357     for (Int_t jj=0; jj< 4; jj++) time[jj] = tdca[(ipmt[jj])] ;
358    
359     if ((time[0]<4095)&&(time[1]<4095)&&(time[2]<4095)&&(time[3]<4095)) {
360     xhelp1 = time[0] + time[1] ;
361     xhelp2 = time[2] + time[3] ;
362     ds = xhelp1-xhelp2;
363     ihelp=48+(tof12_i)*3+tof31_i ;
364     c1 = k1[ihelp] ;
365     beta1[2] = c2*F/(ds-c1);
366     }
367     }
368    
369     //--------- S12 - S32 -------------------------------------------
370    
371     if ((tof12_i>-1)&&(tof32_i>-1)) {
372    
373     dist = zin[1] - zin[5];
374     F = 1./cos(theta);
375     c2 = (2.*0.01*dist)/(3.E08*50.E-12);
376    
377     ipmt[0] = 16+(tof12_i)*2;
378     ipmt[1] = 16+(tof12_i)*2+1;
379     ipmt[2] = 42+(tof32_i)*2;
380     ipmt[3] = 42+(tof32_i)*2+1;
381    
382     for (Int_t jj=0; jj< 4; jj++) time[jj] = tdca[(ipmt[jj])] ;
383    
384     if ((time[0]<4095)&&(time[1]<4095)&&(time[2]<4095)&&(time[3]<4095)) {
385     xhelp1 = time[0] + time[1] ;
386     xhelp2 = time[2] + time[3] ;
387     ds = xhelp1-xhelp2;
388     ihelp=66+(tof12_i)*3+tof32_i ;
389     c1 = k1[ihelp] ;
390     beta1[3] = c2*F/(ds-c1);
391     }
392     }
393    
394     //---------------------- calculate beta mean -----------------
395    
396     sw=0.;
397     sxw=0.;
398     beta_mean_tof=100.;
399    
400     for (Int_t jj=0; jj<4;jj++){
401     if ((beta1[jj]>0.1) && (beta1[jj]<1.5)) {
402     w_i=1./(0.13*0.13);
403     sxw=sxw + beta1[jj]*w_i ;
404     sw =sw + w_i ;
405     }
406     }
407    
408     if (sw>0) beta_mean_tof=sxw/sw;
409    
410     } // if tof11_i > -1 && ...... beta calculation
411    
412    
413     Float_t beta_help = beta_mean_tof ; // pow(beta_mean_tof,1.0) gave best results
414    
415     //----------------------- Select helium --------------------------
416    
417     Int_t icount=0;
418    
419     for (jj=0; jj<24; jj++){
420     a1 = adca[2*jj]*cos(theta);
421     a2 = adca[2*jj+1]*cos(theta);
422    
423     xhelp = 100000.;
424     if ((a1 < 3000) && (a2 < 3000)) xhelp = sqrt(a1*a2); // geometric mean
425     // if geometric mean multiplied by beta_help is inside helium limits, increase counter
426     if ((beta_mean_tof>0.6) && (beta_mean_tof<1.1) &&
427     ((beta_help*xhelp)>A_l[jj]) && ((beta_help*xhelp)<A_h[jj])) icount++ ;
428     }
429    
430     Int_t iz=0;
431     // if (icount > 3) iz=2; // if more than three paddles see helium, then set Z=2
432     if (icount > 4) iz=2;
433    
434     //---------------------- Z=2 fill histograms -----------------------------
435    
436     if (iz==2) {
437    
438     heevent++;
439     for (jj=0; jj<48; jj++) adche[jj]->Fill(adca[jj]);
440    
441     } // iz0==2
442    
443    
444     //===================== end beta and helium part ===========================
445    
446     } // i < nevents
447    
448    
449 campana 1.1 float *X = new float[48];
450     float *means = new float[48];
451     float *entries = new float[48];
452     int *entriestdc = new int[48];
453     int *entriesadc = new int[48];
454 campana 1.3
455 campana 1.1 const char *saveas = format;
456    
457 campana 1.3 int i=0;
458 campana 1.1
459     gStyle->SetStatW(0.4);
460     gStyle->SetStatH(0.4);
461     gStyle->SetOptStat("nmri");
462     gStyle->SetTitleH(0.10);
463     gStyle->SetTitleW(0.96);
464    
465     TCanvas *SCanvas = new TCanvas("SCanvas","SCanvas", 1280, 1024);
466     SCanvas->Divide(4,2);
467 campana 1.3
468 campana 1.1 j = 0;
469     while (j < 12){
470     k = 0;
471     z = 0;
472     if (gROOT->IsBatch()) {
473     SCanvas = new TCanvas("SCanvas","SCanvas", 1280, 1024);
474     SCanvas->Divide(4,2);
475     } else {
476     if (j > 0) SCanvas->DrawClone();
477     }
478    
479    
480     while(k < 4){
481     if (k > 1) z = 2;
482     i = j*4 + k;
483     X[i] = i;
484    
485     SCanvas->cd(k+3+z);
486 campana 1.3 htdc[i] = (TH1F*)hhtdc->At(i);
487     entriestdc[i] = (Int_t)htdc[i]->Integral();
488 campana 1.1 sst.str("");
489     sst << "TDC - " << photoS[i].Data() << " (Nev < 4096 = " << entriestdc[i] << ")";
490 campana 1.3 htdc[i]->SetTitle(sst.str().c_str());
491     htdc[i]->SetTitleSize(10);
492     htdc[i]->SetAxisRange(690,1510);
493     htdc[i]->DrawCopy();
494     htdc[i]->ComputeIntegral();
495     entries[i] = htdc[i]->Integral();
496 campana 1.1
497     SCanvas->cd(k+1+z);
498 campana 1.3 hadc[i] = (TH1F*)hhadc->At(i);
499     entriesadc[i] = (Int_t)hadc[i]->Integral();
500 campana 1.1 sst.str("");
501     sst << "ADC - " << photoS[i].Data() << " (Nev < 4096 = " << entriesadc[i] << ")";
502 campana 1.3 hadc[i]->SetTitle(sst.str().c_str());
503     hadc[i]->SetAxisRange(-10,710);
504     hadc[i]->DrawCopy();
505     means[i] = hadc[i]->GetMean();
506 campana 1.1
507     k++;
508     }
509 campana 1.3
510    
511 campana 1.1 if ( !strcmp(saveas,"ps") ) {
512     sst.str("");
513     sst << outDirectory.Data() << filename.Data() << "TOFScan.ps(";
514     SCanvas->Print(sst.str().c_str());
515     } else {
516     sst.str("");
517     sst << outDirectory.Data() << filename.Data() << "TOFScan" << j+1 << "." << saveas;
518     SCanvas->SaveAs(sst.str().c_str());
519    
520     }
521     j++;
522     }
523 campana 1.3
524 campana 1.1 if (gROOT->IsBatch()) SCanvas->Close();
525    
526     /*
527     * This Canvas will represent a summary of the performances for TOF TDC/ADC channels
528     */
529 pam-de 1.4 // TCanvas *performanceCanvas = new TCanvas("performanceCanvas","performanceCanvas", 1280, 1024);
530     TCanvas *performanceCanvas = new TCanvas("performanceCanvas","performanceCanvas", 1024, 1280);
531     performanceCanvas->Divide(1,3);
532 campana 1.1
533     gStyle->SetTitleW(.9);
534    
535     performanceCanvas->cd(1);
536     TGraph *adcMeans = new TGraph(48, X, means);
537     sst.str("");
538     sst << "ADCMean" << " - Data in " << base.Data() << " - Nevents in the run = " << nevents;
539     adcMeans->SetTitle(sst.str().c_str());
540 campana 1.2 adcMeans->SetFillColor(35);
541 campana 1.1 adcMeans->GetXaxis()->SetTitle("Photomultipliers");
542     adcMeans->GetXaxis()->CenterTitle();
543     adcMeans->GetXaxis()->SetLimits(-0.5, 47.5);
544     adcMeans->GetYaxis()->SetTitle("ADCMean");
545     adcMeans->GetYaxis()->CenterTitle();
546     adcMeans->Draw("AB");
547    
548     performanceCanvas->cd(2);
549     TGraph *tdcEntries = new TGraph(48, X, entries);
550     sst.str("");
551     sst << "TDCEntries" << " - Data in " << base.Data() << " - Nevents in the run = " << nevents;
552     tdcEntries->SetTitle(sst.str().c_str());
553 campana 1.2 tdcEntries->SetFillColor(35);
554 campana 1.1 tdcEntries->GetXaxis()->SetTitle("Photomultipliers");
555     tdcEntries->GetXaxis()->CenterTitle();
556     tdcEntries->GetXaxis()->SetLimits(-0.5, 47.5);
557     tdcEntries->GetYaxis()->SetTitle("TDCIntegral");
558     tdcEntries->GetYaxis()->CenterTitle();
559     tdcEntries->Draw("AB");
560 pam-de 1.4
561     //--------- new part PMT gain check -----------------------------
562    
563     performanceCanvas->cd(3);
564    
565     Float_t xc[48],xmean1[48],xmeana[48];
566     Float_t xmean_arr[12][48];
567    
568     // xmean values from 2-3 april 2007
569    
570     char date_info[]="Reference Data: apr-2007";
571    
572     Float_t xmean[48] = {
573     491.609,509.241,400.786,530.122,699.674,555.747,521.04,486.363,
574     470.173,227.752,611.038,455.889,553.601,520.54,403.527,382.099,
575     349.697,365.113,447.653,377.667,517.815,572.932,338.501,436.681,
576     485.696,450.491,395.375,329.631,751.258,626.681,385.561,578.476,
577     374.454,356.733,641.888,562.767,582.849,521.748,527.043,505.89,
578     489.828,628.408,532.924,506.511,482.872,532.236,554.554,498.849 };
579    
580     // new 01-oct-2007
581     int channelmap[] = {0,7,3,6,2,8,1,5,3,7,3,6,1,7,2,10,
582     10,10,10,5,0,7,0,5,0,6,1,5,
583     2,8,3,8,2,6,1,8,
584     11,9,11,11,9,11,4,4,4,9,9,4};
585    
586    
587     int colormap[] = {46,2,29,4,5,6,7,8,9,11,28,34};
588     //int colormap[] = {417,400,632,617,603,600,434,419,591,625,403,424};
589    
590    
591     for (Int_t j=0; j<48; j++) xmeana[j]=0.;
592     for (Int_t j=0; j<24; j++) xmeana[2*j]=xmean[2*j];
593    
594     for (Int_t i=0; i<12; i++) {
595     for (Int_t j=0; j<48; j++) {
596     xmean_arr[i][j]=0.;
597     }
598     }
599    
600     for (Int_t j=0; j<48; j++) {
601     Int_t ichan = channelmap[j];
602     xmean_arr[ichan][j]=xmean[j];
603     }
604    
605     // get results from ADC histogram
606     for (Int_t j=0; j<48; j++) {
607     xc[j]=j;
608     xmean1[j]=adche[j]->GetMean();
609     }
610    
611    
612     gStyle->SetTitleW(.5);
613     gStyle->SetTitleH(.05);
614    
615     TH2F *hr = new TH2F("frame","2-Dim",2,-0.5,47.5,2,-300.,100.);
616     hr->SetStats(kFALSE);
617     hr->GetXaxis()->CenterTitle();
618     hr->GetXaxis()->SetTitle("Photomultipliers");
619     hr->GetYaxis()->CenterTitle();
620     hr->GetYaxis()->SetTitle("Mean ADC Difference");
621     hr->SetTitle("Difference between Reference and Actual Values");
622     hr->Draw();
623    
624     Int_t npoint=48;
625    
626     for (Int_t j=0; j<12; j++) {
627     for (Int_t i=0; i<48; i++) xmeana[i] = 0.;
628     for (Int_t i=0; i<48; i++) {
629     if (xmean_arr[j][i] != 0) xmeana[i] = xmean1[i] - xmean_arr[j][i];
630     }
631    
632    
633     TGraph *graph1 = new TGraph(npoint,xc,xmeana);
634     graph1->SetFillColor(colormap[j]);
635     graph1->GetXaxis()->SetLimits(-0.5, 47.5);
636     graph1->Draw("BP");
637     }
638    
639     Float_t tp[10];
640     tp[0] = 15.5;
641     tp[1] = 27.5;
642     tp[2] = 31.5;
643     tp[3] = 35.5;
644     tp[4] = 41.5;
645    
646     for (Int_t ii=0; ii<5; ii++) {
647     TLine *l1=new TLine(tp[ii],-300,tp[ii],100);
648     l1->SetLineColor(38);
649     l1->Draw("same");
650     }
651    
652     for (Int_t j=0; j<12; j++) {
653     sprintf(htitle, "HV_%d",j);
654     TText *text1 = new TText(0+j*4,80,htitle);
655     text1->SetTextColor(colormap[j]);
656     //text1->SetTextSize(0.03);
657     text1->SetTextSize(0.05);
658     text1->Draw();
659     }
660    
661    
662     TText *text1 = new TText(0,-185,date_info);
663     text1->SetTextColor(kBlack);
664     text1->SetTextSize(0.023);
665     text1->Draw();
666    
667    
668     sprintf(htitle, "Helium Events: %d",heevent);
669     TText *text2 = new TText(20,-185,htitle);
670     text2->SetTextColor(kBlack);
671     text2->SetTextSize(0.023);
672     text2->Draw();
673    
674    
675     for (Int_t i=0; i<6; i++) {
676     for (Int_t j=0; j<8; j++) {
677     Int_t ihelp = i*8+j;
678     sprintf(htitle, "%d: %.0f/%.0f",(ihelp+1),xmean[ihelp],xmean1[ihelp]);
679     TText *text1 = new TText(0+j*6,-200-i*15,htitle);
680     text1->SetTextColor(kBlack);
681     text1->SetTextSize(0.023);
682     text1->Draw();
683     }
684     }
685    
686     //-------- end new part -------------------------
687    
688    
689 campana 1.1 //------print the ps
690    
691     if ( !strcmp(saveas,"ps") ) {
692     sst.str("");
693     sst << outDirectory.Data() << filename.Data() << "TOFScan.ps)";
694     performanceCanvas->Print(sst.str().c_str());
695    
696     } else {
697     sst.str("");
698     sst << outDirectory.Data() << filename.Data() << "TOFScan13." << saveas;
699     performanceCanvas->SaveAs(sst.str().c_str());
700     }
701     if (gROOT->IsBatch()) {
702     SCanvas->Close();
703     performanceCanvas->Close();
704     }
705    
706 pam-de 1.4
707    
708 campana 1.1 }
709    
710     int main(int argc, char* argv[]){
711     TString path;
712     TString outDir ="./";
713     TString format ="ps";
714    
715     if (argc < 2){
716     printf("You have to insert at least the file to analyze \n");
717     printf("Try '--help' for more information. \n");
718     exit(1);
719     }
720    
721     if (!strcmp(argv[1], "--help")){
722     printf( "Usage: TofScan FILE [OPTION] \n");
723     printf( "\t --help Print this help and exit \n");
724     printf( "\t -outDir[path] Path where to put the output [default ./] \n");
725     printf( "\t -format[ps] Format for output files [default 'ps'] \n");
726     exit(1);
727     }
728    
729    
730     path=argv[1];
731    
732     for (int i = 2; i < argc; i++){
733    
734     if (!strcmp(argv[i], "-outDir")){
735     if (++i >= argc){
736     printf( "-outDir needs arguments. \n");
737     printf( "Try '--help' for more information. \n");
738     exit(1);
739     }
740     else{
741     outDir = argv[i];
742     continue;
743     }
744     }
745    
746    
747    
748     if (!strcmp(argv[i], "-format")){
749     if (++i >= argc){
750     printf( "-format needs arguments. \n");
751     printf( "Try '--help' for more information. \n");
752     exit(1);
753     }
754     else{
755     format = argv[i];
756     continue;
757     }
758     }
759     }
760    
761     TofScan(argv[1], outDir, format);
762    
763     }
764 pam-de 1.4

  ViewVC Help
Powered by ViewVC 1.1.23