/[PAMELA software]/DarthVader/OrbitalInfo/src/OrbitalInfoCore.cpp
ViewVC logotype

Annotation of /DarthVader/OrbitalInfo/src/OrbitalInfoCore.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.13 - (hide annotations) (download)
Tue Mar 20 20:56:30 2007 UTC (17 years, 8 months ago) by mocchiut
Branch: MAIN
Changes since 1.12: +1 -0 lines
Problems with classes streamerinfo fixed

1 mocchiut 1.1 //
2     // C/C++ headers
3     //
4     #include <fstream>
5     #include <string.h>
6     #include <iostream>
7     #include <cstring>
8     #include <stdio.h>
9     //
10     // ROOT headers
11     //
12     #include <TTree.h>
13     #include <TClassEdit.h>
14     #include <TObject.h>
15     #include <TList.h>
16 mocchiut 1.5 #include <TArrayI.h>
17 mocchiut 1.1 #include <TSystem.h>
18     #include <TSystemDirectory.h>
19     #include <TString.h>
20     #include <TFile.h>
21     #include <TClass.h>
22     #include <TSQLServer.h>
23     #include <TSQLRow.h>
24     #include <TSQLResult.h>
25     //
26 mocchiut 1.8 // RunInfo header
27     //
28     #include <RunInfo.h>
29     #include <GLTables.h>
30     //
31 mocchiut 1.1 // YODA headers
32     //
33     #include <PamelaRun.h>
34     #include <PscuHeader.h>
35     #include <PscuEvent.h>
36     #include <EventHeader.h>
37     //
38     // This program headers
39     //
40     #include <OrbitalInfo.h>
41 mocchiut 1.7 #include <OrbitalInfoVerl2.h>
42 mocchiut 1.1 #include <OrbitalInfoCore.h>
43    
44     using namespace std;
45    
46     //
47     // CORE ROUTINE
48     //
49     //
50 mocchiut 1.5 int OrbitalInfoCore(UInt_t run, TFile *file, TSQLServer *dbc, Int_t OrbitalInfoargc, char *OrbitalInfoargv[]){
51 mocchiut 1.7
52     // // Temporary check to use igrf routines. We need dat files in the
53     // // current directory.
54     // fstream igrfdat1("igrf05.dat");
55     // fstream igrfdat2("igrf05s.dat");
56     // if( (!igrfdat1) && (!igrfdat2)) {
57     // cerr<<"\n**************************************\n"
58     // <<"igrf05.dat or igrf05s.dat not in the current directory. Exiting.\n"
59     // <<"**************************************\n";
60     // exit(EXIT_FAILURE);
61     // }
62     // igrfdat1.close();
63     // igrfdat2.close();
64     // // end of temporary code
65     // //
66    
67 mocchiut 1.1 Int_t i = 0;
68     //
69     TString processFolder = "OrbitalInfoFolder";
70     //
71     // Set these to true to have a very verbose output.
72     //
73     Bool_t debug = false;
74     //
75     Bool_t verbose = false;
76    
77     if ( OrbitalInfoargc > 0 ){
78     i = 0;
79     while ( i < OrbitalInfoargc ){
80     if ( !strcmp(OrbitalInfoargv[i],"-processFolder") ) {
81     if ( OrbitalInfoargc < i+1 ){
82     throw -3;
83     };
84     processFolder = (TString)OrbitalInfoargv[i+1];
85     i++;
86     };
87     if ( (!strcmp(OrbitalInfoargv[i],"--debug")) || (!strcmp(OrbitalInfoargv[i],"-g")) ) {
88     verbose = true;
89     };
90     if ( (!strcmp(OrbitalInfoargv[i],"--verbose")) || (!strcmp(OrbitalInfoargv[i],"-v")) ) {
91     verbose = true;
92     };
93     i++;
94     };
95     };
96     //
97     const char* outDir = gSystem->DirName(gSystem->DirName(file->GetPath()));
98     //
99     TTree *OrbitalInfotr = 0;
100 mocchiut 1.5 UInt_t nevents = 0;
101 mocchiut 1.1 //
102     // variables needed to reprocess data
103     //
104 mocchiut 1.6 Long64_t maxsize = 10000000000LL;
105     TTree::SetMaxTreeSize(maxsize);
106     //
107 mocchiut 1.1 TString OrbitalInfoversion;
108     ItoRunInfo *runinfo = 0;
109 mocchiut 1.5 TArrayI *runlist = 0;
110 mocchiut 1.1 TTree *OrbitalInfotrclone = 0;
111     Bool_t reproc = false;
112     Bool_t reprocall = false;
113     UInt_t nobefrun = 0;
114     UInt_t noaftrun = 0;
115     UInt_t numbofrun = 0;
116     stringstream ftmpname;
117     TString fname;
118 mocchiut 1.5 UInt_t totfileentries = 0;
119     UInt_t idRun = 0;
120 mocchiut 1.1 //
121     // variables needed to handle error signals
122     //
123     Int_t code = 0;
124     Int_t sgnl;
125     //
126     // OrbitalInfo classes
127     //
128     OrbitalInfo *orbitalinfo = new OrbitalInfo();
129     OrbitalInfo *orbitalinfoclone = new OrbitalInfo();
130     //
131     // define variables for opening and reading level0 file
132     //
133     TFile *l0File = 0;
134     TTree *l0tr = 0;
135 mocchiut 1.2 // EM: open also header branch
136     TBranch *l0head = 0;
137     pamela::EventHeader *eh = 0;
138     pamela::PscuHeader *ph = 0;
139     // end EM
140 mocchiut 1.1 //
141     // Define other basic variables
142     //
143     UInt_t procev = 0;
144     stringstream file2;
145     stringstream file3;
146     stringstream qy;
147     Int_t totevent = 0;
148 mocchiut 1.5 UInt_t atime = 0;
149     UInt_t re = 0;
150 mocchiut 1.7
151     // Position
152     Float_t lon, lat, alt;
153    
154     //
155     // IGRF stuff
156     //
157     float dimo = 0.0; // dipole moment (computed from dat files)
158     float bnorth, beast, bdown, babs;
159     float xl; // L value
160     float icode; // code value for L accuracy (see fortran code)
161     float bab1; // What's the difference with babs?
162     float stps = 0.005; // step size for field line tracing
163     float bdel = 0.01; // required accuracy
164     float bequ; // equatorial b value (also called b_0)
165     bool value = 0; // false if bequ is not the minimum b value
166     float rr0; // equatorial radius normalized to earth radius
167    
168 mocchiut 1.1 //
169     // Working filename
170     //
171     TString outputfile;
172     stringstream name;
173     name.str("");
174     name << outDir << "/";
175     //
176     // temporary file and folder
177     //
178     TFile *tempfile = 0;
179     TTree *tempOrbitalInfo = 0;
180     stringstream tempname;
181     stringstream OrbitalInfofolder;
182     tempname.str("");
183     tempname << outDir;
184     tempname << "/" << processFolder.Data();
185     OrbitalInfofolder.str("");
186     OrbitalInfofolder << tempname.str().c_str();
187     gSystem->MakeDirectory(OrbitalInfofolder.str().c_str());
188     tempname << "/OrbitalInfotree_run";
189     tempname << run << ".root";
190     //
191     // DB classes
192     //
193     GL_ROOT *glroot = new GL_ROOT();
194 mocchiut 1.5 GL_TIMESYNC *dbtime = 0;
195 mocchiut 1.7 GL_TLE *gltle = new GL_TLE();
196     // Initialize fortran routines!!!
197     Int_t ltp2 = 0;
198     Int_t ltp3 = 0;
199     Int_t uno = 1;
200     char *niente = " ";
201     GL_PARAM *glparam = new GL_PARAM();
202 mocchiut 1.10 GL_PARAM *glparam2 = new GL_PARAM();
203 mocchiut 1.7 Int_t parerror=glparam->Query_GL_PARAM(1,301,dbc); // parameters stored in DB in GL_PRAM table
204     if ( parerror<0 ) {
205     code = parerror;
206     goto closeandexit;
207     };
208     ltp2 = (Int_t)(glparam->PATH+glparam->NAME).Length();
209     if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
210     //
211     parerror=glparam2->Query_GL_PARAM(1,302,dbc); // parameters stored in DB in GL_PRAM table
212     if ( parerror<0 ) {
213     code = parerror;
214     goto closeandexit;
215     };
216     ltp3 = (Int_t)(glparam2->PATH+glparam2->NAME).Length();
217     if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam2->PATH+glparam2->NAME).Data());
218     //
219     initize_((char *)niente,&uno,(char *)(glparam->PATH+glparam->NAME).Data(),&ltp2,(char *)(glparam2->PATH+glparam2->NAME).Data(),&ltp3);
220    
221     //
222     // End IGRF stuff//
223     //
224    
225 mocchiut 1.1 //
226     // Let's start!
227     //
228     // As a first thing we must check what we have to do: if run = 0 we must process all events in the file has been passed
229     // if run != 0 we must process only that run but first we have to check if the tree MyDetector2 already exist in the file
230     // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
231     //
232 mocchiut 1.5 if ( run == 0 ) reproc = true;
233 mocchiut 1.1 //
234     //
235     // Output file is "outputfile"
236     //
237     if ( !file->IsOpen() ){
238     //printf(" OrbitalInfo - ERROR: cannot open file for writing\n");
239     throw -901;
240     };
241     //
242     // Retrieve GL_RUN variables from the level2 file
243     //
244     OrbitalInfoversion = OrbitalInfoInfo(false); // we should decide how to handle versioning system
245     //
246     // create an interface to RunInfo called "runinfo"
247     //
248     runinfo = new ItoRunInfo(file);
249     //
250     // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
251     //
252     sgnl = 0;
253     sgnl = runinfo->Update(run, "ORB", OrbitalInfoversion);
254     //sgnl = runinfo->Read(run);
255    
256     if ( sgnl ){
257     //printf("OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
258     code = sgnl;
259     goto closeandexit;
260     } else {
261     sgnl = 0;
262     };
263     //
264     // number of events in the file BEFORE the first event of our run
265     //
266     nobefrun = runinfo->GetFirstEntry();
267     //
268     // total number of events in the file
269     //
270     totfileentries = runinfo->GetFileEntries();
271     //
272     // first file entry AFTER the last event of our run
273     //
274     noaftrun = runinfo->GetLastEntry() + 1;
275     //
276     // number of run to be processed
277     //
278     numbofrun = runinfo->GetNoRun();
279     //
280     // Try to access the OrbitalInfo tree in the file, if it exists we are reprocessing data if not we are processing a new run
281     //
282     OrbitalInfotrclone = (TTree*)file->Get("OrbitalInfo");
283     //
284     if ( !OrbitalInfotrclone ){
285     //
286     // tree does not exist, we are not reprocessing
287     //
288     reproc = false;
289 mocchiut 1.5 if ( run == 0 ){
290 mocchiut 1.1 if (verbose) printf(" OrbitalInfo - WARNING: you are reprocessing data but OrbitalInfo tree does not exist!\n");
291     }
292 mocchiut 1.5 if ( runinfo->IsReprocessing() && run != 0 ) {
293 mocchiut 1.1 if (verbose) printf(" OrbitalInfo - WARNING: it seems you are not reprocessing data but OrbitalInfo\n versioning information already exists in RunInfo.\n");
294     }
295     } else {
296     //
297     // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
298     //
299 mocchiut 1.6 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
300 mocchiut 1.1 reproc = true;
301     //
302     //
303     if (verbose) printf("\n Preparing the pre-processing...\n");
304     //
305 mocchiut 1.5 if ( run == 0 ){
306 mocchiut 1.1 //
307     // we are reprocessing all the file
308     // if we are reprocessing everything we don't need to copy any old event and we can just work with the new tree and delete the old one immediately
309     //
310     reprocall = true;
311     //
312     if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing all runs\n");
313     //
314     } else {
315     //
316     // we are reprocessing a single run, we must copy to the new tree the events in the file which preceed the first event of the run
317     //
318     reprocall = false;
319     //
320 mocchiut 1.5 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing run number %u \n",run);
321 mocchiut 1.1 //
322     // copying old tree to a new file
323     //
324     tempfile = new TFile(tempname.str().c_str(),"RECREATE");
325     tempOrbitalInfo = OrbitalInfotrclone->CloneTree(-1,"fast");
326     tempOrbitalInfo->SetName("OrbitalInfo-old");
327     tempfile->Write();
328     tempfile->Close();
329     }
330     //
331     // Delete the old tree from old file and memory
332     //
333     OrbitalInfotrclone->Delete("all");
334     //
335     if (verbose) printf(" ...done!\n");
336     //
337     };
338     //
339     // create mydetector tree mydect
340     //
341     file->cd();
342     OrbitalInfotr = new TTree("OrbitalInfo-new","PAMELA OrbitalInfo data");
343 mocchiut 1.6 OrbitalInfotr->SetAutoSave(900000000000000LL);
344 mocchiut 1.1 OrbitalInfotr->Branch("OrbitalInfo","OrbitalInfo",&orbitalinfo);
345     //
346     if ( reproc && !reprocall ){
347     //
348     // open new file and retrieve also tree informations
349     //
350     tempfile = new TFile(tempname.str().c_str(),"READ");
351     OrbitalInfotrclone = (TTree*)tempfile->Get("OrbitalInfo-old");
352 mocchiut 1.6 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
353 mocchiut 1.1 OrbitalInfotrclone->SetBranchAddress("OrbitalInfo",&orbitalinfoclone);
354     //
355     if ( nobefrun > 0 ){
356     if (verbose){
357 mocchiut 1.7 printf("\n Pre-processing: copying events from the old tree before the processed run\n");
358     printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
359     printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
360 mocchiut 1.1 }
361     for (UInt_t j = 0; j < nobefrun; j++){
362     //
363     OrbitalInfotrclone->GetEntry(j);
364     //
365     // copy orbitalinfoclone to mydec
366     //
367 mocchiut 1.3 orbitalinfo->Clear();
368 mocchiut 1.5 //
369 mocchiut 1.1 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
370     //
371     // Fill entry in the new tree
372     //
373     OrbitalInfotr->Fill();
374     //
375     };
376     if (verbose) printf(" Finished successful copying!\n");
377     };
378     };
379     //
380     // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
381     //
382     runlist = runinfo->GetRunList();
383     //
384     // Loop over the run to be processed
385     //
386     for (UInt_t irun=0; irun < numbofrun; irun++){
387     //
388     // retrieve the first run ID to be processed using the RunInfo list
389     //
390     idRun = runlist->At(irun);
391     if (verbose){
392     printf("\n\n\n ####################################################################### \n");
393     printf(" PROCESSING RUN NUMBER %i \n",(int)idRun);
394     printf(" ####################################################################### \n\n\n");
395     }
396     //
397 mocchiut 1.5 runinfo->ID_ROOT_L0 = 0;
398 mocchiut 1.1 //
399     // store in the runinfo class the GL_RUN variables for our run
400     //
401     sgnl = 0;
402     sgnl = runinfo->GetRunInfo(idRun);
403     if ( sgnl ){
404 mocchiut 1.5 if ( debug ) printf("\n OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
405 mocchiut 1.1 code = sgnl;
406     goto closeandexit;
407     } else {
408     sgnl = 0;
409     };
410     //
411     // now you can access that variables using the RunInfo class this way runinfo->ID_REG_RUN
412     //
413 mocchiut 1.5 if ( runinfo->ID_ROOT_L0 == 0 ){
414     if ( debug ) printf("\n OrbitalInfo - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
415 mocchiut 1.1 code = -5;
416     goto closeandexit;
417     };
418     //
419 mocchiut 1.5 // prepare the timesync for the db
420     //
421     dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
422     //
423 mocchiut 1.1 // Search in the DB the path and name of the LEVEL0 file to be processed.
424     //
425 mocchiut 1.5 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
426 mocchiut 1.1 //
427     ftmpname.str("");
428     ftmpname << glroot->PATH.Data() << "/";
429     ftmpname << glroot->NAME.Data();
430     fname = ftmpname.str().c_str();
431     //
432     // print out informations
433     //
434 mocchiut 1.5 totevent = runinfo->NEVENTS;
435 mocchiut 1.1 if (verbose){
436     printf("\n LEVEL0 data file: %s \n",fname.Data());
437 mocchiut 1.5 printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
438     printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
439     printf(" %i events to be processed for run %u: from %i to %i \n\n",totevent,idRun,runinfo->EV_FROM,runinfo->EV_FROM+totevent);
440 mocchiut 1.1 }//
441     // Open Level0 file
442     l0File = new TFile(fname.Data());
443     if ( !l0File ) {
444 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: problems opening Level0 file\n");
445 mocchiut 1.1 code = -6;
446     goto closeandexit;
447     };
448     l0tr = (TTree*)l0File->Get("Physics");
449     if ( !l0tr ) {
450 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: no Physics tree in Level0 file\n");
451 mocchiut 1.1 l0File->Close();
452     code = -7;
453     goto closeandexit;
454     };
455 mocchiut 1.2 // EM: open header branch as well
456     l0head = l0tr->GetBranch("Header");
457     if ( !l0head ) {
458 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: no Header branch in Level0 tree\n");
459 mocchiut 1.2 l0File->Close();
460     code = -8;
461     goto closeandexit;
462     };
463     l0tr->SetBranchAddress("Header", &eh);
464     // end EM
465 mocchiut 1.5 nevents = l0head->GetEntries();
466 mocchiut 1.1 //
467     if ( nevents < 1 ) {
468 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: Level0 file is empty\n\n");
469 mocchiut 1.1 l0File->Close();
470     code = -11;
471     goto closeandexit;
472     };
473     //
474 mocchiut 1.5 if ( runinfo->EV_TO > nevents-1 ) {
475     if ( debug ) printf(" OrbitalInfo - ERROR: too few entries in the registry tree\n");
476 mocchiut 1.1 l0File->Close();
477     code = -12;
478     goto closeandexit;
479     };
480     //
481     // run over all the events of the run
482     //
483     if (verbose) printf("\n Ready to start! \n\n Processed events: \n\n");
484     //
485 mocchiut 1.5 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
486 mocchiut 1.1 //
487 mocchiut 1.5 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
488 mocchiut 1.1 //
489 mocchiut 1.5 l0head->GetEntry(re);
490 mocchiut 1.1 //
491     // absolute time of this event
492     //
493 mocchiut 1.5 ph = eh->GetPscuHeader();
494     atime = dbtime->DBabsTime(ph->GetOrbitalTime());
495 mocchiut 1.1 //
496     // paranoid check
497     //
498     if ( (atime > runinfo->RUNTRAILER_TIME) || (atime < runinfo->RUNHEADER_TIME) ) {
499     if (verbose) printf(" OrbitalInfo - WARNING: event at time outside the run time window, skipping it\n");
500 mocchiut 1.7 debug = true;
501     continue;
502     }
503 mocchiut 1.1 //
504     procev++;
505     //
506     // start processing
507     //
508 mocchiut 1.3 orbitalinfo->Clear();
509 mocchiut 1.5 //
510     // CHANGE HERE!!!!
511     //
512     orbitalinfo->absTime = atime;
513 mocchiut 1.2 // EM: add OBT and plt_num infos from the header
514     ph = eh->GetPscuHeader();
515     orbitalinfo->pkt_num = ph->GetCounter();
516     orbitalinfo->OBT = ph->GetOrbitalTime();
517 mocchiut 1.7
518     // If the absolute time of the event overpass the time of the
519     // tle, get a new tle. GL_TLE::GetToTime() default to zero.
520     //
521     // I also use this condition to compute the dipole moment dimo.
522     // It's really redundant to compute it so often because
523     // probably it will not change at all. But the overhead is
524     // minimum.
525     float jyear=0;
526    
527     if(atime >= gltle->GetToTime()) {
528 mocchiut 1.9 if ( !gltle->Query(atime, dbc) ){
529    
530     // Compute the magnetic dipole moment.
531     UInt_t year, month, day, hour, min, sec;
532    
533     TTimeStamp t = TTimeStamp(atime, kTRUE);
534     t.GetDate(kTRUE, 0, &year, &month, &day);
535     t.GetTime(kTRUE, 0, &hour, &min, &sec);
536     jyear = (float) year
537     + (month*31.+ (float) day)/365.
538     + (hour*3600.+min*60.+(float)sec)/(24*3600*365.);
539    
540     feldcof_(&jyear, &dimo); // get dipole moment for year
541     } else {
542     code = -56;
543     goto closeandexit;
544     };
545 mocchiut 1.7 }
546    
547     // Propagate the orbit from the tle time to atime, using SGP(D)4.
548     cCoordGeo coo = getCoo(atime, gltle->GetFromTime(), gltle->GetTle());
549    
550     // Build coordinates in the right range. We want to convert,
551     // longitude from (0, 2*pi) to (-180deg, 180deg). Altitude is
552     // in meters.
553     lon = (coo.m_Lon > M_PI) ? rad2deg(coo.m_Lon - 2*M_PI) : rad2deg(coo.m_Lon);
554     lat = rad2deg(coo.m_Lat);
555     alt = coo.m_Alt;
556    
557    
558     // if((lon>180) || (lon<-180) || (lat>90) || (lat<-90) || (alt<0))
559     // continue;
560     if( lon<180 && lon>-180 && lat<90 && lat>-90 && alt>0 ){
561    
562     orbitalinfo->lon = lon;
563     orbitalinfo->lat = lat;
564     orbitalinfo->alt = alt ;
565    
566     // compute mag field components and L shell.
567     feldg_(&lat, &lon, &alt, &bnorth, &beast, &bdown, &babs);
568     shellg_(&lat, &lon, &alt, &dimo, &xl, &icode, &bab1);
569     findb0_(&stps, &bdel, &value, &bequ, &rr0);
570    
571     orbitalinfo->Bnorth = bnorth;
572     orbitalinfo->Beast = beast;
573     orbitalinfo->Bdown = bdown;
574     orbitalinfo->Babs = babs;
575     orbitalinfo->BB0 = babs/bequ;
576     orbitalinfo->L = xl;
577    
578     // Set Stormer vertical cutoff using L shell.
579     orbitalinfo->cutoff[0] = 14.9/(xl*xl);
580     };
581 pamelaprod 1.11
582     //
583     // inclination
584     //
585    
586    
587 mocchiut 1.7
588 mocchiut 1.12 // orbitalinfo->yaw = incl->;
589 mocchiut 1.2 // end EM
590 mocchiut 1.1 OrbitalInfotr->Fill();
591 mocchiut 1.7 }
592    
593 mocchiut 1.1 //
594     // Here you may want to clear some variables before processing another run
595     //
596 mocchiut 1.5 delete dbtime;
597 mocchiut 1.13 // l0File->Close();
598 mocchiut 1.1 }; // process all the runs
599     //
600     if (verbose) printf("\n Finished processing data \n");
601     //
602     closeandexit:
603     //
604     // we have finished processing the run(s). If we processed a single run now we must copy all the events after our run from the old tree to the new one and delete the old tree.
605     //
606     if ( !reprocall && reproc && code >= 0 ){
607     if ( totfileentries > noaftrun ){
608     if (verbose){
609     printf("\n Post-processing: copying events from the old tree after the processed run\n");
610     printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
611     printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
612     }
613     for (UInt_t j = noaftrun; j < totfileentries; j++ ){
614     //
615     // Get entry from old tree
616     //
617     OrbitalInfotrclone->GetEntry(j);
618     //
619     // copy orbitalinfoclone to OrbitalInfo
620     //
621 mocchiut 1.3 orbitalinfo->Clear();
622 mocchiut 1.5 //
623 mocchiut 1.1 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
624     //
625     // Fill entry in the new tree
626     //
627     OrbitalInfotr->Fill();
628     };
629     if (verbose) printf(" Finished successful copying!\n");
630     };
631     };
632     //
633     // Close files, delete old tree(s), write and close level2 file
634     //
635     if ( l0File ) l0File->Close();
636     if ( tempfile ) tempfile->Close();
637     gSystem->Unlink(tempname.str().c_str());
638 mocchiut 1.5 //
639 mocchiut 1.1 if ( runinfo ) runinfo->Close();
640     if ( OrbitalInfotr ) OrbitalInfotr->SetName("OrbitalInfo");
641     if ( file ){
642     file->cd();
643     file->Write();
644     };
645     //
646     gSystem->Unlink(OrbitalInfofolder.str().c_str());
647     //
648     // the end
649     //
650     if (verbose) printf("\n Exiting...\n");
651     if(OrbitalInfotr)OrbitalInfotr->Delete();
652 mocchiut 1.4 //
653     if ( orbitalinfo ) delete orbitalinfo;
654     if ( orbitalinfoclone ) delete orbitalinfoclone;
655     if ( glroot ) delete glroot;
656     if ( runinfo ) delete runinfo;
657     //
658 mocchiut 1.1 if(code < 0) throw code;
659     return(code);
660     }
661    
662 mocchiut 1.7
663     //
664     // Returns the cCoordGeo structure holding the geographical
665     // coordinates for the event (see sgp4.h).
666     //
667     // atime is the abstime of the event in UTC unix time.
668     // tletime is the time of the tle in UTC unix time.
669     // tle is the previous and nearest tle (compared to atime).
670     cCoordGeo getCoo(UInt_t atime, UInt_t tletime, cTle *tle)
671     {
672     cEci eci;
673     cOrbit orbit(*tle);
674    
675     orbit.getPosition((double) (atime - tletime)/60., &eci);
676    
677     return eci.toGeo();
678     }

  ViewVC Help
Powered by ViewVC 1.1.23