/[PAMELA software]/quicklook/QLflightTmtc_Header/HeaderScan.cpp
ViewVC logotype

Contents of /quicklook/QLflightTmtc_Header/HeaderScan.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.9 - (show annotations) (download)
Fri Sep 8 08:26:47 2006 UTC (18 years, 2 months ago) by pam-rm2
Branch: MAIN
CVS Tags: v3r00
Changes since 1.8: +136 -153 lines
risolto il crash dell'istogramma e della pad di testo; aggiornato il build.xml

1 /**
2 * Header_Scan
3 * Author Nagni
4 * version 1.0
5 *
6 * Version 1.1 - 28 December 2004
7 * If outList does not exist the function exit to prompt
8 * If nevents = 0 the function exit to promt
9 *
10 * Version 1.2 - 3 January 2005
11 * Two canvases are created to see the graphs better
12 *
13 * Version 1.25 - 13 January 2005
14 * Changed Int_t to Float because of variable range size
15 * (UInt_t has been excluded beacuse of uncompatibility with TGraph)
16 *
17 * version 1.3 - 22 February 2005 - Nagni
18 * For compatibility with batch mode excution:
19 * 1) Added "include <iostream>" and "using namespace std"
20 * 2) Removed gROOT->Reset()
21 *
22 * Version 1.4
23 * Date 08 March 2005
24 * Added "format" parameter for saving the produced image in various formats
25 * (for a complete list of types refer to TPad::SaveAs method)
26 *
27 * Version 1.5
28 * Date 09 February 2006 - Marcelli
29 * Update to work with new Yoda output
30 *
31 * Version 1.6
32 * Date 27 February 2006 - Marcelli
33 * For compilation:
34 * Added function "int main(int argc, char* argv[])"
35 *
36 *
37 * Description: This script creates two canvases with five pads. The first pad shows packetID variable (for all packets) vs. OBT.
38 * The second pad shows the number of physic packets vs. OBT. The third pad shows the lenght of Physic packets (byte) vs. OBT.
39 * The fourth pad shows the packetcounter of physic packets vs. OBT. The fifth pad shows PacketCounter vs. File Offset.
40 *
41 * Parameters:
42 * TSTring base - the path to the root directory for the specific Pamela unpack session
43 * There is no default value, without this input the program will not run
44 * TString outDir - the path where to save the output image (Default = ./)
45 * TString format - the format which will be used fo rsave the produced images (Default = "jpg")
46 *
47 *
48 * Version 1.7
49 * Date 16 June 2006 - Malvezzi
50 *
51 * Description of changes:
52 * Implementation of the case: numebr of events <= 0.
53 * Remove graph "grPcktId1"; see PacketScan for the same information.
54 * Fixed bugs: for a large namber of events is not possible to have vectors, so I have subsituted graphs with histograms
55 * or divided the graphs in more than one canvas.
56 *
57 * Version 1.8
58 * Date 8 August 2006 - Malvezzi
59 *
60 * Description: changed the scale in the second and third graph of the first canvas; added a pad of text in the second canvas
61 *
62 */
63
64
65 #include <fstream>
66 #include <math.h>
67 #include "TLatex.h"
68 #include "TF1.h"
69 #include "TPaveText.h"
70 #include "TMultiGraph.h"
71 #include <sstream>
72 #include <iostream>
73 #include "TString.h"
74 #include "TStyle.h"
75 #include "TFile.h"
76 #include "TList.h"
77 #include "TTree.h"
78 #include "TObjString.h"
79 #include "TCanvas.h"
80 #include "TGraph.h"
81 #include "TH1F.h"
82 #include "TGaxis.h"
83 #include "EventHeader.h"
84 #include "PscuHeader.h"
85 #include "RunHeaderEvent.h"
86 #include "TPaveText.h"
87
88 using namespace std;
89
90 void HeaderScan(TString base, TString outDir, TString format){
91
92 //------------------- Variables initilization -------------------------//
93 Long64_t nevents=0, runnevents;
94 ULong_t lastime, firstime, primotempo, ultimotempo, primoffset=500000000, ultimoffset;
95 double obmin=0.;
96 double obmax=0.;
97 stringstream oss, oss1, oss2, oss3, noentries, stringa;
98
99 //-------------- Load root file, tree and branch -------------------//
100 TFile *file = new TFile(base.Data());
101 if (!file){
102 printf("No such file in the directory has been found");
103 return;
104 }
105
106 TTree *PhysicsTr = (TTree*)file->Get("Physics");
107 TBranch *headBr = PhysicsTr->GetBranch("Header");
108
109 pamela::EventHeader *eh = 0;
110 pamela::PscuHeader *ph = 0;
111
112 PhysicsTr->SetBranchAddress("Header", &eh);
113 nevents = PhysicsTr->GetEntries();
114 const Int_t size = nevents;
115
116 TTree *RunHeadTr = (TTree*)file->Get("RunHeader"); ///run header tree
117 pamela::EventHeader *eH=0;
118 pamela::RunHeaderEvent *reh=0;
119
120 RunHeadTr->SetBranchAddress("Header",&eH);
121 RunHeadTr->SetBranchAddress("RunHeader",&reh);
122 runnevents = RunHeadTr->GetEntries();
123
124 TString filename = ((TObjString*)base.Tokenize('/')->Last())->GetString();
125 filename = ((TObjString*)filename.Tokenize('.')->First())->GetString();
126
127 //----------- If nevents < = 0 ---------------------------------/
128 if (nevents<=0) {
129 printf("nevents = %i \n", nevents);
130 printf(" \n");
131 TCanvas *canv = new TCanvas("No entries", "No entries ", 400, 200);
132 canv->SetFillColor(10);
133 canv->cd();
134 TLatex *l = new TLatex();
135 l->SetTextAlign(12);
136 l->SetTextSize(0.15);
137 l->SetTextColor(2);
138 noentries.str("");
139 noentries<< "HeaderScan_QL:";
140 l->DrawLatex(0.05, 0.7, noentries.str().c_str());
141 noentries.str("");
142 noentries<< "No Physics entries for this files";
143 l->DrawLatex(0.05, 0.5, noentries.str().c_str());
144 if (outDir == "./") {
145 oss.str("");
146 oss << filename.Data() << "_HeaderScan_QL." << format.Data();
147 } else {
148 oss.str("");
149 oss << outDir.Data() << filename.Data() << "_HeaderScan_QL." << format.Data();
150 }
151 canv->Update();
152 canv->SaveAs(oss.str().c_str());
153 return;
154 }
155 //-------------- to know the max and min OBT ----------------------------//
156 headBr->GetEntry(0);
157 ph = eh->GetPscuHeader();
158 firstime = ph->GetOrbitalTime();
159 headBr->GetEntry(nevents);
160 ph = eh->GetPscuHeader();
161 lastime = ph->GetOrbitalTime();
162 int i =0;
163 while(i < nevents){
164 headBr->GetEntry(i);
165 ph = eh->GetPscuHeader();
166 if((ph->GetOrbitalTime()) <= firstime) firstime=ph->GetOrbitalTime();
167 if((ph->GetOrbitalTime()) >= lastime) lastime=ph->GetOrbitalTime();
168 i++;
169 }
170
171 //------------------------ First histogram -----------------------------------//
172 obmin=firstime;
173 obmax=lastime;
174 Int_t nbin = (lastime-firstime)/60000;
175 TH1F *h1 = new TH1F ("histo1", "" , nbin, obmin, obmax);
176
177 //------------------- fill vectors and histogram -----------------------------//
178 Double_t *PscuCounter = new Double_t[size];
179 Double_t *PacketLenght = new Double_t[size];
180 Double_t *OBTime = new Double_t[size];
181 Double_t *Eventsperminute= new Double_t[nbin];
182 Double_t *Minute= new Double_t[nbin];
183 Int_t max=0;
184 for (Int_t k = 0; k < nevents; k++){
185 headBr->GetEntry(k);
186 ph = eh->GetPscuHeader();
187 h1->Fill(ph->GetOrbitalTime());
188 PscuCounter[k]= ph->GetCounter();
189 PacketLenght[k]=ph->GetPacketLenght();
190 OBTime[k]=ph->GetOrbitalTime();
191 }
192 int l=0;
193 while(l<nbin){
194 Eventsperminute[l]=h1->GetBinContent(l);
195 Minute[l]=firstime+l*60000;
196 if(h1->GetBinContent(l) >= max)max =(Int_t)h1->GetBinContent(l);
197 l++;
198 }
199
200 //----------- Graph and MultiGraph -----------------------------------------------//
201 TMultiGraph *rate = new TMultiGraph();
202 TMultiGraph *packetLength = new TMultiGraph();
203 TMultiGraph *packeCounter = new TMultiGraph();
204
205 oss1.str("");
206 oss1 << filename.Data() <<": Physics Packet per minute. Start time = " << obmin << ", End time = "<< obmax <<" ms";
207 TGraph *rate1= new TGraph(nbin, (const Double_t*)Minute, (const Double_t*)Eventsperminute);
208 rate1->SetMarkerColor(kBlack);
209 rate1->SetMarkerSize(.1);
210 rate1->SetMarkerStyle(21);
211 rate->Add(rate1);
212
213 TGraph *packetLength1= new TGraph(nevents, (const Double_t*)OBTime, (const Double_t*)PacketLenght);
214 oss2.str("");
215 oss2 << filename.Data() <<": Lenght of Physic packets";
216 packetLength1->SetMarkerColor(2);
217 packetLength1->SetMarkerSize(.5);
218 packetLength1->SetMarkerStyle(21);
219 packetLength->Add(packetLength1);
220
221 TGraph *packeCounter1= new TGraph(nevents, (const Double_t*)OBTime, (const Double_t*)PscuCounter);
222 oss3.str("");
223 oss3 << filename.Data() <<": Physics Counter vs. OBT";
224 packeCounter1->SetMarkerColor(4);
225 packeCounter1->SetMarkerSize(.2);
226 packeCounter1->SetMarkerStyle(21);
227 packeCounter->Add(packeCounter1);
228
229 //------------ Create and Draw Canvas ---------------------//
230 TCanvas *finalCanv = new TCanvas("Header", base, 1200, 1600);
231 finalCanv->Divide(1,5);
232 finalCanv->SetFillColor(10);
233
234 TPad *all2= new TPad ("","", 0, 0, 1, 1);
235 all2->SetFillColor(10);
236 TPad *all3= new TPad ("","", 0, 0, 1, 1);
237 all3->SetFillColor(10);
238 TPad *all4= new TPad ("","", 0, 0, 1, 1);
239 all4->SetFillColor(10);
240 TPad *all= new TPad ("","", 0, 0, 1, 1);
241 all->SetFillColor(10);
242 TPad *all1= new TPad ("","", 0, 0, 1, 1);
243 all1->SetFillColor(10);
244 TPad *pad = new TPad("pad","pad", .80,.45,.90,.75);
245 pad->SetFillColor(10);
246
247 TLine li;
248 li.SetLineStyle(4);
249 li.SetLineWidth(2);
250
251 //----------- First PAD -------------------------------//
252 finalCanv->cd(1);
253 all2->Draw();
254 all2->cd();
255 rate->Draw("ALP");
256 rate->GetXaxis()->SetTitle("OBT (ms)");
257 rate->GetXaxis()->SetTitleSize(0.05);
258 rate->GetXaxis()->CenterTitle();
259 rate->GetXaxis()->SetLabelSize(0.05);
260 rate->GetYaxis()->SetTitle("Number of events ");
261 rate->GetYaxis()->CenterTitle();
262 rate->GetYaxis()->SetLabelSize(0.05);
263 rate->GetYaxis()->SetTitleSize(0.06);
264 rate->GetYaxis()->SetTitleOffset(0.6);
265 for (Int_t l = 0; l < runnevents; l++){
266 RunHeadTr->GetEntry(l);
267 ph = eH->GetPscuHeader();
268 int ws= reh->RM_ACQ_SETTING_MODE;
269 int id = ph->GetPacketId1();
270 Int_t obt = ph->GetOrbitalTime();
271 if (ws==1){
272 li.SetLineColor(3);
273 li.DrawLine(obt,0,obt,max);
274 }
275 else if (ws==2){
276 li.SetLineColor(4);
277 li.DrawLine(obt,0,obt,max);
278 }
279 }
280
281 finalCanv->cd(1);
282 stringstream ws1, ws2;
283 ws1.str("");
284 ws2.str("");
285 ws1<<"ACQ_SETTING_MODE = 1";
286 ws2<<"ACQ_SETTING_MODE = 2";
287 TPaveText *pt=0;
288 pt = new TPaveText (.60,.92,.76,.98);
289 pt->AddText(ws1.str().c_str());
290 pt->SetTextColor(3);
291 pt->SetFillColor(10);
292 pt->SetBorderSize(0);
293 pt->Draw();
294 TPaveText *pt1=0;
295 pt1 = new TPaveText (.76,.92,.92,.98);
296 pt1->AddText(ws2.str().c_str());
297 pt1->SetTextColor(4);
298 pt1->SetFillColor(10);
299 pt1->SetBorderSize(0);
300 pt1->Draw();
301 pt1 = new TPaveText (.05,.91,.6,1);
302 pt1->AddText(oss1.str().c_str());
303 pt1->SetTextColor(1);
304 pt1->SetFillColor(10);
305 pt1->SetBorderSize(0);
306 pt1->Draw();
307
308 //----------- Second PAD -------------------------------//
309 finalCanv->cd(2);
310 all3->Draw();
311 all3->cd();
312 packetLength->Draw("AP");
313 packetLength->GetXaxis()->SetTitle("OBT (ms)");
314 packetLength->GetXaxis()->SetTitleSize(0.05);
315 packetLength->GetXaxis()->CenterTitle();
316 packetLength->GetXaxis()->SetLabelSize(0.05);
317 packetLength->GetYaxis()->SetTitle("Lenght (byte)");
318 packetLength->GetYaxis()->CenterTitle();
319 packetLength->GetYaxis()->SetLabelSize(0.05);
320 packetLength->GetYaxis()->SetTitleSize(0.06);
321 packetLength->GetYaxis()->SetTitleOffset(0.6);
322
323
324 finalCanv->cd(2);
325 pt = new TPaveText (.6,.91,.90,1);
326 pt->AddText(oss2.str().c_str());
327 pt->SetTextColor(2);
328 pt->SetFillColor(10);
329 pt->SetBorderSize(0);
330 pt->Draw();
331
332 //----------- Third PAD -------------------------------//
333 finalCanv->cd(3);
334 all4->Draw();
335 all4->cd();
336 packeCounter->Draw("AP");
337 packeCounter->GetXaxis()->SetTitle("OBT (ms)");
338 packeCounter->GetXaxis()->SetTitleSize(0.05);
339 packeCounter->GetXaxis()->CenterTitle();
340 packeCounter->GetXaxis()->SetLabelSize(0.05);
341 packeCounter->GetYaxis()->SetTitle("Counter");
342 packeCounter->GetYaxis()->SetTitleSize(0.05);
343 packeCounter->GetYaxis()->CenterTitle();
344 packeCounter->GetYaxis()->SetLabelSize(0.05);
345 packeCounter->GetYaxis()->SetTitleSize(0.06);
346 packeCounter->GetYaxis()->SetTitleOffset(0.6);
347
348
349 finalCanv->cd(3);
350 TPaveText *pt2=0;
351 pt2 = new TPaveText (.6,.91,.90,1);
352 pt2->AddText(oss3.str().c_str());
353 pt2->SetTextColor(4);
354 pt2->SetFillColor(10);
355 pt2->SetBorderSize(0);
356 pt2->Draw();
357
358 /**********************************************************************************************/
359
360 TMultiGraph *mg1 = new TMultiGraph();
361 TMultiGraph *mg2 = new TMultiGraph();
362 //---------------------- fill vectors and histogram --------------------------------------------------//
363 TList *list = new TList;
364 Int_t numkey;
365 TObject *key = new TObject;
366 const char *name;
367 TTree* tr = new TTree;
368 Long64_t nevntskey=0;
369 list = file->GetListOfKeys();
370 numkey = file->GetNkeys();
371 Double_t salto;
372 for (Int_t m=0; m<numkey; m++){
373 key = list->At(m);
374 name=(char *)(key->GetName());
375 tr = (TTree*)file->Get(name);
376 if (tr->IsZombie()) continue;
377 tr->SetBranchAddress("Header", &eh);
378 TBranch *Br = tr->GetBranch("Header");
379 nevntskey = tr->GetEntries();
380
381 if(nevntskey !=0){
382 Int_t size1=nevntskey;
383 Double_t *PscuCounter1 = new Double_t[size1];
384 Double_t *FileOffset1 = new Double_t[size1];
385 Double_t *tempo1 = new Double_t[size1];
386
387 int n=0;
388 while(n<nevntskey){
389 Br->GetEntry(n);
390 ph = eh->GetPscuHeader();
391 PscuCounter1[n]= ph->GetCounter();
392 FileOffset1[n]=ph->GetFileOffset();
393 tempo1[n]=ph->GetOrbitalTime();
394 if((m==0) && (n==0)){
395 primotempo=ph->GetOrbitalTime();
396 salto=ph->GetOrbitalTime();
397 }
398 if(salto > ph->GetOrbitalTime())salto=ph->GetOrbitalTime();
399 if(ph->GetFileOffset()<= primoffset){
400 primoffset=ph->GetFileOffset();
401 primotempo=ph->GetOrbitalTime();
402 }
403 if(ph->GetFileOffset()>=ultimoffset){
404 ultimotempo=ph->GetOrbitalTime();
405 ultimoffset=ph->GetFileOffset();
406 }
407 n++;
408 }
409 TGraph *graph3= new TGraph(nevntskey, (const Double_t*)FileOffset1, (const Double_t*)PscuCounter1);
410 graph3->SetMarkerColor(3);
411 graph3->SetMarkerSize(.2);
412 graph3->SetMarkerStyle(21);
413 mg1->Add(graph3);
414
415 TGraph *graph4= new TGraph(nevntskey, (const Double_t*)FileOffset1, (const Double_t*)tempo1);
416 graph4->SetMarkerColor(kBlue);
417 graph4->SetMarkerSize(.2);
418 graph4->SetMarkerStyle(21);
419 mg2->Add(graph4);
420 }
421 }
422
423 TLatex *lat = new TLatex();
424 lat->SetTextAlign(12);
425 lat->SetTextSize(0.15);
426 lat->SetTextColor(kBlue);
427
428 //------------ Fourth PAD ---------------------//
429 finalCanv->cd(4);
430 all1->Draw();
431 all1->cd();
432
433 oss1.str("");
434 oss1 << filename.Data() <<": PscuCounter vs FileOffset.";
435 mg1->Draw("AP");
436 mg1->GetXaxis()->SetTitle("File Offset");
437 mg1->GetXaxis()->CenterTitle();
438 mg1->GetXaxis()->SetTitleOffset(0.8);
439 mg1->GetXaxis()->SetTitleSize(0.05);
440 mg1->GetXaxis()->SetLabelSize(0.05);
441 mg1->GetYaxis()->SetTitle("Counter");
442 mg1->GetYaxis()->CenterTitle();
443 mg1->GetYaxis()->SetTitleSize(0.06);
444 mg1->GetYaxis()->SetLabelSize(0.06);
445 mg1->GetYaxis()->SetTitleOffset(0.6);
446 finalCanv->cd(4);
447 TPaveText *pt3=0;
448 pt3 = new TPaveText (.60,.91,.90,1);
449 pt3->AddText(oss1.str().c_str());
450 pt3->SetTextColor(3);
451 pt3->SetFillColor(10);
452 pt3->SetBorderSize(0);
453 pt3->Draw();
454
455 //------------ Fifth PAD ---------------------//
456 finalCanv->cd(5);
457 all->Draw();
458 all->cd();
459
460 oss3.str("");
461 oss3 << filename.Data() <<" OBT vs FileOffset. First packet at "<<primotempo <<" ms, last packet at "<<ultimotempo<<" ms.";
462 mg2->Draw("AP");
463 mg2->GetXaxis()->SetTitle("File Offset");
464 mg2->GetXaxis()->CenterTitle();
465 mg2->GetXaxis()->SetTitleSize(0.05);
466 mg2->GetXaxis()->SetLabelSize(0.05);
467 mg2->GetYaxis()->SetTitle("OBT");
468 mg2->GetYaxis()->CenterTitle();
469 mg2->GetYaxis()->SetTitleSize(0.06);
470 mg2->GetYaxis()->SetLabelSize(0.05);
471 mg2->GetYaxis()->SetTitleOffset(0.6);
472 pad->Draw();
473 pad->cd();
474 stringa.str("");
475 stringa << "jump at: "<<salto<<" ms \n";
476 if(salto != primotempo) lat->DrawLatex(0.08, 0.8,stringa.str().c_str());
477
478 finalCanv->cd(5);
479 TPaveText *pt4=0;
480 pt4 = new TPaveText (.38 ,.91,.92,1);
481 pt4->AddText(oss3.str().c_str());
482 pt4->SetTextColor(kBlue);
483 pt4->SetFillColor(10);
484 pt4->SetBorderSize(0);
485 pt4->Draw();
486
487 finalCanv->Update();
488
489 oss1.str("");
490 oss1 << outDir.Data() << filename.Data();
491 oss1 << "_HeaderScan"<<"." << format.Data();
492
493 finalCanv->SaveAs(oss1.str().c_str());
494
495 file->Close();
496
497 }
498
499 int main(int argc, char* argv[]){
500 TString path;
501 TString outDir = "./";
502 TString format = "jpg";
503 if (argc < 2){
504 printf("You have to insert at least the file to analyze \n");
505 printf("Try '--help' for more information. \n");
506 exit(1);
507 }
508 if (!strcmp(argv[1], "--help")){
509 printf( "Usage: HeaderScan FILE [OPTION] \n");
510 printf( "\t --help Print this help and exit \n");
511 printf( "\t -outDir[path] Path where to put the output [default ./] \n");
512 printf( "\t -format[jpg|ps|gif] Format for output files [default 'jpg'] \n");
513 exit(1);
514 }
515 path=argv[1];
516 for (int i = 2; i < argc; i++){
517 if (!strcmp(argv[i], "-outDir")){
518 if (++i >= argc){
519 printf( "-outDir needs arguments. \n");
520 printf( "Try '--help' for more information. \n");
521 exit(1);
522 }
523 else{
524 outDir = argv[i];
525 continue;
526 }
527 }
528 if (!strcmp(argv[i], "-format")){
529 if (++i >= argc){
530 printf( "-format needs arguments. \n");
531 printf( "Try '--help' for more information. \n");
532 exit(1);
533 }
534 else{
535 format = argv[i];
536 continue;
537 }
538 }
539 }
540 HeaderScan(argv[1], outDir, format);
541 }
542

  ViewVC Help
Powered by ViewVC 1.1.23