/[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.8 - (show annotations) (download)
Fri Nov 17 10:08:09 2006 UTC (18 years, 3 months ago) by mocchiut
Branch: MAIN
CVS Tags: v2r01
Changes since 1.7: +2 -0 lines
Added some checks on MySQL connection

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

  ViewVC Help
Powered by ViewVC 1.1.23