/[PAMELA software]/DarthVader/NDLevel2/src/NDCore.cpp
ViewVC logotype

Contents of /DarthVader/NDLevel2/src/NDCore.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.10 - (show annotations) (download)
Tue Jan 23 15:36:43 2007 UTC (17 years, 10 months ago) by mocchiut
Branch: MAIN
CVS Tags: v3r00
Changes since 1.9: +7 -2 lines
Compiling bug under fedora4 fixed, Syntax error inc/TrkLevelX 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 <TArrayI.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 #include <TClonesArray.h>
26 #include <stdlib.h>
27 #include <math.h>
28 //
29 // RunInfo header
30 //
31 #include <RunInfo.h>
32 #include <GLTables.h>
33 //
34 // YODA headers
35 //
36 #include <PamelaRun.h>
37 #include <PscuHeader.h>
38 #include <PscuEvent.h>
39 #include <EventHeader.h>
40 #include <physics/neutronDetector/NeutronEvent.h>
41 #include <physics/neutronDetector/NeutronRecord.h>
42 //#include <yodaUtility.h>
43 //
44 // RunInfo header
45 //
46 //#include <RunInfo.h>
47 //#include <GLTables.h>
48 //
49 // This program headers
50 //
51 #include <NDLevel2.h>
52 #include <NDCore.h>
53 #include <NDVerl2.h>
54
55 using namespace std;
56
57 //
58 // CORE ROUTINE
59 //
60 //
61 int NDCore(UInt_t run, TFile *file, TSQLServer *dbc, Int_t NDargc, char *NDargv[]){
62 //
63 Int_t i = 0;
64 //
65 TString processFolder = "NDFolder";
66 //
67 // Set these to true to have a very verbose output.
68 //
69 Bool_t debug = false;
70 //
71 Bool_t verbose = false;
72 //
73 if ( NDargc > 0 ){
74 i = 0;
75 while ( i < NDargc ){
76 if ( !strcmp(NDargv[i],"-processFolder") ) {
77 if ( NDargc < i+1 ){
78 throw -3;
79 };
80 processFolder = (TString)NDargv[i+1];
81 i++;
82 };
83 if ( (!strcmp(NDargv[i],"--verbose")) || (!strcmp(NDargv[i],"-v")) ) {
84 verbose = true;
85 };
86 if ( (!strcmp(NDargv[i],"--debug")) || (!strcmp(NDargv[i],"-g")) ) {
87 verbose= true;
88 };
89 i++;
90 };
91 };
92 //
93 const char* outDir = gSystem->DirName(gSystem->DirName(file->GetPath()));
94 // Variables for level2
95 //
96 TTree *NDtr = 0;
97 UInt_t nevents = 0;
98 //
99 // variables needed to reprocess data
100 //
101 TString NDversion;
102 ItoRunInfo *runinfo = 0;
103 TArrayI *runlist = 0;
104 TTree *NDtrclone = 0;
105 Bool_t reproc = false;
106 Bool_t reprocall = false;
107 UInt_t nobefrun = 0;
108 UInt_t noaftrun = 0;
109 UInt_t numbofrun = 0;
110 stringstream ftmpname;
111 TString fname;
112 UInt_t totfileentries = 0;
113 UInt_t idRun = 0;
114 Int_t tmpSize =0;
115 Float_t yTrig =0;
116 Float_t yUpperBackground =0;
117 Float_t yBottomBackground =0;
118 //
119 // variables needed to handle error signals
120 //
121 Int_t code = 0;
122 Int_t sgnl;
123 //
124 Long64_t maxsize = 10000000000LL;
125 TTree::SetMaxTreeSize(maxsize);
126 //
127 //
128 // ND level2 classes
129 //
130 NDLevel2 *nd = new NDLevel2();
131 NDLevel2 *ndclone = new NDLevel2();
132 //
133 // define variables for opening and reading level0 file
134 //
135 TFile *l0File = 0;
136 TTree *l0tr = 0;
137 TBranch *l0head = 0;
138 TBranch *l0ND =0;
139 pamela::EventHeader *eh = 0;
140 pamela::PscuHeader *ph = 0;
141 pamela::neutron::NeutronRecord *l0nr = 0;
142 pamela::neutron::NeutronEvent *l0ne = 0;
143 //
144 // Define other basic variables
145 //
146 UInt_t procev = 0;
147 stringstream file2;
148 stringstream file3;
149 stringstream qy;
150 Int_t totevent = 0;
151 UInt_t atime = 0;
152 UInt_t re = 0;
153 //
154 // Working filename
155 //
156 TString outputfile;
157 stringstream name;
158 name.str("");
159 name << outDir << "/";
160 //
161 // temporary file and folder
162 //
163 TFile *tempfile = 0;
164 TTree *tempND = 0;
165 stringstream tempname;
166 stringstream NDfolder;
167 tempname.str("");
168 tempname << outDir;
169 tempname << "/" << processFolder.Data();
170 NDfolder.str("");
171 NDfolder << tempname.str().c_str();
172 gSystem->MakeDirectory(NDfolder.str().c_str());
173 tempname << "/NDtree_run";
174 tempname << run << ".root";
175 //
176 // DB classes
177 //
178 GL_ROOT *glroot = new GL_ROOT();
179 GL_TIMESYNC *dbtime = 0;
180 //
181 // Let's start!
182 //
183 // 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
184 // if run != 0 we must process only that run but first we have to check if the tree MyDetector2 already exist in the file
185 // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
186 //
187 if ( run == 0 ) reproc = true;
188 //
189 //
190 // Output file is "outputfile"
191 //
192 if ( !file->IsOpen() ){
193 //printf(" ND - ERROR: cannot open file for writing\n");
194 throw -601 ;
195 };
196 //
197 // Retrieve GL_RUN variables from the level2 file
198 //
199 NDversion = NDInfo(false); // we should decide how to handle versioning system
200 //
201 // create an interface to RunInfo called "runinfo"
202 //
203 runinfo = new ItoRunInfo(file);
204 //
205 // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
206 //
207 sgnl = 0;
208 sgnl = runinfo->Update(run, "ND",NDversion);
209 if ( sgnl ){
210 if ( debug ) printf(" ND - ERROR: RunInfo exited with non-zero status\n");
211 code = sgnl;
212 goto closeandexit;
213 } else {
214 sgnl = 0;
215 };
216 //
217 // number of events in the file BEFORE the first event of our run
218 //
219 nobefrun = runinfo->GetFirstEntry();
220 //
221 // total number of events in the file
222 //
223 totfileentries = runinfo->GetFileEntries();
224 //
225 // first file entry AFTER the last event of our run
226 //
227 noaftrun = runinfo->GetLastEntry() + 1;
228 //
229 // number of run to be processed
230 //
231 numbofrun = runinfo->GetNoRun();
232 //
233 // Try to access the S4 tree in the file, if it exists we are reprocessing data if not we are processing a new run
234 //
235 NDtrclone = (TTree*)file->Get("NeutronD");
236 //
237 if ( !NDtrclone ){
238 //
239 // tree does not exist, we are not reprocessing
240 //
241 reproc = false;
242 if ( run == 0 ) {
243 if (verbose) printf(" ND - WARNING: you are reprocessing data but ND tree does not exist!\n");
244 }
245 if ( runinfo->IsReprocessing() && run != 0 ) {
246 if (verbose) printf(" ND - WARNING: it seems you are not reprocessing data but ND\n versioning information already exists in RunInfo.\n");
247 }
248 } else {
249 //
250 NDtrclone->SetAutoSave(900000000000000LL);
251 //
252 // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
253 //
254 reproc = true;
255 //
256 // update versioning information
257 //
258 if (verbose) printf("\n Preparing the pre-processing...\n");
259 //
260 if ( run == 0 ){
261 //
262 // we are reprocessing all the file
263 // 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
264 //
265 reprocall = true;
266 //
267 if (verbose) printf("\n ND - WARNING: Reprocessing all runs\n");
268 //
269 } else {
270 //
271 // 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
272 //
273 reprocall = false;
274 //
275 if (verbose) printf("\n ND - WARNING: Reprocessing run number %u \n",run);
276 //
277 // copying old tree to a new file
278 //
279 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
280 tempND = NDtrclone->CloneTree(-1,"fast");
281 tempND->SetName("NeutronD-old");
282 tempfile->Write();
283 tempfile->Close();
284 }
285 //
286 // Delete the old tree from old file and memory
287 //
288 NDtrclone->Delete("all");
289 //
290 if (verbose) printf(" ...done!\n");
291 //
292 };
293 //
294 // create mydetector tree mydect
295 //
296 file->cd();
297 NDtr = new TTree("NeutronD-new","PAMELA Level2 NeutronD data");
298 NDtr->SetAutoSave(900000000000000LL);
299 NDtr->Branch("NDLevel2","NDLevel2",&nd);
300 //
301 if ( reproc && !reprocall ){
302 //
303 // open new file and retrieve also tree informations
304 //
305 tempfile = new TFile(tempname.str().c_str(),"READ");
306 NDtrclone = (TTree*)tempfile->Get("NeutronD-old");
307 NDtrclone->SetAutoSave(900000000000000LL);
308 NDtrclone->SetBranchAddress("NDLevel2",&ndclone);
309 //
310 if ( nobefrun > 0 ){
311 if (verbose) {
312 printf("\n Pre-processing: copying events from the old tree before the processed run\n");
313 printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
314 printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
315 }
316 for (UInt_t j = 0; j < nobefrun; j++){
317 //
318 NDtrclone->GetEntry(j);
319 //
320 // copy ndclone to mydec
321 //
322 // nd = new NDLevel2();
323 nd->Clear();
324 memcpy(&nd,&ndclone,sizeof(ndclone));
325 //
326 // Fill entry in the new tree
327 //
328 NDtr->Fill();
329 //
330 };
331 if (verbose) printf(" Finished successful copying!\n");
332 };
333 };
334 //
335 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
336 //
337 runlist = runinfo->GetRunList();
338 //
339 // Loop over the run to be processed
340 //
341 for (UInt_t irun=0; irun < numbofrun; irun++){
342 //
343 // retrieve the first run ID to be processed using the RunInfo list
344 //
345 idRun = runlist->At(irun);
346 if (verbose) {
347 printf("\n\n\n ####################################################################### \n");
348 printf(" PROCESSING RUN NUMBER %u \n",idRun);
349 printf(" ####################################################################### \n\n\n");
350 }
351 //
352 runinfo->ID_ROOT_L0 = 0;
353 //
354 // store in the runinfo class the GL_RUN variables for our run
355 //
356 sgnl = 0;
357 sgnl = runinfo->GetRunInfo(idRun);
358 if ( sgnl ){
359 if ( debug ) printf(" ND - ERROR: RunInfo exited with non-zero status\n");
360 code = sgnl;
361 goto closeandexit;
362 } else {
363 sgnl = 0;
364 };
365 //
366 // now you can access that variables using the RunInfo class this way runinfo->ID_ROOT_L0
367 //
368 if ( runinfo->ID_ROOT_L0 == 0 ){
369 if ( debug ) printf("\n ND - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
370 code = -5;
371 goto closeandexit;
372 };
373 //
374 // prepare the timesync for the db
375 //
376 if ( !dbc->IsConnected() ) throw -604;
377 dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
378 //
379 // Search in the DB the path and name of the LEVEL0 file to be processed.
380 //
381 if ( !dbc->IsConnected() ) throw -604;
382 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
383 //
384 ftmpname.str("");
385 ftmpname << glroot->PATH.Data() << "/";
386 ftmpname << glroot->NAME.Data();
387 fname = ftmpname.str().c_str();
388 //
389 // print out informations
390 //
391 totevent = runinfo->NEVENTS;
392 if (verbose) {
393 printf("\n LEVEL0 data file: %s \n",fname.Data());
394 printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
395 printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
396 printf(" %i events to be processed for run %u: from %i to %i \n\n",totevent,idRun,runinfo->EV_FROM,runinfo->EV_FROM+totevent);
397 }
398 //
399 // Open Level0 file
400 //
401 l0File = new TFile(fname.Data());
402 if ( !l0File ) {
403 if ( debug ) printf(" ND - ERROR: problems opening Level0 file\n");
404 code = -6;
405 goto closeandexit;
406 };
407 l0tr = (TTree*)l0File->Get("Physics");
408 if ( !l0tr ) {
409 if ( debug ) printf(" ND - ERROR: no Physics tree in Level0 file\n");
410 l0File->Close();
411 code = -7;
412 goto closeandexit;
413 };
414 l0head = l0tr->GetBranch("Header");
415 if ( !l0head ) {
416 if ( debug ) printf(" ND - ERROR: no Header branch in Level0 tree\n");
417 l0File->Close();
418 code = -8;
419 goto closeandexit;
420 };
421 l0ND = l0tr->GetBranch("Neutron");
422 if ( !l0ND ) {
423 if ( debug ) printf(" ND - ERROR: no ND branch in Level0 tree\n");
424 l0File->Close();
425 code = -603;
426 goto closeandexit;
427 };
428
429 //
430 l0tr->SetBranchAddress("Neutron", &l0ne);
431 l0tr->SetBranchAddress("Header", &eh);
432 //
433 nevents = l0ND->GetEntries();
434 //
435 if ( nevents < 1 ) {
436 if ( debug ) printf(" ND - ERROR: Level0 file is empty\n\n");
437 l0File->Close();
438 code = -11;
439 goto closeandexit;
440 };
441 //
442 if ( runinfo->EV_TO > nevents-1 ) {
443 if ( debug ) printf(" ND - ERROR: too few entries in the tree\n");
444 l0File->Close();
445 code = -12;
446 goto closeandexit;
447 };
448 //
449 // run over all the events of the run
450 //
451 if (verbose) printf("\n Ready to start! \n\n Processed events: \n\n");
452 //
453 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
454 //
455 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
456 //
457 l0head->GetEntry(re);
458 //
459 // absolute time of this event
460 //
461 ph = eh->GetPscuHeader();
462 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
463 //
464 // paranoid check
465 //
466 if ( (atime > runinfo->RUNTRAILER_TIME) || (atime < runinfo->RUNHEADER_TIME) ) {
467 if (verbose) printf(" ND - WARNING: event at time outside the run time window, skipping it\n");
468 goto jumpev;
469 };
470 //
471 procev++;
472 //
473 // start processing
474 //
475 nd->Clear();
476 //
477 l0ND->GetEntry(re);
478 tmpSize = l0ne->Records->GetEntries();
479 if ( tmpSize && l0ne->unpackError == 0 ){
480 for (Int_t j = 0; j < tmpSize; j++){
481 l0nr = (pamela::neutron::NeutronRecord*)l0ne->Records->At(j);
482 yTrig = yTrig + (float)l0nr->trigPhysics;
483 yUpperBackground = yUpperBackground + (float)l0nr->upperBack;
484 yBottomBackground = yBottomBackground + (float)l0nr->bottomBack;
485 }
486 nd->upperBack = yUpperBackground;
487 nd->bottomBack = yBottomBackground;
488 nd->trigPhysics = yTrig;
489 };
490 //
491 nd->unpackError = l0ne->unpackError;
492 //
493 NDtr->Fill(); // Fill tree NDLevel2 in the root file
494 //
495 yUpperBackground =0;
496 yBottomBackground =0;
497 yTrig=0;
498 //
499 //
500 jumpev:
501 debug = false;
502 //
503 };
504 //
505 // Here you may want to clear some variables before processing another run
506 //
507 delete dbtime;
508 }; // process all the runs
509 //
510 if (verbose) printf("\n Finished processing data \n");
511 //
512 closeandexit:
513 //
514 // 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.
515 //
516 if ( !reprocall && reproc && code >= 0 ){
517 if ( totfileentries > noaftrun ){
518 if (verbose) {
519 printf("\n Post-processing: copying events from the old tree after the processed run\n");
520 printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
521 printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
522 }
523 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
524 //
525 // Get entry from old tree
526 //
527 NDtrclone->GetEntry(j);
528 //
529 // copy ndclone to nd
530 //
531 nd->Clear();
532 memcpy(&nd,&ndclone,sizeof(ndclone));
533 //
534 // Fill entry in the new tree
535 //
536 NDtr->Fill();
537 };
538 if (verbose) printf(" Finished successful copying!\n");
539 };
540 };
541 //
542 // Close files, delete old tree(s), write and close level2 file
543 //
544 if ( l0File ) l0File->Close();
545 if ( tempfile ) tempfile->Close();
546 gSystem->Unlink(tempname.str().c_str());
547 //
548 if ( runinfo ) runinfo->Close();
549 if ( NDtr ) NDtr->SetName("NeutronD");
550 if ( file ){
551 file->cd();
552 file->Write();
553 };
554 //
555 gSystem->Unlink(NDfolder.str().c_str());
556 //
557 // the end
558 //
559 if ( verbose ) printf("\n Exiting...\n");
560 if ( NDtr ) NDtr->Delete();
561 //
562 if ( nd ) delete nd;
563 if ( ndclone ) delete ndclone;
564 if ( glroot ) delete glroot;
565 if ( runinfo ) delete runinfo;
566 //
567 if(code < 0) throw code;
568 return(code);
569 }
570
571
572

  ViewVC Help
Powered by ViewVC 1.1.23