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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.52 - (hide annotations) (download)
Fri Jun 6 14:18:25 2014 UTC (10 years, 7 months ago) by mocchiut
Branch: MAIN
Changes since 1.51: +361 -147 lines
New extended algorithm track-related added to CALO, TOF and ORB

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

  ViewVC Help
Powered by ViewVC 1.1.23