/[PAMELA software]/DarthVader/TrackerLevel2/src/TrkCore.cpp
ViewVC logotype

Annotation of /DarthVader/TrackerLevel2/src/TrkCore.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.42 - (hide annotations) (download)
Tue Nov 5 13:58:45 2013 UTC (11 years, 3 months ago) by mocchiut
Branch: MAIN
Changes since 1.41: +119 -0 lines
New code for retracking

1 mocchiut 1.1 /**
2     * \file TrkCore.cpp
3     * \author Elena Vannuccini
4     */
5     // .........................................................
6     // ROOT header files
7     // .........................................................
8     #include <TTree.h>
9     #include <TClassEdit.h>
10     #include <TObject.h>
11     #include <TList.h>
12     #include <TSystem.h>
13     #include <TSystemDirectory.h>
14     #include <TString.h>
15     #include <TFile.h>
16     #include <TClass.h>
17     #include <TCanvas.h>
18     #include <TH1.h>
19     #include <TH1F.h>
20     #include <TH2D.h>
21     #include <TLatex.h>
22     #include <TPad.h>
23     #include <TSQLServer.h>
24     #include <TSQLRow.h>
25     #include <TSQLResult.h>
26 pam-fi 1.11 #include <TBenchmark.h>
27    
28 mocchiut 1.1 // .........................................................
29     // other general header files
30     // .........................................................
31     #include <fstream>
32     #include <iostream>
33     // .........................................................
34     // files in the inc directory
35     // .........................................................
36     #include <RunInfo.h>
37     #include <GLTables.h>
38     #include <TrkVerl2.h>
39     #include <TrkProcess.h>
40 pam-fi 1.2 //#include <TrkStruct.h>
41 mocchiut 1.1 #include <TrkLevel2.h>
42 pam-fi 1.6 #include <TrkLevel1.h>
43 mocchiut 1.1 #include <TrkLevel0.h>
44 pam-fi 1.15 #include <TrkHough.h>
45 mocchiut 1.1 // .........................................................
46     // YODA header files
47     // .........................................................
48     #include <PamelaRun.h>
49 pam-fi 1.7 //#include <RegistryEvent.h>
50 mocchiut 1.1 #include <physics/tracker/TrackerEvent.h>
51     #include <CalibTrk1Event.h>
52     #include <CalibTrk2Event.h>
53 mocchiut 1.42 //
54     // Calorimeter level1 classes headers and definitions
55     //
56     #include <CaloLevel1.h> //EM
57 pam-fi 1.2
58 mocchiut 1.1 // .........................................................
59     // namespaces
60     // .........................................................
61     using namespace std;
62     using namespace pamela;
63     // ================================================================================
64     //
65     //
66     // ================================================================================
67     /**
68     * \brief Tracker data reduction routine.
69     *
70     * It performs data reduction LEVEL0->LEVEL1->LEVEL2, producing LEVEL1 and or LEVEL2 output, in ROOT or HBOOK format.
71     * Input parameters:
72     * @param run id of the run to be processed (if run=0 a whole file is reprocessed)
73     * @param dbc pointer to BD server
74     * @param file1 output LEVEL1 file name (null if no LEVEL1 output is required)
75     * @param file2 output LEVEL2 file name (null if no LEVEL2 output is required)
76     * @param standalone (flag to run the program in standalone mode, that is without reading RunInfo)
77     */
78     //short int TrkCore(Int_t run, TSQLServer *dbc, TString file1, TString file2, Bool_t standalone)
79 mocchiut 1.31 //int TrkCore(UInt_t run, TFile *f2, TSQLServer *dbc, int ncustom, char*vcustom[]) // EMILIANO
80     int TrkCore(UInt_t run, TFile *f2, GL_TABLES *glt, int ncustom, char*vcustom[]) // EMILIANO
81 mocchiut 1.1 {
82     // ---------------------------
83     // Define some basic variables
84     // ---------------------------
85    
86 pizzolot 1.36 TString filename = "";
87 pam-fi 1.7 Long64_t nentries = 0LL;
88 pam-fi 1.2
89 mocchiut 1.1 // LEVEL0 input
90     TFile *f0 = 0;
91     TTree *physicsTree = 0;
92     TBranch *b_trk = 0;
93     TBranch *b_header = 0;
94     EventHeader *header = 0;
95     PscuHeader *pscu = 0;
96 pam-fi 1.22 TrkLevel0 *l0_event = new TrkLevel0();
97 mocchiut 1.1 // RunInfo
98     ItoRunInfo *runinfo = 0;
99 pam-fi 1.7 TArrayI *runlist = 0;
100     UInt_t from_run;
101     UInt_t to_run;
102 mocchiut 1.1 Bool_t reprocessing = false;
103     //LEVEL2 output - ROOT
104     TTree *t_level2 = 0;
105     TTree *t_clone = 0;
106     TrkLevel2 *l2_event = new TrkLevel2();
107     TrkLevel2 *l2_clone = new TrkLevel2();
108 pam-fi 1.13 //LEVEL1 output - ROOT
109     TrkLevel1 *l1_event = new TrkLevel1();
110 pam-fi 1.15 TrkHough *lh_event = new TrkHough();
111 mocchiut 1.42
112     // EM: how to use CaloLevel1 data inside tracker core routine, example:
113     TTree *T = (TTree*)f2->Get("Calorimeter");
114     if ( T ){
115     CaloStrip *cs = new CaloStrip(false);
116     CaloLevel1 *cl1 = new CaloLevel1();
117     T->SetBranchAddress("CaloLevel1",&cl1);
118     //
119     for (Int_t i=0; i < T->GetEntries(); i++){
120     T->GetEntry(i);
121     Bool_t OK = false;
122     Int_t k = 0;
123     Int_t view = 0;
124     Int_t plane = 0;
125     Int_t strip = 0;
126     Float_t mip = 0.;
127     Float_t ene[2][22][96];
128     memset (ene,0,sizeof(Float_t)*2*22*96);
129     //
130     while ( k < cl1->istrip ){
131     mip = cl1->DecodeEstrip(k,view,plane,strip);
132     // if ( view==0 && plane==18 ) mip = 0.; // 9thRED - set to zero signals from plane 18X
133     //
134     ene[view][plane][strip] = mip;
135     k++;
136     };
137     // check 5
138     Float_t totX = 0.;
139     Float_t totY = 0.;
140     Float_t posX = 0.;
141     Float_t posY = 0.;
142     Int_t minYstrip = 100;
143     Int_t maxYstrip = -1;
144     // calculate baricenter strip on first x and y planes
145     for (Int_t j = 0; j<2; j++){
146     for (Int_t i = 0; i<96; i++){
147     totX += ene[0][j][i];
148     totY += ene[1][j][i];
149     if ( j==0 && ene[1][0][i]>0. && i > maxYstrip ) maxYstrip = i;
150     if ( j==0 && ene[1][0][i]>0. && i < minYstrip ) minYstrip = i;
151     posX += (float)i*ene[0][j][i];
152     posY += (float)i*ene[1][j][i];
153     }
154     }
155     //
156     // printf(" totX %f totY %f maxYstrip %i minYstrip %i \n",totX,totY,maxYstrip,minYstrip);
157     Int_t bariX = 0;
158     Int_t bariY = 0;
159     if ( totX>0. && totY > 0.){
160     // printf(" totX %f totY %f \n",totX,totY);
161     posX /= totX;
162     posY /= totY;
163     bariX = (int)(posX);
164     bariY = (int)(posY);
165     Int_t minX = bariX - 10;
166     Int_t maxX = bariX + 10;
167     Int_t minY = bariY - 10;
168     Int_t maxY = bariY + 10;
169     if ( minX < 0 ) minX = 0;
170     if ( minY < 0 ) minY = 0;
171     if ( maxX > 95 ) maxX = 95;
172     if ( maxY > 95 ) maxY = 95;
173     totX = 0.;
174     totY = 0.;
175     Float_t parX = 0.;
176     Float_t parY = 0.;
177    
178     for (Int_t j = 0; j<3; j++){
179     for (Int_t i = 0; i<96; i++){
180     totX += ene[0][j][i];
181    
182    
183     totY += ene[1][j][i];
184     if ( i>=minX && i<=maxX){
185     parX += ene[0][j][i];
186     }
187     if ( i>=minY && i<=maxY){
188     parY += ene[1][j][i];
189     }
190     }
191     }
192     // printf(" bariX %i bariY %i minX %i maxX %i minY %i maxY %i fractionX %f fractionY %f \n",bariX,bariY,minX,maxX,minY,maxY,parX/totX,parY/totY);
193    
194     if ( parX/totX > 0.65 && parY/totY > 0.75 ) OK = true; // less strict cut on the X due to possible bremsstrhalung double showers
195     // if ( OK ) printf(" CHECK5 PASSED\n");
196     }
197    
198     // check 1
199     if ( totY > 0. && abs(maxYstrip-minYstrip)<10 ){
200     posY /= totY;
201     bariY = (int)(posY);
202     Float_t ntotX = 0.;
203     Float_t ntotY = 0.;
204     for (Int_t j = 0; j<2; j++){
205     for (Int_t i = 0; i<96; i++){
206     ntotX += ene[0][j][i];
207     if ( j > 0 ) ntotY += ene[1][j][i];
208     }
209     }
210     if ( ntotX < 0.3 && ntotY < 0.3 ){
211     OK = true;
212     // printf(" CHECK1 PASSED\n");
213     }
214     }
215    
216     if ( OK ){
217     printf("%i: event is good! \n",i);
218     cs->Set(1,0,bariY);
219     printf(" Y0 hit strip coordinates: X %f Y %f Z %f \n",cs->GetX(),cs->GetY(),cs->GetZ());
220     }
221     }
222     T->Delete();
223     }
224     // end EM
225    
226 mocchiut 1.1 // -----------------------
227     // -----------------------
228     // Handle input parameters
229     // (data)
230     //
231     // - open/create output files, determining the environment root/hbook from the estension
232     // - create the level1/level2 tree+branch/nt-uple
233     // -----------------------
234     // -----------------------
235     TrkProcess *p = new TrkProcess(run,f2);
236 pam-fi 1.26 if( p->HandleCustomPar(ncustom,vcustom) )return(p->ostatus);;
237 pam-fi 1.25 if( TrkParams::VerboseMode() )p->Dump();
238 mocchiut 1.1
239     // ===================================
240     // Open/Create level2 output file
241     // ===================================
242     if(p->get2){
243 pam-fi 1.13 //-------------------------------------------
244     // read RunInfo
245     //-------------------------------------------
246     if(!(p->standalone)){
247     // Open "RunInfo" tree and get list of run
248     runinfo = new ItoRunInfo(f2);
249     char *trkversion = TrkInfo(false);
250     int runinfo_error = runinfo->Update(run,"TRK",trkversion);
251     if( runinfo_error ) throw runinfo_error;
252     runlist = runinfo->GetRunList();
253     reprocessing = runinfo->IsReprocessing();//????
254 pam-fi 1.25 if(TrkParams::VerboseMode()){
255 pam-fi 1.13 cout << "#events "<< runinfo->GetFileEntries() << endl;// #events in the file
256     cout << "#runs "<< runinfo->GetRunEntries() << endl;// #runs in the file
257     cout << "1-st run "<< runlist->At(0) << endl;
258     cout << "last run "<< runlist->At(runinfo->GetRunEntries()-1) << endl;
259     cout << "1-st event "<< runinfo->GetFirstEntry() << endl;// first event of our run
260     cout << "last event+1 "<< runinfo->GetLastEntry() << endl;// first event AFTER the last event of our run
261     cout << "reprocessing "<< runinfo->IsReprocessing() << endl << endl;
262     };
263     };
264     //-------------------------------------------
265     //
266     //-------------------------------------------
267     // Take (if present) the old Tracker tree
268     t_clone = (TTree*)f2->Get("Tracker");
269     if( t_clone != NULL ) t_clone->SetBranchAddress("TrkLevel2",&l2_clone);
270     // Create NEW Tracker tree
271     t_level2 = new TTree("Tracker","PAMELA tracker level2 data ");
272 pam-fi 1.21 l2_event->Set(); // ****NBNBNBN*****
273 pam-fi 1.13 t_level2->Branch("TrkLevel2","TrkLevel2",&l2_event);
274     if(p->get1){
275 pam-fi 1.25 if(TrkParams::VerboseMode())cout << endl << "Requested LEVEL1 output" << endl;
276 pam-fi 1.21 l1_event->Set(); // ****NBNBNBN*****
277 pam-fi 1.13 t_level2->Branch("TrkLevel1","TrkLevel1",&l1_event);
278 pam-fi 1.16 }
279     if(p->geth){
280 pam-fi 1.25 if(TrkParams::VerboseMode())cout << endl << "Requested Hough-Transform output" << endl;
281 pam-fi 1.16 t_level2->Branch("TrkHough","TrkHough",&lh_event);
282 pam-fi 1.7 };
283 pam-fi 1.13 };
284 mocchiut 1.1
285    
286     // -------------------------------------------
287     // define runs to be processed/reprocessed
288     // -------------------------------------------
289     if(run == 0){
290 pam-fi 1.13 // reprocessing ALL runs
291     if(p->standalone)throw -298; // reprocessing not implemented
292     from_run = runlist->At(0);
293     to_run = runlist->At(runinfo->GetRunEntries()-1);
294 mocchiut 1.1 }else{
295     // processing/reprocessing ONE single run
296 pam-fi 1.13 from_run = run;
297     to_run = run;
298 mocchiut 1.1 };
299    
300     //
301     // init event counter
302     //
303     Int_t ev_count =0;
304     //
305     // create query-results objects
306     //
307 pam-fi 1.7 GL_RUN q1 = GL_RUN();
308 mocchiut 1.1 GL_TRK_CALIB q2 = GL_TRK_CALIB();
309     GL_ROOT q3 = GL_ROOT();
310     GL_PARAM q4 = GL_PARAM();
311    
312     // ------------------------------------------------------------
313     // if reprocessing one run, copy all the events BEFORE the run
314     // ------------------------------------------------------------
315     if( !(p->standalone) ){
316 pam-fi 1.13 for(UInt_t i=0; i<runinfo->GetFirstEntry(); i++){
317 mocchiut 1.37 if ( t_clone->GetEntry(i) <= 0 ) throw -36;// EM
318 pam-fi 1.13 *l2_event = *l2_clone;
319     t_level2->Fill();
320     l2_event->Clear();
321     // COPY COPY COPY
322     };
323 mocchiut 1.1 };
324     // ------------------------------------------------------------
325     // ------------------------------------------------------------
326 pam-fi 1.25 // START LOOP OVER RUNS TO PROCESS/REPROCESS
327 mocchiut 1.1 // ------------------------------------------------------------
328     // ------------------------------------------------------------
329 pam-fi 1.7 for(UInt_t idRun = from_run; idRun <= to_run; idRun++){
330 mocchiut 1.1
331 pam-fi 1.25 if(TrkParams::VerboseMode()) cout << endl<<" ========================= Run: "<< idRun << endl;
332 pam-fi 1.7 UInt_t runheadtime = 0;
333     UInt_t runtrailtime = 0;
334 pam-fi 1.32 UInt_t runheadobt = 0;
335     UInt_t runtrailobt = 0;
336 mocchiut 1.41 UInt_t abstime = 0;
337 mocchiut 1.1 UInt_t evfrom = 0;
338     UInt_t evto = 0;
339 pam-fi 1.7 UInt_t nevents = 0;
340     Int_t id_root_l0 =-1;
341 mocchiut 1.1 Int_t trk_calib_used = 0;
342 pam-fi 1.7
343 mocchiut 1.31 TString host = glt->CGetHost(); // EMILIANO
344     TString user = glt->CGetUser(); // EMILIANO
345     TString psw = glt->CGetPsw(); // EMILIANO
346     // TString host = "mysql://localhost/pamelaprod";
347     // TString user = "anonymous";
348     // TString psw = "";
349     // const char *pamdbhost=gSystem->Getenv("PAM_DBHOST");
350     // const char *pamdbuser=gSystem->Getenv("PAM_DBUSER");
351     // const char *pamdbpsw=gSystem->Getenv("PAM_DBPSW");
352     // if ( !pamdbhost ) pamdbhost = "";
353     // if ( !pamdbuser ) pamdbuser = "";
354     // if ( !pamdbpsw ) pamdbpsw = "";
355     // if ( strcmp(pamdbhost,"") ) host = pamdbhost;
356     // if ( strcmp(pamdbuser,"") ) user = pamdbuser;
357     // if ( strcmp(pamdbpsw,"") ) psw = pamdbpsw;
358 mocchiut 1.29 TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
359 mocchiut 1.31 stringstream myquery; // EMILIANO
360     myquery.str(""); // EMILIANO
361     myquery << "SET time_zone='+0:00'"; // EMILIANO
362 mocchiut 1.39 delete dbc->Query(myquery.str().c_str()); // EMILIANO
363 mocchiut 1.1 if(p->standalone){
364     // ==============================
365     // first query: retrieve run info
366     // ==============================
367 pam-fi 1.7 if (q1.Query_GL_RUN(idRun,dbc) )throw -50;
368     id_root_l0 = q1.ID_ROOT_L0;
369 mocchiut 1.1 runheadtime = q1.RUNHEADER_TIME;
370     runtrailtime = q1.RUNTRAILER_TIME;
371 pam-fi 1.32 runheadobt = q1.RUNHEADER_OBT;
372     runtrailobt = q1.RUNTRAILER_OBT;
373 pam-fi 1.7 evfrom = q1.EV_FROM;
374     evto = q1.EV_TO;
375 pam-fi 1.13 nevents = q1.NEVENTS;
376     trk_calib_used = q1.TRK_CALIB_USED;
377 mocchiut 1.1 }else{
378     // ==============================
379     // get run info from RunInfo tree
380     // ==============================
381 pam-fi 1.3 int runinfo_error = runinfo->GetRunInfo(idRun);
382 pam-fi 1.7 if( runinfo_error ) throw runinfo_error;
383     id_root_l0 = runinfo->ID_ROOT_L0;
384 mocchiut 1.1 runheadtime = runinfo->RUNHEADER_TIME;
385     runtrailtime = runinfo->RUNTRAILER_TIME;
386 pam-fi 1.32 runheadobt = runinfo->RUNHEADER_OBT;
387     runtrailobt = runinfo->RUNTRAILER_OBT;
388 pam-fi 1.7 evfrom = runinfo->EV_FROM;
389     evto = runinfo->EV_TO;
390 pam-fi 1.13 nevents = runinfo->NEVENTS;
391     trk_calib_used = runinfo->TRK_CALIB_USED;
392 pam-fi 1.24
393 mocchiut 1.1 };
394 pam-fi 1.24
395 mocchiut 1.1 //
396 pam-fi 1.25 if(TrkParams::VerboseMode()){
397 pam-fi 1.7 cout << "ROOT file ID "<< id_root_l0 << endl;
398 mocchiut 1.1 cout << "RunHeader time "<< runheadtime << endl;
399     cout << "RunTrailer time "<< runtrailtime << endl;
400 pam-fi 1.7 cout << " from event "<< evfrom << endl;
401     cout << " to event "<< evto << endl;
402 pam-fi 1.24 cout << " num. events "<< nevents << endl;
403     cout << "trk-calibration used "<< trk_calib_used << endl;
404 mocchiut 1.1 };
405     // ========================================================
406     // second query: search the LEVEL0 file that contains idRun
407     // ========================================================
408     TString lastfilename = filename;
409 pam-fi 1.7 if( q3.Query_GL_ROOT(id_root_l0,dbc) )throw -51;
410 mocchiut 1.1 filename = q3.PATH + q3.NAME;
411     // ========================================================
412     // Open the input LEVEL0 data file
413     // ========================================================
414     if(filename.CompareTo(lastfilename)){
415     if(!lastfilename.IsNull())f0->Close();
416     //if( debug ) cout << "Opening LEVEL0 file: "<< filename << endl;
417 pam-fi 1.25 if(TrkParams::VerboseMode()) cout << "Opening LEVEL0 file: "<< filename << endl;
418 pam-fi 1.2 FileStat_t t;
419     if( gSystem->GetPathInfo(filename.Data(),t) )throw -6;
420 mocchiut 1.39 if ( f0 ) f0->Close();
421 mocchiut 1.1 f0 = new TFile(filename);
422     if ( !f0 ) throw -6;
423     physicsTree = (TTree*)f0->Get("Physics"); if(!physicsTree) throw -7;
424     b_header = physicsTree ->GetBranch("Header"); if(!b_header ) throw -8;
425     b_trk = physicsTree ->GetBranch("Tracker"); if(!b_trk ) throw -203;
426 pam-fi 1.22 l0_event->Set();
427     physicsTree->SetBranchAddress("Tracker" ,l0_event->GetPointerToTrackerEvent());
428 mocchiut 1.1 physicsTree->SetBranchAddress("Header",&header);
429    
430 pam-fi 1.7 nentries = physicsTree->GetEntries();
431 mocchiut 1.34 if ( nentries < 1 && nevents ) throw -11; // EMILIANO: go on if the file is empty, why not? in principle we should not have any event to process but if the case the next line will take care of exiting; exit only if the file is empty but nevents is not zero
432 mocchiut 1.33 if ( nentries < (evto+1) && nevents > 0 ) throw -12; // EMILIANO: if NEVENTS = 0 and the file is empty everything is ok but due to a mistake in the
433     // variable EV_TO type (UInt_t instead of Int_t) that we don't want to change to avoid changing a lot of code, evto becomes inf hence
434     // the condition ( nentries < (evto+1) ) is satisfied and DV exit with error even if the error is only in the DB.
435 pam-fi 1.7
436 mocchiut 1.1 };
437    
438 pam-fi 1.7 GL_TIMESYNC *dbtime = new GL_TIMESYNC(id_root_l0,"ID",dbc);
439    
440 pam-fi 1.24
441     // =============================================================
442     // retrieve information about parameters to process LEVEL2
443     // =============================================================
444    
445 mocchiut 1.41 TrkParams::Set(runinfo->GetGL_RUN(),dbc);
446     for(int i=0; i<p->npar; i++){
447     if(TrkParams::VerboseMode())cout<<" ((( force parameters from input path )))"<<endl;
448     TrkParams::Set(p->parpath[i],p->partype[i]);
449     }
450 pam-fi 1.25
451 mocchiut 1.41 TrkParams::Load();
452     if( !TrkParams::IsLoaded() )throw -52;
453 pam-fi 1.24
454     // =============================================================
455     // retrieve calibration file needed to reduce data
456     // =============================================================
457    
458     TrkParams::SetCalib(runinfo->GetGL_RUN(),dbc);
459     TrkParams::LoadCalib( );
460     if( !TrkParams::CalibIsLoaded() )throw -53;
461    
462 pam-fi 1.11 TBenchmark *reduction = new TBenchmark();
463 pam-fi 1.25 if(TrkParams::VerboseMode())reduction->Start("reduction");
464 pam-fi 1.11 Int_t ntrk=0;
465 mocchiut 1.1 // ====================================================
466     // start looping on events cointained in the data file
467     // ====================================================
468 mocchiut 1.41 if ( TrkParams::GetSimuFlag() ){
469     abstime = runheadtime;
470     } else {
471     if(dbc){
472     dbc->Close();
473     delete dbc;
474     }
475     }
476    
477 pam-fi 1.32 for (UInt_t re = evfrom+min(p->nskip,nevents); re < evfrom+nevents; re++){
478 mocchiut 1.1
479     ev_count++;
480 pam-fi 1.32
481    
482 pam-fi 1.25 // if ( TrkParams::DebugMode() && re%100 == 0 && re > 0 ) cout << ".";
483 mocchiut 1.1
484 mocchiut 1.37 if ( b_trk->GetEntry(re) <= 0 ) throw -36;//EM
485     if ( b_header->GetEntry(re) <= 0 ) throw -36;//EM
486 pam-fi 1.13 pscu = header->GetPscuHeader();
487 mocchiut 1.41
488     // =============================================================
489     // The following 6 lines have been moved here by VALERIO.
490     if(TrkParams::GetSimuFlag()){
491     abstime = runheadtime + (int) floor(0.03*(re-evfrom)); //If simulated data we need to assign a fake abstime. 30ms between each event.
492     if(TrkParams::VerboseMode())cout << "Event: " << re-evfrom << " - Attempting to retrieve Mask Info for abstime=" << abstime << endl;
493     if(!TrkParams::Set(runinfo->GetGL_RUN(),dbc,6,abstime))throw -52; // Setting to load mask (VALERIO)
494     TrkParams::Load(6);
495     if( !TrkParams::IsLoaded() )throw -52;
496     if(TrkParams::VerboseMode())cout << "Mask Info for abstime=" << abstime << " retrieved" << endl;
497     }
498    
499 pam-fi 1.25 if( TrkParams::DebugMode() )cout << ">>> "<<ev_count-1<<" @ OBT "<<pscu->GetOrbitalTime()<<endl;
500 pam-fi 1.23
501 mocchiut 1.35 if ( dbtime->DBabsTime(pscu->GetOrbitalTime()) > (runtrailtime+1) || dbtime->DBabsTime(pscu->GetOrbitalTime()) < (runheadtime-1)) {
502 pam-fi 1.7
503 pam-fi 1.25 if (TrkParams::VerboseMode()){
504 pam-fi 1.14 printf(" TrkCore - WARNING: event outside the run time window, skipping it\n");
505 pam-fi 1.32 cout << " OBT "<<pscu->GetOrbitalTime()<<" RUN "<<runheadobt<<"-"<<runtrailobt<<" ABS-time "<<dbtime->DBabsTime(pscu->GetOrbitalTime())<<" RUN "<<runheadtime<<"-"<<runtrailtime<<endl;
506 pam-fi 1.14 };
507 pam-fi 1.13 }else{
508 pam-fi 1.25 if ( TrkParams::DebugMode() )
509 pam-fi 1.23 printf("\n-------------------------------\nevent %d\n",re-evfrom);
510 pam-fi 1.16
511 pam-fi 1.13 p->ProcessEvent(l0_event);
512 pam-fi 1.7
513 pam-fi 1.13 // ----------------
514     // LEVEL1 output
515     // ----------------
516     if(p->get1){
517     if(p->ifroot1){ // root
518     l1_event->Clear();
519 pam-fi 1.22 // l1_event->SetFromLevel1Struct(&level1event_,p->full1);
520     l1_event->SetFromLevel1Struct(p->full1);
521 pam-fi 1.13 // t_level1->Fill();
522 pam-fi 1.23 }else{ // hbook
523 pam-fi 1.13 throw -299;
524     };
525     };
526     // ----------------
527     // HOUGH output
528     // ----------------
529     if(p->geth){
530     if(p->ifrooth){ // root
531 pam-fi 1.23 lh_event->Delete();
532     lh_event->SetFromHoughStruct(&houghevent_);
533 pam-fi 1.13 }else{ // hbook
534     throw -299;
535     };
536     };
537     // ----------------
538     // LEVEL2 output
539     // ----------------
540     if(p->get2){
541 pam-fi 1.23 l2_event->Clear();
542     if(p->get1){
543     l2_event->SetFromLevel2Struct(&level2event_,l1_event);//set references to level1
544     }else{
545     l2_event->SetFromLevel2Struct(&level2event_);
546     }
547 pam-fi 1.10 // l2_event->Dump();
548 pam-fi 1.23 t_level2->Fill();
549     if( l2_event->ntrk()>0 )ntrk++;
550 pam-fi 1.27 if(TrkParams::VerboseMode())l2_event->Dump();
551 mocchiut 1.1 };
552 pam-fi 1.13 };
553 mocchiut 1.1 }; // end loop on events
554 pam-fi 1.25 if(TrkParams::VerboseMode()){
555 pam-fi 1.13 cout << " >>> processed "<< ev_count <<" events"<< endl;
556     if(ev_count)cout << ntrk << " events with at least one track ("<<(Int_t)(100*ntrk/ev_count)<<"%)\n";
557     reduction->Show("reduction");
558 pam-fi 1.11 }
559     delete reduction;
560 pam-fi 1.7
561     delete dbtime;
562    
563 pam-fi 1.13 }; // end loop on runs
564    
565    
566 mocchiut 1.1 // ------------------------------------------------------------
567     // if reprocessing one run, copy all the events AFTER the run
568     // ------------------------------------------------------------
569     if( !(p->standalone) ){
570     for(UInt_t i=runinfo->GetLastEntry()+1; i<runinfo->GetFileEntries(); i++){
571 mocchiut 1.37 if ( t_clone->GetEntry(i) <= 0 ) throw -36;//EM
572 mocchiut 1.1 *l2_event = *l2_clone;
573     t_level2->Fill();
574     l2_event->Clear();
575     // COPY COPY COPY
576     };
577     };
578     // ---------------
579     // close the files
580     // ---------------
581     if(p->get2){
582 mocchiut 1.40 if( t_clone )t_clone->Delete("all");//delete old tree from file
583 pam-fi 1.13 if( !(p->standalone) )runinfo->Close();
584 mocchiut 1.40 // t_level2->Write("Tracker", TObject::kOverwrite);
585     f2->Write("Tracker", TObject::kOverwrite);
586     if( t_level2 )t_level2->Delete(); //delete new tree from memory
587 pam-fi 1.13
588     };
589 pam-fi 1.6
590 pam-fi 1.24 if(f0) if(f0->IsOpen()) f0->Close();
591     // if( f0_c->IsOpen() )f0_c->Close();
592 pam-fi 1.13
593     // lh_event->Delete();
594     l1_event->Delete();
595     l2_event->Delete();
596     l2_clone->Delete();
597 pam-fi 1.2
598     return(p->ostatus);
599 mocchiut 1.1 }
600    

  ViewVC Help
Powered by ViewVC 1.1.23