/[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.41 - (hide annotations) (download)
Wed Mar 6 14:18:32 2013 UTC (11 years, 11 months ago) by mocchiut
Branch: MAIN
Changes since 1.40: +29 -12 lines
Simulation flag (Valerio's style) added for Tracker

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 mocchiut 1.31 //int TrkCore(UInt_t run, TFile *f2, TSQLServer *dbc, int ncustom, char*vcustom[]) // EMILIANO
76     int TrkCore(UInt_t run, TFile *f2, GL_TABLES *glt, int ncustom, char*vcustom[]) // EMILIANO
77 mocchiut 1.1 {
78     // ---------------------------
79     // Define some basic variables
80     // ---------------------------
81    
82 pizzolot 1.36 TString filename = "";
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     TTree *physicsTree = 0;
88     TBranch *b_trk = 0;
89     TBranch *b_header = 0;
90     EventHeader *header = 0;
91     PscuHeader *pscu = 0;
92 pam-fi 1.22 TrkLevel0 *l0_event = new TrkLevel0();
93 mocchiut 1.1 // RunInfo
94     ItoRunInfo *runinfo = 0;
95 pam-fi 1.7 TArrayI *runlist = 0;
96     UInt_t from_run;
97     UInt_t to_run;
98 mocchiut 1.1 Bool_t reprocessing = false;
99     //LEVEL2 output - ROOT
100     TTree *t_level2 = 0;
101     TTree *t_clone = 0;
102     TrkLevel2 *l2_event = new TrkLevel2();
103     TrkLevel2 *l2_clone = new TrkLevel2();
104 pam-fi 1.13 //LEVEL1 output - ROOT
105     TrkLevel1 *l1_event = new TrkLevel1();
106 pam-fi 1.15 TrkHough *lh_event = new TrkHough();
107 mocchiut 1.1 // -----------------------
108     // -----------------------
109     // Handle input parameters
110     // (data)
111     //
112     // - open/create output files, determining the environment root/hbook from the estension
113     // - create the level1/level2 tree+branch/nt-uple
114     // -----------------------
115     // -----------------------
116     TrkProcess *p = new TrkProcess(run,f2);
117 pam-fi 1.26 if( p->HandleCustomPar(ncustom,vcustom) )return(p->ostatus);;
118 pam-fi 1.25 if( TrkParams::VerboseMode() )p->Dump();
119 mocchiut 1.1
120     // ===================================
121     // Open/Create level2 output file
122     // ===================================
123     if(p->get2){
124 pam-fi 1.13 //-------------------------------------------
125     // read RunInfo
126     //-------------------------------------------
127     if(!(p->standalone)){
128     // Open "RunInfo" tree and get list of run
129     runinfo = new ItoRunInfo(f2);
130     char *trkversion = TrkInfo(false);
131     int runinfo_error = runinfo->Update(run,"TRK",trkversion);
132     if( runinfo_error ) throw runinfo_error;
133     runlist = runinfo->GetRunList();
134     reprocessing = runinfo->IsReprocessing();//????
135 pam-fi 1.25 if(TrkParams::VerboseMode()){
136 pam-fi 1.13 cout << "#events "<< runinfo->GetFileEntries() << endl;// #events in the file
137     cout << "#runs "<< runinfo->GetRunEntries() << endl;// #runs in the file
138     cout << "1-st run "<< runlist->At(0) << endl;
139     cout << "last run "<< runlist->At(runinfo->GetRunEntries()-1) << endl;
140     cout << "1-st event "<< runinfo->GetFirstEntry() << endl;// first event of our run
141     cout << "last event+1 "<< runinfo->GetLastEntry() << endl;// first event AFTER the last event of our run
142     cout << "reprocessing "<< runinfo->IsReprocessing() << endl << endl;
143     };
144     };
145     //-------------------------------------------
146     //
147     //-------------------------------------------
148     // Take (if present) the old Tracker tree
149     t_clone = (TTree*)f2->Get("Tracker");
150     if( t_clone != NULL ) t_clone->SetBranchAddress("TrkLevel2",&l2_clone);
151     // Create NEW Tracker tree
152     t_level2 = new TTree("Tracker","PAMELA tracker level2 data ");
153 pam-fi 1.21 l2_event->Set(); // ****NBNBNBN*****
154 pam-fi 1.13 t_level2->Branch("TrkLevel2","TrkLevel2",&l2_event);
155     if(p->get1){
156 pam-fi 1.25 if(TrkParams::VerboseMode())cout << endl << "Requested LEVEL1 output" << endl;
157 pam-fi 1.21 l1_event->Set(); // ****NBNBNBN*****
158 pam-fi 1.13 t_level2->Branch("TrkLevel1","TrkLevel1",&l1_event);
159 pam-fi 1.16 }
160     if(p->geth){
161 pam-fi 1.25 if(TrkParams::VerboseMode())cout << endl << "Requested Hough-Transform output" << endl;
162 pam-fi 1.16 t_level2->Branch("TrkHough","TrkHough",&lh_event);
163 pam-fi 1.7 };
164 pam-fi 1.13 };
165 mocchiut 1.1
166    
167     // -------------------------------------------
168     // define runs to be processed/reprocessed
169     // -------------------------------------------
170     if(run == 0){
171 pam-fi 1.13 // reprocessing ALL runs
172     if(p->standalone)throw -298; // reprocessing not implemented
173     from_run = runlist->At(0);
174     to_run = runlist->At(runinfo->GetRunEntries()-1);
175 mocchiut 1.1 }else{
176     // processing/reprocessing ONE single run
177 pam-fi 1.13 from_run = run;
178     to_run = run;
179 mocchiut 1.1 };
180    
181     //
182     // init event counter
183     //
184     Int_t ev_count =0;
185     //
186     // create query-results objects
187     //
188 pam-fi 1.7 GL_RUN q1 = GL_RUN();
189 mocchiut 1.1 GL_TRK_CALIB q2 = GL_TRK_CALIB();
190     GL_ROOT q3 = GL_ROOT();
191     GL_PARAM q4 = GL_PARAM();
192    
193     // ------------------------------------------------------------
194     // if reprocessing one run, copy all the events BEFORE the run
195     // ------------------------------------------------------------
196     if( !(p->standalone) ){
197 pam-fi 1.13 for(UInt_t i=0; i<runinfo->GetFirstEntry(); i++){
198 mocchiut 1.37 if ( t_clone->GetEntry(i) <= 0 ) throw -36;// EM
199 pam-fi 1.13 *l2_event = *l2_clone;
200     t_level2->Fill();
201     l2_event->Clear();
202     // COPY COPY COPY
203     };
204 mocchiut 1.1 };
205     // ------------------------------------------------------------
206     // ------------------------------------------------------------
207 pam-fi 1.25 // START LOOP OVER RUNS TO PROCESS/REPROCESS
208 mocchiut 1.1 // ------------------------------------------------------------
209     // ------------------------------------------------------------
210 pam-fi 1.7 for(UInt_t idRun = from_run; idRun <= to_run; idRun++){
211 mocchiut 1.1
212 pam-fi 1.25 if(TrkParams::VerboseMode()) cout << endl<<" ========================= Run: "<< idRun << endl;
213 pam-fi 1.7 UInt_t runheadtime = 0;
214     UInt_t runtrailtime = 0;
215 pam-fi 1.32 UInt_t runheadobt = 0;
216     UInt_t runtrailobt = 0;
217 mocchiut 1.41 UInt_t abstime = 0;
218 mocchiut 1.1 UInt_t evfrom = 0;
219     UInt_t evto = 0;
220 pam-fi 1.7 UInt_t nevents = 0;
221     Int_t id_root_l0 =-1;
222 mocchiut 1.1 Int_t trk_calib_used = 0;
223 pam-fi 1.7
224 mocchiut 1.31 TString host = glt->CGetHost(); // EMILIANO
225     TString user = glt->CGetUser(); // EMILIANO
226     TString psw = glt->CGetPsw(); // EMILIANO
227     // TString host = "mysql://localhost/pamelaprod";
228     // TString user = "anonymous";
229     // TString psw = "";
230     // const char *pamdbhost=gSystem->Getenv("PAM_DBHOST");
231     // const char *pamdbuser=gSystem->Getenv("PAM_DBUSER");
232     // const char *pamdbpsw=gSystem->Getenv("PAM_DBPSW");
233     // if ( !pamdbhost ) pamdbhost = "";
234     // if ( !pamdbuser ) pamdbuser = "";
235     // if ( !pamdbpsw ) pamdbpsw = "";
236     // if ( strcmp(pamdbhost,"") ) host = pamdbhost;
237     // if ( strcmp(pamdbuser,"") ) user = pamdbuser;
238     // if ( strcmp(pamdbpsw,"") ) psw = pamdbpsw;
239 mocchiut 1.29 TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
240 mocchiut 1.31 stringstream myquery; // EMILIANO
241     myquery.str(""); // EMILIANO
242     myquery << "SET time_zone='+0:00'"; // EMILIANO
243 mocchiut 1.39 delete dbc->Query(myquery.str().c_str()); // EMILIANO
244 mocchiut 1.1 if(p->standalone){
245     // ==============================
246     // first query: retrieve run info
247     // ==============================
248 pam-fi 1.7 if (q1.Query_GL_RUN(idRun,dbc) )throw -50;
249     id_root_l0 = q1.ID_ROOT_L0;
250 mocchiut 1.1 runheadtime = q1.RUNHEADER_TIME;
251     runtrailtime = q1.RUNTRAILER_TIME;
252 pam-fi 1.32 runheadobt = q1.RUNHEADER_OBT;
253     runtrailobt = q1.RUNTRAILER_OBT;
254 pam-fi 1.7 evfrom = q1.EV_FROM;
255     evto = q1.EV_TO;
256 pam-fi 1.13 nevents = q1.NEVENTS;
257     trk_calib_used = q1.TRK_CALIB_USED;
258 mocchiut 1.1 }else{
259     // ==============================
260     // get run info from RunInfo tree
261     // ==============================
262 pam-fi 1.3 int runinfo_error = runinfo->GetRunInfo(idRun);
263 pam-fi 1.7 if( runinfo_error ) throw runinfo_error;
264     id_root_l0 = runinfo->ID_ROOT_L0;
265 mocchiut 1.1 runheadtime = runinfo->RUNHEADER_TIME;
266     runtrailtime = runinfo->RUNTRAILER_TIME;
267 pam-fi 1.32 runheadobt = runinfo->RUNHEADER_OBT;
268     runtrailobt = runinfo->RUNTRAILER_OBT;
269 pam-fi 1.7 evfrom = runinfo->EV_FROM;
270     evto = runinfo->EV_TO;
271 pam-fi 1.13 nevents = runinfo->NEVENTS;
272     trk_calib_used = runinfo->TRK_CALIB_USED;
273 pam-fi 1.24
274 mocchiut 1.1 };
275 pam-fi 1.24
276 mocchiut 1.1 //
277 pam-fi 1.25 if(TrkParams::VerboseMode()){
278 pam-fi 1.7 cout << "ROOT file ID "<< id_root_l0 << endl;
279 mocchiut 1.1 cout << "RunHeader time "<< runheadtime << endl;
280     cout << "RunTrailer time "<< runtrailtime << endl;
281 pam-fi 1.7 cout << " from event "<< evfrom << endl;
282     cout << " to event "<< evto << endl;
283 pam-fi 1.24 cout << " num. events "<< nevents << endl;
284     cout << "trk-calibration used "<< trk_calib_used << endl;
285 mocchiut 1.1 };
286     // ========================================================
287     // second query: search the LEVEL0 file that contains idRun
288     // ========================================================
289     TString lastfilename = filename;
290 pam-fi 1.7 if( q3.Query_GL_ROOT(id_root_l0,dbc) )throw -51;
291 mocchiut 1.1 filename = q3.PATH + q3.NAME;
292     // ========================================================
293     // Open the input LEVEL0 data file
294     // ========================================================
295     if(filename.CompareTo(lastfilename)){
296     if(!lastfilename.IsNull())f0->Close();
297     //if( debug ) cout << "Opening LEVEL0 file: "<< filename << endl;
298 pam-fi 1.25 if(TrkParams::VerboseMode()) cout << "Opening LEVEL0 file: "<< filename << endl;
299 pam-fi 1.2 FileStat_t t;
300     if( gSystem->GetPathInfo(filename.Data(),t) )throw -6;
301 mocchiut 1.39 if ( f0 ) f0->Close();
302 mocchiut 1.1 f0 = new TFile(filename);
303     if ( !f0 ) throw -6;
304     physicsTree = (TTree*)f0->Get("Physics"); if(!physicsTree) throw -7;
305     b_header = physicsTree ->GetBranch("Header"); if(!b_header ) throw -8;
306     b_trk = physicsTree ->GetBranch("Tracker"); if(!b_trk ) throw -203;
307 pam-fi 1.22 l0_event->Set();
308     physicsTree->SetBranchAddress("Tracker" ,l0_event->GetPointerToTrackerEvent());
309 mocchiut 1.1 physicsTree->SetBranchAddress("Header",&header);
310    
311 pam-fi 1.7 nentries = physicsTree->GetEntries();
312 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
313 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
314     // 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
315     // the condition ( nentries < (evto+1) ) is satisfied and DV exit with error even if the error is only in the DB.
316 pam-fi 1.7
317 mocchiut 1.1 };
318    
319 pam-fi 1.7 GL_TIMESYNC *dbtime = new GL_TIMESYNC(id_root_l0,"ID",dbc);
320    
321 pam-fi 1.24
322     // =============================================================
323     // retrieve information about parameters to process LEVEL2
324     // =============================================================
325    
326 mocchiut 1.41 TrkParams::Set(runinfo->GetGL_RUN(),dbc);
327     for(int i=0; i<p->npar; i++){
328     if(TrkParams::VerboseMode())cout<<" ((( force parameters from input path )))"<<endl;
329     TrkParams::Set(p->parpath[i],p->partype[i]);
330     }
331 pam-fi 1.25
332 mocchiut 1.41 TrkParams::Load();
333     if( !TrkParams::IsLoaded() )throw -52;
334 pam-fi 1.24
335     // =============================================================
336     // retrieve calibration file needed to reduce data
337     // =============================================================
338    
339     TrkParams::SetCalib(runinfo->GetGL_RUN(),dbc);
340     TrkParams::LoadCalib( );
341     if( !TrkParams::CalibIsLoaded() )throw -53;
342    
343 pam-fi 1.11 TBenchmark *reduction = new TBenchmark();
344 pam-fi 1.25 if(TrkParams::VerboseMode())reduction->Start("reduction");
345 pam-fi 1.11 Int_t ntrk=0;
346 mocchiut 1.1 // ====================================================
347     // start looping on events cointained in the data file
348     // ====================================================
349 mocchiut 1.41 if ( TrkParams::GetSimuFlag() ){
350     abstime = runheadtime;
351     } else {
352     if(dbc){
353     dbc->Close();
354     delete dbc;
355     }
356     }
357    
358 pam-fi 1.32 for (UInt_t re = evfrom+min(p->nskip,nevents); re < evfrom+nevents; re++){
359 mocchiut 1.1
360     ev_count++;
361 pam-fi 1.32
362    
363 pam-fi 1.25 // if ( TrkParams::DebugMode() && re%100 == 0 && re > 0 ) cout << ".";
364 mocchiut 1.1
365 mocchiut 1.37 if ( b_trk->GetEntry(re) <= 0 ) throw -36;//EM
366     if ( b_header->GetEntry(re) <= 0 ) throw -36;//EM
367 pam-fi 1.13 pscu = header->GetPscuHeader();
368 mocchiut 1.41
369     // =============================================================
370     // The following 6 lines have been moved here by VALERIO.
371     if(TrkParams::GetSimuFlag()){
372     abstime = runheadtime + (int) floor(0.03*(re-evfrom)); //If simulated data we need to assign a fake abstime. 30ms between each event.
373     if(TrkParams::VerboseMode())cout << "Event: " << re-evfrom << " - Attempting to retrieve Mask Info for abstime=" << abstime << endl;
374     if(!TrkParams::Set(runinfo->GetGL_RUN(),dbc,6,abstime))throw -52; // Setting to load mask (VALERIO)
375     TrkParams::Load(6);
376     if( !TrkParams::IsLoaded() )throw -52;
377     if(TrkParams::VerboseMode())cout << "Mask Info for abstime=" << abstime << " retrieved" << endl;
378     }
379    
380 pam-fi 1.25 if( TrkParams::DebugMode() )cout << ">>> "<<ev_count-1<<" @ OBT "<<pscu->GetOrbitalTime()<<endl;
381 pam-fi 1.23
382 mocchiut 1.35 if ( dbtime->DBabsTime(pscu->GetOrbitalTime()) > (runtrailtime+1) || dbtime->DBabsTime(pscu->GetOrbitalTime()) < (runheadtime-1)) {
383 pam-fi 1.7
384 pam-fi 1.25 if (TrkParams::VerboseMode()){
385 pam-fi 1.14 printf(" TrkCore - WARNING: event outside the run time window, skipping it\n");
386 pam-fi 1.32 cout << " OBT "<<pscu->GetOrbitalTime()<<" RUN "<<runheadobt<<"-"<<runtrailobt<<" ABS-time "<<dbtime->DBabsTime(pscu->GetOrbitalTime())<<" RUN "<<runheadtime<<"-"<<runtrailtime<<endl;
387 pam-fi 1.14 };
388 pam-fi 1.13 }else{
389 pam-fi 1.25 if ( TrkParams::DebugMode() )
390 pam-fi 1.23 printf("\n-------------------------------\nevent %d\n",re-evfrom);
391 pam-fi 1.16
392 pam-fi 1.13 p->ProcessEvent(l0_event);
393 pam-fi 1.7
394 pam-fi 1.13 // ----------------
395     // LEVEL1 output
396     // ----------------
397     if(p->get1){
398     if(p->ifroot1){ // root
399     l1_event->Clear();
400 pam-fi 1.22 // l1_event->SetFromLevel1Struct(&level1event_,p->full1);
401     l1_event->SetFromLevel1Struct(p->full1);
402 pam-fi 1.13 // t_level1->Fill();
403 pam-fi 1.23 }else{ // hbook
404 pam-fi 1.13 throw -299;
405     };
406     };
407     // ----------------
408     // HOUGH output
409     // ----------------
410     if(p->geth){
411     if(p->ifrooth){ // root
412 pam-fi 1.23 lh_event->Delete();
413     lh_event->SetFromHoughStruct(&houghevent_);
414 pam-fi 1.13 }else{ // hbook
415     throw -299;
416     };
417     };
418     // ----------------
419     // LEVEL2 output
420     // ----------------
421     if(p->get2){
422 pam-fi 1.23 l2_event->Clear();
423     if(p->get1){
424     l2_event->SetFromLevel2Struct(&level2event_,l1_event);//set references to level1
425     }else{
426     l2_event->SetFromLevel2Struct(&level2event_);
427     }
428 pam-fi 1.10 // l2_event->Dump();
429 pam-fi 1.23 t_level2->Fill();
430     if( l2_event->ntrk()>0 )ntrk++;
431 pam-fi 1.27 if(TrkParams::VerboseMode())l2_event->Dump();
432 mocchiut 1.1 };
433 pam-fi 1.13 };
434 mocchiut 1.1 }; // end loop on events
435 pam-fi 1.25 if(TrkParams::VerboseMode()){
436 pam-fi 1.13 cout << " >>> processed "<< ev_count <<" events"<< endl;
437     if(ev_count)cout << ntrk << " events with at least one track ("<<(Int_t)(100*ntrk/ev_count)<<"%)\n";
438     reduction->Show("reduction");
439 pam-fi 1.11 }
440     delete reduction;
441 pam-fi 1.7
442     delete dbtime;
443    
444 pam-fi 1.13 }; // end loop on runs
445    
446    
447 mocchiut 1.1 // ------------------------------------------------------------
448     // if reprocessing one run, copy all the events AFTER the run
449     // ------------------------------------------------------------
450     if( !(p->standalone) ){
451     for(UInt_t i=runinfo->GetLastEntry()+1; i<runinfo->GetFileEntries(); i++){
452 mocchiut 1.37 if ( t_clone->GetEntry(i) <= 0 ) throw -36;//EM
453 mocchiut 1.1 *l2_event = *l2_clone;
454     t_level2->Fill();
455     l2_event->Clear();
456     // COPY COPY COPY
457     };
458     };
459     // ---------------
460     // close the files
461     // ---------------
462     if(p->get2){
463 mocchiut 1.40 if( t_clone )t_clone->Delete("all");//delete old tree from file
464 pam-fi 1.13 if( !(p->standalone) )runinfo->Close();
465 mocchiut 1.40 // t_level2->Write("Tracker", TObject::kOverwrite);
466     f2->Write("Tracker", TObject::kOverwrite);
467     if( t_level2 )t_level2->Delete(); //delete new tree from memory
468 pam-fi 1.13
469     };
470 pam-fi 1.6
471 pam-fi 1.24 if(f0) if(f0->IsOpen()) f0->Close();
472     // if( f0_c->IsOpen() )f0_c->Close();
473 pam-fi 1.13
474     // lh_event->Delete();
475     l1_event->Delete();
476     l2_event->Delete();
477     l2_clone->Delete();
478 pam-fi 1.2
479     return(p->ostatus);
480 mocchiut 1.1 }
481    

  ViewVC Help
Powered by ViewVC 1.1.23