/[PAMELA software]/DarthVader/OrbitalInfo/src/OrbitalInfoCore.cpp
ViewVC logotype

Contents of /DarthVader/OrbitalInfo/src/OrbitalInfoCore.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.3 - (show annotations) (download)
Fri Jun 30 09:21:58 2006 UTC (18 years, 8 months ago) by mocchiut
Branch: MAIN
CVS Tags: v1r01beta, v1r01
Changes since 1.2: +6 -3 lines
Memory leak bug fixed

1 //
2 // C/C++ headers
3 //
4 #include <fstream>
5 #include <string.h>
6 #include <iostream>
7 #include <cstring>
8 #include <stdio.h>
9 //
10 // ROOT headers
11 //
12 #include <TTree.h>
13 #include <TClassEdit.h>
14 #include <TObject.h>
15 #include <TList.h>
16 #include <TArrayL.h>
17 #include <TSystem.h>
18 #include <TSystemDirectory.h>
19 #include <TString.h>
20 #include <TFile.h>
21 #include <TClass.h>
22 #include <TSQLServer.h>
23 #include <TSQLRow.h>
24 #include <TSQLResult.h>
25 //
26 // YODA headers
27 //
28 #include <PamelaRun.h>
29 #include <RegistryEvent.h>
30 #include <PscuHeader.h>
31 #include <PscuEvent.h>
32 #include <EventHeader.h>
33 #include <RegistryEvent.h>
34 //
35 // RunInfo header
36 //
37 #include <RunInfo.h>
38 #include <GLTables.h>
39 //
40 // This program headers
41 //
42 #include <OrbitalInfo.h>
43 #include <OrbitalInfoCore.h>
44 #include <OrbitalInfoVerl2.h>
45
46 using namespace std;
47
48 //
49 // CORE ROUTINE
50 //
51 //
52
53 int OrbitalInfoCore(ULong64_t run, TFile *file, TSQLServer *dbc, Int_t OrbitalInfoargc, char *OrbitalInfoargv[]){
54 Int_t i = 0;
55 //
56 TString processFolder = "OrbitalInfoFolder";
57 //
58 // Set these to true to have a very verbose output.
59 //
60 Bool_t debug = false;
61 //
62 Bool_t verbose = false;
63
64 if ( OrbitalInfoargc > 0 ){
65 i = 0;
66 while ( i < OrbitalInfoargc ){
67 if ( !strcmp(OrbitalInfoargv[i],"-processFolder") ) {
68 if ( OrbitalInfoargc < i+1 ){
69 throw -3;
70 };
71 processFolder = (TString)OrbitalInfoargv[i+1];
72 i++;
73 };
74 if ( (!strcmp(OrbitalInfoargv[i],"--debug")) || (!strcmp(OrbitalInfoargv[i],"-g")) ) {
75 verbose = true;
76 };
77 if ( (!strcmp(OrbitalInfoargv[i],"--verbose")) || (!strcmp(OrbitalInfoargv[i],"-v")) ) {
78 verbose = true;
79 };
80 i++;
81 };
82 };
83 //
84 const char* outDir = gSystem->DirName(gSystem->DirName(file->GetPath()));
85 //
86 TTree *OrbitalInfotr = 0;
87 Long64_t nevents = 0LL;
88 //
89 // variables needed to reprocess data
90 //
91 TString OrbitalInfoversion;
92 ItoRunInfo *runinfo = 0;
93 TArrayL *runlist = 0;
94 TTree *OrbitalInfotrclone = 0;
95 Bool_t reproc = false;
96 Bool_t reprocall = false;
97 UInt_t nobefrun = 0;
98 UInt_t noaftrun = 0;
99 UInt_t numbofrun = 0;
100 stringstream ftmpname;
101 TString fname;
102 Long64_t totfileentries = 0ULL;
103 Long64_t idRun = 0LL;
104 //
105 // variables needed to handle error signals
106 //
107 Int_t code = 0;
108 Int_t sgnl;
109 //
110 // OrbitalInfo classes
111 //
112 OrbitalInfo *orbitalinfo = new OrbitalInfo();
113 OrbitalInfo *orbitalinfoclone = new OrbitalInfo();
114 //
115 // define variables for opening and reading level0 file
116 //
117 TFile *l0File = 0;
118 TTree *l0tr = 0;
119 TBranch *l0registry = 0;
120 pamela::RegistryEvent *l0reg=0;
121 // EM: open also header branch
122 TBranch *l0head = 0;
123 pamela::EventHeader *eh = 0;
124 pamela::PscuHeader *ph = 0;
125 // end EM
126 //
127 // Define other basic variables
128 //
129 UInt_t procev = 0;
130 stringstream file2;
131 stringstream file3;
132 stringstream qy;
133 Int_t totevent = 0;
134 ULong64_t atime = 0ULL;
135 Int_t ei = 0;
136 Int_t re = 0;
137 //
138 // Working filename
139 //
140 TString outputfile;
141 stringstream name;
142 name.str("");
143 name << outDir << "/";
144 //
145 // temporary file and folder
146 //
147 TFile *tempfile = 0;
148 TTree *tempOrbitalInfo = 0;
149 stringstream tempname;
150 stringstream OrbitalInfofolder;
151 tempname.str("");
152 tempname << outDir;
153 tempname << "/" << processFolder.Data();
154 OrbitalInfofolder.str("");
155 OrbitalInfofolder << tempname.str().c_str();
156 gSystem->MakeDirectory(OrbitalInfofolder.str().c_str());
157 tempname << "/OrbitalInfotree_run";
158 tempname << run << ".root";
159 //
160 // DB classes
161 //
162 GL_ROOT *glroot = new GL_ROOT();
163 //
164 // Let's start!
165 //
166 // As a first thing we must check what we have to do: if run = 0 we must process all events in the file has been passed
167 // if run != 0 we must process only that run but first we have to check if the tree MyDetector2 already exist in the file
168 // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
169 //
170 if ( run == 0ULL ) reproc = true;
171 //
172 //
173 // Output file is "outputfile"
174 //
175 if ( !file->IsOpen() ){
176 //printf(" OrbitalInfo - ERROR: cannot open file for writing\n");
177 throw -901;
178 };
179 //
180 // Retrieve GL_RUN variables from the level2 file
181 //
182 OrbitalInfoversion = OrbitalInfoInfo(false); // we should decide how to handle versioning system
183 //
184 // create an interface to RunInfo called "runinfo"
185 //
186 runinfo = new ItoRunInfo(file);
187 //
188 // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
189 //
190 sgnl = 0;
191 sgnl = runinfo->Update(run, "ORB", OrbitalInfoversion);
192 //sgnl = runinfo->Read(run);
193
194 if ( sgnl ){
195 //printf("OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
196 code = sgnl;
197 goto closeandexit;
198 } else {
199 sgnl = 0;
200 };
201 //
202 // number of events in the file BEFORE the first event of our run
203 //
204 nobefrun = runinfo->GetFirstEntry();
205 //
206 // total number of events in the file
207 //
208 totfileentries = runinfo->GetFileEntries();
209 //
210 // first file entry AFTER the last event of our run
211 //
212 noaftrun = runinfo->GetLastEntry() + 1;
213 //
214 // number of run to be processed
215 //
216 numbofrun = runinfo->GetNoRun();
217 //
218 // Try to access the OrbitalInfo tree in the file, if it exists we are reprocessing data if not we are processing a new run
219 //
220 OrbitalInfotrclone = (TTree*)file->Get("OrbitalInfo");
221 //
222 if ( !OrbitalInfotrclone ){
223 //
224 // tree does not exist, we are not reprocessing
225 //
226 reproc = false;
227 if ( run == 0ULL ){
228 if (verbose) printf(" OrbitalInfo - WARNING: you are reprocessing data but OrbitalInfo tree does not exist!\n");
229 }
230 if ( runinfo->IsReprocessing() && run != 0ULL ) {
231 if (verbose) printf(" OrbitalInfo - WARNING: it seems you are not reprocessing data but OrbitalInfo\n versioning information already exists in RunInfo.\n");
232 }
233 } else {
234 //
235 // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
236 //
237 reproc = true;
238 //
239 //
240 if (verbose) printf("\n Preparing the pre-processing...\n");
241 //
242 if ( run == 0ULL ){
243 //
244 // we are reprocessing all the file
245 // if we are reprocessing everything we don't need to copy any old event and we can just work with the new tree and delete the old one immediately
246 //
247 reprocall = true;
248 //
249 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing all runs\n");
250 //
251 } else {
252 //
253 // we are reprocessing a single run, we must copy to the new tree the events in the file which preceed the first event of the run
254 //
255 reprocall = false;
256 //
257 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing run number %llu \n",run);
258 //
259 // copying old tree to a new file
260 //
261 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
262 tempOrbitalInfo = OrbitalInfotrclone->CloneTree(-1,"fast");
263 tempOrbitalInfo->SetName("OrbitalInfo-old");
264 tempfile->Write();
265 tempfile->Close();
266 }
267 //
268 // Delete the old tree from old file and memory
269 //
270 OrbitalInfotrclone->Delete("all");
271 //
272 if (verbose) printf(" ...done!\n");
273 //
274 };
275 //
276 // create mydetector tree mydect
277 //
278 file->cd();
279 OrbitalInfotr = new TTree("OrbitalInfo-new","PAMELA OrbitalInfo data");
280 OrbitalInfotr->Branch("OrbitalInfo","OrbitalInfo",&orbitalinfo);
281 //
282 if ( reproc && !reprocall ){
283 //
284 // open new file and retrieve also tree informations
285 //
286 tempfile = new TFile(tempname.str().c_str(),"READ");
287 OrbitalInfotrclone = (TTree*)tempfile->Get("OrbitalInfo-old");
288 OrbitalInfotrclone->SetBranchAddress("OrbitalInfo",&orbitalinfoclone);
289 //
290 if ( nobefrun > 0 ){
291 if (verbose){
292 printf("\n Pre-processing: copying events from the old tree before the processed run\n");
293 printf(" Copying %u events in the file which are before the beginning of the run %llu \n",nobefrun,run);
294 printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
295 }
296 for (UInt_t j = 0; j < nobefrun; j++){
297 //
298 OrbitalInfotrclone->GetEntry(j);
299 //
300 // copy orbitalinfoclone to mydec
301 //
302 // orbitalinfo = new OrbitalInfo();
303 orbitalinfo->Clear();
304 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
305 //
306 // Fill entry in the new tree
307 //
308 OrbitalInfotr->Fill();
309 //
310 };
311 if (verbose) printf(" Finished successful copying!\n");
312 };
313 };
314 //
315 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
316 //
317 runlist = runinfo->GetRunList();
318 //
319 // Loop over the run to be processed
320 //
321 for (UInt_t irun=0; irun < numbofrun; irun++){
322 //
323 // retrieve the first run ID to be processed using the RunInfo list
324 //
325 idRun = runlist->At(irun);
326 if (verbose){
327 printf("\n\n\n ####################################################################### \n");
328 printf(" PROCESSING RUN NUMBER %i \n",(int)idRun);
329 printf(" ####################################################################### \n\n\n");
330 }
331 //
332 runinfo->ID_REG_RUN = 0ULL;
333 //
334 // store in the runinfo class the GL_RUN variables for our run
335 //
336 sgnl = 0;
337 sgnl = runinfo->GetRunInfo(idRun);
338 if ( sgnl ){
339 //printf("\n OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
340 code = sgnl;
341 goto closeandexit;
342 } else {
343 sgnl = 0;
344 };
345 //
346 // now you can access that variables using the RunInfo class this way runinfo->ID_REG_RUN
347 //
348 if ( runinfo->ID_REG_RUN == 0 ){
349 //printf("\n OrbitalInfo - ERROR: no run with ID_RUN = %i \n\n Exiting... \n\n",(int)idRun);
350 code = -5;
351 goto closeandexit;
352 };
353 //
354 // Search in the DB the path and name of the LEVEL0 file to be processed.
355 //
356 glroot->Query_GL_ROOT(runinfo->ID_REG_RUN,dbc);
357 //
358 ftmpname.str("");
359 ftmpname << glroot->PATH.Data() << "/";
360 ftmpname << glroot->NAME.Data();
361 fname = ftmpname.str().c_str();
362 //
363 // print out informations
364 //
365 totevent = runinfo->EV_REG_PHYS_TO - runinfo->EV_REG_PHYS_FROM + 1;
366 if (verbose){
367 printf("\n LEVEL0 data file: %s \n",fname.Data());
368 printf(" RUN HEADER absolute time is: %llu \n",runinfo->RUNHEADER_TIME);
369 printf(" RUN TRAILER absolute time is: %llu \n",runinfo->RUNTRAILER_TIME);
370 printf(" %i events to be processed for run %llu: from %i to %i (reg entries)\n\n",totevent,idRun,runinfo->EV_REG_PHYS_FROM,runinfo->EV_REG_PHYS_TO);
371 }//
372 // Open Level0 file
373 l0File = new TFile(fname.Data());
374 if ( !l0File ) {
375 //printf(" OrbitalInfo - ERROR: problems opening Level0 file\n");
376 code = -6;
377 goto closeandexit;
378 };
379 l0tr = (TTree*)l0File->Get("Physics");
380 if ( !l0tr ) {
381 //printf(" OrbitalInfo - ERROR: no Physics tree in Level0 file\n");
382 l0File->Close();
383 code = -7;
384 goto closeandexit;
385 };
386 // EM: open header branch as well
387 l0head = l0tr->GetBranch("Header");
388 if ( !l0head ) {
389 //if ( verbose ) printf(" CALORIMETER - ERROR: no Header branch in Level0 tree\n");
390 l0File->Close();
391 code = -8;
392 goto closeandexit;
393 };
394 l0tr->SetBranchAddress("Header", &eh);
395 // end EM
396 l0registry = l0tr->GetBranch("Registry");
397 if ( !l0registry ) {
398 //printf(" OrbitalInfo - ERROR: no Registry branch in Level0 tree\n");
399 l0File->Close();
400 code = -9;
401 goto closeandexit;
402 };
403 //
404 l0tr->SetBranchAddress("Registry", &l0reg);
405 //
406 nevents = l0registry->GetEntries();
407 //
408 if ( nevents < 1 ) {
409 //printf(" OrbitalInfo - ERROR: Level0 file is empty\n\n");
410 l0File->Close();
411 code = -11;
412 goto closeandexit;
413 };
414 //
415 if ( runinfo->EV_REG_PHYS_TO > nevents-1 ) {
416 //printf(" OrbitalInfo - ERROR: too few entries in the registry tree\n");
417 l0File->Close();
418 code = -12;
419 goto closeandexit;
420 };
421 //
422 // run over all the events of the run
423 //
424 if (verbose) printf("\n Ready to start! \n\n Processed events: \n\n");
425 //
426 for ( re = runinfo->EV_REG_PHYS_FROM; re <= runinfo->EV_REG_PHYS_TO; re++){
427 //
428 if ( procev%1000 == 0 && procev > 0 && verbose) printf(" %iK \n",procev/1000);
429 //
430 l0registry->GetEntry(re);
431 //
432 // absolute time of this event
433 //
434 atime = l0reg->absTime;
435 //
436 // physics events is at entry number ei where
437 //
438 ei = l0reg->event;
439 //
440 // paranoid check
441 //
442 if ( (atime > runinfo->RUNTRAILER_TIME) || (atime < runinfo->RUNHEADER_TIME) ) {
443 if (verbose) printf(" OrbitalInfo - WARNING: event at time outside the run time window, skipping it\n");
444 goto jumpev;
445 };
446 //
447 procev++;
448 //
449 // start processing
450 //
451 orbitalinfo->Clear();
452 //orbitalinfo = new OrbitalInfo();
453 orbitalinfo->absTime = l0reg->absTime;
454 // EM: add OBT and plt_num infos from the header
455 l0head->GetEntry(ei);
456 ph = eh->GetPscuHeader();
457 orbitalinfo->pkt_num = ph->GetCounter();
458 orbitalinfo->OBT = ph->GetOrbitalTime();
459 // end EM
460 OrbitalInfotr->Fill();
461 //
462 //
463 jumpev:
464 debug = false;
465 //
466 };
467 //
468 // Here you may want to clear some variables before processing another run
469 //
470 ei = 0;
471 }; // process all the runs
472 //
473 if (verbose) printf("\n Finished processing data \n");
474 //
475 closeandexit:
476 //
477 // we have finished processing the run(s). If we processed a single run now we must copy all the events after our run from the old tree to the new one and delete the old tree.
478 //
479 if ( !reprocall && reproc && code >= 0 ){
480 if ( totfileentries > noaftrun ){
481 if (verbose){
482 printf("\n Post-processing: copying events from the old tree after the processed run\n");
483 printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
484 printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
485 }
486 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
487 //
488 // Get entry from old tree
489 //
490 OrbitalInfotrclone->GetEntry(j);
491 //
492 // copy orbitalinfoclone to OrbitalInfo
493 //
494 orbitalinfo->Clear();
495 // orbitalinfo = new OrbitalInfo();
496 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
497 //
498 // Fill entry in the new tree
499 //
500 OrbitalInfotr->Fill();
501 };
502 if (verbose) printf(" Finished successful copying!\n");
503 };
504 };
505 //
506 // Close files, delete old tree(s), write and close level2 file
507 //
508 if ( l0File ) l0File->Close();
509 if ( tempfile ) tempfile->Close();
510 gSystem->Unlink(tempname.str().c_str());
511
512 //if ( code < 0 ) printf("\n OrbitalInfo - ERROR: an error occurred, try to save anyway...\n");
513 //printf("\n Writing and closing rootple\n");
514 if ( runinfo ) runinfo->Close();
515 if ( OrbitalInfotr ) OrbitalInfotr->SetName("OrbitalInfo");
516 if ( file ){
517 file->cd();
518 file->Write();
519 };
520 //
521 gSystem->Unlink(OrbitalInfofolder.str().c_str());
522 //
523 // the end
524 //
525 if (verbose) printf("\n Exiting...\n");
526 if(OrbitalInfotr)OrbitalInfotr->Delete();
527 if(code < 0) throw code;
528 return(code);
529 }
530

  ViewVC Help
Powered by ViewVC 1.1.23