/[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.8 - (show annotations) (download)
Thu Sep 28 13:31:53 2006 UTC (18 years, 5 months ago) by mocchiut
Branch: MAIN
CVS Tags: v2r00BETA
Changes since 1.7: +3 -0 lines
Call trk->Clear() before reading tracker data in Calo/ToF-Core.cpp

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 []);
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 dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
404 //
405 // Search in the DB the path and name of the LEVEL0 file to be processed.
406 //
407 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
408 //
409 ftmpname.str("");
410 ftmpname << glroot->PATH.Data() << "/";
411 ftmpname << glroot->NAME.Data();
412 fname = ftmpname.str().c_str();
413 //
414 // print out informations
415 //
416 totevent = runinfo->NEVENTS;
417 if ( verbose ) printf("\n LEVEL0 data file: %s \n",fname.Data());
418 if ( verbose ) printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
419 if ( verbose ) printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
420 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);
421 //
422 // Open Level0 file
423 //
424 l0File = new TFile(fname.Data());
425 if ( !l0File ) {
426 if ( verbose ) printf(" TOF - ERROR: problems opening Level0 file\n");
427 code = -6;
428 goto closeandexit;
429 };
430 l0tr = (TTree*)l0File->Get("Physics");
431 if ( !l0tr ) {
432 if ( verbose ) printf(" TOF - ERROR: no Physics tree in Level0 file\n");
433 l0File->Close();
434 code = -7;
435 goto closeandexit;
436 };
437 l0head = l0tr->GetBranch("Header");
438 if ( !l0head ) {
439 if ( verbose ) printf(" TOF - ERROR: no Header branch in Level0 tree\n");
440 l0File->Close();
441 code = -8;
442 goto closeandexit;
443 };
444 l0trig = l0tr->GetBranch("Trigger");
445 if ( !l0trig ) {
446 if ( verbose ) printf(" TOF - ERROR: no Trigger branch in Level0 tree\n");
447 l0File->Close();
448 code = -300;
449 goto closeandexit;
450 };
451 l0tof = l0tr->GetBranch("Tof");
452 if ( !l0tof ) {
453 if ( verbose ) printf(" TOF - ERROR: no ToF branch in Level0 tree\n");
454 l0File->Close();
455 code = -303;
456 goto closeandexit;
457 };
458 //
459 l0tr->SetBranchAddress("Trigger", &trig);
460 l0tr->SetBranchAddress("Tof", &tofEvent);
461 l0tr->SetBranchAddress("Header", &eh);
462 //
463 nevents = l0tof->GetEntries();
464 //
465 if ( nevents < 1 ) {
466 if ( verbose ) printf(" TOF - ERROR: Level0 file is empty\n\n");
467 l0File->Close();
468 code = -11;
469 goto closeandexit;
470 };
471 //
472 if ( runinfo->EV_TO > nevents-1 ) {
473 if ( verbose ) printf(" TOF - ERROR: too few entries in the registry tree\n");
474 l0File->Close();
475 code = -12;
476 goto closeandexit;
477 };
478 //
479 // Check if we have to load parameter files (or calibration associated to runs and not to events)
480 //
481 // for example let's assume that we could have different magnetic field maps for different runs:
482 //
483 if ( trkpar1 || ( tttrkpar1 != 0 && tttrkpar1 < runinfo->RUNHEADER_TIME ) ){
484 trkpar1 = false;
485 // read from DB infos about Magnetic filed maps
486 glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,1,dbc); // parameters stored in DB in GL_PRAM table
487 tttrkpar1 = glparam->TO_TIME;
488 // ----------------------------
489 // Read the magnetic field
490 // ----------------------------
491 if ( verbose ) printf(" Reading magnetic field maps: \n");
492 trk->LoadField(glparam->PATH+glparam->NAME);
493 if ( verbose ) printf("\n");
494 };
495 //
496 if ( tofpar1 || ( tttofpar1 != 0 && tttofpar1 < runinfo->RUNHEADER_TIME ) ){
497 tofpar1 = false;
498 //
499 Int_t error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,201,dbc); // parameters stored in DB in GL_PRAM table
500 if ( error<0 ) {
501 code = error;
502 goto closeandexit;
503 };
504 //
505 tttofpar1 = glparam->TO_TIME;
506 rdtofcal((char *)(glparam->PATH+glparam->NAME).Data());
507 };
508 //
509 // run over all the events of the run
510 //
511 if ( verbose ) printf("\n Ready to start! \n\n Processed events: \n\n");
512 //
513 jumped = 0;
514 //
515 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
516 //
517 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
518 //
519 l0head->GetEntry(re);
520 //
521 // absolute time of this event
522 //
523 ph = eh->GetPscuHeader();
524 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
525 //
526 tof->Clear();
527 Int_t pmt_id = 0;
528 ToFPMT *t_pmt = new ToFPMT();
529 TClonesArray &tpmt = *tof->PMT;
530 ToFTrkVar *t_tof = new ToFTrkVar();
531 TClonesArray &t = *tof->ToFTrk;
532 //
533 // paranoid check
534 //
535 if ( atime > runinfo->RUNTRAILER_TIME || atime < runinfo->RUNHEADER_TIME ) {
536 if ( verbose ) printf(" TOF - WARNING: event at time outside the run time window, skipping it\n");
537 jumped++;
538 goto jumpev;
539 };
540 //
541 // retrieve tracker informations, the LEVEL2 entry which correspond to our event will be "itr"
542 //
543 if ( !reprocall ){
544 itr = nobefrun + (re - runinfo->EV_FROM -jumped);
545 } else {
546 itr = runinfo->GetFirstEntry() + (re - runinfo->EV_FROM -jumped);
547 };
548 if ( itr > nevtrkl2 ){ // nevtrkl2 tracker entry number
549 if ( verbose ) printf(" TOF - ERROR: no tracker events with entry = %i in Level2 file\n",itr);
550 l0File->Close();
551 code = -313;
552 goto closeandexit;
553 };
554 //
555 trk->Clear();
556 //
557 tracker->GetEntry(itr);
558 ///
559 //
560 l0tof->GetEntry(re);
561 l0trig->GetEntry(re);
562 ///
563 //
564 procev++;
565 //
566 // start processing
567 //
568 //
569 // Here we will use some procedure to calibrate our data and put some kind of informations in the cinput structure
570
571 for (Int_t gg=0; gg<4;gg++){
572 for (Int_t hh=0; hh<12;hh++){
573 tofinput_.tdc[hh][gg]=tofEvent->tdc[gg][hh];
574 tofinput_.adc[hh][gg]=tofEvent->adc[gg][hh];
575 };
576 };
577
578 for (Int_t hh=0; hh<5;hh++){
579 tofinput_.patterntrig[hh]=trig->patterntrig[hh];
580 };
581 //
582 // Here we have calibrated data, ready to be passed to the FORTRAN routine which will extract common and track-related variables.
583 //
584 npmtentry = 0;
585 //
586 ntrkentry = 0;
587 //
588 // Calculate tracks informations from ToF alone
589 //
590 tofl2com();
591 //
592 memcpy(tof->tof_j_flag,tofoutput_.tof_j_flag,6*sizeof(Int_t));
593 //
594 t_tof->trkseqno = -1;
595 //
596 // and now we must copy from the output structure to the level2 class:
597 //
598 t_tof->npmttdc = 0;
599 //
600 for (Int_t hh=0; hh<12;hh++){
601 for (Int_t kk=0; kk<4;kk++){
602 if ( tofoutput_.tofmask[hh][kk] != 0 ){
603 pmt_id = tof->GetPMTid(kk,hh);
604 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
605 t_tof->npmttdc++;
606 };
607 };
608 };
609 for (Int_t kk=0; kk<13;kk++){
610 t_tof->beta[kk] = tofoutput_.betatof_a[kk];
611 }
612 //
613 t_tof->npmtadc = 0;
614 for (Int_t hh=0; hh<12;hh++){
615 for (Int_t kk=0; kk<4;kk++){
616 if ( tofoutput_.adctof_c[hh][kk] < 1000 ){
617 t_tof->dedx.AddAt(tofoutput_.adctof_c[hh][kk],t_tof->npmtadc);
618 pmt_id = tof->GetPMTid(kk,hh);
619 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
620 t_tof->npmtadc++;
621 };
622 };
623 };
624 //
625 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
626 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
627 //
628 new(t[ntrkentry]) ToFTrkVar(*t_tof);
629 ntrkentry++;
630 t_tof->Clear();
631 //
632 //
633 //
634 t_pmt->Clear();
635 //
636 for (Int_t hh=0; hh<12;hh++){
637 for (Int_t kk=0; kk<4;kk++){
638 if ( tofoutput_.tdc_c[hh][kk] < 4095 || tofEvent->adc[kk][hh] < 4095 ){
639 //
640 t_pmt->pmt_id = tof->GetPMTid(kk,hh);
641 t_pmt->tdc_tw = tofoutput_.tdc_c[hh][kk];
642 t_pmt->adc = tofEvent->adc[kk][hh];
643 //
644 new(tpmt[npmtentry]) ToFPMT(*t_pmt);
645 npmtentry++;
646 t_pmt->Clear();
647 };
648 };
649 };
650 //
651 // Calculate track-related variables
652 //
653 if ( trk->ntrk() > 0 ){
654 //
655 // We have at least one track
656 //
657 //
658 // Run over tracks
659 //
660 for(Int_t nt=0; nt < trk->ntrk(); nt++){
661 //
662 TrkTrack *ptt = trk->GetStoredTrack(nt);
663 //
664 // Copy the alpha vector in the input structure
665 //
666 for (Int_t e = 0; e < 5 ; e++){
667 tofinput_.al_pp[e] = ptt->al[e];
668 };
669 //
670 // Get tracker related variables for this track
671 //
672 toftrk();
673 //
674 // Copy values in the class from the structure (we need to use a temporary class to store variables).
675 //
676 t_tof->npmttdc = 0;
677 for (Int_t hh=0; hh<12;hh++){
678 for (Int_t kk=0; kk<4;kk++){
679 if ( tofoutput_.tofmask[hh][kk] != 0 ){
680 pmt_id = tof->GetPMTid(kk,hh);
681 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
682 t_tof->npmttdc++;
683 };
684 };
685 };
686 for (Int_t kk=0; kk<13;kk++){
687 t_tof->beta[kk] = tofoutput_.beta_a[kk];
688 };
689 //
690 t_tof->npmtadc = 0;
691 for (Int_t hh=0; hh<12;hh++){
692 for (Int_t kk=0; kk<4;kk++){
693 if ( tofoutput_.adc_c[hh][kk] < 1000 ){
694 t_tof->dedx.AddAt(tofoutput_.adc_c[hh][kk],t_tof->npmtadc);
695 pmt_id = tof->GetPMTid(kk,hh);
696 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
697 t_tof->npmtadc++;
698 };
699 };
700 };
701 //
702 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
703 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
704 //
705 // Store the tracker track number in order to be sure to have shyncronized data during analysis
706 //
707 t_tof->trkseqno = nt;
708 //
709 // create a new object for this event with track-related variables
710 //
711 new(t[ntrkentry]) ToFTrkVar(*t_tof);
712 ntrkentry++;
713 t_tof->Clear();
714 //
715 }; // loop on all the tracks
716 };
717 //
718 // Fill the rootple
719 //
720 toft->Fill();
721 //
722 //
723 //
724 delete t_tof;
725 //
726 //
727 //
728 jumpev:
729 debug = false;
730 //
731 };
732 //
733 // Here you may want to clear some variables before processing another run
734 //
735 delete dbtime;
736 }; // process all the runs
737 //
738 if ( verbose ) printf("\n Finished processing data \n");
739 //
740 closeandexit:
741 //
742 // 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.
743 //
744 if ( !reprocall && reproc && code >= 0 ){
745 if ( totfileentries > noaftrun ){
746 if ( verbose ) printf("\n Post-processing: copying events from the old tree after the processed run\n");
747 if ( verbose ) printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
748 if ( verbose ) printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
749 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
750 //
751 // Get entry from old tree
752 //
753 toftclone->GetEntry(j);
754 //
755 // copy tofclone to tof
756 //
757 tof->Clear();
758 memcpy(&tof,&tofclone,sizeof(tofclone));
759 //
760 // Fill entry in the new tree
761 //
762 toft->Fill();
763 };
764 if ( verbose ) printf(" Finished successful copying!\n");
765 };
766 };
767 //
768 // Close files, delete old tree(s), write and close level2 file
769 //
770 if ( l0File ) l0File->Close();
771 if ( tempfile ) tempfile->Close();
772 gSystem->Unlink(tempname.str().c_str());
773 if ( tracker ) tracker->Delete(); // delete tracker tree from memory only to avoid writing a copy to file!
774 //
775 if ( code < 0 && verbose ) printf("\n TOF - ERROR: an error occurred, try to save anyway...\n");
776 if ( verbose ) printf("\n Writing and closing rootple\n");
777 if ( runinfo ) runinfo->Close();
778 if ( toft ) toft->SetName("ToF");
779 if ( file ){
780 file->cd();
781 file->Write("ToF");
782 };
783 //
784 gSystem->Unlink(toffolder.str().c_str());
785 //
786 // the end
787 //
788 if ( verbose ) printf("\n Exiting...\n");
789 if(toft)toft->Delete();
790 //
791 if ( tof ) delete tof;
792 if ( tofclone ) delete tofclone;
793 if ( glroot ) delete glroot;
794 if ( runinfo ) delete runinfo;
795 //
796 if ( code < 0 ) throw code;
797 return(code);
798 }

  ViewVC Help
Powered by ViewVC 1.1.23