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

  ViewVC Help
Powered by ViewVC 1.1.23