/[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.12 - (hide annotations) (download)
Mon Oct 16 12:56:57 2006 UTC (18 years, 4 months ago) by pam-fi
Branch: MAIN
Changes since 1.11: +1 -1 lines
small bug fixed

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

  ViewVC Help
Powered by ViewVC 1.1.23