/[PAMELA software]/quicklook/tracker/flight/macros/FTrkCalibQLook_BASIC.cxx
ViewVC logotype

Annotation of /quicklook/tracker/flight/macros/FTrkCalibQLook_BASIC.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (hide annotations) (download)
Thu Jun 1 15:23:12 2006 UTC (18 years, 6 months ago) by pam-fi
Branch: MAIN
Changes since 1.1: +60 -50 lines
*** empty log message ***

1 pam-fi 1.1 /**
2     * FTrkCalibQLook_BASIC.cxx
3     *
4     * autor: D.Fedele
5     * version 2.0
6     * Parameters:
7     * file - the data file to analyze
8     * fromevent - first event to analyze
9     * toevent - last event to analyze
10     * outdir - total path of output file
11     * outfile - extension of output file (pdf,ps,gif,jpg)
12     *
13     */
14     //
15     #include <TLatex.h>
16     #include <TCanvas.h>
17     #include <TGraph.h>
18     #include <TTree.h>
19     #include <TStyle.h>
20     #include <TString.h>
21     //
22     #include <PscuHeader.h>
23     #include <EventHeader.h>
24     #include <CalibTrk1Event.h>
25     #include <CalibTrk2Event.h>
26     //
27    
28     void stringcopy(TString& s1, const TString& s2, Int_t from=0, Int_t to=0){
29     if ( to == 0 ){
30     Int_t t2length = s2.Length();
31     s1 = "";
32     to = t2length;
33     };
34     for (Int_t i = from; i<to; i++){
35     s1.Append(s2[i],1);
36     };
37     };
38    
39     typedef struct caltrk_def{
40     Int_t good0[2];
41     Int_t daqmode[12];
42     Int_t dspnum[12];
43     Int_t calibnum[12];
44     Int_t ncalev[12];
45     Int_t calfl[12];
46     Int_t ped1[12];
47     Int_t ped2[12];
48     Int_t ped3[12];
49     Int_t sig1[12];
50     Int_t sig2[12];
51     Int_t sig3[12];
52     Int_t nbad1[12];
53     Int_t nbad2[12];
54     Int_t nbad3[12];
55     Float_t dspped[12][3702];
56     Float_t dspsig[12][3702];
57     Float_t dspbad[12][3702];
58     Int_t crc_c[12][3];
59     Int_t crc_hc[12];
60     };
61    
62    
63     void FTrkCalibQLook_BASIC(TString file,Int_t fromevent,Int_t toevent,TString outdir,TString outfile)
64     {
65     //
66     // obtain information about the data file and select the output dir
67     const string filepath=file.Data();
68 pam-fi 1.2 Int_t dwpos = filepath.rfind("/");
69     Int_t dwpos1 = filepath.rfind(".root");
70 pam-fi 1.1 TString fpath=(filepath.c_str());
71     TString base,ffile;
72     stringcopy(base,fpath,0,dwpos);
73 pam-fi 1.2 stringcopy(ffile,fpath,dwpos+1,dwpos1);
74     if(dwpos>0) base+="/";
75 pam-fi 1.1
76     TString out;
77     if(outdir.Length()==0){
78     out = base;
79     }else{
80     out = outdir;
81     };
82    
83     //
84     // inizialise the variables and open the file
85     struct caltrk_def ctrk;
86     Int_t nevents=0;
87     Int_t minevent = 0;
88     Int_t maxevent = 0;
89     ULong64_t OBT[2];
90    
91     OBT[0]=0;
92     OBT[1]=0;
93     ctrk.good0[0]=0;
94     ctrk.good0[1]=0;
95     for(Int_t i=0;i<12;i++){
96     ctrk.daqmode[i]=0;
97     ctrk.dspnum[i]=0;
98     ctrk.calibnum[i]=0;
99     ctrk.ncalev[i]=0;
100     ctrk.calfl[i]=0;
101     ctrk.ped1[i]=0;
102     ctrk.ped2[i]=0;
103     ctrk.ped3[i]=0;
104     ctrk.sig1[i]=0;
105     ctrk.sig2[i]=0;
106     ctrk.sig3[i]=0;
107     ctrk.nbad1[i]=0;
108     ctrk.nbad2[i]=0;
109     ctrk.nbad3[i]=0;
110     ctrk.crc_hc[i]=0;
111     ctrk.crc_c[i][0]=0;
112     ctrk.crc_c[i][1]=0;
113     ctrk.crc_c[i][2]=0;
114     for(Int_t iii=0;iii<3072;iii++){
115     ctrk.dspped[i][iii]=0;
116     ctrk.dspsig[i][iii]=0;
117     ctrk.dspbad[i][iii]=0;
118     }
119     }
120    
121     pamela::EventHeader *eh1=0,*eh2=0;
122     pamela::PscuHeader *ph1=0,*ph2=0;
123     pamela::CalibTrk1Event *trk1 = 0;
124     pamela::CalibTrk2Event *trk2 = 0;
125    
126     TFile *datafile = new TFile(file);
127     if ( !datafile ){
128     printf("No data file, exiting...\n");
129     return;
130     };
131     TTree *otr1,*otr2;
132    
133     otr1 = (TTree*)datafile->Get("CalibTrk1");
134     otr1->SetBranchAddress("CalibTrk1", &trk1);
135     otr1->SetBranchAddress("Header",&eh1);
136     otr2 = (TTree*)datafile->Get("CalibTrk2");
137     otr2->SetBranchAddress("CalibTrk2", &trk2);
138     otr2->SetBranchAddress("Header",&eh2);
139    
140    
141     if(otr1->GetEntries()==otr2->GetEntries())
142     nevents = otr1->GetEntries();
143     else{
144     printf("WARNING: CalibTrk1 entries is different from CalibTrk2 entries");
145     return;}
146    
147     if (nevents<=0) {
148     datafile->Close();
149     printf("No calibration packets found, exiting...\n");
150     return;
151     };
152     printf("Number of calibration packets: %i\n",nevents);
153    
154     if ( fromevent > toevent && toevent ){
155     printf("It must be fromevent < toevent \n");
156     return;
157     };
158    
159     if ( fromevent > nevents || fromevent < 0 ) {
160     printf("You can choose fromevent between 0 (all) and %i \n",nevents);
161     return;
162     };
163    
164     if ( toevent > nevents || toevent < 0 ) {
165     printf("You can choose toevent between 0 (all) and %i \n",nevents);
166     return;
167     };
168     if ( fromevent == 0 ) {
169     minevent = 0;
170     maxevent = nevents;
171     } else {
172     minevent = fromevent - 1;
173     if ( toevent > 0 ){
174     maxevent = toevent;
175     } else if (toevent > nevents) {
176     maxevent = nevents;
177     } else {
178     maxevent = fromevent;
179     };
180     };
181    
182     //
183     // other variables definitions
184     stringstream fromfile,rep,tit;
185     fromfile<<"FTrkCalibQLook_BASIC File: "<<ffile;
186    
187     gStyle->SetLabelSize(0.08,"x");
188     gStyle->SetLabelSize(0.08,"y");
189     gStyle->SetFillColor(0);
190     gStyle->SetTitleFillColor(0);
191     gStyle->SetTitleFontSize(0.1);
192     gStyle->SetTitleOffset(0.8,"y");
193     gStyle->SetTitleOffset(1.,"x");
194     gStyle->SetTitleSize(0.06,"y");
195     gStyle->SetTitleSize(0.06,"x");
196     gStyle->SetOptStat(0);
197    
198     TLatex *tzz=new TLatex();
199     tzz->SetTextFont(32);
200     tzz->SetTextColor(1);
201     tzz->SetTextAlign(12);
202     tzz->SetTextSize(0.02);
203    
204     Int_t canvasx=900;
205     Int_t canvasy=1200;
206    
207     Int_t ndsp =0,alarm=0;
208     Float_t pedav[12][12],pedavtemp[12][12],sigav[12][12],sigavtemp[12][12];
209     Int_t flpedav[12][12],flsigav[12][12];
210     Float_t siglimsup[12][12],sigliminf[12][12],pedlimsup[12][12],pedliminf[12][12];
211    
212     //
213     // inizialize the limits for simga and pedestall
214     for(Int_t i=0;i<12;i++){
215     for(Int_t ii=0;ii<12;ii++){
216 pam-fi 1.2 siglimsup[i][ii]=30.;
217     sigliminf[i][ii]=1.5;
218     if(!(i%2)){
219     pedlimsup[i][ii]=3700.;
220     pedliminf[i][ii]=1700.;
221     }
222     else{
223     pedlimsup[i][ii]=2200.;
224     pedliminf[i][ii]=200.;
225     }
226 pam-fi 1.1 }
227     }
228    
229     //
230     // count of possible alarm to set the number of output pages
231     alarm+=(maxevent-minevent)/9 +1;
232     for (Int_t i = minevent; i < maxevent; i++){
233     otr1->GetEntry(i);
234     otr2->GetEntry(i);
235     for(Int_t m=0;m<6;m++){
236     for(Int_t mm=0;mm<3;mm++){
237     if(trk1->crc_cal[m][mm]!=0) alarm+=1;
238     if(trk2->crc_cal[m][mm]!=0 ) alarm+=1;
239     }
240     if(trk1->crc_hcal[m]!=0 ) alarm+=1;
241     if(trk1->cal_flag[m]!=0 ) alarm+=1;
242     if(trk2->crc_hcal[m]!=0) alarm+=1;
243     if(trk2->cal_flag[m]!=0) alarm+=1;
244     }
245     }
246 pam-fi 1.2 const Int_t cnum=alarm/30 + 100;
247 pam-fi 1.1 Int_t flcanvas=1;
248     Float_t spacep=1.5,space[cnum];
249 pam-fi 1.2 TCanvas *c[cnum];
250 pam-fi 1.1 //
251     // create output canvas
252     for(Int_t i=0;i<cnum;i++){
253     space[i]=96.0;
254     rep.str("");
255     tit.str("");
256     tit<<"c"<<i;
257     rep<<"FTrkCalibQLook_BASIC_pag"<<i+1;
258     c[i]=new TCanvas(tit.str().c_str(),rep.str().c_str(),canvasx,canvasy);
259     c[i]->Range(0,0,100,100);
260 pam-fi 1.2 c[i]->SetFillColor(10);
261 pam-fi 1.1 tzz->DrawLatex(1,98.5,fromfile.str().c_str());
262     rep.str("");
263     rep<<"CALIBRATION REPORT pag"<<i+1;
264     tzz->DrawLatex(70,98.5,rep.str().c_str());
265     rep.str("");
266     }
267    
268    
269     //**********************************************************************
270     //
271     // LOOP OVER EVENTS
272     //
273     //**********************************************************************
274    
275     Int_t wc=0;
276     printf("\n Scan of calibration packets from %i to %i ... \n",minevent+1,maxevent);
277     for (Int_t i = minevent; i < maxevent; i++){
278    
279     otr1->GetEntry(i);
280     otr2->GetEntry(i);
281    
282     ctrk.good0[0]=trk1->good0;
283     ctrk.good0[1]=trk2->good0;
284     for (Int_t m = 0; m < 6; m++){
285     ph1 = eh1->GetPscuHeader();
286     OBT[0]= ph1->GetOrbitalTime();
287     ctrk.daqmode[trk1->DSPnumber[m]-1]=trk1->DAQmode[m];
288     ctrk.dspnum[trk1->DSPnumber[m]-1]=trk1->DSPnumber[m];
289     ctrk.calibnum[trk1->DSPnumber[m]-1]=trk1->calibnumber[m];
290     ctrk.ncalev[trk1->DSPnumber[m]-1]=trk1->ncalib_event[m];
291     ctrk.ped1[trk1->DSPnumber[m]-1]=trk1->ped_l1[m];
292     ctrk.ped2[trk1->DSPnumber[m]-1]=trk1->ped_l2[m];
293     ctrk.ped3[trk1->DSPnumber[m]-1]=trk1->ped_l3[m];
294     ctrk.sig1[trk1->DSPnumber[m]-1]=trk1->sig_l1[m];
295     ctrk.sig2[trk1->DSPnumber[m]-1]=trk1->sig_l2[m];
296     ctrk.sig3[trk1->DSPnumber[m]-1]=trk1->sig_l3[m];
297     ctrk.nbad1[trk1->DSPnumber[m]-1]=trk1->nbad_l1[m];
298     ctrk.nbad2[trk1->DSPnumber[m]-1]=trk1->nbad_l2[m];
299     ctrk.nbad3[trk1->DSPnumber[m]-1]=trk1->nbad_l3[m];
300     ctrk.calfl[trk1->DSPnumber[m]-1]=trk1->cal_flag[m];
301     ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk1->crc_cal[m][0];
302     ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk1->crc_cal[m][1];
303     ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk1->crc_cal[m][2];
304     ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk1->crc_hcal[m];
305     for (Int_t j = 0; j < 3072; j++){
306     ctrk.dspped[trk1->DSPnumber[m]-1][j]=trk1->DSPped_par[m][j];
307     ctrk.dspsig[trk1->DSPnumber[m]-1][j]=trk1->DSPsig_par[m][j];
308     ctrk.dspbad[trk1->DSPnumber[m]-1][j]=trk1->DSPbad_par[m][j];
309     }
310     ph2 = eh2->GetPscuHeader();
311     OBT[1]= ph2->GetOrbitalTime();
312     ctrk.daqmode[trk2->DSPnumber[m]-1]=trk2->DAQmode[m];
313     ctrk.dspnum[trk2->DSPnumber[m]-1]=trk2->DSPnumber[m];
314     ctrk.calibnum[trk2->DSPnumber[m]-1]=trk2->calibnumber[m];
315     ctrk.ncalev[trk2->DSPnumber[m]-1]=trk2->ncalib_event[m];
316     ctrk.ped1[trk2->DSPnumber[m]-1]=trk2->ped_l1[m];
317     ctrk.ped2[trk2->DSPnumber[m]-1]=trk2->ped_l2[m];
318     ctrk.ped3[trk2->DSPnumber[m]-1]=trk2->ped_l3[m];
319     ctrk.sig1[trk2->DSPnumber[m]-1]=trk2->sig_l1[m];
320     ctrk.sig2[trk2->DSPnumber[m]-1]=trk2->sig_l2[m];
321     ctrk.sig3[trk2->DSPnumber[m]-1]=trk2->sig_l3[m];
322     ctrk.nbad1[trk2->DSPnumber[m]-1]=trk2->nbad_l1[m];
323     ctrk.nbad2[trk2->DSPnumber[m]-1]=trk2->nbad_l2[m];
324     ctrk.nbad3[trk2->DSPnumber[m]-1]=trk2->nbad_l3[m];
325     ctrk.calfl[trk2->DSPnumber[m]-1]=trk2->cal_flag[m];
326     ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk2->crc_cal[m][0];
327     ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk2->crc_cal[m][1];
328     ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk2->crc_cal[m][2];
329     ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk2->crc_hcal[m];
330     for (Int_t j = 0; j < 3072; j++){
331     ctrk.dspped[trk2->DSPnumber[m]-1][j]=trk2->DSPped_par[m][j];
332     ctrk.dspsig[trk2->DSPnumber[m]-1][j]=trk2->DSPsig_par[m][j];
333     ctrk.dspbad[trk2->DSPnumber[m]-1][j]=trk2->DSPbad_par[m][j];
334     }
335     }
336    
337     for(Int_t n = 0; n<12; n++){
338     for(Int_t nm = 0; nm<12; nm++){
339     pedav[n][nm]=0;
340     pedavtemp[n][nm]=0;
341     sigav[n][nm]=0;
342     sigavtemp[n][nm]=0;
343     flpedav[n][nm]=0;
344     flsigav[n][nm]=0;
345     }
346     }
347    
348     Int_t nn,ok=0;
349    
350     //
351     // write warning if it occur
352     for(Int_t n = 0; n<12; n++){
353    
354     ndsp = ctrk.dspnum[n];
355     nn = ndsp-1;
356    
357     for(Int_t iii=0;iii<3;iii++){
358     if(ctrk.crc_c[nn][iii]!=0){
359     ok=1;
360     if(space[wc]<=3){
361     wc+=1;
362     flcanvas+=1;
363     }
364     c[wc]->cd();
365     tzz->SetTextFont(40);
366     tzz->SetTextSize(0.02);
367     tzz->SetTextAlign(13);
368     tzz->SetTextColor(2);
369     rep<<"***************************************************************************************************************************";
370     tzz->DrawLatex(2,space[wc],rep.str().c_str());
371     rep.str("");
372     space[wc]-=spacep;
373     rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> crc_cal["<<nn+1<<"]["<<iii+1<<"]= "<<ctrk.crc_c[nn][iii];
374     tzz->DrawLatex(2,space[wc],rep.str().c_str());
375     rep.str("");
376     space[wc]-=spacep;
377     }
378     }
379    
380     if(ctrk.crc_hc[nn]!=0){
381     ok=1;
382     if(space[wc]<=3){
383     wc+=1;
384     flcanvas+=1;
385     }
386     c[wc]->cd();
387     tzz->SetTextFont(40);
388     tzz->SetTextSize(0.02);
389     tzz->SetTextAlign(13);
390     tzz->SetTextColor(2);
391     rep<<"***************************************************************************************************************************";
392     tzz->DrawLatex(2,space[wc],rep.str().c_str());
393     rep.str("");
394     space[wc]-=spacep;
395     rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> crc_hcal["<<nn+1<<"]= "<<ctrk.crc_hc[nn];
396     tzz->DrawLatex(2,space[wc],rep.str().c_str());
397     rep.str("");
398     space[wc]-=spacep;
399     }
400    
401     if(ctrk.calfl[nn]!=0){
402     ok=1;
403     if(space[wc]<=3){
404     wc+=1;
405     flcanvas+=1;
406     }
407     c[wc]->cd();
408     tzz->SetTextFont(40);
409     tzz->SetTextSize(0.02);
410     tzz->SetTextAlign(13);
411     tzz->SetTextColor(2);
412     rep<<"***************************************************************************************************************************";
413     tzz->DrawLatex(2,space[wc],rep.str().c_str());
414     rep.str("");
415     space[wc]-=spacep;
416     rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> cal_flag["<<nn+1<<"]= "<<ctrk.calfl[nn];
417     tzz->DrawLatex(2,space[wc],rep.str().c_str());
418     rep.str("");
419     space[wc]-=spacep;
420     }
421    
422     if(ctrk.good0[0]==1 && ctrk.good0[1]==1){
423     //
424     // evaluate the mean value of the sigma and pedestal
425     for(Int_t j = 0; j < 3072; j++){
426     sigavtemp[nn][j/256]+=ctrk.dspsig[nn][j];
427     pedavtemp[nn][j/256]+=ctrk.dspped[nn][j];
428     }
429    
430     for(Int_t ii=0;ii<12;ii++){
431     pedav[nn][ii]=pedavtemp[nn][ii]/256;
432     sigav[nn][ii]=sigavtemp[nn][ii]/256;
433    
434     if(pedav[nn][ii]>pedlimsup[nn][ii] || pedav[nn][ii]<pedliminf[nn][ii]) flpedav[nn][ii]=1;
435     if(sigav[nn][ii]>siglimsup[nn][ii] || sigav[nn][ii]<sigliminf[nn][ii]) flsigav[nn][ii]=1;
436    
437 pam-fi 1.2 if((nn==1 && ii==11)||(nn==6 && ii==2)||(nn==6 && ii==4)||(nn==6 && ii==5)||(nn==6 && ii==6)||(nn==6 && ii==10)||(nn==11 && ii==3))
438     continue;
439     else{
440     if(flpedav[nn][ii]==1){
441     ok=1;
442     if(space[wc]<=3){
443     wc+=1;
444     flcanvas+=1;
445     }
446     c[wc]->cd();
447     tzz->SetTextFont(40);
448     tzz->SetTextSize(0.02);
449     tzz->SetTextAlign(13);
450     tzz->SetTextColor(50);
451     rep<<"********************************************************************************************************************************";
452     tzz->DrawLatex(2,space[wc],rep.str().c_str());
453     rep.str("");
454     space[wc]-=spacep;
455     rep<<"WARNING >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<"-->DSP "<<nn+1<<" -VA1 "<<2*ii+1<<"-"<<2*ii+2<<" --> <PED>= "<<pedav[nn][ii];
456     tzz->DrawLatex(2,space[wc],rep.str().c_str());
457     rep.str("");
458     space[wc]-=spacep;
459 pam-fi 1.1 }
460 pam-fi 1.2
461     if(flsigav[nn][ii]==1){
462     ok=1;
463     if(space[wc]<=3){
464     wc+=1;
465     flcanvas+=1;
466     }
467     c[wc]->cd();
468     tzz->SetTextFont(40);
469     tzz->SetTextSize(0.02);
470     tzz->SetTextAlign(13);
471     tzz->SetTextColor(50);
472     rep<<"********************************************************************************************************************************";
473     tzz->DrawLatex(2,space[wc],rep.str().c_str());
474     rep.str("");
475     space[wc]-=spacep;
476     rep<<"WARNING >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<"-->DSP "<<nn+1<<" -VA1 "<<2*ii+1<<"-"<<2*ii+2<<" --> <SIG>= "<<sigav[nn][ii];
477     tzz->DrawLatex(2,space[wc],rep.str().c_str());
478     rep.str("");
479     space[wc]-=spacep;
480 pam-fi 1.1 }
481     }
482     }
483     }
484     }
485 pam-fi 1.2 if(ok==0 && ctrk.good0[0]==1 && ctrk.good0[1]==1){
486 pam-fi 1.1 if(space[wc]<=10){
487     wc+=1;
488     flcanvas+=1;
489     }
490     c[wc]->cd();
491     tzz->SetTextFont(40);
492     tzz->SetTextSize(0.03);
493     tzz->SetTextAlign(13);
494     tzz->SetTextColor(1);
495     rep<<"*********************************** CALIBRATION pkt "<<i+1<<"********************************";
496     tzz->DrawLatex(2,space[wc],rep.str().c_str());
497     rep.str("");
498     space[wc]-=spacep+1;
499     rep<<">>>>>>> CalibTrk1 at OBT: "<<OBT[0]<<" ---------> OK ";
500     tzz->DrawLatex(2,space[wc],rep.str().c_str());
501     rep.str("");
502     space[wc]-=spacep+1;
503     rep<<">>>>>>> CalibTrk2 at OBT: "<<OBT[1]<<" ---------> OK ";
504     tzz->DrawLatex(2,space[wc],rep.str().c_str());
505     rep.str("");
506     space[wc]-=spacep+1;
507     rep<<"*******************************************************************************************";
508     tzz->DrawLatex(2,space[wc],rep.str().c_str());
509     rep.str("");
510     space[wc]-=spacep+1;
511     }
512     };//end loop on events
513    
514     printf("... end of packets. \n");
515     //
516     // Save output Files
517     stringstream nom1,nom2,nom3;
518    
519     for(Int_t i=0;i<cnum;i++)
520     c[i]->Update();//draw pads in canvas
521    
522     for(Int_t fl=0;fl<flcanvas;fl++){
523     if(flcanvas==1){
524     nom1.str("");
525     nom1<<ffile<<"_FTrkCalibQLook_BASIC."<<outfile.Data();
526     c[fl]->Print(out+nom1.str().c_str());
527     nom1.str("");
528     }
529    
530     if(flcanvas>=2){
531     if(!strcmp(outfile.Data(),"ps")||!strcmp(outfile.Data(),"pdf")){
532     nom1.str("");
533     nom2.str("");
534     nom3.str("");
535     nom1<<ffile<<"_FTrkCalibQLook_BASIC.ps(";
536     nom2<<ffile<<"_FTrkCalibQLook_BASIC.ps";
537     nom3<<ffile<<"_FTrkCalibQLook_BASIC.ps)";
538     if(fl==0) c[fl]->Print(out+nom1.str().c_str(),"portrait");
539     else if(fl==flcanvas-1) c[fl]->Print(out+nom3.str().c_str(),"portrait");
540     else c[fl]->Print(out+nom2.str().c_str(),"portrait");
541    
542     }
543     else{
544     nom1.str("");
545     nom1<<ffile<<"_FTrkCalibQLook_BASIC-pag"<<fl+1<<"."<<outfile.Data();
546     c[fl]->Print(out+nom1.str().c_str());
547     }
548     }
549     }
550    
551     if(!strcmp(outfile.Data(),"pdf")&&flcanvas>=2){
552     stringstream com;
553     com<<"ps2pdf13 "<<out<<ffile<<"_FTrkCalibQLook_BASIC.ps "<<out<<ffile<<"_FTrkCalibQLook_BASIC.pdf";
554     system(com.str().c_str());
555     printf("\n---> ps file converted in pdf format!\n");
556     com.str("");
557     com<<"rm -f "<<out<<ffile<<"_FTrkCalibQLook_BASIC.ps";
558     system(com.str().c_str());
559     printf("---> ps file removed!\n\n");
560     com.str("");
561     }
562    
563    
564     gROOT->Reset();
565     return;
566    
567     }

  ViewVC Help
Powered by ViewVC 1.1.23