/[PAMELA software]/DarthVader/ToFLevel2/src/ToFCore.cpp
ViewVC logotype

Contents of /DarthVader/ToFLevel2/src/ToFCore.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.14 - (show annotations) (download)
Mon Jan 22 10:45:24 2007 UTC (18 years, 1 month ago) by mocchiut
Branch: MAIN
Changes since 1.13: +4 -9 lines
ToF routines updated

1 //
2 // C/C++ headers
3 //
4 #include <fstream>
5 #include <string.h>
6 //
7 // ROOT headers
8 //
9 #include <TTree.h>
10 #include <TClassEdit.h>
11 #include <TObject.h>
12 #include <TList.h>
13 #include <TArrayI.h>
14 #include <TSystem.h>
15 #include <TSystemDirectory.h>
16 #include <TString.h>
17 #include <TFile.h>
18 #include <TClass.h>
19 #include <TCanvas.h>
20 #include <TH1.h>
21 #include <TH1F.h>
22 #include <TH2D.h>
23 #include <TLatex.h>
24 #include <TPad.h>
25 #include <TSQLServer.h>
26 #include <TSQLRow.h>
27 #include <TSQLResult.h>
28 #include <TClonesArray.h>
29 //
30 // YODA headers
31 //
32 #include <PamelaRun.h>
33 #include <physics/trigger/TriggerEvent.h>
34 #include <physics/tof/TofEvent.h>
35 //
36 // RunInfo header
37 //
38 #include <RunInfo.h>
39 //
40 // This program headers
41 //
42 #include <ToFCore.h>
43 #include <ToFLevel2.h>
44 #include <ToFVerl2.h>
45 //
46 //
47 // Declaration of the core fortran routines
48 //
49 #define tofl2com tofl2com_
50 extern "C" int tofl2com();
51 #define toftrk toftrk_
52 extern "C" int toftrk();
53 #define rdtofcal rdtofcal_
54 extern "C" int rdtofcal(char [], int *);
55
56 //
57 // Tracker classes headers and definitions
58 //
59 #include <TrkLevel2.h>
60 //
61 using namespace std;
62 //
63 //
64 // CORE ROUTINE
65 //
66 //
67 int ToFCore(UInt_t run, TFile *file, TSQLServer *dbc, Int_t ToFargc, char *ToFargv[]){
68 //
69 //
70 // Set these to true to have a very verbose output.
71 //
72 Bool_t verbose = false;
73 Bool_t debug = false;
74 //
75 TString processFolder = "ToFFolder";
76 if ( ToFargc > 0 ){
77 Int_t i = 0;
78 while ( i < ToFargc ){
79 if ( !strcmp(ToFargv[i],"-processFolder") ) {
80 if ( ToFargc < i+1 ){
81 throw -3;
82 };
83 processFolder = (TString)ToFargv[i+1];
84 i++;
85 };
86 if ( !strcmp(ToFargv[i],"-v") || !strcmp(ToFargv[i],"--verbose") ) {
87 verbose = true;
88 };
89 if ( !strcmp(ToFargv[i],"-g") || !strcmp(ToFargv[i],"--debug") ) {
90 debug = true;
91 };
92 i++;
93 };
94 };
95 //
96 //
97 // Output directory is the working directoy.
98 //
99 const char* outdir = gSystem->DirName(gSystem->DirName(file->GetPath()));
100 //
101 // Variables for level2
102 //
103 TTree *tracker = 0;
104 TTree *toft = 0;
105 UInt_t nevents = 0;
106 //
107 // variables needed to reprocess data
108 //
109 TString tofversion;
110 ItoRunInfo *runinfo = 0;
111 TArrayI *runlist = 0;
112 TTree *toftclone = 0;
113 Bool_t reproc = false;
114 Bool_t reprocall = false;
115 UInt_t nobefrun = 0;
116 UInt_t noaftrun = 0;
117 UInt_t numbofrun = 0;
118 stringstream ftmpname;
119 TString fname;
120 UInt_t totfileentries = 0;
121 UInt_t idRun = 0;
122 //
123 // variables needed to handle error signals
124 //
125 Int_t code = 0;
126 Int_t sgnl;
127 //
128 // tof level2 classes
129 //
130 ToFLevel2 *tof = new ToFLevel2();
131 ToFLevel2 *tofclone = new ToFLevel2();
132 //
133 // tracker level2 variables
134 //
135 TrkLevel2 *trk = new TrkLevel2();
136 Int_t nevtrkl2 = 0;
137 //
138 // define variables for opening and reading level0 file
139 //
140 TFile *l0File = 0;
141 TTree *l0tr = 0;
142 TBranch *l0head = 0;
143 TBranch *l0trig = 0;
144 TBranch *l0tof = 0;
145 pamela::EventHeader *eh = 0;
146 pamela::PscuHeader *ph = 0;
147 pamela::trigger::TriggerEvent *trig = 0;
148 pamela::tof::TofEvent *tofEvent = 0;
149 //
150 // Define other basic variables
151 //
152 UInt_t procev = 0;
153 stringstream file2;
154 stringstream file3;
155 stringstream qy;
156 Int_t itr = -1;
157 Int_t totevent = 0;
158 UInt_t atime = 0;
159 UInt_t re = 0;
160 UInt_t jumped = 0;
161 //
162 // Working filename
163 //
164 TString outputfile;
165 stringstream name;
166 name.str("");
167 name << outdir << "/";
168 //
169 // temporary file and folder
170 //
171 TFile *tempfile = 0;
172 TTree *temptof = 0;
173 stringstream tempname;
174 stringstream toffolder;
175 tempname.str("");
176 tempname << outdir;
177 tempname << "/" << processFolder.Data();
178 toffolder.str("");
179 toffolder << tempname.str().c_str();
180 gSystem->MakeDirectory(toffolder.str().c_str());
181 tempname << "/toftree_run";
182 tempname << run << ".root";
183 //
184 // variables needed to load magnetic field maps
185 //
186 Int_t ntrkentry = 0;
187 Int_t npmtentry = 0;
188 UInt_t tttrkpar1 = 0;
189 Bool_t trkpar1 = true;
190 UInt_t tttofpar1 = 0;
191 Bool_t tofpar1 = true;
192 //
193 // DB classes
194 //
195 GL_ROOT *glroot = new GL_ROOT();
196 GL_PARAM *glparam = new GL_PARAM();
197 GL_TIMESYNC *dbtime = 0;
198 //
199 // declaring external output and input structures
200 //
201 extern struct ToFInput tofinput_;
202 extern struct ToFOutput tofoutput_;
203 //
204 // Let's start!
205 //
206 //
207 // 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
208 // if run != 0 we must process only that run but first we have to check if the tree ToF already exist in the file
209 // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
210 //
211 if ( run == 0 ) reproc = true;
212 //
213 //
214 // Output file is "outputfile"
215 //
216 if ( !file->IsOpen() ){
217 if ( verbose ) printf(" ToF - ERROR: cannot open file for writing\n");
218 throw -301;
219 };
220 //
221 // Does it contain the Tracker tree?
222 //
223 tracker = (TTree*)file->Get("Tracker");
224 if ( !tracker ) {
225 if ( verbose ) printf(" TOF - ERROR: no tracker tree\n");
226 code = -302;
227 goto closeandexit;
228 };
229 //
230 // get tracker level2 data pointer
231 //
232 tracker->SetBranchAddress("TrkLevel2",&trk);
233 nevtrkl2 = tracker->GetEntries();
234 //
235 // Retrieve GL_RUN variables from the level2 file
236 //
237 tofversion = ToFInfo(false); // we should decide how to handle versioning system
238 //
239 // create an interface to RunInfo called "runinfo"
240 //
241 // ItoRunInfo= interface with RunInfo and GL_RUN
242 runinfo = new ItoRunInfo(file);
243 //
244 // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
245 //
246 sgnl = 0;
247 sgnl = runinfo->Update(run, "TOF",tofversion);
248 if ( sgnl ){
249 if ( verbose ) printf(" TOF - ERROR: RunInfo exited with non-zero status\n");
250 code = sgnl;
251 goto closeandexit;
252 } else {
253 sgnl = 0;
254 };
255 //
256 // number of events in the file BEFORE the first event of our run
257 //
258 nobefrun = runinfo->GetFirstEntry();
259 //
260 // total number of events in the file
261 //
262 totfileentries = runinfo->GetFileEntries();
263 //
264 // first file entry AFTER the last event of our run
265 //
266 noaftrun = runinfo->GetLastEntry() + 1;
267 //
268 // number of run to be processed
269 //
270 numbofrun = runinfo->GetNoRun();
271 //
272 // Try to access the ToF tree in the file, if it exists we are reprocessing data if not we are processing a new run
273 //
274 toftclone = (TTree*)file->Get("ToF");
275 //
276 if ( !toftclone ){
277 //
278 // tree does not exist, we are not reprocessing
279 //
280 reproc = false;
281 if ( run == 0 && verbose ) printf(" ToF - WARNING: you are reprocessing data but ToF tree does not exist!\n");
282 if ( runinfo->IsReprocessing() && run != 0 && verbose ) printf(" ToF - WARNING: it seems you are not reprocessing data but ToF\n versioning information already exists in RunInfo.\n");
283
284 } else {
285 //
286 // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
287 //
288 reproc = true;
289 //
290 // update versioning information
291 //
292 if ( verbose ) printf("\n Preparing the pre-processing...\n");
293 //
294 if ( run == 0 ){
295 //
296 // we are reprocessing all the file
297 // 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
298 //
299 reprocall = true;
300 //
301 if ( verbose ) printf("\n ToF - WARNING: Reprocessing all runs\n");
302 //
303 } else {
304 //
305 // 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
306 //
307 reprocall = false;
308 //
309 if ( verbose ) printf("\n ToF - WARNING: Reprocessing run number %u \n",run);
310 //
311 // copying old tree to a new file
312 //
313 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
314 temptof = toftclone->CloneTree(-1,"fast");
315 temptof->SetName("ToF-old");
316 tempfile->Write();
317 tempfile->Close();
318 }
319 //
320 // Delete the old tree from old file and memory
321 //
322 toftclone->Delete("all");
323 //
324 if ( verbose ) printf(" ...done!\n");
325 //
326 };
327 //
328 // create ToF detector tree toft
329 //
330 file->cd();
331 toft = new TTree("ToF-new","PAMELA Level2 ToF data");
332 toft->Branch("ToFLevel2","ToFLevel2",&tof);
333 //
334 if ( reproc && !reprocall ){
335 //
336 // open new file and retrieve all tree informations
337 //
338 tempfile = new TFile(tempname.str().c_str(),"READ");
339 toftclone = (TTree*)tempfile->Get("ToF-old");
340 toftclone->SetBranchAddress("ToFLevel2",&tofclone);
341 //
342 if ( nobefrun > 0 ){
343 if ( verbose ) printf("\n Pre-processing: copying events from the old tree before the processed run\n");
344 if ( verbose ) printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
345 if ( verbose ) printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
346 for (UInt_t j = 0; j < nobefrun; j++){
347 //
348 toftclone->GetEntry(j);
349 //
350 // copy tofclone to tof
351 //
352 tof->Clear();
353 memcpy(&tof,&tofclone,sizeof(tofclone));
354 //
355 // Fill entry in the new tree
356 //
357 toft->Fill();
358 //
359 };
360 if ( verbose ) printf(" Finished successful copying!\n");
361 };
362 };
363 //
364 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
365 //
366 runlist = runinfo->GetRunList();
367 //
368 // Loop over the run to be processed
369 //
370 for (UInt_t irun=0; irun < numbofrun; irun++){
371 //
372 // retrieve the first run ID to be processed using the RunInfo list
373 //
374 idRun = runlist->At(irun);
375 if ( verbose ) printf("\n\n\n ####################################################################### \n");
376 if ( verbose ) printf(" PROCESSING RUN NUMBER %u \n",idRun);
377 if ( verbose ) printf(" ####################################################################### \n\n\n");
378 //
379 runinfo->ID_ROOT_L0 = 0;
380 //
381 // store in the runinfo class the GL_RUN variables for our run
382 //
383 sgnl = 0;
384 sgnl = runinfo->GetRunInfo(idRun);
385 if ( sgnl ){
386 if ( verbose ) printf(" TOF - ERROR: RunInfo exited with non-zero status\n");
387 code = sgnl;
388 goto closeandexit;
389 } else {
390 sgnl = 0;
391 };
392 //
393 // now you can access that variables using the RunInfo class this way runinfo->ID_ROOT_L0
394 //
395 if ( runinfo->ID_ROOT_L0 == 0 ){
396 if ( verbose ) printf("\n TOF - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
397 code = -5;
398 goto closeandexit;
399 };
400 //
401 // prepare the timesync for the db
402 //
403 if ( !dbc->IsConnected() ) throw -314;
404 dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
405 //
406 // Search in the DB the path and name of the LEVEL0 file to be processed.
407 //
408 if ( !dbc->IsConnected() ) throw -314;
409 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
410 //
411 ftmpname.str("");
412 ftmpname << glroot->PATH.Data() << "/";
413 ftmpname << glroot->NAME.Data();
414 fname = ftmpname.str().c_str();
415 //
416 // print out informations
417 //
418 totevent = runinfo->NEVENTS;
419 if ( verbose ) printf("\n LEVEL0 data file: %s \n",fname.Data());
420 if ( verbose ) printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
421 if ( verbose ) printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
422 if ( verbose ) printf(" %i events to be processed for run %u: from %i to %i \n\n",totevent,idRun,runinfo->EV_FROM,runinfo->EV_FROM+totevent);
423 //
424 // Open Level0 file
425 //
426 l0File = new TFile(fname.Data());
427 if ( !l0File ) {
428 if ( verbose ) printf(" TOF - ERROR: problems opening Level0 file\n");
429 code = -6;
430 goto closeandexit;
431 };
432 l0tr = (TTree*)l0File->Get("Physics");
433 if ( !l0tr ) {
434 if ( verbose ) printf(" TOF - ERROR: no Physics tree in Level0 file\n");
435 l0File->Close();
436 code = -7;
437 goto closeandexit;
438 };
439 l0head = l0tr->GetBranch("Header");
440 if ( !l0head ) {
441 if ( verbose ) printf(" TOF - ERROR: no Header branch in Level0 tree\n");
442 l0File->Close();
443 code = -8;
444 goto closeandexit;
445 };
446 l0trig = l0tr->GetBranch("Trigger");
447 if ( !l0trig ) {
448 if ( verbose ) printf(" TOF - ERROR: no Trigger branch in Level0 tree\n");
449 l0File->Close();
450 code = -300;
451 goto closeandexit;
452 };
453 l0tof = l0tr->GetBranch("Tof");
454 if ( !l0tof ) {
455 if ( verbose ) printf(" TOF - ERROR: no ToF branch in Level0 tree\n");
456 l0File->Close();
457 code = -303;
458 goto closeandexit;
459 };
460 //
461 l0tr->SetBranchAddress("Trigger", &trig);
462 l0tr->SetBranchAddress("Tof", &tofEvent);
463 l0tr->SetBranchAddress("Header", &eh);
464 //
465 nevents = l0tof->GetEntries();
466 //
467 if ( nevents < 1 ) {
468 if ( verbose ) printf(" TOF - ERROR: Level0 file is empty\n\n");
469 l0File->Close();
470 code = -11;
471 goto closeandexit;
472 };
473 //
474 if ( runinfo->EV_TO > nevents-1 ) {
475 if ( verbose ) printf(" TOF - ERROR: too few entries in the registry tree\n");
476 l0File->Close();
477 code = -12;
478 goto closeandexit;
479 };
480 //
481 // Check if we have to load parameter files (or calibration associated to runs and not to events)
482 //
483 // for example let's assume that we could have different magnetic field maps for different runs:
484 //
485 if ( trkpar1 || ( tttrkpar1 != 0 && tttrkpar1 < runinfo->RUNHEADER_TIME ) ){
486 trkpar1 = false;
487 // read from DB infos about Magnetic filed maps
488 if ( !dbc->IsConnected() ) throw -314;
489 glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,1,dbc); // parameters stored in DB in GL_PRAM table
490 tttrkpar1 = glparam->TO_TIME;
491 // ----------------------------
492 // Read the magnetic field
493 // ----------------------------
494 if ( verbose ) printf(" Reading magnetic field maps: \n");
495 trk->LoadField(glparam->PATH+glparam->NAME);
496 if ( verbose ) printf("\n");
497 };
498 //
499 if ( tofpar1 || ( tttofpar1 != 0 && tttofpar1 < runinfo->RUNHEADER_TIME ) ){
500 tofpar1 = false;
501 //
502 if ( !dbc->IsConnected() ) throw -314;
503 Int_t error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,201,dbc); // parameters stored in DB in GL_PRAM table
504 if ( error<0 ) {
505 code = error;
506 goto closeandexit;
507 };
508 //
509 if ( verbose ) printf(" Reading ToF parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
510 //
511 tttofpar1 = glparam->TO_TIME;
512 Int_t nlen = (Int_t)(glparam->PATH+glparam->NAME).Length();
513 rdtofcal((char *)(glparam->PATH+glparam->NAME).Data(),&nlen);
514 };
515 //
516 // run over all the events of the run
517 //
518 if ( verbose ) printf("\n Ready to start! \n\n Processed events: \n\n");
519 //
520 jumped = 0;
521 //
522 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
523 //
524 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
525 //
526 l0head->GetEntry(re);
527 //
528 // absolute time of this event
529 //
530 ph = eh->GetPscuHeader();
531 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
532 //
533 tof->Clear();
534 Int_t pmt_id = 0;
535 ToFPMT *t_pmt = new ToFPMT();
536 TClonesArray &tpmt = *tof->PMT;
537 ToFTrkVar *t_tof = new ToFTrkVar();
538 TClonesArray &t = *tof->ToFTrk;
539 //
540 // paranoid check
541 //
542 if ( atime > runinfo->RUNTRAILER_TIME || atime < runinfo->RUNHEADER_TIME ) {
543 if ( verbose ) printf(" TOF - WARNING: event at time outside the run time window, skipping it\n");
544 jumped++;
545 goto jumpev;
546 };
547 //
548 // retrieve tracker informations, the LEVEL2 entry which correspond to our event will be "itr"
549 //
550 if ( !reprocall ){
551 itr = nobefrun + (re - runinfo->EV_FROM -jumped);
552 } else {
553 itr = runinfo->GetFirstEntry() + (re - runinfo->EV_FROM -jumped);
554 };
555 if ( itr > nevtrkl2 ){ // nevtrkl2 tracker entry number
556 if ( verbose ) printf(" TOF - ERROR: no tracker events with entry = %i in Level2 file\n",itr);
557 l0File->Close();
558 code = -313;
559 goto closeandexit;
560 };
561 //
562 trk->Clear();
563 //
564 tracker->GetEntry(itr);
565 ///
566 //
567 l0tof->GetEntry(re);
568 l0trig->GetEntry(re);
569 ///
570 //
571 procev++;
572 //
573 // start processing
574 //
575 //
576 // Here we will use some procedure to calibrate our data and put some kind of informations in the cinput structure
577
578 for (Int_t gg=0; gg<4;gg++){
579 for (Int_t hh=0; hh<12;hh++){
580 tofinput_.tdc[hh][gg]=tofEvent->tdc[gg][hh];
581 tofinput_.adc[hh][gg]=tofEvent->adc[gg][hh];
582 };
583 };
584
585 for (Int_t hh=0; hh<5;hh++){
586 tofinput_.patterntrig[hh]=trig->patterntrig[hh];
587 };
588 //
589 // Here we have calibrated data, ready to be passed to the FORTRAN routine which will extract common and track-related variables.
590 //
591 npmtentry = 0;
592 //
593 ntrkentry = 0;
594 //
595 // Calculate tracks informations from ToF alone
596 //
597 tofl2com();
598 //
599 memcpy(tof->tof_j_flag,tofoutput_.tof_j_flag,6*sizeof(Int_t));
600 //
601 t_tof->trkseqno = -1;
602 //
603 // and now we must copy from the output structure to the level2 class:
604 //
605 t_tof->npmttdc = 0;
606 //
607 for (Int_t hh=0; hh<12;hh++){
608 for (Int_t kk=0; kk<4;kk++){
609 if ( tofoutput_.tofmask[hh][kk] != 0 ){
610 pmt_id = tof->GetPMTid(kk,hh);
611 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
612 t_tof->tdcflag.AddAt(tofoutput_.tdcflagtof[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
613 t_tof->npmttdc++;
614 };
615 };
616 };
617 for (Int_t kk=0; kk<13;kk++){
618 t_tof->beta[kk] = tofoutput_.betatof_a[kk];
619 }
620 //
621 t_tof->npmtadc = 0;
622 for (Int_t hh=0; hh<12;hh++){
623 for (Int_t kk=0; kk<4;kk++){
624 if ( tofoutput_.adctof_c[hh][kk] < 1000 ){
625 t_tof->dedx.AddAt(tofoutput_.adctof_c[hh][kk],t_tof->npmtadc);
626 pmt_id = tof->GetPMTid(kk,hh);
627 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
628 t_tof->adcflag.AddAt(tofoutput_.adcflagtof[hh][kk],t_tof->npmtadc); // gf: Jan 09/07
629 t_tof->npmtadc++;
630 };
631 };
632 };
633 //
634 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
635 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
636 //
637 new(t[ntrkentry]) ToFTrkVar(*t_tof);
638 ntrkentry++;
639 t_tof->Clear();
640 //
641 //
642 //
643 t_pmt->Clear();
644 //
645 for (Int_t hh=0; hh<12;hh++){
646 for (Int_t kk=0; kk<4;kk++){
647 if ( tofoutput_.tdc_c[hh][kk] < 4095 || tofEvent->adc[kk][hh] < 4095 ){
648 //
649 t_pmt->pmt_id = tof->GetPMTid(kk,hh);
650 t_pmt->tdc_tw = tofoutput_.tdc_c[hh][kk];
651 t_pmt->adc = tofEvent->adc[kk][hh];
652 //
653 new(tpmt[npmtentry]) ToFPMT(*t_pmt);
654 npmtentry++;
655 t_pmt->Clear();
656 };
657 };
658 };
659 //
660 // Calculate track-related variables
661 //
662 if ( trk->ntrk() > 0 ){
663 //
664 // We have at least one track
665 //
666 //
667 // Run over tracks
668 //
669 for(Int_t nt=0; nt < trk->ntrk(); nt++){
670 //
671 TrkTrack *ptt = trk->GetStoredTrack(nt);
672 //
673 // Copy the alpha vector in the input structure
674 //
675 for (Int_t e = 0; e < 5 ; e++){
676 tofinput_.al_pp[e] = ptt->al[e];
677 };
678 //
679 // Get tracker related variables for this track
680 //
681 toftrk();
682 //
683 // Copy values in the class from the structure (we need to use a temporary class to store variables).
684 //
685 t_tof->npmttdc = 0;
686 for (Int_t hh=0; hh<12;hh++){
687 for (Int_t kk=0; kk<4;kk++){
688 if ( tofoutput_.tofmask[hh][kk] != 0 ){
689 pmt_id = tof->GetPMTid(kk,hh);
690 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
691 t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
692 t_tof->npmttdc++;
693 };
694 };
695 };
696 for (Int_t kk=0; kk<13;kk++){
697 t_tof->beta[kk] = tofoutput_.beta_a[kk];
698 };
699 //
700 t_tof->npmtadc = 0;
701 for (Int_t hh=0; hh<12;hh++){
702 for (Int_t kk=0; kk<4;kk++){
703 if ( tofoutput_.adc_c[hh][kk] < 1000 ){
704 t_tof->dedx.AddAt(tofoutput_.adc_c[hh][kk],t_tof->npmtadc);
705 pmt_id = tof->GetPMTid(kk,hh);
706 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
707 t_tof->adcflag.AddAt(tofoutput_.adcflag[hh][kk],t_tof->npmtadc); // gf: Jan 09/07
708 t_tof->npmtadc++;
709 };
710 };
711 };
712 //
713 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
714 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
715 //
716 // Store the tracker track number in order to be sure to have shyncronized data during analysis
717 //
718 t_tof->trkseqno = nt;
719 //
720 // create a new object for this event with track-related variables
721 //
722 new(t[ntrkentry]) ToFTrkVar(*t_tof);
723 ntrkentry++;
724 t_tof->Clear();
725 //
726 }; // loop on all the tracks
727 };
728 //
729 tof->unpackError = tofEvent->unpackError;
730 //
731 // Fill the rootple
732 //
733 toft->Fill();
734 //
735 //
736 //
737 delete t_tof;
738 //
739 //
740 //
741 jumpev:
742 debug = false;
743 //
744 };
745 //
746 // Here you may want to clear some variables before processing another run
747 //
748 delete dbtime;
749 }; // process all the runs
750 //
751 if ( verbose ) printf("\n Finished processing data \n");
752 //
753 closeandexit:
754 //
755 // 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.
756 //
757 if ( !reprocall && reproc && code >= 0 ){
758 if ( totfileentries > noaftrun ){
759 if ( verbose ) printf("\n Post-processing: copying events from the old tree after the processed run\n");
760 if ( verbose ) printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
761 if ( verbose ) printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
762 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
763 //
764 // Get entry from old tree
765 //
766 toftclone->GetEntry(j);
767 //
768 // copy tofclone to tof
769 //
770 tof->Clear();
771 memcpy(&tof,&tofclone,sizeof(tofclone));
772 //
773 // Fill entry in the new tree
774 //
775 toft->Fill();
776 };
777 if ( verbose ) printf(" Finished successful copying!\n");
778 };
779 };
780 //
781 // Close files, delete old tree(s), write and close level2 file
782 //
783 if ( l0File ) l0File->Close();
784 if ( tempfile ) tempfile->Close();
785 gSystem->Unlink(tempname.str().c_str());
786 if ( tracker ) tracker->Delete(); // delete tracker tree from memory only to avoid writing a copy to file!
787 //
788 if ( code < 0 && verbose ) printf("\n TOF - ERROR: an error occurred, try to save anyway...\n");
789 if ( verbose ) printf("\n Writing and closing rootple\n");
790 if ( runinfo ) runinfo->Close();
791 if ( toft ) toft->SetName("ToF");
792 if ( file ){
793 file->cd();
794 file->Write("ToF");
795 };
796 //
797 gSystem->Unlink(toffolder.str().c_str());
798 //
799 // the end
800 //
801 if ( verbose ) printf("\n Exiting...\n");
802 if(toft)toft->Delete();
803 //
804 if ( tof ) delete tof;
805 if ( tofclone ) delete tofclone;
806 if ( glroot ) delete glroot;
807 if ( runinfo ) delete runinfo;
808 //
809 if ( code < 0 ) throw code;
810 return(code);
811 }

  ViewVC Help
Powered by ViewVC 1.1.23