/[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.59 - (show annotations) (download)
Tue Oct 14 14:07:34 2014 UTC (10 years, 2 months ago) by mocchiut
Branch: MAIN
CVS Tags: v10RED, v10REDr01, HEAD
Changes since 1.58: +2 -1 lines
10RED: lost sync bug fixed

1 // // C/C++ headers //
2 #include <fstream>
3 #include <string.h>
4 #include <iomanip>
5 //
6 // ROOT headers
7 //
8 #include <TTree.h>
9 #include <TClassEdit.h>
10 #include <TObject.h>
11 #include <TList.h>
12 #include <TArrayI.h>
13 #include <TSystem.h>
14 #include <TSystemDirectory.h>
15 #include <TString.h>
16 #include <TFile.h>
17 #include <TClass.h>
18 #include <TCanvas.h>
19 #include <TH1.h>
20 #include <TH1F.h>
21 #include <TH2D.h>
22 #include <TLatex.h>
23 #include <TPad.h>
24 #include <TSQLServer.h>
25 #include <TSQLRow.h>
26 #include <TSQLResult.h>
27 #include <TClonesArray.h>
28 //
29 // RunInfo header
30 //
31 #include <RunInfo.h>
32 //
33 // YODA headers
34 //
35 #include <PamelaRun.h>
36 //#include <physics/trigger/TriggerEvent.h>
37 #include <physics/tof/TofEvent.h>
38 //
39 // This program headers
40 //
41 #include <ToFCore.h>
42 #include <ToFLevel2.h>
43 #include <ToFVerl2.h>
44 //
45 //
46 // Declaration of the core fortran routines
47 //
48 #define tofl2com tofl2com_
49 extern "C" int tofl2com();
50 #define toftrk toftrk_
51 extern "C" int toftrk();
52 #define rdtofcal rdtofcal_
53 //extern "C" int rdtofcal(char [], int *);
54 extern "C" int rdtofcal(const char *, int *);
55
56 //
57 // Tracker classes headers and definitions
58 //
59 #include <TrkLevel2.h>
60 #include <ExtTrack.h> // new tracking code
61 //
62 using namespace std;
63 //
64 //
65 // CORE ROUTINE
66 //
67 //
68 int ToFCore(UInt_t run, TFile *file, GL_TABLES *glt, Int_t ToFargc, char *ToFargv[]){
69 //
70 //
71 // Set these to true to have a very verbose output.
72 //
73 Bool_t verbose = false;
74 Bool_t debug = false;
75 //
76 Bool_t l1only = false;
77 //
78 Bool_t deltree = false;
79 //
80 //
81 TString processFolder = Form("ToFFolder_%u",run);
82 if ( ToFargc > 0 ){
83 Int_t i = 0;
84 while ( i < ToFargc ){
85 if ( !strcmp(ToFargv[i],"-processFolder") ) {
86 if ( ToFargc < i+1 ){
87 throw -3;
88 };
89 processFolder = (TString)ToFargv[i+1];
90 i++;
91 };
92 if ( !strcmp(ToFargv[i],"-v") || !strcmp(ToFargv[i],"--verbose") ) {
93 verbose = true;
94 };
95 if ( !strcmp(ToFargv[i],"-g") || !strcmp(ToFargv[i],"--debug") ) {
96 verbose = true;
97 debug = true;
98 };
99 if ( !strcmp(ToFargv[i],"--level1-only") ) {
100 l1only = true;
101 }
102 if ( !strcmp(ToFargv[i],"--delete-tree") ) {
103 deltree = true;
104 }
105 i++;
106 };
107 };
108 //
109 //
110 // Output directory is the working directoy.
111 //
112 const char* outdir = gSystem->DirName(gSystem->DirName(file->GetPath()));
113 //
114 // Variables for level2
115 //
116 TTree *tracker = 0;
117 TTree *trigger = 0;
118 TTree *toft = 0;
119 UInt_t nevents = 0;
120 Long64_t maxsize = 10000000000LL;
121 TTree::SetMaxTreeSize(maxsize);
122 //
123 // variables needed to reprocess data
124 //
125 TString tofversion;
126 ItoRunInfo *runinfo = 0;
127 TArrayI *runlist = 0;
128 TTree *toftclone = 0;
129 Bool_t reproc = false;
130 Bool_t reprocall = false;
131 UInt_t nobefrun = 0;
132 UInt_t noaftrun = 0;
133 UInt_t numbofrun = 0;
134 stringstream ftmpname;
135 TString fname;
136 UInt_t totfileentries = 0;
137 UInt_t idRun = 0;
138 //
139 // variables needed to handle error signals
140 //
141 Int_t code = 0;
142 Int_t sgnl;
143 //
144 // tof level2 classes
145 //
146 ToFLevel2 *tof = new ToFLevel2();
147 ToFLevel2 *tofclone = new ToFLevel2();
148 ToFdEdx *tofdedx = new ToFdEdx();
149 //
150 // tracker level2 variables
151 //
152 TrkLevel2 *trk = new TrkLevel2();
153 Int_t nevtrkl2 = 0;
154 //
155 // trigger level2 variables
156 //
157 TrigLevel2 *trg = new TrigLevel2();
158 Int_t nevtrgl2 = 0;
159 //
160 // define variables for opening and reading level0 file
161 //
162 TFile *l0File = 0;
163 TTree *l0tr = 0;
164 TBranch *l0head = 0;
165 // TBranch *l0trig = 0;
166 TBranch *l0tof = 0;
167 pamela::EventHeader *eh = 0;
168 pamela::PscuHeader *ph = 0;
169 // pamela::trigger::TriggerEvent *trig = 0;
170 pamela::tof::TofEvent *tofEvent = 0;
171 //
172 // Define other basic variables
173 //
174 UInt_t procev = 0;
175 stringstream file2;
176 stringstream file3;
177 stringstream qy;
178 Int_t itr = -1;
179 Int_t totevent = 0;
180 UInt_t atime = 0;
181 UInt_t re = 0;
182 UInt_t jumped = 0;
183 //
184 // Working filename
185 //
186 TString outputfile;
187 stringstream name;
188 name.str("");
189 name << outdir << "/";
190 //
191 // temporary file and folder
192 //
193 TFile *tempfile = 0;
194 TTree *temptof = 0;
195 stringstream tempname;
196 stringstream toffolder;
197 Bool_t myfold = false;
198 tempname.str("");
199 tempname << outdir;
200 tempname << "/" << processFolder.Data();
201 toffolder.str("");
202 toffolder << tempname.str().c_str();
203 tempname << "/toftree_run";
204 tempname << run << ".root";
205 UInt_t totnorun = 0;
206 //
207 // variables needed to load magnetic field maps
208 //
209 Int_t ntrkentry = 0;
210 Int_t npmtentry = 0;
211 UInt_t tttrkpar1 = 0;
212 Bool_t trkpar1 = true;
213 UInt_t tttofpar1 = 0;
214 Bool_t tofpar1 = true;
215 //
216 // DB classes
217 //
218 GL_ROOT *glroot = new GL_ROOT();
219 GL_PARAM *glparam = new GL_PARAM();
220 GL_TIMESYNC *dbtime = 0;
221 //
222 // declaring external output and input structures
223 //
224 extern struct ToFInput tofinput_;
225 extern struct ToFOutput tofoutput_;
226 //
227 // WM variables perform dE/dx II order corrections
228 //
229 //Float_t dedx_corr_m[100][48],dedx_corr[48];
230 Double_t mtime[100],t1,t2,tm;
231 //Float_t yhelp1,yhelp2,slope,inter,thelp1,thelp2;
232
233 //RC variables for new dEdx II order correction (10th reduction)
234 Float_t Heyhelp1,Heyhelp2,Heslope,Heinter,thelp1,thelp2;
235 Float_t pyhelp1,pyhelp2,pslope,pinter;
236 Float_t dedx_Hepeak[48],dedx_ppeak[48];
237 Float_t dedx_Hepeak_m[100][48],dedx_ppeak_m[100][48];
238 Float_t inter_dedx[48],slope_dedx[48];
239
240 Float_t xmean1,xwidth1;
241 Int_t ical,ii,wj,jj;
242 Float_t xleft=0;
243 Float_t xright=0;
244 Float_t yleft=0;
245 Float_t yright=0;
246
247 Int_t warning = 0;
248 int a=0, b=0;
249
250 //
251 // Let's start!
252 //
253 //
254 // 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
255 // if run != 0 we must process only that run but first we have to check if the tree ToF already exist in the file
256 // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
257 //
258 if ( run == 0 ) reproc = true;
259 //
260 //
261 // Output file is "outputfile"
262 //
263 if ( !file->IsOpen() ){
264 if ( verbose ) printf(" ToF - ERROR: cannot open file for writing\n");
265 throw -301;
266 }
267
268 if ( debug ) file->ls();
269
270 //
271 // Delete tree if requested
272 //
273 if ( deltree ){
274 TTree *T = (TTree*)file->Get("ToF");
275 if ( T ){
276 if ( verbose ) printf(" ToF - REMOVING ToF TTree \n");
277 T->Delete("all");
278 }
279 }
280
281 //
282 // Does it contain the Tracker tree?
283 //
284 //
285 TClonesArray *tcNucleiTrk = NULL;
286 TClonesArray *tcExtNucleiTrk = NULL;
287 TClonesArray *tcExtTrk = NULL;
288 TClonesArray *ttofNucleiTrk = NULL;
289 TClonesArray *ttofExtNucleiTrk = NULL;
290 TClonesArray *ttofExtTrk = NULL;
291 Bool_t hasNucleiTrk = false;
292 Bool_t hasExtNucleiTrk = false;
293 Bool_t hasExtTrk = false;
294 if ( !l1only ){
295 tracker = (TTree*)file->Get("Tracker");
296 if ( !tracker ) {
297 if ( verbose ) printf(" TOF - ERROR: no tracker tree\n");
298 code = -302;
299 goto closeandexit;
300 }
301 //
302 // get tracker level2 data pointer
303 //
304 // tracker->SetMaxVirtualSize(2500000000LL); // EM residual Tracker-new tree in level2 files when NEVENTS is big
305 tracker->SetBranchAddress("TrkLevel2",&trk);
306 nevtrkl2 = tracker->GetEntries();
307 //
308 // Look for extended tracking algorithm
309 //
310 if ( verbose ) printf("Look for extended and nuclei tracking algorithms\n");
311 // Nuclei tracking algorithm
312 Int_t checkAlgo = 0;
313 tcNucleiTrk = new TClonesArray("TrkTrack");
314 checkAlgo = tracker->SetBranchAddress("TrackNuclei",&tcNucleiTrk);
315 if ( !checkAlgo ){
316 if ( verbose ) printf(" Nuclei tracking algorithm branch found! :D \n");
317 hasNucleiTrk = true;
318 } else {
319 if ( verbose ) printf(" Nuclei tracking algorithm branch not found :( !\n");
320 printf(" ok, this is not a problem (it depends on tracker settings) \n");
321 delete tcNucleiTrk;
322 tcNucleiTrk = NULL; // 10RED reprocessing bug when missing extratracks
323 }
324 // Nuclei tracking algorithm using calorimeter points
325 tcExtNucleiTrk = new TClonesArray("ExtTrack");
326 checkAlgo = tracker->SetBranchAddress("RecoveredTrackNuclei",&tcExtNucleiTrk);
327 if ( !checkAlgo ){
328 if ( verbose ) printf(" Recovered nuclei tracking algorithm branch found! :D \n");
329 hasExtNucleiTrk = true;
330 } else {
331 if ( verbose ) printf(" Recovered nuclei tracking algorithm branch not found :( !\n");
332 printf(" ok, this is not a problem (it depends on tracker settings) \n");
333 delete tcExtNucleiTrk;
334 tcExtNucleiTrk = NULL; // 10RED reprocessing bug when missing extratracks
335 }
336 // Tracking algorithm using calorimeter points
337 tcExtTrk = new TClonesArray("ExtTrack");
338 checkAlgo = tracker->SetBranchAddress("RecoveredTrack",&tcExtTrk);
339 if ( !checkAlgo ){
340 if ( verbose ) printf(" Recovered track algorithm branch found! :D \n");
341 hasExtTrk = true;
342 } else {
343 if ( verbose ) printf(" Recovered track algorithm branch not found :( !\n");
344 printf(" ok, this is not a problem (it depends on tracker settings) \n");
345 delete tcExtTrk;
346 tcExtTrk = NULL; // 10RED reprocessing bug when missing extratracks
347 }
348 }
349 //
350 // Does it contain the Trigger tree?
351 //
352 trigger = (TTree*)file->Get("Trigger");
353 if ( !trigger ) {
354 if ( verbose ) printf(" TOF - ERROR: no trigger tree\n");
355 code = -302;
356 goto closeandexit;
357 };
358 //
359 // get trigger level2 data pointer
360 //
361 // trigger->SetMaxVirtualSize(2500000000LL); // EM residual Tracker-new tree in level2 files when NEVENTS is big
362 trigger->SetBranchAddress("TrigLevel2",&trg);
363 nevtrgl2 = trigger->GetEntries();
364
365 //
366 // Retrieve GL_RUN variables from the level2 file
367 //
368 tofversion = ToFInfo(false); // we should decide how to handle versioning system
369 //
370 // create an interface to RunInfo called "runinfo"
371 //
372 // ItoRunInfo= interface with RunInfo and GL_RUN
373 runinfo = new ItoRunInfo(file);
374 //
375 // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
376 //
377 sgnl = 0;
378 sgnl = runinfo->Update(run, "TOF",tofversion);
379 if ( sgnl ){
380 if ( verbose ) printf(" TOF - ERROR: RunInfo exited with non-zero status\n");
381 code = sgnl;
382 goto closeandexit;
383 } else {
384 sgnl = 0;
385 };
386 //
387 // number of events in the file BEFORE the first event of our run
388 //
389 nobefrun = runinfo->GetFirstEntry();
390 //
391 // total number of events in the file
392 //
393 totfileentries = runinfo->GetFileEntries();
394 //
395 // first file entry AFTER the last event of our run
396 //
397 noaftrun = runinfo->GetLastEntry() + 1;
398 //
399 // number of run to be processed
400 //
401 numbofrun = runinfo->GetNoRun();
402 totnorun = runinfo->GetRunEntries();
403 //
404 // Try to access the ToF tree in the file, if it exists we are reprocessing data if not we are processing a new run
405 //
406 toftclone = (TTree*)file->Get("ToF");
407 //
408 if ( !toftclone ){
409 //
410 // tree does not exist, we are not reprocessing
411 //
412 reproc = false;
413 if ( run == 0 && verbose ) printf(" ToF - WARNING: you are reprocessing data but ToF tree does not exist!\n");
414 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");
415
416 } else {
417 //
418 // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
419 //
420 // toftclone->SetMaxVirtualSize(2500000000LL); // EM residual Tracker-new tree in level2 files when NEVENTS is big
421 toftclone->SetAutoSave(900000000000000LL);
422 reproc = true;
423 //
424 // update versioning information
425 //
426 if ( verbose ) printf("\n Preparing the pre-processing...\n");
427 //
428 if ( run == 0 || totnorun == 1 ){
429 //
430 // we are reprocessing all the file
431 // 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
432 //
433 reprocall = true;
434 //
435 if ( verbose ) printf("\n ToF - WARNING: Reprocessing all runs\n");
436 //
437 } else {
438 //
439 // 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
440 //
441 reprocall = false;
442 //
443 if ( verbose ) printf("\n ToF - WARNING: Reprocessing run number %u \n",run);
444 //
445 // copying old tree to a new file
446 //
447 gSystem->MakeDirectory(toffolder.str().c_str());
448 myfold = true;
449 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
450 temptof = toftclone->CloneTree(-1,"fast");
451 temptof->SetName("ToF-old");
452 tempfile->Write();
453 tempfile->Close();
454 }
455 //
456 // Delete the old tree from old file and memory
457 //
458 toftclone->Delete("all");
459 //
460 if ( verbose ) printf(" ...done!\n");
461 //
462 };
463 //
464 // create ToF detector tree toft
465 //
466 file->cd();
467 toft = new TTree("ToF-new","PAMELA Level2 ToF data");
468 // toft->SetMaxVirtualSize(2500000000LL); // EM residual Tracker-new tree in level2 files when NEVENTS is big
469 toft->SetAutoSave(900000000000000LL);
470 tof->Set();//ELENA **TEMPORANEO?**
471 toft->Branch("ToFLevel2","ToFLevel2",&tof);
472 //
473 // create new branches for new tracking algorithms
474 //
475 if ( hasNucleiTrk ){
476 ttofNucleiTrk = new TClonesArray("ToFTrkVar",1);
477 toft->Branch("TrackNuclei",&ttofNucleiTrk);
478 }
479 if ( hasExtNucleiTrk ){
480 ttofExtNucleiTrk = new TClonesArray("ToFTrkVar",1);
481 toft->Branch("RecoveredTrackNuclei",&ttofExtNucleiTrk);
482 }
483 if ( hasExtTrk ){
484 ttofExtTrk = new TClonesArray("ToFTrkVar",1);
485 toft->Branch("RecoveredTrack",&ttofExtTrk);
486 }
487
488 //
489 if ( reproc && !reprocall ){
490 //
491 // open new file and retrieve all tree informations
492 //
493 tempfile = new TFile(tempname.str().c_str(),"READ");
494 toftclone = (TTree*)tempfile->Get("ToF-old");
495 // toftclone->SetMaxVirtualSize(2500000000LL); // EM residual Tracker-new tree in level2 files when NEVENTS is big
496 toftclone->SetAutoSave(900000000000000LL);
497 toftclone->SetBranchAddress("ToFLevel2",&tofclone); // EM reprocessing bug fixed
498 //
499 if ( nobefrun > 0 ){
500 if ( verbose ) printf("\n Pre-processing: copying events from the old tree before the processed run\n");
501 if ( verbose ) printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
502 if ( verbose ) printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
503 for (UInt_t j = 0; j < nobefrun; j++){
504 //
505 if ( toftclone->GetEntry(j) <= 0 ) throw -36;
506 //
507 // copy tofclone to tof
508 //
509 memcpy(&tof,&tofclone,sizeof(tofclone));
510 //
511 // Fill entry in the new tree
512 //
513 toft->Fill();
514 //
515 }
516 if ( verbose ) printf(" Finished successful copying!\n");
517 }
518 }
519 //
520 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
521 //
522 runlist = runinfo->GetRunList();
523 //
524 // Loop over the run to be processed
525 //
526 for (UInt_t irun=0; irun < numbofrun; irun++){
527 //
528 // retrieve the first run ID to be processed using the RunInfo list
529 //
530 idRun = runlist->At(irun);
531 if ( verbose ) printf("\n\n\n ####################################################################### \n");
532 if ( verbose ) printf(" PROCESSING RUN NUMBER %u \n",idRun);
533 if ( verbose ) printf(" ####################################################################### \n\n\n");
534 //
535 runinfo->ID_ROOT_L0 = 0;
536 //
537 // store in the runinfo class the GL_RUN variables for our run
538 //
539 sgnl = 0;
540 sgnl = runinfo->GetRunInfo(idRun);
541 if ( sgnl ){
542 if ( verbose ) printf(" TOF - ERROR: RunInfo exited with non-zero status\n");
543 code = sgnl;
544 goto closeandexit;
545 } else {
546 sgnl = 0;
547 };
548 //
549 // now you can access that variables using the RunInfo class this way runinfo->ID_ROOT_L0
550 //
551 if ( runinfo->ID_ROOT_L0 == 0 ){
552 if ( verbose ) printf("\n TOF - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
553 code = -5;
554 goto closeandexit;
555 };
556 //
557 // prepare the timesync for the db
558 //
559 TString host = glt->CGetHost();
560 TString user = glt->CGetUser();
561 TString psw = glt->CGetPsw();
562 TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
563 if ( !dbc->IsConnected() ) throw -314;
564 stringstream myquery;
565 myquery.str("");
566 myquery << "SET time_zone='+0:00';";
567 delete dbc->Query(myquery.str().c_str());
568 delete dbc->Query("SET sql_mode = 'NO_UNSIGNED_SUBTRACTION';");
569 dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
570 //
571 // Search in the DB the path and name of the LEVEL0 file to be processed.
572 //
573 // if ( !dbc->IsConnected() ) throw -314;
574 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
575 //
576 ftmpname.str("");
577 ftmpname << glroot->PATH.Data() << "/";
578 ftmpname << glroot->NAME.Data();
579 fname = ftmpname.str().c_str();
580 //
581 // print out informations
582 //
583 totevent = runinfo->NEVENTS;
584 if ( verbose ) printf("\n LEVEL0 data file: %s \n",fname.Data());
585 if ( verbose ) printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
586 if ( verbose ) printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
587 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);
588 //
589 // if ( !totevent ) goto closeandexit;
590 //
591 // Open Level0 file
592 //
593 if ( l0File ) l0File->Close();
594 l0File = new TFile(fname.Data());
595 if ( !l0File ) {
596 if ( verbose ) printf(" TOF - ERROR: problems opening Level0 file\n");
597 code = -6;
598 goto closeandexit;
599 };
600 l0tr = (TTree*)l0File->Get("Physics");
601 if ( !l0tr ) {
602 if ( verbose ) printf(" TOF - ERROR: no Physics tree in Level0 file\n");
603 l0File->Close();
604 code = -7;
605 goto closeandexit;
606 };
607 l0head = l0tr->GetBranch("Header");
608 if ( !l0head ) {
609 if ( verbose ) printf(" TOF - ERROR: no Header branch in Level0 tree\n");
610 l0File->Close();
611 code = -8;
612 goto closeandexit;
613 };
614 // l0trig = l0tr->GetBranch("Trigger");
615 // if ( !l0trig ) {
616 // if ( verbose ) printf(" TOF - ERROR: no Trigger branch in Level0 tree\n");
617 // l0File->Close();
618 // code = -300;
619 // goto closeandexit;
620 // };
621 l0tof = l0tr->GetBranch("Tof");
622 if ( !l0tof ) {
623 if ( verbose ) printf(" TOF - ERROR: no ToF branch in Level0 tree\n");
624 l0File->Close();
625 code = -303;
626 goto closeandexit;
627 };
628 //
629 // l0tr->SetBranchAddress("Trigger", &trig);
630 l0tr->SetBranchAddress("Tof", &tofEvent);
631 l0tr->SetBranchAddress("Header", &eh);
632 //
633 nevents = l0tof->GetEntries();
634 //
635 if ( nevents < 1 && totevent ) {
636 if ( verbose ) printf(" TOF - ERROR: Level0 file is empty\n\n");
637 l0File->Close();
638 code = -11;
639 goto closeandexit;
640 };
641 //
642 if ( runinfo->EV_TO > nevents-1 && totevent ) {
643 if ( verbose ) printf(" TOF - ERROR: too few entries in the registry tree\n");
644 l0File->Close();
645 code = -12;
646 goto closeandexit;
647 };
648 //
649 // Check if we have to load parameter files (or calibration associated to runs and not to events)
650 //
651 // for example let's assume that we could have different magnetic field maps for different runs:
652 //
653 if ( !l1only ){
654 if ( trkpar1 || ( tttrkpar1 != 0 && tttrkpar1 < runinfo->RUNHEADER_TIME ) ){
655 trkpar1 = false;
656 // read from DB infos about Magnetic filed maps
657 // if ( !dbc->IsConnected() ) throw -314;
658 glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,1,dbc); // parameters stored in DB in GL_PRAM table
659 tttrkpar1 = glparam->TO_TIME;
660 // ----------------------------
661 // Read the magnetic field
662 // ----------------------------
663 if ( verbose ) printf(" Reading magnetic field maps: \n");
664 trk->LoadField(glparam->PATH+glparam->NAME);
665 if ( verbose ) printf("\n");
666 }
667 }
668 //
669 // variable to save information about the tof calibration used
670 //
671 Bool_t defcal = true;
672 //
673 if ( tofpar1 || ( tttofpar1 != 0 && tttofpar1 < runinfo->RUNHEADER_TIME ) ){
674 tofpar1 = false;
675 //
676 // if ( !dbc->IsConnected() ) throw -314;
677 Int_t error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,201,dbc); // parameters stored in DB in GL_PRAM table
678 if ( error<0 ) {
679 code = error;
680 goto closeandexit;
681 };
682 //
683 if ( verbose ) printf(" Reading ToF parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
684 //
685 if ( (UInt_t)glparam->TO_TIME != (UInt_t)4294967295UL ) defcal = false;
686 //
687 tttofpar1 = glparam->TO_TIME;
688 Int_t nlen = (Int_t)(glparam->PATH+glparam->NAME).Length();
689 // rdtofcal((char *)(glparam->PATH+glparam->NAME).Data(),&nlen);
690 rdtofcal((const char *)(glparam->PATH+glparam->NAME).Data(),&nlen);
691 //
692 };
693 //
694 Int_t error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,204,dbc); // parameters stored in DB in GL_PRAM table
695 if ( error<0 ) {
696 code = error;
697 goto closeandexit;
698 };
699 //
700 if ( verbose ) printf(" Reading ToF attenuation parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
701 tofdedx->ReadParAtt((glparam->PATH+glparam->NAME).Data());
702
703 //
704 error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,205,dbc); // parameters stored in DB in GL_PRAM table
705 if ( error<0 ) {
706 code = error;
707 goto closeandexit;
708 };
709 //
710 if ( verbose ) printf(" Reading ToF desaturation on position parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
711 tofdedx->ReadParPos((glparam->PATH+glparam->NAME).Data());
712
713 //
714 error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,206,dbc); // parameters stored in DB in GL_PRAM table
715 if ( error<0 ) {
716 code = error;
717 goto closeandexit;
718 };
719 //
720 if ( verbose ) printf(" Reading ToF BetheBloch parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
721 tofdedx->ReadParBBneg((glparam->PATH+glparam->NAME).Data());
722
723 //
724 error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,207,dbc); // parameters stored in DB in GL_PRAM table
725 if ( error<0 ) {
726 code = error;
727 goto closeandexit;
728 };
729 //
730 if ( verbose ) printf(" Reading ToF Bethe-Bloch parameter file for beta gt1: %s \n",(glparam->PATH+glparam->NAME).Data());
731 tofdedx->ReadParBBpos((glparam->PATH+glparam->NAME).Data());
732
733 //
734 error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,208,dbc); // parameters stored in DB in GL_PRAM table
735 if ( error<0 ) {
736 code = error;
737 goto closeandexit;
738 };
739 //
740 if ( verbose ) printf(" Reading ToF desaturation on beta parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
741 tofdedx->ReadParDesatBB((glparam->PATH+glparam->NAME).Data());
742
743
744 tofdedx->CheckConnectors(runinfo->RUNHEADER_TIME,glparam,dbc);
745
746 //
747 // WM reading parameter file for dE/dx II order corrections
748 //
749 //memset(dedx_corr_m,0,100*48*sizeof(Float_t));
750 //memset(dedx_corr,0,48*sizeof(Float_t));
751 //memset(mtime,0,100*sizeof(Double_t));
752
753 //
754 // RC reading parameter file for new dE/dx II order correction (10th red)
755 //
756 memset(dedx_Hepeak_m,0,100*48*sizeof(Float_t));
757 memset(dedx_ppeak_m,0,100*48*sizeof(Float_t));
758 memset(dedx_Hepeak,0,48*sizeof(Float_t));
759 memset(dedx_ppeak,0,48*sizeof(Float_t));
760 memset(mtime,0,100*sizeof(Double_t));
761
762 //
763 // Query the DB to get the file
764 //
765 error=glparam->Query_GL_PARAM(runinfo->RUNHEADER_TIME,203,dbc); // parameters stored in DB in GL_PRAM table
766 if ( error<0 ) {
767 code = error;
768 goto closeandexit;
769 };
770 //
771 if ( verbose ) printf(" Reading ToF dE/dx II order correction parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
772 //
773 ical=0; // counter set to zero if first-time reading
774 //-----------------------------------------------------------
775 // Here I read the dEdx_korr parameters
776 //-----------------------------------------------------------
777 jj=0;
778 ifstream fin((glparam->PATH+glparam->NAME).Data());
779 UInt_t window = 200000;
780 Bool_t first = true;
781 Bool_t last = true;
782 //Float_t sdedx_corr_m[48];
783 //memset(sdedx_corr_m,0,48*sizeof(Float_t));
784
785 Float_t sdedx_Hepeak_m[48];
786 memset(sdedx_Hepeak_m,0,48*sizeof(Float_t));
787 Float_t sdedx_ppeak_m[48];
788 memset(sdedx_ppeak_m,0,48*sizeof(Float_t));
789
790 Double_t stm = 0;
791 while ( !fin.eof() ){
792 stm = tm;
793 // if ( jj > 0 ) memcpy(sdedx_corr_m,dedx_corr_m[jj-1],48*sizeof(Float_t)); // BUG sdedx should be the previous in time not the previous saved [absurd dE/dx for 8th reduction March and > March 2008 data - fixed on 2009/02/04
794 fin>>t1>>tm>>t2;
795 if ( debug ) cout << setiosflags(ios::fixed) << setw(10) << setprecision(3) << tm << endl;
796 if ( (tm >= (runinfo->RUNHEADER_TIME-window) && tm <= (runinfo->RUNTRAILER_TIME+window)) || (tm > (runinfo->RUNTRAILER_TIME+window) && last) ){
797 if ( first ){
798 mtime[jj]=stm;
799 jj++;
800 if ( jj >= 100 ){
801 code = -318;
802 goto closeandexit;
803 };
804 };
805 mtime[jj]=tm;
806 };
807 for (ii=0; ii<48;ii++){
808 fin>>wj>>xmean1>>xwidth1;
809 if ( (tm >= (runinfo->RUNHEADER_TIME-window) && tm <= (runinfo->RUNTRAILER_TIME+window)) || (tm > (runinfo->RUNTRAILER_TIME+window) && last) ){
810 if ( first ){
811 //dedx_corr_m[jj-1][ii]=sdedx_corr_m[ii];
812
813 dedx_Hepeak_m[jj-1][ii]=sdedx_Hepeak_m[ii];
814 dedx_ppeak_m[jj-1][ii]=sdedx_ppeak_m[ii];
815 };
816 //dedx_corr_m[jj][ii]=xmean1;
817
818 dedx_Hepeak_m[jj][ii]=xmean1;
819 dedx_ppeak_m[jj][ii]=xwidth1;
820 };
821 //sdedx_corr_m[ii]=xmean1; // BUG sdedx should be the previous in time not the previous saved [absurd dE/dx for 8th reduction March and > March 2008 data - fixed on 2009/02/04
822 sdedx_Hepeak_m[ii]=xmean1;
823 sdedx_ppeak_m[ii]=xwidth1;
824
825 };
826 if ( (tm >= (runinfo->RUNHEADER_TIME-window) && tm <= (runinfo->RUNTRAILER_TIME+window)) || (tm > (runinfo->RUNTRAILER_TIME+window) && last)){
827 if ( first ) first = false;
828 if ( tm > (runinfo->RUNTRAILER_TIME+window) ) last = false;
829 jj++;
830 };
831 if ( jj >= 100 ){
832 code = -318;
833 goto closeandexit;
834 };
835 };
836 //
837 fin.close();
838 // this is a possible bug...
839 // Bool_t ff = false;
840 // while ( runinfo->RUNHEADER_TIME > mtime[ical] && ical < 100 ) {
841 // ical = ical+1;
842 // ff = true;
843 // };
844 while ( (mtime[ical] > runinfo->RUNHEADER_TIME || runinfo->RUNHEADER_TIME > mtime[ical+1] ) && ical < 99 ) {
845 ical = ical+1;
846 // ff = true;
847 };
848 // if ( ff ) ical = ical-1;
849 if ( verbose ) cout<<"rh time "<<runinfo->RUNHEADER_TIME<<" rt time "<<runinfo->RUNTRAILER_TIME<<" limit low "<<mtime[ical]<<" limit up "<<mtime[ical+1]<<" ical "<<ical<< " jj " << jj<< endl;
850 if ( ical < 0 || ical >= 98 ){
851 code = -315;
852 goto closeandexit;
853 };
854 //
855 // run over all the events of the run
856 //
857 if ( verbose ) printf("\n Ready to start! \n\n Processed events: \n\n");
858 //
859 if ( dbc ){
860 dbc->Close();
861 delete dbc;
862 dbc = 0;
863 };
864 //
865 jumped = 0;
866 //
867 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
868 // for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+100); re++){ // QUIIIIIII
869 //
870 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
871 //
872 if ( l0head->GetEntry(re) <= 0 ) throw -36;
873 //
874 // absolute time of this event
875 //
876 ph = eh->GetPscuHeader();
877 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
878 //
879 tof->Clear();
880 Int_t pmt_id = 0;
881 ToFPMT *t_pmt = new ToFPMT();
882 if(!(tof->PMT))tof->PMT = new TClonesArray("ToFPMT",12); //ELENA
883 TClonesArray &tpmt = *tof->PMT;
884 ToFTrkVar *t_tof = new ToFTrkVar();
885 if(!(tof->ToFTrk))tof->ToFTrk = new TClonesArray("ToFTrkVar",2); //ELENA
886 TClonesArray &t = *tof->ToFTrk;
887 //
888 // paranoid check
889 //
890 if ( atime > (runinfo->RUNTRAILER_TIME+1) || atime < (runinfo->RUNHEADER_TIME-1) ) {
891 if ( verbose ) printf(" TOF - WARNING: event at time outside the run time window, skipping it\n");
892 jumped++;
893 goto jumpev;
894 };
895 //
896 // retrieve tracker informations, the LEVEL2 entry which correspond to our event will be "itr"
897 //
898 if ( !reprocall ){
899 itr = nobefrun + (re - runinfo->EV_FROM -jumped);
900 } else {
901 itr = runinfo->GetFirstEntry() + (re - runinfo->EV_FROM -jumped);
902 };
903 if ( !l1only ){
904 if ( itr > nevtrkl2 ){ // nevtrkl2 tracker entry number
905 if ( verbose ) printf(" TOF - ERROR: no tracker events with entry = %i in Level2 file\n",itr);
906 l0File->Close();
907 code = -313;
908 goto closeandexit;
909 }
910 }
911 if ( itr > nevtrgl2 ){ // nevtrgl2 trigger entry number
912 if ( verbose ) printf(" TOF - ERROR: no trigger events with entry = %i in Level2 file\n",itr);
913 l0File->Close();
914 code = -319;
915 goto closeandexit;
916 }
917 //
918 if ( !l1only ){
919 trk->Clear();
920 //
921 // Clones array must be cleared before going on
922 //
923 if ( hasNucleiTrk ){
924 tcNucleiTrk->Delete();
925 ttofNucleiTrk->Delete();
926 }
927 if ( hasExtNucleiTrk ){
928 tcExtNucleiTrk->Delete();
929 ttofExtNucleiTrk->Delete();
930 }
931 if ( hasExtTrk ){
932 tcExtTrk->Delete();
933 ttofExtTrk->Delete();
934 }
935 }
936 trg->Clear();
937 //
938 if ( !l1only && tracker->GetEntry(itr) <= 0 ) throw -36;
939 if ( trigger->GetEntry(itr) <= 0 ) throw -36;
940 ///
941 //
942 if ( l0tof->GetEntry(re) <= 0 ) throw -36;
943 // if ( l0trig->GetEntry(re) <= 0 ) throw -36;
944 ///
945 //
946 procev++;
947 //
948 // start processing
949 //
950 // dE/dx II order correction: check time limits and interpolate time correction
951 //==================================================================
952 //== if time is outside time limits:
953 //==================================================================
954 if ( atime<mtime[ical] || atime>mtime[ical+1] ){
955 if ( verbose ) cout<<"Checking Time Limits!"<<endl;
956 ical=0;
957 while ( ( mtime[ical] > atime || atime > mtime[ical+1]) && ical < 99 ){
958 ical = ical+1;
959 }
960 //
961 if ( ical < 0 || ical >= 98 ){
962 code = -317;
963 goto closeandexit;
964 };
965 if ( verbose ) cout<<"abs time "<<atime<<" limit low "<<mtime[ical]<<" limit up "<<mtime[ical+1]<<" ical "<<ical<<endl;
966 };
967 //==================================================================
968 //== interpolate betwen time limits
969 //==================================================================
970 thelp1 = mtime[ical];
971 thelp2 = mtime[ical+1];
972 for (ii=0; ii<48;ii++) {
973
974 Heyhelp1 = fabs(dedx_Hepeak_m[ical][ii]);
975 if ( Heyhelp1 < 0.1 ) Heyhelp1 = 4.;
976 Heyhelp2 = fabs(dedx_Hepeak_m[ical+1][ii]);
977 if ( Heyhelp2 < 0.1 ) Heyhelp2 = 4.;
978 Heslope = (Heyhelp2-Heyhelp1)/(thelp2-thelp1);
979 Heinter = Heyhelp1 - Heslope*thelp1;
980 dedx_Hepeak[ii] = Heslope*atime + Heinter;
981
982 pyhelp1 = fabs(dedx_ppeak_m[ical][ii]);
983 if ( pyhelp1 < 0.1 ) pyhelp1 = 1.;
984 pyhelp2 = fabs(dedx_ppeak_m[ical+1][ii]);
985 if ( pyhelp2 < 0.1 ) pyhelp2 = 1.;
986 pslope = (pyhelp2-pyhelp1)/(thelp2-thelp1);
987 pinter = pyhelp1 - pslope*thelp1;
988 dedx_ppeak[ii] = pslope*atime + pinter;
989
990 if(dedx_Hepeak[ii]>dedx_ppeak[ii])slope_dedx[ii]=3./(dedx_Hepeak[ii]-dedx_ppeak[ii]);
991 else slope_dedx[ii]=4.;
992 if(dedx_Hepeak[ii]>dedx_ppeak[ii])inter_dedx[ii]=1.-(slope_dedx[ii]*dedx_ppeak[ii]);
993 else inter_dedx[ii]=0.;
994
995 if ( fabs(dedx_ppeak[ii]) <= 1e-15 ){
996 if ( verbose ) printf("ii %i pslope %f atime %u pinter %f dedx_ppeak %f \n",ii,pslope,atime,pinter,dedx_ppeak[ii]);
997 if ( verbose ) printf("ical %i pyhelp2 %f pyhelp1 %f thelp2 %f thelp1 %f \n",ical,pyhelp2,pyhelp1,thelp2,thelp1);
998 code = -316;
999 goto closeandexit;
1000 }
1001 }
1002 //================================================================
1003 //================================================================
1004
1005 //
1006 // Here we will use some procedure to calibrate our data and put some kind of informations in the cinput structure
1007 //
1008 for (Int_t gg=0; gg<4;gg++){
1009 for (Int_t hh=0; hh<12;hh++){
1010 tofinput_.tdc[hh][gg] = (0xFFF & tofEvent->tdc[gg][hh]); // exclude warning bits
1011 tofinput_.adc[hh][gg] = (0xFFF & tofEvent->adc[gg][hh]); // exclude warning bits
1012 }
1013 }
1014 //
1015 tofdedx->Init(tofEvent);
1016 warning = 0;
1017 //
1018 for (Int_t hh=0; hh<5;hh++){
1019 tofinput_.patterntrig[hh]=trg->patterntrig[hh];
1020 }
1021 //
1022 // Here we have calibrated data, ready to be passed to the FORTRAN routine which will extract common and track-related variables.
1023 //
1024 npmtentry = 0;
1025 //
1026 ntrkentry = 0;
1027 //
1028 // Calculate tracks informations from ToF alone
1029 //
1030 tofl2com();
1031 //
1032 memcpy(tof->tof_j_flag,tofoutput_.tof_j_flag,6*sizeof(Int_t));
1033 //
1034 if ( !l1only ){
1035 //
1036 t_tof->trkseqno = -1;
1037 //
1038 // and now we must copy from the output structure to the level2 class:
1039 //
1040 t_tof->npmttdc = 0;
1041 //
1042 for (Int_t hh=0; hh<12;hh++){
1043 for (Int_t kk=0; kk<4;kk++){
1044 if ( tofoutput_.tofmask[hh][kk] != 0 ){
1045 pmt_id = tof->GetPMTid(kk,hh);
1046 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
1047 t_tof->tdcflag.AddAt(tofoutput_.tdcflagtof[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
1048 t_tof->npmttdc++;
1049 };
1050 };
1051 };
1052 for (Int_t kk=0; kk<13;kk++){
1053 t_tof->beta[kk] = tofoutput_.betatof_a[kk];
1054 }
1055 //
1056 //
1057 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
1058 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
1059 memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
1060 memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
1061 //
1062 {
1063 Float_t xtof_temp[6]={100.,100.,100.,100.,100.,100.};
1064 Float_t ytof_temp[6]={100.,100.,100.,100.,100.,100.};
1065
1066 if(t_tof->xtofpos[0]<100. && t_tof->ytofpos[0]<100.){
1067 xtof_temp[1]=t_tof->xtofpos[0];
1068 ytof_temp[0]=t_tof->ytofpos[0];
1069 }else if(t_tof->xtofpos[0]>=100. && t_tof->ytofpos[0]<100.){
1070 ytof_temp[0]=t_tof->ytofpos[0];
1071 tof->GetPaddleGeometry(0,(Int_t)log2(tof->tof_j_flag[0]),xleft, xright, yleft, yright);
1072 xtof_temp[1]=xleft+2.55;
1073 }else if(t_tof->ytofpos[0]>=100. && t_tof->xtofpos[0]<100.){
1074 xtof_temp[1]=t_tof->xtofpos[0];
1075 tof->GetPaddleGeometry(1,(Int_t)log2(tof->tof_j_flag[1]),xleft, xright, yleft, yright);
1076 ytof_temp[0]=yleft+2.75;
1077 }
1078
1079 if(t_tof->xtofpos[1]<100. && t_tof->ytofpos[1]<100.){
1080 xtof_temp[2]=t_tof->xtofpos[1];
1081 ytof_temp[3]=t_tof->ytofpos[1];
1082 }else if(t_tof->xtofpos[1]>=100. && t_tof->ytofpos[1]<100.){
1083 ytof_temp[3]=t_tof->ytofpos[1];
1084 tof->GetPaddleGeometry(3,(Int_t)log2(tof->tof_j_flag[3]),xleft, xright, yleft, yright);
1085 xtof_temp[2]=xleft+4.5;
1086 }else if(t_tof->ytofpos[1]>=100. && t_tof->xtofpos[1]<100.){
1087 xtof_temp[2]=t_tof->xtofpos[1];
1088 tof->GetPaddleGeometry(2,(Int_t)log2(tof->tof_j_flag[2]),xleft, xright, yleft, yright);
1089 ytof_temp[3]=yleft+3.75;
1090 }
1091
1092 if(t_tof->xtofpos[2]<100. && t_tof->ytofpos[2]<100.){
1093 xtof_temp[5]=t_tof->xtofpos[2];
1094 ytof_temp[4]=t_tof->ytofpos[2];
1095 }else if(t_tof->xtofpos[2]>=100. && t_tof->ytofpos[2]<100.){
1096 ytof_temp[4]=t_tof->ytofpos[2];
1097 tof->GetPaddleGeometry(4,(Int_t)log2(tof->tof_j_flag[4]),xleft, xright, yleft, yright);
1098 xtof_temp[5]=xleft+3;
1099 }else if(t_tof->ytofpos[2]>=100. && t_tof->xtofpos[2]<100.){
1100 xtof_temp[5]=t_tof->xtofpos[2];
1101 tof->GetPaddleGeometry(5,(Int_t)log2(tof->tof_j_flag[5]),xleft, xright, yleft, yright);
1102 ytof_temp[4]=yleft+2.5;
1103 }
1104 //
1105 tofdedx->Process(atime,t_tof->beta[12], (Float_t *)xtof_temp,(Float_t *)ytof_temp);
1106 t_tof->npmtadc = 0;
1107 for (Int_t hh=0; hh<12;hh++){
1108 for (Int_t kk=0; kk<4;kk++){
1109 pmt_id = tof->GetPMTid(kk,hh);
1110 Int_t Iplane=-1;
1111 Int_t Ipaddle=-1;
1112 tof->GetPMTPaddle(pmt_id, Iplane, Ipaddle);
1113 tof->GetPaddleGeometry(Iplane,Ipaddle,xleft,xright,yleft,yright);
1114 if (tofEvent->tdc[kk][hh] < 4095 || tofEvent->adc[kk][hh] < 4095 || tofinput_.tdc[hh][kk] < 4095 || tofinput_.adc[hh][kk] < 4095 ) {
1115 if ( tofdedx->GetdEdx_pmt(pmt_id)>-1. && (inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)) > 0. &&((xtof_temp[Iplane]>=xleft&&xtof_temp[Iplane]<=xright) || (ytof_temp[Iplane]>=yleft&&ytof_temp[Iplane]<=yright)) ){
1116
1117 t_tof->dedx.AddAt((inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),t_tof->npmtadc);// RC new dE/dx II order correction
1118
1119 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
1120 t_tof->adcflag.AddAt(0,t_tof->npmtadc); // gf: Jan 09/07
1121 t_tof->npmtadc++;
1122 }
1123 }
1124 }
1125 }
1126 }
1127 new(t[ntrkentry]) ToFTrkVar(*t_tof);
1128 ntrkentry++;
1129 t_tof->Clear();
1130 t_pmt->Clear();
1131 //
1132 for (Int_t gg=0; gg<4;gg++){
1133 for (Int_t hh=0; hh<12;hh++){
1134 // new WM
1135 if ( tofoutput_.tdc_c[hh][gg] < 4095 || (0xFFF & tofEvent->adc[gg][hh]) < 4095 || (0xFFF & tofEvent->tdc[gg][hh]) < 4095 ){
1136 t_pmt->pmt_id = tof->GetPMTid(gg,hh);
1137 t_pmt->tdc_tw = tofoutput_.tdc_c[hh][gg];
1138 t_pmt->adc = (Float_t)(0xFFF & tofEvent->adc[gg][hh]);
1139 t_pmt->tdc = (Float_t)(0xFFF & tofEvent->tdc[gg][hh]);
1140 t_pmt->l0flag_adc = (Float_t)(tofEvent->adc[gg][hh]>>12);
1141 t_pmt->l0flag_tdc = (Float_t)(tofEvent->tdc[gg][hh]>>12);
1142 if ( t_pmt->l0flag_adc || t_pmt->l0flag_tdc ) warning |= 1 << 0;
1143 //
1144 new(tpmt[npmtentry]) ToFPMT(*t_pmt);
1145 npmtentry++;
1146 t_pmt->Clear();
1147 }
1148 }
1149 }
1150 //
1151 if ( debug ) printf(" ATIME %u re %u \n",atime,(UInt_t)re);
1152 //
1153 // Calculate track-related variables
1154 //
1155
1156 //
1157 // Run over tracks - standard algorithm
1158 //
1159 for(Int_t nt=0; nt < trk->ntrk(); nt++){
1160 //
1161 TrkTrack *ptt = trk->GetStoredTrack(nt);
1162 //
1163 // Copy the alpha vector in the input structure
1164 //
1165 for (Int_t e = 0; e < 5 ; e++){
1166 tofinput_.al_pp[e] = ptt->al[e];
1167 }
1168 // new input for 9th reduction: tracker dEdx
1169 tofinput_.trkmip = ptt->GetDEDX();
1170 //
1171 // Get tracker related variables for this track
1172 //
1173 toftrk();
1174 //
1175 // Copy values in the class from the structure (we need to use a temporary class to store variables).
1176 //
1177 t_tof->npmttdc = 0;
1178 for (Int_t hh=0; hh<12;hh++){
1179 for (Int_t kk=0; kk<4;kk++){
1180 if ( tofoutput_.tofmask[hh][kk] != 0 ){
1181 pmt_id = tof->GetPMTid(kk,hh);
1182 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
1183 t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
1184 t_tof->npmttdc++;
1185 }
1186 }
1187 }
1188 for (Int_t kk=0; kk<13;kk++){
1189 t_tof->beta[kk] = tofoutput_.beta_a[kk];
1190 }
1191 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
1192 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
1193 memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
1194 memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
1195 //
1196 tofdedx->Process(atime,t_tof->beta[12], (Float_t *)t_tof->xtr_tof,(Float_t *)t_tof->ytr_tof);
1197 t_tof->npmtadc = 0;
1198 for (Int_t hh=0; hh<12;hh++){
1199 for (Int_t kk=0; kk<4;kk++){
1200 pmt_id = tof->GetPMTid(kk,hh);
1201 Int_t Iplane=-1;
1202 Int_t Ipaddle=-1;
1203 Int_t IpaddleT=-1;
1204 tof->GetPMTPaddle(pmt_id, Iplane, Ipaddle);
1205 IpaddleT=tof->GetPaddleIdOfTrack(t_tof->xtr_tof[Iplane],t_tof->ytr_tof[Iplane], Iplane,0.0);
1206 if ( debug ) printf(" 1nt %i pmt_id %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,pmt_id,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1207 if (tofEvent->tdc[kk][hh] < 4095 || tofEvent->adc[kk][hh] < 4095 || tofinput_.tdc[hh][kk] < 4095 || tofinput_.adc[hh][kk] < 4095 ) {
1208 if ( tofdedx->GetdEdx_pmt(pmt_id) > -1. && (inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)) > 0. && Ipaddle==IpaddleT ){
1209 t_tof->dedx.AddAt((inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),t_tof->npmtadc);// RC new dE/dx II order correction
1210 if ( debug ) printf(" 2nt %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1211 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
1212 t_tof->adcflag.AddAt(0,t_tof->npmtadc); // gf: Jan 09/07
1213 t_tof->npmtadc++;
1214 }
1215 }
1216 }
1217 }
1218 //
1219 // Store the tracker track number in order to be sure to have shyncronized data during analysis
1220 //
1221 t_tof->trkseqno = nt;
1222 //
1223 // create a new object for this event with track-related variables
1224 //
1225 new(t[ntrkentry]) ToFTrkVar(*t_tof);
1226 ntrkentry++;
1227 t_tof->Clear();
1228 //
1229 } // loop on all the tracks
1230 //
1231 // Code for extended tracking algorithm:
1232 //
1233 //
1234 // Run over tracks - nuclei algorithm
1235 //
1236 if ( hasNucleiTrk ){
1237 Int_t ttentry = 0;
1238 for(Int_t nt=0; nt < tcNucleiTrk->GetEntries(); nt++){
1239 //
1240 TrkTrack *ptt = (TrkTrack*)(tcNucleiTrk->At(nt));
1241 //
1242 // Copy the alpha vector in the input structure
1243 //
1244 for (Int_t e = 0; e < 5 ; e++){
1245 tofinput_.al_pp[e] = ptt->al[e];
1246 }
1247 // new input for 9th reduction: tracker dEdx
1248 tofinput_.trkmip = ptt->GetDEDX();
1249 //
1250 // Get tracker related variables for this track
1251 //
1252 toftrk();
1253 //
1254 // Copy values in the class from the structure (we need to use a temporary class to store variables).
1255 //
1256 t_tof->npmttdc = 0;
1257 for (Int_t hh=0; hh<12;hh++){
1258 for (Int_t kk=0; kk<4;kk++){
1259 if ( tofoutput_.tofmask[hh][kk] != 0 ){
1260 pmt_id = tof->GetPMTid(kk,hh);
1261 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
1262 t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
1263 t_tof->npmttdc++;
1264 }
1265 }
1266 }
1267 for (Int_t kk=0; kk<13;kk++){
1268 t_tof->beta[kk] = tofoutput_.beta_a[kk];
1269 }
1270 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
1271 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
1272 memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
1273 memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
1274 //
1275 tofdedx->Process(atime,t_tof->beta[12], (Float_t *)t_tof->xtr_tof,(Float_t *)t_tof->ytr_tof);
1276 t_tof->npmtadc = 0;
1277 for (Int_t hh=0; hh<12;hh++){
1278 for (Int_t kk=0; kk<4;kk++){
1279 pmt_id = tof->GetPMTid(kk,hh);
1280 Int_t Iplane=-1;
1281 Int_t Ipaddle=-1;
1282 Int_t IpaddleT=-1;
1283 tof->GetPMTPaddle(pmt_id, Iplane, Ipaddle);
1284 IpaddleT=tof->GetPaddleIdOfTrack(t_tof->xtr_tof[Iplane],t_tof->ytr_tof[Iplane], Iplane,0.0);
1285 if ( debug ) printf(" 1nt %i pmt_id %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,pmt_id,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1286 if (tofEvent->tdc[kk][hh] < 4095 || tofEvent->adc[kk][hh] < 4095 || tofinput_.tdc[hh][kk] < 4095 || tofinput_.adc[hh][kk] < 4095 ) {
1287 if ( tofdedx->GetdEdx_pmt(pmt_id) > -1. && (inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)) > 0. && Ipaddle==IpaddleT ){
1288 t_tof->dedx.AddAt((inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),t_tof->npmtadc);// RC new dE/dx II order correction
1289 if ( debug ) printf(" 2nt %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1290 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
1291 t_tof->adcflag.AddAt(0,t_tof->npmtadc); // gf: Jan 09/07
1292 t_tof->npmtadc++;
1293 }
1294 }
1295 }
1296 }
1297 //
1298 // Store the tracker track number in order to be sure to have shyncronized data during analysis
1299 //
1300 t_tof->trkseqno = nt;
1301 //
1302 // create a new object for this event with track-related variables
1303 //
1304 TClonesArray &tt1 = *ttofNucleiTrk;
1305 new(tt1[ttentry]) ToFTrkVar(*t_tof);
1306 ttentry++;
1307 t_tof->Clear();
1308 //
1309 } // loop on all the tracks, nuclei algorithm
1310 }
1311 //
1312 // Run over tracks - extended nuclei algorithm
1313 //
1314 if ( hasExtNucleiTrk ){
1315 Int_t ttentry = 0;
1316 for(Int_t nt=0; nt < tcExtNucleiTrk->GetEntries(); nt++){
1317 //
1318 ExtTrack *ptt = (ExtTrack*)(tcExtNucleiTrk->At(nt));
1319 //
1320 // Copy the alpha vector in the input structure
1321 //
1322 for (Int_t e = 0; e < 5 ; e++){
1323 tofinput_.al_pp[e] = ptt->al[e];
1324 }
1325 // new input for 9th reduction: tracker dEdx
1326 tofinput_.trkmip = ptt->GetDEDX();
1327 //
1328 // Get tracker related variables for this track
1329 //
1330 toftrk();
1331 //
1332 // Copy values in the class from the structure (we need to use a temporary class to store variables).
1333 //
1334 t_tof->npmttdc = 0;
1335 for (Int_t hh=0; hh<12;hh++){
1336 for (Int_t kk=0; kk<4;kk++){
1337 if ( tofoutput_.tofmask[hh][kk] != 0 ){
1338 pmt_id = tof->GetPMTid(kk,hh);
1339 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
1340 t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
1341 t_tof->npmttdc++;
1342 }
1343 }
1344 }
1345 for (Int_t kk=0; kk<13;kk++){
1346 t_tof->beta[kk] = tofoutput_.beta_a[kk];
1347 }
1348 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
1349 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
1350 memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
1351 memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
1352 //
1353 tofdedx->Process(atime,t_tof->beta[12], (Float_t *)t_tof->xtr_tof,(Float_t *)t_tof->ytr_tof);
1354 t_tof->npmtadc = 0;
1355 for (Int_t hh=0; hh<12;hh++){
1356 for (Int_t kk=0; kk<4;kk++){
1357 pmt_id = tof->GetPMTid(kk,hh);
1358 Int_t Iplane=-1;
1359 Int_t Ipaddle=-1;
1360 Int_t IpaddleT=-1;
1361 tof->GetPMTPaddle(pmt_id, Iplane, Ipaddle);
1362 IpaddleT=tof->GetPaddleIdOfTrack(t_tof->xtr_tof[Iplane],t_tof->ytr_tof[Iplane], Iplane,0.0);
1363 if ( debug ) printf(" 1nt %i pmt_id %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,pmt_id,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1364 if (tofEvent->tdc[kk][hh] < 4095 || tofEvent->adc[kk][hh] < 4095 || tofinput_.tdc[hh][kk] < 4095 || tofinput_.adc[hh][kk] < 4095 ) {
1365 if ( tofdedx->GetdEdx_pmt(pmt_id) > -1. && (inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)) > 0. && Ipaddle==IpaddleT ){
1366 t_tof->dedx.AddAt((inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),t_tof->npmtadc);// RC new dE/dx II order correction
1367 if ( debug ) printf(" 2nt %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1368 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
1369 t_tof->adcflag.AddAt(0,t_tof->npmtadc); // gf: Jan 09/07
1370 t_tof->npmtadc++;
1371 }
1372 }
1373 }
1374 }
1375 //
1376 // Store the tracker track number in order to be sure to have shyncronized data during analysis
1377 //
1378 t_tof->trkseqno = nt;
1379 //
1380 // create a new object for this event with track-related variables
1381 //
1382 TClonesArray &tt2 = *ttofExtNucleiTrk;
1383 new(tt2[ttentry]) ToFTrkVar(*t_tof);
1384 ttentry++;
1385 t_tof->Clear();
1386 //
1387 } // loop on all the tracks, extended nuclei algorithm
1388 }
1389 //
1390 // Run over tracks - extended algorithm
1391 //
1392 if ( hasExtTrk ){
1393 Int_t ttentry = 0;
1394 for(Int_t nt=0; nt < tcExtTrk->GetEntries(); nt++){
1395 //
1396 ExtTrack *ptt = (ExtTrack*)(tcExtTrk->At(nt));
1397 //
1398 // Copy the alpha vector in the input structure
1399 //
1400 for (Int_t e = 0; e < 5 ; e++){
1401 tofinput_.al_pp[e] = ptt->al[e];
1402 }
1403 // new input for 9th reduction: tracker dEdx
1404 tofinput_.trkmip = ptt->GetDEDX();
1405 //
1406 // Get tracker related variables for this track
1407 //
1408 toftrk();
1409 //
1410 // Copy values in the class from the structure (we need to use a temporary class to store variables).
1411 //
1412 t_tof->npmttdc = 0;
1413 for (Int_t hh=0; hh<12;hh++){
1414 for (Int_t kk=0; kk<4;kk++){
1415 if ( tofoutput_.tofmask[hh][kk] != 0 ){
1416 pmt_id = tof->GetPMTid(kk,hh);
1417 t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
1418 t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
1419 t_tof->npmttdc++;
1420 }
1421 }
1422 }
1423 for (Int_t kk=0; kk<13;kk++){
1424 t_tof->beta[kk] = tofoutput_.beta_a[kk];
1425 }
1426 memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
1427 memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
1428 memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
1429 memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
1430 //
1431 tofdedx->Process(atime,t_tof->beta[12], (Float_t *)t_tof->xtr_tof,(Float_t *)t_tof->ytr_tof);
1432 t_tof->npmtadc = 0;
1433 for (Int_t hh=0; hh<12;hh++){
1434 for (Int_t kk=0; kk<4;kk++){
1435 pmt_id = tof->GetPMTid(kk,hh);
1436 Int_t Iplane=-1;
1437 Int_t Ipaddle=-1;
1438 Int_t IpaddleT=-1;
1439 tof->GetPMTPaddle(pmt_id, Iplane, Ipaddle);
1440 IpaddleT=tof->GetPaddleIdOfTrack(t_tof->xtr_tof[Iplane],t_tof->ytr_tof[Iplane], Iplane,0.0);
1441 if ( debug ) printf(" 1nt %i pmt_id %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,pmt_id,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1442 if (tofEvent->tdc[kk][hh] < 4095 || tofEvent->adc[kk][hh] < 4095 || tofinput_.tdc[hh][kk] < 4095 || tofinput_.adc[hh][kk] < 4095 ) {
1443 if ( tofdedx->GetdEdx_pmt(pmt_id) > -1. && (inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)) > 0. && Ipaddle==IpaddleT ){
1444 t_tof->dedx.AddAt((inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),t_tof->npmtadc);// RC new dE/dx II order correction
1445 if ( debug ) printf(" 2nt %i npmtadc %i dedx %f dedx slope %f dedx inter %f\n",nt,t_tof->npmtadc,(inter_dedx[pmt_id]+slope_dedx[pmt_id]*tofdedx->GetdEdx_pmt(pmt_id)),inter_dedx[pmt_id],slope_dedx[pmt_id]);
1446 t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
1447 t_tof->adcflag.AddAt(0,t_tof->npmtadc); // gf: Jan 09/07
1448 t_tof->npmtadc++;
1449 }
1450 }
1451 }
1452 }
1453 //
1454 // Store the tracker track number in order to be sure to have shyncronized data during analysis
1455 //
1456 t_tof->trkseqno = nt;
1457 //
1458 // create a new object for this event with track-related variables
1459 //
1460 TClonesArray &tt3 = *ttofExtTrk;
1461 new(tt3[ttentry]) ToFTrkVar(*t_tof);
1462 ttentry++;
1463 t_tof->Clear();
1464 //
1465 } // loop on all the tracks, extended algorithm
1466 }
1467
1468 } // if !l1only
1469 //
1470 tof->unpackError = tofEvent->unpackError;
1471
1472 a = 0;
1473 b = 0;
1474 if ( !tof->checkPMTpatternPMThit(trg, a, b) ) warning |= 1 << 1;
1475 if ( !tof->checkPMTpmttrig(trg) ) warning |= 1 << 2;
1476 if ( !trg->checkPMTpatterntrig() ) warning |= 1 << 3;
1477 tof->unpackWarning = warning;
1478
1479 if ( defcal ){
1480 tof->default_calib = 1;
1481 } else {
1482 tof->default_calib = 0;
1483 }
1484 //
1485 // Fill the rootple
1486 //
1487 toft->Fill();
1488 //
1489 //
1490 //
1491 delete t_tof;
1492 //
1493 //
1494 //
1495 jumpev:
1496 if ( !debug ) debug = false;
1497 //
1498 }
1499 //
1500 // Here you may want to clear some variables before processing another run
1501 //
1502 delete dbtime;
1503 } // process all the runs
1504 //
1505 if ( verbose ) printf("\n Finished processing data \n");
1506 //
1507 closeandexit:
1508 //
1509 // 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.
1510 //
1511 if ( !reprocall && reproc && code >= 0 ){
1512 if ( totfileentries > noaftrun ){
1513 if ( verbose ) printf("\n Post-processing: copying events from the old tree after the processed run\n");
1514 if ( verbose ) printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
1515 if ( verbose ) printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
1516 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
1517 //
1518 // Get entry from old tree
1519 //
1520 if ( toftclone->GetEntry(j) <= 0 ) throw -36;
1521 //
1522 // copy tofclone to tof
1523 //
1524 memcpy(&tof,&tofclone,sizeof(tofclone));// EM reprocessing bug fixed
1525 //
1526 // Fill entry in the new tree
1527 //
1528 toft->Fill();// EM reprocessing bug fixed
1529 };
1530 if ( verbose ) printf(" Finished successful copying!\n");
1531 };
1532 };
1533 //
1534 // Close files, delete old tree(s), write and close level2 file
1535 //
1536 if ( l0File ) l0File->Close();
1537 if ( tempfile ) tempfile->Close();
1538 if ( myfold ) gSystem->Unlink(tempname.str().c_str());
1539 //
1540 if ( code < 0 && verbose ) printf("\n TOF - ERROR: an error occurred, try to save anyway...\n");
1541 if ( verbose ) printf("\n Writing and closing rootple\n");
1542 if ( toft ) toft->SetName("ToF");
1543 if ( file ){
1544 file->cd();
1545 if ( toft ) toft->Write(0, TObject::kOverwrite); // 10RED bug fixed
1546 };
1547 //
1548 if ( myfold ) gSystem->Unlink(toffolder.str().c_str());
1549 //
1550 // the end
1551 //
1552 if ( tcNucleiTrk ){
1553 tcNucleiTrk->Delete();
1554 delete tcNucleiTrk;
1555 tcNucleiTrk = NULL;
1556 }
1557 if ( tcExtNucleiTrk ){
1558 tcExtNucleiTrk->Delete();
1559 delete tcExtNucleiTrk;
1560 tcExtNucleiTrk = NULL;
1561 }
1562 if ( tcExtTrk ){
1563 tcExtTrk->Delete();
1564 delete tcExtTrk;
1565 tcExtTrk = NULL;
1566 }
1567 //
1568 if ( verbose ) printf("\n Exiting...\n");
1569 //
1570 if ( tofdedx ) delete tofdedx;
1571 if ( glroot ) delete glroot;
1572 if ( glparam ) delete glparam;
1573 if ( runinfo ) runinfo->Close();
1574 if ( runinfo ) delete runinfo;
1575 //
1576 if ( code < 0 ) throw code;
1577 if ( debug ) file->ls();
1578 return(code);
1579 }

  ViewVC Help
Powered by ViewVC 1.1.23