/[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.15 - (hide annotations) (download)
Fri Nov 10 11:38:44 2006 UTC (18 years ago) by pam-fi
Branch: MAIN
Changes since 1.14: +8 -7 lines
TrkHough class added for track recognition debug

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 pam-fi 1.2
54 mocchiut 1.1 // .........................................................
55     // namespaces
56     // .........................................................
57     using namespace std;
58     using namespace pamela;
59     // ================================================================================
60     //
61     //
62     // ================================================================================
63     /**
64     * \brief Tracker data reduction routine.
65     *
66     * It performs data reduction LEVEL0->LEVEL1->LEVEL2, producing LEVEL1 and or LEVEL2 output, in ROOT or HBOOK format.
67     * Input parameters:
68     * @param run id of the run to be processed (if run=0 a whole file is reprocessed)
69     * @param dbc pointer to BD server
70     * @param file1 output LEVEL1 file name (null if no LEVEL1 output is required)
71     * @param file2 output LEVEL2 file name (null if no LEVEL2 output is required)
72     * @param standalone (flag to run the program in standalone mode, that is without reading RunInfo)
73     */
74     //short int TrkCore(Int_t run, TSQLServer *dbc, TString file1, TString file2, Bool_t standalone)
75 pam-fi 1.7 int TrkCore(UInt_t run, TFile *f2, TSQLServer *dbc, int ncustom, char*vcustom[])
76 mocchiut 1.1 {
77     // ---------------------------
78     // Define some basic variables
79     // ---------------------------
80    
81     TString filename = 0;
82     TString trkcalibfile = 0;
83     Int_t last_trk_calib_used = 0;
84 pam-fi 1.7 Long64_t nentries = 0LL;
85 pam-fi 1.2
86 mocchiut 1.1 // LEVEL0 input
87     TFile *f0 = 0;
88 pam-fi 1.2 TFile *f0_c = 0;
89 mocchiut 1.1 TTree *physicsTree = 0;
90     TBranch *b_trk = 0;
91     TBranch *b_header = 0;
92     EventHeader *header = 0;
93     PscuHeader *pscu = 0;
94     TrkLevel0 *l0_event = 0;
95     // RunInfo
96     ItoRunInfo *runinfo = 0;
97 pam-fi 1.7 TArrayI *runlist = 0;
98     UInt_t from_run;
99     UInt_t to_run;
100 mocchiut 1.1 Bool_t reprocessing = false;
101     //LEVEL2 output - ROOT
102     TTree *t_level2 = 0;
103     TTree *t_clone = 0;
104     TrkLevel2 *l2_event = new TrkLevel2();
105     TrkLevel2 *l2_clone = new TrkLevel2();
106 pam-fi 1.13 //LEVEL1 output - ROOT
107     TrkLevel1 *l1_event = new TrkLevel1();
108 pam-fi 1.15 TrkHough *lh_event = new TrkHough();
109 mocchiut 1.1 // -----------------------
110     // -----------------------
111     // Handle input parameters
112     // (data)
113     //
114     // - open/create output files, determining the environment root/hbook from the estension
115     // - create the level1/level2 tree+branch/nt-uple
116     // -----------------------
117     // -----------------------
118     TrkProcess *p = new TrkProcess(run,f2);
119     p->HandleCustomPar(ncustom,vcustom);
120 pam-fi 1.2 if( p->DebugMode() )p->Dump();
121 mocchiut 1.1
122     // ===================================
123     // Open/Create level2 output file
124     // ===================================
125     if(p->get2){
126 pam-fi 1.13 //-------------------------------------------
127     // read RunInfo
128     //-------------------------------------------
129     if(!(p->standalone)){
130     // Open "RunInfo" tree and get list of run
131     runinfo = new ItoRunInfo(f2);
132     char *trkversion = TrkInfo(false);
133     int runinfo_error = runinfo->Update(run,"TRK",trkversion);
134     if( runinfo_error ) throw runinfo_error;
135     runlist = runinfo->GetRunList();
136     reprocessing = runinfo->IsReprocessing();//????
137     if(p->VerboseMode()){
138     cout << "#events "<< runinfo->GetFileEntries() << endl;// #events in the file
139     cout << "#runs "<< runinfo->GetRunEntries() << endl;// #runs in the file
140     cout << "1-st run "<< runlist->At(0) << endl;
141     cout << "last run "<< runlist->At(runinfo->GetRunEntries()-1) << endl;
142     cout << "1-st event "<< runinfo->GetFirstEntry() << endl;// first event of our run
143     cout << "last event+1 "<< runinfo->GetLastEntry() << endl;// first event AFTER the last event of our run
144     cout << "reprocessing "<< runinfo->IsReprocessing() << endl << endl;
145     };
146     };
147     //-------------------------------------------
148     //
149     //-------------------------------------------
150     // Take (if present) the old Tracker tree
151     t_clone = (TTree*)f2->Get("Tracker");
152     if( t_clone != NULL ) t_clone->SetBranchAddress("TrkLevel2",&l2_clone);
153     // Create NEW Tracker tree
154     t_level2 = new TTree("Tracker","PAMELA tracker level2 data ");
155     t_level2->Branch("TrkLevel2","TrkLevel2",&l2_event);
156     if(p->get1){
157     if(p->VerboseMode())cout << endl << "Requested LEVEL1 output" << endl;
158     t_level2->Branch("TrkLevel1","TrkLevel1",&l1_event);
159     t_level2->BranchRef();
160     }else if(p->geth){
161     if(p->VerboseMode())cout << endl << "Requested Hough-Transform output" << endl;
162 pam-fi 1.15 t_level2->Branch("TrkHough","TrkHough",&lh_event);
163 pam-fi 1.13 // t_level2->BranchRef();
164 pam-fi 1.7 };
165 pam-fi 1.13 };
166 mocchiut 1.1
167    
168     // -------------------------------------------
169     // define runs to be processed/reprocessed
170     // -------------------------------------------
171     if(run == 0){
172 pam-fi 1.13 // reprocessing ALL runs
173     if(p->standalone)throw -298; // reprocessing not implemented
174     from_run = runlist->At(0);
175     to_run = runlist->At(runinfo->GetRunEntries()-1);
176 mocchiut 1.1 }else{
177     // processing/reprocessing ONE single run
178 pam-fi 1.13 from_run = run;
179     to_run = run;
180 mocchiut 1.1 };
181    
182     //
183     // init "expiration date" of calibrations and parameters
184     //
185 pam-fi 1.7 UInt_t pedsig_time =0;
186 pam-fi 1.13 UInt_t B_time =0;
187     UInt_t mip_time =0;
188     UInt_t charge_time =0;
189     UInt_t eta_time =0;
190     UInt_t align_time =0;
191     UInt_t mask_time =0;
192 mocchiut 1.1 //
193     // init event counter
194     //
195     Int_t ev_count =0;
196     //
197     // create query-results objects
198     //
199 pam-fi 1.7 GL_RUN q1 = GL_RUN();
200 mocchiut 1.1 GL_TRK_CALIB q2 = GL_TRK_CALIB();
201     GL_ROOT q3 = GL_ROOT();
202     GL_PARAM q4 = GL_PARAM();
203    
204     // ------------------------------------------------------------
205     // if reprocessing one run, copy all the events BEFORE the run
206     // ------------------------------------------------------------
207     if( !(p->standalone) ){
208 pam-fi 1.13 for(UInt_t i=0; i<runinfo->GetFirstEntry(); i++){
209     t_clone->GetEntry(i);
210     *l2_event = *l2_clone;
211     t_level2->Fill();
212     l2_event->Clear();
213     // COPY COPY COPY
214     };
215 mocchiut 1.1 };
216     // ------------------------------------------------------------
217     // ------------------------------------------------------------
218     // START LOOP OVER RUNS TO PROCESSED/REPROCESSED
219     // ------------------------------------------------------------
220     // ------------------------------------------------------------
221 pam-fi 1.7 for(UInt_t idRun = from_run; idRun <= to_run; idRun++){
222 mocchiut 1.1
223 pam-fi 1.2 if(p->VerboseMode()) cout << endl<<" ========================= Run: "<< idRun << endl;
224 pam-fi 1.7 UInt_t runheadtime = 0;
225     UInt_t runtrailtime = 0;
226 mocchiut 1.1 UInt_t evfrom = 0;
227     UInt_t evto = 0;
228 pam-fi 1.7 UInt_t nevents = 0;
229     Int_t id_root_l0 =-1;
230 mocchiut 1.1 Int_t trk_calib_used = 0;
231 pam-fi 1.7
232 mocchiut 1.1 if(p->standalone){
233     // ==============================
234     // first query: retrieve run info
235     // ==============================
236 pam-fi 1.7 if (q1.Query_GL_RUN(idRun,dbc) )throw -50;
237     id_root_l0 = q1.ID_ROOT_L0;
238 mocchiut 1.1 runheadtime = q1.RUNHEADER_TIME;
239     runtrailtime = q1.RUNTRAILER_TIME;
240 pam-fi 1.7 evfrom = q1.EV_FROM;
241     evto = q1.EV_TO;
242 pam-fi 1.13 nevents = q1.NEVENTS;
243     trk_calib_used = q1.TRK_CALIB_USED;
244 mocchiut 1.1 }else{
245     // ==============================
246     // get run info from RunInfo tree
247     // ==============================
248 pam-fi 1.3 int runinfo_error = runinfo->GetRunInfo(idRun);
249 pam-fi 1.7 if( runinfo_error ) throw runinfo_error;
250     id_root_l0 = runinfo->ID_ROOT_L0;
251 mocchiut 1.1 runheadtime = runinfo->RUNHEADER_TIME;
252     runtrailtime = runinfo->RUNTRAILER_TIME;
253 pam-fi 1.7 evfrom = runinfo->EV_FROM;
254     evto = runinfo->EV_TO;
255 pam-fi 1.13 nevents = runinfo->NEVENTS;
256     trk_calib_used = runinfo->TRK_CALIB_USED;
257 mocchiut 1.1 };
258     //
259 pam-fi 1.2 if(p->VerboseMode()){
260 pam-fi 1.7 cout << "ROOT file ID "<< id_root_l0 << endl;
261 mocchiut 1.1 cout << "RunHeader time "<< runheadtime << endl;
262     cout << "RunTrailer time "<< runtrailtime << endl;
263 pam-fi 1.7 cout << " from event "<< evfrom << endl;
264     cout << " to event "<< evto << endl;
265     cout << " num. events "<< nevents << endl;
266     cout << "trk-calibration used "<< trk_calib_used << endl;
267 mocchiut 1.1 };
268     // ========================================================
269     // second query: search the LEVEL0 file that contains idRun
270     // ========================================================
271     TString lastfilename = filename;
272 pam-fi 1.7 if( q3.Query_GL_ROOT(id_root_l0,dbc) )throw -51;
273 mocchiut 1.1 filename = q3.PATH + q3.NAME;
274     // ========================================================
275     // Open the input LEVEL0 data file
276     // ========================================================
277     if(filename.CompareTo(lastfilename)){
278     if(!lastfilename.IsNull())f0->Close();
279     //if( debug ) cout << "Opening LEVEL0 file: "<< filename << endl;
280 pam-fi 1.2 if(p->VerboseMode()) cout << "Opening LEVEL0 file: "<< filename << endl;
281     FileStat_t t;
282     if( gSystem->GetPathInfo(filename.Data(),t) )throw -6;
283 mocchiut 1.1 f0 = new TFile(filename);
284     if ( !f0 ) throw -6;
285     physicsTree = (TTree*)f0->Get("Physics"); if(!physicsTree) throw -7;
286     b_header = physicsTree ->GetBranch("Header"); if(!b_header ) throw -8;
287 pam-fi 1.7 // b_registry = physicsTree ->GetBranch("Registry"); if(!b_registry ) throw -9;
288 mocchiut 1.1 b_trk = physicsTree ->GetBranch("Tracker"); if(!b_trk ) throw -203;
289     physicsTree->SetBranchAddress("Tracker" ,&l0_event);
290 pam-fi 1.7 // physicsTree->SetBranchAddress("Registry",&reg);
291 mocchiut 1.1 physicsTree->SetBranchAddress("Header",&header);
292    
293 pam-fi 1.7 nentries = physicsTree->GetEntries();
294     if ( nentries < 1 ) throw -11;
295     if ( nentries < (evto+1)) throw -12;
296    
297 mocchiut 1.1 };
298    
299 pam-fi 1.7 GL_TIMESYNC *dbtime = new GL_TIMESYNC(id_root_l0,"ID",dbc);
300    
301 mocchiut 1.1 // =============================================================
302 pam-fi 1.2 // retrieve calibration file needed to reduce data
303 mocchiut 1.1 // =============================================================
304 pam-fi 1.2 // if run OBT is > last calibration "expiration date"
305     // - search for new calibration packet
306     // - load calibration parameters (full + truncated)
307 pam-fi 1.5 if(p->VerboseMode() )cout << "Full pedestals for cluster finding:";
308 pam-fi 1.2 if(runheadtime > pedsig_time){
309     if( q2.Query_GL_TRK_CALIB(runheadtime,dbc) )throw -53;
310     pedsig_time = q2.TO_TIME;
311 pam-fi 1.7 if(q2.EV_ROOT_CALIBTRK1 != q2.EV_ROOT_CALIBTRK2)
312     printf("WARNING!! ---> EV_ROOT_CALIBTRK1=%d it's different from EV_ROOT_CALIBTRK2=%d \n\n",q2.EV_ROOT_CALIBTRK1,q2.EV_ROOT_CALIBTRK2);
313     if( q3.Query_GL_ROOT(q2.ID_ROOT_L0,dbc) )throw -51;
314 pam-fi 1.5 trkcalibfile = q3.PATH + q3.NAME;
315     if(p->VerboseMode()) cout << " >> Loading from LEVEL0 file: "<< trkcalibfile << endl;
316 pam-fi 1.2 if(trkcalibfile.CompareTo(filename)){
317 pam-fi 1.5 FileStat_t t;
318     if( gSystem->GetPathInfo(trkcalibfile.Data(),t) )throw -6;
319     f0_c=new TFile(trkcalibfile);
320     if ( !f0_c ) throw -6;
321 pam-fi 1.2 }else{
322 pam-fi 1.5 f0_c = f0;
323 pam-fi 1.2 };
324     if(p->VerboseMode()) {
325 pam-fi 1.7 cout << " calibration entries "<< q2.EV_ROOT_CALIBTRK1 << " " << q2.EV_ROOT_CALIBTRK2;
326     cout << " (from time "<< q2.FROM_TIME <<" to time "<< q2.TO_TIME <<")"<<endl;
327 pam-fi 1.2 };
328 pam-fi 1.7 pedsigbad_.FillACalibFrom(f0_c,q2.EV_ROOT_CALIBTRK1,q2.EV_ROOT_CALIBTRK2);
329 pam-fi 1.5 }else{
330     if(p->VerboseMode() )cout<< " already loaded. "<< endl;
331 pam-fi 1.2 };
332     if( p->DebugMode() ) for(int i=0; i<12; i++) cout << " DSP "<< i << " "<< pedsigbad_.pedestal[64][12][i] << endl;
333     // =============================================================
334     // retrieve calibration file needed to uncompress data
335     // =============================================================
336     // if the run was compressed using default calib
337     // load truncated pedestals from default
338     // otherwise reload them from on-line calibration
339     if(p->VerboseMode() )cout << "Truncated pedestals for uncompression:";
340 mocchiut 1.1 if(trk_calib_used==104 && last_trk_calib_used!=104){
341    
342 pam-fi 1.5 if(p->VerboseMode() )cout << " >> Loading default calibration " << endl;
343 pam-fi 1.7 if (q4.Query_GL_PARAM(runheadtime,7,dbc) )throw -52;
344 pam-fi 1.2 pedsigbad_.FillTCalibFrom(q4.PATH+q4.NAME);
345    
346     }else if( (trk_calib_used==1 | trk_calib_used==2) && last_trk_calib_used==104 ){
347    
348     if(p->VerboseMode()) cout << ">> Re-loading on-line calibration " << endl;
349 pam-fi 1.7 pedsigbad_.FillTCalibFrom(f0_c,q2.EV_ROOT_CALIBTRK1,q2.EV_ROOT_CALIBTRK2);
350 pam-fi 1.2 }else{
351     if(p->VerboseMode()) cout << ">> Using on-line calibration " << endl;
352     };
353     if( p->DebugMode() ) for(int i=0; i<12; i++) cout << " DSP "<< i << " "<< pedsigbad_.pedestal_t[64][12][i] << endl;
354     last_trk_calib_used = trk_calib_used;
355    
356 mocchiut 1.1 // =============================================================
357     // retrieve information about parameters to process LEVEL2
358     // =============================================================
359     //
360     // magnetic field
361     //
362 pam-fi 1.2 if(p->VerboseMode()) cout << "Loading parameter files : "<< endl;
363    
364 mocchiut 1.1 if(runheadtime > B_time){
365 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,1,dbc) )throw -52;
366 mocchiut 1.1 B_time = q4.TO_TIME;
367     l2_event->LoadField(q4.PATH+q4.NAME);
368 pam-fi 1.2 // l2_event->LoadField(q4.PATH+q4.NAME);
369 mocchiut 1.1 if(path_.error) throw -215;
370     };
371     //
372     // mip conversion parameters
373     //
374     if(runheadtime > mip_time){
375 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,2,dbc) )throw -52;
376 mocchiut 1.1 mip_time = q4.TO_TIME;
377     path_.FillWith(q4.PATH+q4.NAME);
378     readmipparam_();
379     if(path_.error) throw -212;
380     };
381     //
382     // charge correlation parameters
383     //
384     if(runheadtime > charge_time){
385 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,3,dbc) )throw -52;
386 mocchiut 1.1 charge_time = q4.TO_TIME;
387     path_.FillWith(q4.PATH+q4.NAME);
388     readchargeparam_();
389     if(path_.error) throw -213;
390     };
391     //
392     // eta p.f.a. parameters
393     //
394     if(runheadtime > eta_time){
395 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,4,dbc) )throw -52;
396 mocchiut 1.1 eta_time = q4.TO_TIME;
397     path_.FillWith(q4.PATH+q4.NAME);
398     readetaparam_();
399     if(path_.error) throw -214;
400     };
401     //
402     // alignment parameters
403     //
404     if(runheadtime > align_time){
405 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,5,dbc) )throw -52;
406 mocchiut 1.1 align_time = q4.TO_TIME;
407     path_.FillWith(q4.PATH+q4.NAME);
408     readalignparam_();
409     if(path_.error) throw -211;
410     };
411     //
412     // viking mask
413     //
414     if(runheadtime > mask_time){
415 pam-fi 1.7 if( q4.Query_GL_PARAM(runheadtime,6,dbc) )throw -52;
416 mocchiut 1.1 mask_time = q4.TO_TIME;
417     path_.FillWith(q4.PATH+q4.NAME);
418     readvkmask_();
419     if(path_.error) throw -210;
420     };
421    
422 pam-fi 1.11 TBenchmark *reduction = new TBenchmark();
423     if(p->VerboseMode())reduction->Start("reduction");
424     Int_t ntrk=0;
425 mocchiut 1.1 // ====================================================
426     // start looping on events cointained in the data file
427     // ====================================================
428 pam-fi 1.7 for (UInt_t re = evfrom; re < evfrom+nevents; re++){
429 mocchiut 1.1
430     ev_count++;
431    
432 pam-fi 1.13 if ( p->DebugMode() && re%100 == 0 && re > 0 ) cout << ".";
433 mocchiut 1.1
434 pam-fi 1.13 b_trk->GetEntry(re);
435     b_header->GetEntry(re);
436     pscu = header->GetPscuHeader();
437 pam-fi 1.7
438 pam-fi 1.13 if ( dbtime->DBabsTime(pscu->GetOrbitalTime()) > runtrailtime || dbtime->DBabsTime(pscu->GetOrbitalTime()) < runheadtime) {
439 pam-fi 1.7
440 pam-fi 1.14 if (p->VerboseMode()){
441     printf(" TrkCore - WARNING: event outside the run time window, skipping it\n");
442     cout << " OBT "<<pscu->GetOrbitalTime()<<" ABS-time "<<dbtime->DBabsTime(pscu->GetOrbitalTime())<<" RUN "<<runheadtime<<"-"<<runtrailtime<<endl;
443     };
444 pam-fi 1.13 }else{
445 pam-fi 1.7
446 pam-fi 1.13 p->ProcessEvent(l0_event);
447 pam-fi 1.7
448 pam-fi 1.13 // ----------------
449     // LEVEL1 output
450     // ----------------
451     if(p->get1){
452     if(p->ifroot1){ // root
453     l1_event->Clear();
454     l1_event->SetFromLevel1Struct(&level1event_);
455     // t_level1->Fill();
456     }else{ // hbook
457     throw -299;
458     };
459     };
460     // ----------------
461     // HOUGH output
462     // ----------------
463     if(p->geth){
464     if(p->ifrooth){ // root
465 pam-fi 1.15 lh_event->Delete();
466     lh_event->SetFromHoughStruct(&houghevent_);
467 pam-fi 1.13 }else{ // hbook
468     throw -299;
469     };
470     };
471     // ----------------
472     // LEVEL2 output
473     // ----------------
474     if(p->get2){
475 pam-fi 1.15 l2_event->Clear();
476     if(p->get1) l2_event->SetFromLevel2Struct(&level2event_,l1_event);//set references to level1
477     else l2_event->SetFromLevel2Struct(&level2event_);
478 pam-fi 1.10 // l2_event->Dump();
479 pam-fi 1.13 t_level2->Fill();
480     if( l2_event->ntrk()>0 )ntrk++;
481 mocchiut 1.1 };
482 pam-fi 1.13 };
483 mocchiut 1.1 }; // end loop on events
484 pam-fi 1.11 if(p->VerboseMode()){
485 pam-fi 1.13 cout << " >>> processed "<< ev_count <<" events"<< endl;
486     if(ev_count)cout << ntrk << " events with at least one track ("<<(Int_t)(100*ntrk/ev_count)<<"%)\n";
487     reduction->Show("reduction");
488 pam-fi 1.11 }
489     delete reduction;
490 pam-fi 1.7
491     delete dbtime;
492    
493 pam-fi 1.13 }; // end loop on runs
494    
495    
496 mocchiut 1.1 // ------------------------------------------------------------
497     // if reprocessing one run, copy all the events AFTER the run
498     // ------------------------------------------------------------
499     if( !(p->standalone) ){
500     for(UInt_t i=runinfo->GetLastEntry()+1; i<runinfo->GetFileEntries(); i++){
501     t_clone->GetEntry(i);
502     *l2_event = *l2_clone;
503     t_level2->Fill();
504     l2_event->Clear();
505     // COPY COPY COPY
506     };
507     };
508     // ---------------
509     // close the files
510     // ---------------
511     if(p->get2){
512 pam-fi 1.13 if( t_clone )t_clone->Delete("all");//delete old tree from file
513     if( !(p->standalone) )runinfo->Close();
514     f2->Write("Tracker");
515     if( t_level2 )t_level2->Delete(); //delete new tree from memory
516    
517     };
518 pam-fi 1.6
519 pam-fi 1.13 if( f0->IsOpen() ) f0->Close();
520 pam-fi 1.2 if( f0_c->IsOpen() )f0_c->Close();
521 pam-fi 1.13
522     // lh_event->Delete();
523     l1_event->Delete();
524     l2_event->Delete();
525     l2_clone->Delete();
526 pam-fi 1.2
527     return(p->ostatus);
528 mocchiut 1.1 }
529    

  ViewVC Help
Powered by ViewVC 1.1.23