/[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.73 - (hide annotations) (download)
Fri Jun 6 14:18:20 2014 UTC (10 years, 10 months ago) by mocchiut
Branch: MAIN
Changes since 1.72: +584 -134 lines
New extended algorithm track-related added to CALO, TOF and ORB

1 mocchiut 1.1 // C/C++ headers
2     //
3     #include <fstream>
4     #include <string.h>
5     #include <iostream>
6     #include <cstring>
7     #include <stdio.h>
8     //
9     // ROOT headers
10     //
11 mocchiut 1.44 //#include <TCanvas.h>
12 pam-mep 1.60 #include <TH2F.h> //for test only. Vitaly.
13     #include <TVector3.h>
14 mocchiut 1.44 //#include <TF1.h>
15    
16 mocchiut 1.1 #include <TTree.h>
17     #include <TClassEdit.h>
18     #include <TObject.h>
19     #include <TList.h>
20 mocchiut 1.5 #include <TArrayI.h>
21 mocchiut 1.1 #include <TSystem.h>
22     #include <TSystemDirectory.h>
23     #include <TString.h>
24     #include <TFile.h>
25     #include <TClass.h>
26     #include <TSQLServer.h>
27     #include <TSQLRow.h>
28     #include <TSQLResult.h>
29     //
30 mocchiut 1.8 // RunInfo header
31     //
32     #include <RunInfo.h>
33     #include <GLTables.h>
34     //
35 mocchiut 1.1 // YODA headers
36     //
37     #include <PamelaRun.h>
38     #include <PscuHeader.h>
39     #include <PscuEvent.h>
40     #include <EventHeader.h>
41 mocchiut 1.15 #include <mcmd/McmdEvent.h>
42     #include <mcmd/McmdRecord.h>
43 mocchiut 1.1 //
44     // This program headers
45     //
46     #include <OrbitalInfo.h>
47 mocchiut 1.7 #include <OrbitalInfoVerl2.h>
48 mocchiut 1.1 #include <OrbitalInfoCore.h>
49 mocchiut 1.15 #include <InclinationInfo.h>
50 mocchiut 1.1
51 pam-mep 1.69 //
52     // Tracker and ToF classes headers and definitions
53     //
54     #include <ToFLevel2.h>
55     #include <TrkLevel2.h>
56 mocchiut 1.73 #include <ExtTrack.h> // new tracking code
57 mocchiut 1.40
58 mocchiut 1.1 using namespace std;
59    
60     //
61     // CORE ROUTINE
62     //
63     //
64 mocchiut 1.25 int OrbitalInfoCore(UInt_t run, TFile *file, GL_TABLES *glt, Int_t OrbitalInfoargc, char *OrbitalInfoargv[]){
65 mocchiut 1.15 //
66 mocchiut 1.1 Int_t i = 0;
67 mocchiut 1.25 TString host = glt->CGetHost();
68     TString user = glt->CGetUser();
69     TString psw = glt->CGetPsw();
70     TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
71 mocchiut 1.1 //
72 mocchiut 1.26 stringstream myquery;
73     myquery.str("");
74     myquery << "SET time_zone='+0:00'";
75 mocchiut 1.59 delete dbc->Query(myquery.str().c_str());
76 mocchiut 1.26 //
77 mocchiut 1.16 TString processFolder = Form("OrbitalInfoFolder_%u",run);
78 mocchiut 1.1 //
79     // Set these to true to have a very verbose output.
80     //
81     Bool_t debug = false;
82     //
83     Bool_t verbose = false;
84 mocchiut 1.30 //
85 mocchiut 1.31 Bool_t standalone = false;
86 mocchiut 1.30 //
87 mocchiut 1.1 if ( OrbitalInfoargc > 0 ){
88     i = 0;
89     while ( i < OrbitalInfoargc ){
90     if ( !strcmp(OrbitalInfoargv[i],"-processFolder") ) {
91     if ( OrbitalInfoargc < i+1 ){
92     throw -3;
93     };
94     processFolder = (TString)OrbitalInfoargv[i+1];
95     i++;
96     };
97     if ( (!strcmp(OrbitalInfoargv[i],"--debug")) || (!strcmp(OrbitalInfoargv[i],"-g")) ) {
98     verbose = true;
99 mocchiut 1.15 debug = true;
100 mocchiut 1.1 };
101     if ( (!strcmp(OrbitalInfoargv[i],"--verbose")) || (!strcmp(OrbitalInfoargv[i],"-v")) ) {
102     verbose = true;
103     };
104 mocchiut 1.30 if ( (!strcmp(OrbitalInfoargv[i],"--standalone")) ) {
105     standalone = true;
106     };
107     if ( (!strcmp(OrbitalInfoargv[i],"--calculate-pitch")) ) {
108     standalone = false;
109     };
110 mocchiut 1.1 i++;
111     };
112     };
113     //
114     const char* outDir = gSystem->DirName(gSystem->DirName(file->GetPath()));
115     //
116     TTree *OrbitalInfotr = 0;
117 mocchiut 1.5 UInt_t nevents = 0;
118 mocchiut 1.15 UInt_t neventsm = 0;
119 mocchiut 1.1 //
120     // variables needed to reprocess data
121     //
122 mocchiut 1.6 Long64_t maxsize = 10000000000LL;
123     TTree::SetMaxTreeSize(maxsize);
124     //
125 mocchiut 1.1 TString OrbitalInfoversion;
126     ItoRunInfo *runinfo = 0;
127 mocchiut 1.5 TArrayI *runlist = 0;
128 mocchiut 1.1 TTree *OrbitalInfotrclone = 0;
129     Bool_t reproc = false;
130     Bool_t reprocall = false;
131 mocchiut 1.66 Bool_t igrfloaded = false;
132 mocchiut 1.1 UInt_t nobefrun = 0;
133     UInt_t noaftrun = 0;
134     UInt_t numbofrun = 0;
135     stringstream ftmpname;
136     TString fname;
137 mocchiut 1.5 UInt_t totfileentries = 0;
138 mocchiut 1.15 UInt_t idRun = 0;
139 pam-fi 1.53 UInt_t anni5 = 60 * 60 * 24 * 365 * 5 ;//1576800
140 mocchiut 1.15 //
141     // My variables. Vitaly.
142     //
143 mocchiut 1.44 // UInt_t oi = 0;
144 mocchiut 1.15 Int_t tmpSize = 0;
145 mocchiut 1.1 //
146     // variables needed to handle error signals
147     //
148     Int_t code = 0;
149     Int_t sgnl;
150     //
151     // OrbitalInfo classes
152     //
153     OrbitalInfo *orbitalinfo = new OrbitalInfo();
154     OrbitalInfo *orbitalinfoclone = new OrbitalInfo();
155 mocchiut 1.44
156 mocchiut 1.1 //
157     // define variables for opening and reading level0 file
158     //
159     TFile *l0File = 0;
160     TTree *l0tr = 0;
161 mocchiut 1.37 // TTree *l0trm = 0;
162     TChain *ch = 0;
163 mocchiut 1.2 // EM: open also header branch
164     TBranch *l0head = 0;
165     pamela::EventHeader *eh = 0;
166     pamela::PscuHeader *ph = 0;
167 mocchiut 1.15 pamela::McmdEvent *mcmdev = 0;
168     pamela::McmdRecord *mcmdrc = 0;
169 mocchiut 1.2 // end EM
170 mocchiut 1.15
171     // pamela::RunHeaderEvent *reh = new pamela::RunHeaderEvent;
172     // pamela::EventHeader *eH = new pamela::EventHeader;
173    
174 mocchiut 1.1 //
175     // Define other basic variables
176     //
177     UInt_t procev = 0;
178     stringstream file2;
179     stringstream file3;
180     stringstream qy;
181     Int_t totevent = 0;
182 mocchiut 1.5 UInt_t atime = 0;
183     UInt_t re = 0;
184 mocchiut 1.15 UInt_t ik = 0;
185 mocchiut 1.7
186     // Position
187     Float_t lon, lat, alt;
188    
189     //
190     // IGRF stuff
191     //
192 mocchiut 1.64 Double_t dimo = 0.0; // dipole moment (computed from dat files) // EM GCC 4.7
193 mocchiut 1.40 Float_t bnorth, beast, bdown, babs;
194     Float_t xl; // L value
195     Float_t icode; // code value for L accuracy (see fortran code)
196     Float_t bab1; // What's the difference with babs?
197     Float_t stps = 0.005; // step size for field line tracing
198     Float_t bdel = 0.01; // required accuracy
199     Float_t bequ; // equatorial b value (also called b_0)
200     Bool_t value = 0; // false if bequ is not the minimum b value
201     Float_t rr0; // equatorial radius normalized to earth radius
202 mocchiut 1.7
203 mocchiut 1.1 //
204     // Working filename
205     //
206     TString outputfile;
207     stringstream name;
208     name.str("");
209     name << outDir << "/";
210     //
211     // temporary file and folder
212     //
213     TFile *tempfile = 0;
214     TTree *tempOrbitalInfo = 0;
215     stringstream tempname;
216     stringstream OrbitalInfofolder;
217 mocchiut 1.23 Bool_t myfold = false;
218 mocchiut 1.1 tempname.str("");
219     tempname << outDir;
220     tempname << "/" << processFolder.Data();
221     OrbitalInfofolder.str("");
222     OrbitalInfofolder << tempname.str().c_str();
223     tempname << "/OrbitalInfotree_run";
224     tempname << run << ".root";
225 mocchiut 1.40 UInt_t totnorun = 0;
226 mocchiut 1.1 //
227     // DB classes
228     //
229     GL_ROOT *glroot = new GL_ROOT();
230 mocchiut 1.5 GL_TIMESYNC *dbtime = 0;
231 mocchiut 1.7 GL_TLE *gltle = new GL_TLE();
232 mocchiut 1.15 //
233     //Quaternions classes
234     //
235 mocchiut 1.57 Quaternions *L_QQ_Q_l_lower = 0;
236     InclinationInfo *RYPang_lower = 0;
237     Quaternions *L_QQ_Q_l_upper = 0;
238     InclinationInfo *RYPang_upper = 0;
239 mocchiut 1.15
240     cEci eCi;
241    
242 mocchiut 1.7 // Initialize fortran routines!!!
243 mocchiut 1.52 Int_t ltp1 = 0;
244 mocchiut 1.7 Int_t ltp2 = 0;
245     Int_t ltp3 = 0;
246 mocchiut 1.55 // Int_t uno = 1;
247     // const char *niente = " ";
248 pam-mep 1.69 GL_PARAM *glparam0 = new GL_PARAM();
249 mocchiut 1.7 GL_PARAM *glparam = new GL_PARAM();
250 mocchiut 1.10 GL_PARAM *glparam2 = new GL_PARAM();
251 mocchiut 1.52 GL_PARAM *glparam3 = new GL_PARAM();
252 mocchiut 1.44
253 mocchiut 1.30 //
254 mocchiut 1.44 // Orientation variables. Vitaly
255 mocchiut 1.30 //
256     UInt_t evfrom = 0;
257     UInt_t jumped = 0;
258     Int_t itr = -1;
259 mocchiut 1.64 // Double_t A1;
260     // Double_t A2;
261     // Double_t A3;
262 mocchiut 1.30 Double_t Px = 0;
263     Double_t Py = 0;
264     Double_t Pz = 0;
265     TTree *ttof = 0;
266     ToFLevel2 *tof = new ToFLevel2();
267 pam-mep 1.69 TTree *ttrke = 0;
268     TrkLevel2 *trke = new TrkLevel2();
269 mocchiut 1.30 OrientationInfo *PO = new OrientationInfo();
270     Int_t nz = 6;
271     Float_t zin[6];
272     Int_t nevtofl2 = 0;
273 pam-mep 1.69 Int_t nevtrkl2 = 0;
274 mocchiut 1.44 if ( verbose ) cout<<"Reading quaternions external file"<<endl;
275     cout.setf(ios::fixed,ios::floatfield);
276     /******Reading recovered quaternions...*********/
277     vector<Double_t> recqtime;
278     vector<Float_t> recq0;
279     vector<Float_t> recq1;
280     vector<Float_t> recq2;
281     vector<Float_t> recq3;
282     Float_t Norm = 1;
283 mocchiut 1.73
284     vector<UInt_t> RTtime1;
285     vector<UInt_t> RTtime2;
286     vector<Double_t> RTbank1;
287     vector<Double_t> RTbank2;
288     vector<Int_t> RTazim;
289     vector<Int_t> RTdir1;
290     vector<Int_t> RTdir2;
291     vector<Int_t> RTerrq;
292    
293     TClonesArray *tcNucleiTrk = NULL;
294     TClonesArray *tcExtNucleiTrk = NULL;
295     TClonesArray *tcExtTrk = NULL;
296     TClonesArray *tcNucleiTof = NULL;
297     TClonesArray *tcExtNucleiTof = NULL;
298     TClonesArray *tcExtTof = NULL;
299     TClonesArray *torbNucleiTrk = NULL;
300     TClonesArray *torbExtNucleiTrk = NULL;
301     TClonesArray *torbExtTrk = NULL;
302     Bool_t hasNucleiTrk = false;
303     Bool_t hasExtNucleiTrk = false;
304     Bool_t hasExtTrk = false;
305     Bool_t hasNucleiTof = false;
306     Bool_t hasExtNucleiTof = false;
307     Bool_t hasExtTof = false;
308    
309     ifstream in;
310     ifstream an;
311     // ofstream mc;
312     // TString gr;
313     Int_t parerror2=0;
314    
315 pam-mep 1.69 Int_t parerror=glparam0->Query_GL_PARAM(1,303,dbc); // parameters stored in DB in GL_PRAM table
316 mocchiut 1.72 if ( verbose ) cout<<parerror<<"\t"<<(char*)(glparam0->PATH+glparam0->NAME).Data()<<endl;
317 mocchiut 1.44 if ( parerror<0 ) {
318     code = parerror;
319 mocchiut 1.73 goto closeandexit;
320 pam-mep 1.54 }
321 mocchiut 1.73 in.open((char*)(glparam0->PATH+glparam0->NAME).Data(),ios::in);
322 mocchiut 1.44 while(!in.eof()){
323     recqtime.resize(recqtime.size()+1);
324     Int_t sizee = recqtime.size();
325     recq0.resize(sizee);
326     recq1.resize(sizee);
327     recq2.resize(sizee);
328     recq3.resize(sizee);
329     in>>recqtime[sizee-1];
330     in>>recq0[sizee-1];
331     in>>recq1[sizee-1];
332     in>>recq2[sizee-1];
333     in>>recq3[sizee-1];
334     in>>Norm;
335     }
336 pam-mep 1.69 in.close();
337 mocchiut 1.44 if ( verbose ) cout<<"We have read recovered data"<<endl;
338 pam-mep 1.69 if (debug) cout << "size of recovered quaterions data set is " << recqtime.size() << endl;
339    
340 mocchiut 1.73 // 10RED CHECK
341    
342     // TH2F* DIFFX = new TH2F("diffx","",100,0,100,90,0,90);
343     // TH2F* DIFFY = new TH2F("diffy","",100,0,100,90,0,90);
344     // TH2F* DIFFZ = new TH2F("diffz","",100,0,100,90,0,90);
345    
346     // gr = "methodscomparison_";
347     // gr+=run;
348     // gr+=".txt";
349     // mc.open(gr);
350     // mc.setf(ios::fixed,ios::floatfield);
351     // 10RED CHECK END
352 pam-mep 1.69
353     if ( verbose ) cout<<"read Rotation Table"<<endl;
354 mocchiut 1.73
355     parerror2=glparam0->Query_GL_PARAM(1,305,dbc);
356 mocchiut 1.71 if ( verbose ) cout<<parerror2<<"\t"<<(char*)(glparam0->PATH+glparam0->NAME).Data()<<endl;
357 mocchiut 1.73 if ( parerror2<0 ) {
358     code = parerror;
359     goto closeandexit;
360     }
361     an.open((char*)(glparam0->PATH+glparam0->NAME).Data(),ios::in);
362 pam-mep 1.69
363     while(!an.eof()){
364     RTtime1.resize(RTtime1.size()+1);
365     Int_t sizee = RTtime1.size();
366     RTbank1.resize(sizee+1);
367     RTazim.resize(sizee+1);
368     RTdir1.resize(sizee+1);
369     RTerrq.resize(sizee+1);
370     an>>RTtime1[sizee-1];
371     an>>RTbank1[sizee-1];
372     an>>RTazim[sizee-1];
373     an>>RTdir1[sizee-1];
374     an>>RTerrq[sizee-1];
375     if(sizee>1) {
376     RTtime2.resize(sizee+1);
377     RTbank2.resize(sizee+1);
378     RTdir2.resize(sizee+1);
379     RTtime2[sizee-2]=RTtime1[sizee-1];
380     RTbank2[sizee-2]=RTbank1[sizee-1];
381     RTdir2[sizee-2]=RTdir1[sizee-1];
382     }
383     }
384     an.close();
385     //cout<<"put some number here"<<endl;
386     //Int_t yupi;
387     //cin>>yupi;
388    
389     if ( verbose ) cout<<"We have read Rotation Table"<<endl;
390     //Geomagnetic coordinates calculations staff
391    
392     GMtype_CoordGeodetic location;
393 mocchiut 1.70 // GMtype_CoordDipole GMlocation;
394 pam-mep 1.69 GMtype_Ellipsoid Ellip;
395     GMtype_Data G0, G1, H1;
396    
397     // { // this braces is necessary to avoid jump to label 'closeandexit' error // but it is wrong since the variable "igpath" will not exist outside. To overcome the "jump to label 'closeandexit' error" it is necessary to set the "igpath" before line 276
398     // TString igpath="/data03/Malakhov/pam9Malakhov/installed10/calib/orb-param/";
399     // }
400    
401     // GM_ScanIGRF(glparam->PATH, &G0, &G1, &H1);
402     GM_ScanIGRF(dbc, &G0, &G1, &H1);
403    
404     //cout << G0.element[0] << "\t" << G1.element[0] << "\t" << H1.element[0] << endl;
405     //cout << G0.element[5] << "\t" << G1.element[5] << "\t" << H1.element[5] << endl;
406    
407     GM_SetEllipsoid(&Ellip);
408 mocchiut 1.44
409 mocchiut 1.55 // IGRF stuff moved inside run loop!
410    
411 mocchiut 1.30 for (Int_t ip=0;ip<nz;ip++){
412     zin[ip] = tof->GetZTOF(tof->GetToFPlaneID(ip));
413     };
414     //
415     if ( !standalone ){
416     //
417 mocchiut 1.73 // Does it contain the Tracker and ToF trees?
418 mocchiut 1.30 //
419     ttof = (TTree*)file->Get("ToF");
420     if ( !ttof ) {
421     if ( verbose ) printf(" OrbitalInfo - ERROR: no tof tree\n");
422     code = -900;
423     goto closeandexit;
424 pam-mep 1.69 }
425 mocchiut 1.30 ttof->SetBranchAddress("ToFLevel2",&tof);
426     nevtofl2 = ttof->GetEntries();
427 pam-mep 1.69
428 mocchiut 1.73 //
429     // Look for extended tracking algorithm
430     //
431     if ( verbose ) printf("Look for extended and nuclei tracking algorithms in ToF\n");
432     // Nuclei tracking algorithm
433     Int_t checkAlgo = 0;
434     tcNucleiTof = new TClonesArray("ToFTrkVar");
435     checkAlgo = ttof->SetBranchAddress("TrackNuclei",&tcNucleiTof);
436     if ( !checkAlgo ){
437     if ( verbose ) printf(" Nuclei tracking algorithm ToF branch found! :D \n");
438     hasNucleiTof = true;
439     } else {
440     if ( verbose ) printf(" Nuclei tracking algorithm ToF branch not found :( !\n");
441     printf(" ok, this is not a problem (it depends on tracker settings) \n");
442     delete tcNucleiTof;
443     }
444     // Nuclei tracking algorithm using calorimeter points
445     tcExtNucleiTof = new TClonesArray("ToFTrkVar");
446     checkAlgo = ttof->SetBranchAddress("RecoveredTrackNuclei",&tcExtNucleiTof);
447     if ( !checkAlgo ){
448     if ( verbose ) printf(" Recovered nuclei tracking algorithm ToF branch found! :D \n");
449     hasExtNucleiTof = true;
450     } else {
451     if ( verbose ) printf(" Recovered nuclei tracking algorithm ToF branch not found :( !\n");
452     printf(" ok, this is not a problem (it depends on tracker settings) \n");
453     delete tcExtNucleiTof;
454     }
455     // Tracking algorithm using calorimeter points
456     tcExtTof = new TClonesArray("ToFTrkVar");
457     checkAlgo = ttof->SetBranchAddress("RecoveredTrack",&tcExtTof);
458     if ( !checkAlgo ){
459     if ( verbose ) printf(" Recovered track algorithm ToF branch found! :D \n");
460     hasExtTof = true;
461     } else {
462     if ( verbose ) printf(" Recovered track algorithm ToF branch not found :( !\n");
463     printf(" ok, this is not a problem (it depends on tracker settings) \n");
464     delete tcExtTof;
465     }
466    
467 pam-mep 1.69 ttrke = (TTree*)file->Get("Tracker");
468     if ( !ttrke ) {
469     if ( verbose ) printf(" OrbitalInfo - ERROR: no trk tree\n");
470     code = -903;
471     goto closeandexit;
472     }
473     ttrke->SetBranchAddress("TrkLevel2",&trke);
474     nevtrkl2 = ttrke->GetEntries();
475 mocchiut 1.73
476     //
477     // Look for extended tracking algorithm
478     //
479     if ( verbose ) printf("Look for extended and nuclei tracking algorithms\n");
480     // Nuclei tracking algorithm
481     checkAlgo = 0;
482     tcNucleiTrk = new TClonesArray("TrkTrack");
483     checkAlgo = ttrke->SetBranchAddress("TrackNuclei",&tcNucleiTrk);
484     if ( !checkAlgo ){
485     if ( verbose ) printf(" Nuclei tracking algorithm branch found! :D \n");
486     hasNucleiTrk = true;
487     } else {
488     if ( verbose ) printf(" Nuclei tracking algorithm branch not found :( !\n");
489     printf(" ok, this is not a problem (it depends on tracker settings) \n");
490     delete tcNucleiTrk;
491     }
492     // Nuclei tracking algorithm using calorimeter points
493     tcExtNucleiTrk = new TClonesArray("ExtTrack");
494     checkAlgo = ttrke->SetBranchAddress("RecoveredTrackNuclei",&tcExtNucleiTrk);
495     if ( !checkAlgo ){
496     if ( verbose ) printf(" Recovered nuclei tracking algorithm branch found! :D \n");
497     hasExtNucleiTrk = true;
498     } else {
499     if ( verbose ) printf(" Recovered nuclei tracking algorithm branch not found :( !\n");
500     printf(" ok, this is not a problem (it depends on tracker settings) \n");
501     delete tcExtNucleiTrk;
502     }
503     // Tracking algorithm using calorimeter points
504     tcExtTrk = new TClonesArray("ExtTrack");
505     checkAlgo = ttrke->SetBranchAddress("RecoveredTrack",&tcExtTrk);
506     if ( !checkAlgo ){
507     if ( verbose ) printf(" Recovered track algorithm branch found! :D \n");
508     hasExtTrk = true;
509     } else {
510     if ( verbose ) printf(" Recovered track algorithm branch not found :( !\n");
511     printf(" ok, this is not a problem (it depends on tracker settings) \n");
512     delete tcExtTrk;
513     }
514    
515     if ( (hasNucleiTrk && !hasNucleiTof) || (!hasNucleiTrk && hasNucleiTof) ||
516     (hasExtNucleiTrk && !hasExtNucleiTof) || (!hasExtNucleiTrk && hasExtNucleiTof) ||
517     (hasExtTrk && !hasExtTof) || (!hasExtTrk && hasExtTof)
518     ){
519     if ( verbose ) printf(" ERROR: Mismatch between tracker and tof tree branches concerning extended tracking algorithm(s)\n");
520     if ( debug ) printf("hasNucleiTrk %i hasExtNucleiTrk %i hasExtTrk %i \n",hasNucleiTrk,hasExtNucleiTrk,hasExtTrk);
521     if ( debug ) printf("hasNucleiTof %i hasExtNucleiTof %i hasExtTof %i \n",hasNucleiTof,hasExtNucleiTof,hasExtTof);
522     throw -901;
523     }
524    
525 pam-mep 1.69 }
526 mocchiut 1.1 //
527     // Let's start!
528     //
529     // 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
530     // if run != 0 we must process only that run but first we have to check if the tree MyDetector2 already exist in the file
531     // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
532     //
533 mocchiut 1.5 if ( run == 0 ) reproc = true;
534 mocchiut 1.1 //
535     //
536     // Output file is "outputfile"
537     //
538     if ( !file->IsOpen() ){
539     //printf(" OrbitalInfo - ERROR: cannot open file for writing\n");
540     throw -901;
541     };
542     //
543     // Retrieve GL_RUN variables from the level2 file
544     //
545     OrbitalInfoversion = OrbitalInfoInfo(false); // we should decide how to handle versioning system
546     //
547     // create an interface to RunInfo called "runinfo"
548     //
549     runinfo = new ItoRunInfo(file);
550     //
551     // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
552     //
553     sgnl = 0;
554     sgnl = runinfo->Update(run, "ORB", OrbitalInfoversion);
555     //sgnl = runinfo->Read(run);
556    
557     if ( sgnl ){
558     //printf("OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
559     code = sgnl;
560     goto closeandexit;
561     } else {
562     sgnl = 0;
563     };
564     //
565     // number of events in the file BEFORE the first event of our run
566     //
567     nobefrun = runinfo->GetFirstEntry();
568     //
569     // total number of events in the file
570     //
571     totfileentries = runinfo->GetFileEntries();
572     //
573     // first file entry AFTER the last event of our run
574     //
575     noaftrun = runinfo->GetLastEntry() + 1;
576     //
577     // number of run to be processed
578     //
579     numbofrun = runinfo->GetNoRun();
580 mocchiut 1.40 totnorun = runinfo->GetRunEntries();
581 mocchiut 1.1 //
582     // Try to access the OrbitalInfo tree in the file, if it exists we are reprocessing data if not we are processing a new run
583     //
584     OrbitalInfotrclone = (TTree*)file->Get("OrbitalInfo");
585     //
586     if ( !OrbitalInfotrclone ){
587     //
588     // tree does not exist, we are not reprocessing
589     //
590     reproc = false;
591 mocchiut 1.5 if ( run == 0 ){
592 mocchiut 1.1 if (verbose) printf(" OrbitalInfo - WARNING: you are reprocessing data but OrbitalInfo tree does not exist!\n");
593     }
594 mocchiut 1.5 if ( runinfo->IsReprocessing() && run != 0 ) {
595 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");
596     }
597     } else {
598     //
599     // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
600     //
601 mocchiut 1.6 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
602 mocchiut 1.1 reproc = true;
603     //
604     //
605     if (verbose) printf("\n Preparing the pre-processing...\n");
606     //
607 mocchiut 1.22 if ( run == 0 || totnorun == 1 ){
608 mocchiut 1.1 //
609     // we are reprocessing all the file
610     // 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
611     //
612     reprocall = true;
613     //
614 mocchiut 1.57 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing all runs\n Deleting old tree...\n");
615 mocchiut 1.1 //
616     } else {
617     //
618     // 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
619     //
620     reprocall = false;
621     //
622 mocchiut 1.5 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing run number %u \n",run);
623 mocchiut 1.1 //
624     // copying old tree to a new file
625     //
626 mocchiut 1.23 gSystem->MakeDirectory(OrbitalInfofolder.str().c_str());
627     myfold = true;
628 mocchiut 1.1 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
629     tempOrbitalInfo = OrbitalInfotrclone->CloneTree(-1,"fast");
630     tempOrbitalInfo->SetName("OrbitalInfo-old");
631     tempfile->Write();
632 mocchiut 1.57 tempOrbitalInfo->Delete();
633 mocchiut 1.1 tempfile->Close();
634     }
635     //
636     // Delete the old tree from old file and memory
637     //
638 mocchiut 1.57 OrbitalInfotrclone->Clear();
639 mocchiut 1.1 OrbitalInfotrclone->Delete("all");
640     //
641     if (verbose) printf(" ...done!\n");
642     //
643     };
644     //
645     // create mydetector tree mydect
646     //
647     file->cd();
648     OrbitalInfotr = new TTree("OrbitalInfo-new","PAMELA OrbitalInfo data");
649 mocchiut 1.6 OrbitalInfotr->SetAutoSave(900000000000000LL);
650 mocchiut 1.30 orbitalinfo->Set();//ELENA **TEMPORANEO?**
651 mocchiut 1.1 OrbitalInfotr->Branch("OrbitalInfo","OrbitalInfo",&orbitalinfo);
652     //
653 mocchiut 1.73 // create new branches for new tracking algorithms
654     //
655     if ( hasNucleiTrk ){
656     torbNucleiTrk = new TClonesArray("OrbitalInfoTrkVar",1);
657     OrbitalInfotr->Branch("TrackNuclei",&torbNucleiTrk);
658     }
659     if ( hasExtNucleiTrk ){
660     torbExtNucleiTrk = new TClonesArray("OrbitalInfoTrkVar",1);
661     OrbitalInfotr->Branch("RecoveredTrackNuclei",&torbExtNucleiTrk);
662     }
663     if ( hasExtTrk ){
664     torbExtTrk = new TClonesArray("OrbitalInfoTrkVar",1);
665     OrbitalInfotr->Branch("RecoveredTrack",&torbExtTrk);
666     }
667    
668     //
669 mocchiut 1.1 if ( reproc && !reprocall ){
670     //
671     // open new file and retrieve also tree informations
672     //
673     tempfile = new TFile(tempname.str().c_str(),"READ");
674     OrbitalInfotrclone = (TTree*)tempfile->Get("OrbitalInfo-old");
675 mocchiut 1.6 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
676 mocchiut 1.1 OrbitalInfotrclone->SetBranchAddress("OrbitalInfo",&orbitalinfoclone);
677     //
678     if ( nobefrun > 0 ){
679     if (verbose){
680 mocchiut 1.7 printf("\n Pre-processing: copying events from the old tree before the processed run\n");
681     printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
682     printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
683 mocchiut 1.1 }
684     for (UInt_t j = 0; j < nobefrun; j++){
685     //
686 mocchiut 1.43 if ( OrbitalInfotrclone->GetEntry(j) <= 0 ) throw -36;
687 mocchiut 1.1 //
688     // copy orbitalinfoclone to mydec
689     //
690 mocchiut 1.3 orbitalinfo->Clear();
691 mocchiut 1.5 //
692 mocchiut 1.1 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
693     //
694     // Fill entry in the new tree
695     //
696     OrbitalInfotr->Fill();
697     //
698     };
699     if (verbose) printf(" Finished successful copying!\n");
700 mocchiut 1.57 };
701 mocchiut 1.1 };
702     //
703 mocchiut 1.44 //
704 mocchiut 1.1 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
705     //
706     runlist = runinfo->GetRunList();
707     //
708     // Loop over the run to be processed
709     //
710     for (UInt_t irun=0; irun < numbofrun; irun++){
711 mocchiut 1.57
712     L_QQ_Q_l_lower = new Quaternions();
713     RYPang_lower = new InclinationInfo();
714     L_QQ_Q_l_upper = new Quaternions();
715     RYPang_upper = new InclinationInfo();
716    
717 mocchiut 1.1 //
718     // retrieve the first run ID to be processed using the RunInfo list
719     //
720 mocchiut 1.15
721 mocchiut 1.1 idRun = runlist->At(irun);
722     if (verbose){
723     printf("\n\n\n ####################################################################### \n");
724     printf(" PROCESSING RUN NUMBER %i \n",(int)idRun);
725     printf(" ####################################################################### \n\n\n");
726     }
727     //
728 mocchiut 1.5 runinfo->ID_ROOT_L0 = 0;
729 mocchiut 1.1 //
730     // store in the runinfo class the GL_RUN variables for our run
731     //
732     sgnl = 0;
733     sgnl = runinfo->GetRunInfo(idRun);
734     if ( sgnl ){
735 mocchiut 1.5 if ( debug ) printf("\n OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
736 mocchiut 1.1 code = sgnl;
737     goto closeandexit;
738     } else {
739     sgnl = 0;
740     };
741     //
742     // now you can access that variables using the RunInfo class this way runinfo->ID_REG_RUN
743     //
744 mocchiut 1.5 if ( runinfo->ID_ROOT_L0 == 0 ){
745     if ( debug ) printf("\n OrbitalInfo - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
746 mocchiut 1.1 code = -5;
747     goto closeandexit;
748     };
749     //
750 mocchiut 1.5 // prepare the timesync for the db
751     //
752     dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
753 mocchiut 1.15
754 mocchiut 1.5 //
755 mocchiut 1.1 // Search in the DB the path and name of the LEVEL0 file to be processed.
756     //
757 mocchiut 1.5 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
758 mocchiut 1.1 //
759     ftmpname.str("");
760     ftmpname << glroot->PATH.Data() << "/";
761     ftmpname << glroot->NAME.Data();
762     fname = ftmpname.str().c_str();
763 mocchiut 1.15 ftmpname.str("");
764 mocchiut 1.1 //
765 mocchiut 1.40 // print nout informations
766 mocchiut 1.1 //
767 mocchiut 1.5 totevent = runinfo->NEVENTS;
768 mocchiut 1.30 evfrom = runinfo->EV_FROM;
769 mocchiut 1.15 //cout<<"totevents = "<<totevent<<"\n";
770 mocchiut 1.1 if (verbose){
771     printf("\n LEVEL0 data file: %s \n",fname.Data());
772 mocchiut 1.5 printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
773     printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
774 mocchiut 1.15 printf(" %i events to be processed for run %u: from %i to %i \n\n",totevent,idRun,runinfo->EV_FROM+1,runinfo->EV_FROM+totevent);
775 mocchiut 1.1 }//
776 mocchiut 1.27 //
777 mocchiut 1.28 // if ( !totevent ) goto closeandexit;
778 mocchiut 1.1 // Open Level0 file
779 mocchiut 1.59 if ( l0File ) l0File->Close();
780 mocchiut 1.1 l0File = new TFile(fname.Data());
781     if ( !l0File ) {
782 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: problems opening Level0 file\n");
783 mocchiut 1.1 code = -6;
784     goto closeandexit;
785     };
786     l0tr = (TTree*)l0File->Get("Physics");
787     if ( !l0tr ) {
788 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: no Physics tree in Level0 file\n");
789 mocchiut 1.1 l0File->Close();
790     code = -7;
791     goto closeandexit;
792     };
793 mocchiut 1.2 // EM: open header branch as well
794     l0head = l0tr->GetBranch("Header");
795     if ( !l0head ) {
796 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: no Header branch in Level0 tree\n");
797 mocchiut 1.2 l0File->Close();
798     code = -8;
799     goto closeandexit;
800     };
801     l0tr->SetBranchAddress("Header", &eh);
802     // end EM
803 mocchiut 1.5 nevents = l0head->GetEntries();
804 mocchiut 1.1 //
805 mocchiut 1.28 if ( nevents < 1 && totevent ) {
806 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: Level0 file is empty\n\n");
807 mocchiut 1.1 l0File->Close();
808     code = -11;
809     goto closeandexit;
810     };
811 mocchiut 1.20 //
812 mocchiut 1.28 if ( runinfo->EV_TO > nevents-1 && totevent ) {
813 mocchiut 1.5 if ( debug ) printf(" OrbitalInfo - ERROR: too few entries in the registry tree\n");
814 mocchiut 1.1 l0File->Close();
815     code = -12;
816     goto closeandexit;
817     };
818 mocchiut 1.55
819     //
820     // open IGRF files and do it only once if we are processing a full level2 file
821     //
822 mocchiut 1.68 if ( !igrfloaded ){
823 mocchiut 1.66
824     if ( l0head->GetEntry(runinfo->EV_FROM) > 0 ){
825     igrfloaded = true;
826 mocchiut 1.68 //
827     // absolute time of first event of the run (it should not matter a lot)
828     //
829     ph = eh->GetPscuHeader();
830     atime = dbtime->DBabsTime(ph->GetOrbitalTime());
831    
832     parerror=glparam->Query_GL_PARAM(atime-anni5,301,dbc); // parameters stored in DB in GL_PRAM table
833     if ( parerror<0 ) {
834     code = parerror;
835     goto closeandexit;
836     }
837     ltp1 = (Int_t)(glparam->PATH+glparam->NAME).Length();
838     if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
839     //
840     parerror=glparam2->Query_GL_PARAM(atime,301,dbc); // parameters stored in DB in GL_PRAM table
841     if ( parerror<0 ) {
842     code = parerror;
843     goto closeandexit;
844     }
845     ltp2 = (Int_t)(glparam2->PATH+glparam2->NAME).Length();
846     if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam2->PATH+glparam2->NAME).Data());
847     //
848     parerror=glparam3->Query_GL_PARAM(atime,302,dbc); // parameters stored in DB in GL_PRAM table
849     if ( parerror<0 ) {
850     code = parerror;
851     goto closeandexit;
852     }
853     ltp3 = (Int_t)(glparam3->PATH+glparam3->NAME).Length();
854     if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam3->PATH+glparam3->NAME).Data());
855     //
856     initize_((char *)(glparam->PATH+glparam->NAME).Data(),&ltp1,(char *)(glparam2->PATH+glparam2->NAME).Data(),&ltp2,(char *)(glparam3->PATH+glparam3->NAME).Data(),&ltp3);
857     //
858 mocchiut 1.70 if (debug) cout<<"initize: "<<(char *)(glparam->PATH+glparam->NAME).Data()<<"\t"<<(char *)(glparam2->PATH+glparam2->NAME).Data()<<"\t"<<(char *)(glparam3->PATH+glparam3->NAME).Data()<<endl;
859 mocchiut 1.66 }
860 mocchiut 1.55 }
861     //
862     // End IGRF stuff//
863     //
864    
865 mocchiut 1.15 ULong_t TimeSync = (ULong_t)dbtime->GetTimesync();
866     ULong_t ObtSync = (ULong_t)(dbtime->GetObt0()/1000);
867     ULong_t DeltaOBT = TimeSync - ObtSync;
868    
869     if ( debug ) printf(" 2 TimeSync %lu ObtSync %lu DeltaOBT %lu\n",(ULong_t)(dbtime->GetTimesync()/1000),(ULong_t)dbtime->GetObt0(),TimeSync-ObtSync);
870 mocchiut 1.37 //
871     // Read MCMDs from up to 11 files, 5 before and 5 after the present one in order to have some kind of inclination information
872     //
873     ch = new TChain("Mcmd","Mcmd");
874     //
875     // look in the DB to find the closest files to this run
876     //
877     TSQLResult *pResult = 0;
878     TSQLRow *Row = 0;
879     stringstream myquery;
880     UInt_t l0fid[10];
881     Int_t i = 0;
882     memset(l0fid,0,10*sizeof(Int_t));
883     //
884     myquery.str("");
885     myquery << "select ID_ROOT_L0 from GL_RUN where RUNHEADER_TIME<=" << runinfo->RUNHEADER_TIME << " group by ID_ROOT_L0 order by RUNHEADER_TIME desc limit 5;";
886     //
887     pResult = dbc->Query(myquery.str().c_str());
888     //
889     i = 9;
890     if( pResult ){
891     //
892     Row = pResult->Next();
893     //
894     while ( Row ){
895     //
896     // store infos and exit
897     //
898     l0fid[i] = (UInt_t)atoll(Row->GetField(0));
899     i--;
900     Row = pResult->Next();
901     //
902     };
903     pResult->Delete();
904     };
905     //
906     myquery.str("");
907     myquery << "select ID_ROOT_L0 from GL_RUN where RUNHEADER_TIME>" << runinfo->RUNHEADER_TIME << " group by ID_ROOT_L0 order by RUNHEADER_TIME asc limit 5;";
908     //
909     pResult = dbc->Query(myquery.str().c_str());
910     //
911     i = 0;
912     if( pResult ){
913     //
914     Row = pResult->Next();
915     //
916     while ( Row ){
917     //
918     // store infos and exit
919     //
920     l0fid[i] = (UInt_t)atoll(Row->GetField(0));
921     i++;
922     Row = pResult->Next();
923     //
924     };
925     pResult->Delete();
926     };
927     //
928     i = 0;
929     UInt_t previd = 0;
930     while ( i < 10 ){
931     if ( l0fid[i] && previd != l0fid[i] ){
932     previd = l0fid[i];
933     myquery.str("");
934     myquery << "select PATH,NAME from GL_ROOT where ID=" << l0fid[i] << " ;";
935     //
936     pResult = dbc->Query(myquery.str().c_str());
937     //
938     if( pResult ){
939     //
940     Row = pResult->Next();
941     //
942     if ( debug ) printf(" Using inclination informations from file: %s \n",(((TString)gSystem->ExpandPathName(Row->GetField(0)))+"/"+(TString)Row->GetField(1)).Data());
943     ch->Add(((TString)gSystem->ExpandPathName(Row->GetField(0)))+"/"+(TString)Row->GetField(1));
944     //
945     pResult->Delete();
946     };
947     };
948     i++;
949     };
950     //
951     ch->SetBranchAddress("Mcmd",&mcmdev);
952     neventsm = ch->GetEntries();
953     if ( debug ) printf(" entries %u \n", neventsm);
954 mocchiut 1.15 //
955     if (neventsm == 0){
956 mocchiut 1.18 if ( debug ) printf("InclinationInfo - WARNING: No quaternions in this File");
957     code = 900;
958 mocchiut 1.15 }
959     //
960 mocchiut 1.73 Double_t lowerqtime = 0;
961 mocchiut 1.15 //
962 mocchiut 1.44 // init quaternions information from mcmd-packets
963 mocchiut 1.15 //
964     Bool_t isf = true;
965 mocchiut 1.44
966     vector<Float_t> q0;
967     vector<Float_t> q1;
968     vector<Float_t> q2;
969     vector<Float_t> q3;
970     vector<Double_t> qtime;
971     vector<Float_t> qPitch;
972     vector<Float_t> qRoll;
973     vector<Float_t> qYaw;
974     vector<Int_t> qmode;
975    
976     Int_t nt = 0;
977     UInt_t must = 0;
978    
979 mocchiut 1.15 //
980 mocchiut 1.1 // run over all the events of the run
981     //
982     if (verbose) printf("\n Ready to start! \n\n Processed events: \n\n");
983     //
984 mocchiut 1.44 //
985 mocchiut 1.5 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
986 mocchiut 1.1 //
987 mocchiut 1.5 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
988 mocchiut 1.15 if ( debug ) printf(" %i \n",procev);
989 mocchiut 1.1 //
990 mocchiut 1.43 if ( l0head->GetEntry(re) <= 0 ) throw -36;
991 mocchiut 1.1 //
992     // absolute time of this event
993     //
994 mocchiut 1.5 ph = eh->GetPscuHeader();
995 mocchiut 1.15 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
996 mocchiut 1.40 if ( debug ) printf(" %i absolute time \n",procev);
997 mocchiut 1.1 //
998     // paranoid check
999     //
1000 mocchiut 1.34 if ( (atime > (runinfo->RUNTRAILER_TIME+1)) || (atime < (runinfo->RUNHEADER_TIME-1)) ) {
1001 mocchiut 1.1 if (verbose) printf(" OrbitalInfo - WARNING: event at time outside the run time window, skipping it\n");
1002 mocchiut 1.30 jumped++;
1003 mocchiut 1.44 // debug = true;
1004 mocchiut 1.7 continue;
1005     }
1006 mocchiut 1.30
1007     //
1008     // retrieve tof informations
1009     //
1010     if ( !reprocall ){
1011     itr = nobefrun + (re - evfrom - jumped);
1012     //itr = re-(46438+200241);
1013     } else {
1014     itr = runinfo->GetFirstEntry() + (re - evfrom - jumped);
1015     };
1016     //
1017     if ( !standalone ){
1018     if ( itr > nevtofl2 ){
1019     if ( verbose ) printf(" OrbitalInfo - ERROR: no tof events with entry = %i in Level2 file\n",itr);
1020     if ( debug ) printf(" nobefrun %u re %u evfrom %u jumped %u reprocall %i \n",nobefrun,re,evfrom,jumped,reprocall);
1021     l0File->Close();
1022 pam-mep 1.69 code = -904;
1023 mocchiut 1.30 goto closeandexit;
1024     };
1025     //
1026     tof->Clear();
1027     //
1028 mocchiut 1.73 // Clones array must be cleared before going on
1029     //
1030     if ( hasNucleiTof ){
1031     tcNucleiTof->Delete();
1032     }
1033     if ( hasExtNucleiTof ){
1034     tcExtNucleiTof->Delete();
1035     }
1036     if ( hasExtTof ){
1037     tcExtTof->Delete();
1038     }
1039     //
1040     if ( verbose ) printf(" get tof tree entries... entry = %i in Level2 file\n",itr);
1041 emocchiutti 1.67 if ( ttof->GetEntry(itr) <= 0 ){
1042 mocchiut 1.73 if ( verbose ) printf(" problems with tof tree entries... entry = %i in Level2 file\n",itr);
1043     if ( verbose ) printf(" nobefrun %u re %u evfrom %u jumped %u reprocall %i \n",nobefrun,re,evfrom,jumped,reprocall);
1044     throw -36;
1045 emocchiutti 1.67 }
1046 mocchiut 1.73 if ( verbose ) printf(" gat0\n");
1047 mocchiut 1.30 //
1048 pam-mep 1.69 }
1049     //
1050     // retrieve tracker informations
1051     //
1052     if ( !standalone ){
1053     if ( itr > nevtrkl2 ){
1054     if ( verbose ) printf(" OrbitalInfo - ERROR: no trk events with entry = %i in Level2 file\n",itr);
1055     if ( debug ) printf(" nobefrun %u re %u evfrom %u jumped %u reprocall %i \n",nobefrun,re,evfrom,jumped,reprocall);
1056     l0File->Close();
1057     code = -905;
1058     goto closeandexit;
1059 mocchiut 1.73 }
1060 pam-mep 1.69 //
1061 mocchiut 1.73 if ( verbose ) printf(" gat1\n");
1062 pam-mep 1.69 trke->Clear();
1063     //
1064 mocchiut 1.73 // Clones array must be cleared before going on
1065     //
1066     if ( hasNucleiTrk ){
1067     if ( verbose ) printf(" gat2\n");
1068     tcNucleiTrk->Delete();
1069     if ( verbose ) printf(" gat3\n");
1070     torbNucleiTrk->Delete();
1071     }
1072     if ( hasExtNucleiTrk ){
1073     if ( verbose ) printf(" gat4\n");
1074     tcExtNucleiTrk->Delete();
1075     if ( verbose ) printf(" gat5\n");
1076     torbExtNucleiTrk->Delete();
1077     }
1078     if ( hasExtTrk ){
1079     if ( verbose ) printf(" gat6\n");
1080     tcExtTrk->Delete();
1081     if ( verbose ) printf(" gat7\n");
1082     torbExtTrk->Delete();
1083     }
1084     //
1085     if ( verbose ) printf(" get trk tree entries... entry = %i in Level2 file\n",itr);
1086 pam-mep 1.69 if ( ttrke->GetEntry(itr) <= 0 ) throw -36;
1087     //
1088     }
1089    
1090 mocchiut 1.1 //
1091     procev++;
1092     //
1093     // start processing
1094     //
1095 mocchiut 1.40 if ( debug ) printf(" %i start processing \n",procev);
1096 mocchiut 1.3 orbitalinfo->Clear();
1097 mocchiut 1.5 //
1098 mocchiut 1.30 OrbitalInfoTrkVar *t_orb = new OrbitalInfoTrkVar();
1099     if( !(orbitalinfo->OrbitalInfoTrk) ) orbitalinfo->OrbitalInfoTrk = new TClonesArray("OrbitalInfoTrkVar",2);
1100     TClonesArray &tor = *orbitalinfo->OrbitalInfoTrk;
1101 pam-mep 1.69
1102     // Geomagnetic coordinates calculation variables
1103     GMtype_CoordSpherical CoordSpherical, DipoleSpherical;
1104     GMtype_CoordCartesian CoordCartesian, DipoleCartesian;
1105     GMtype_Model Model;
1106     GMtype_Pole Pole;
1107    
1108 mocchiut 1.30 //
1109 mocchiut 1.15 // Fill OBT, pkt_num and absTime
1110     //
1111 mocchiut 1.2 orbitalinfo->pkt_num = ph->GetCounter();
1112     orbitalinfo->OBT = ph->GetOrbitalTime();
1113 mocchiut 1.15 orbitalinfo->absTime = atime;
1114 mocchiut 1.40 if ( debug ) printf(" %i pktnum obt abstime \n",procev);
1115 mocchiut 1.15 //
1116     // Propagate the orbit from the tle time to atime, using SGP(D)4.
1117     //
1118 mocchiut 1.40 if ( debug ) printf(" %i sgp4 \n",procev);
1119 mocchiut 1.15 cCoordGeo coo;
1120 mocchiut 1.40 Float_t jyear=0.;
1121 mocchiut 1.7 //
1122     if(atime >= gltle->GetToTime()) {
1123 mocchiut 1.9 if ( !gltle->Query(atime, dbc) ){
1124 mocchiut 1.15 //
1125 mocchiut 1.9 // Compute the magnetic dipole moment.
1126 mocchiut 1.15 //
1127 mocchiut 1.40 if ( debug ) printf(" %i compute magnetic dipole moment \n",procev);
1128 mocchiut 1.9 UInt_t year, month, day, hour, min, sec;
1129 mocchiut 1.15 //
1130 mocchiut 1.9 TTimeStamp t = TTimeStamp(atime, kTRUE);
1131     t.GetDate(kTRUE, 0, &year, &month, &day);
1132     t.GetTime(kTRUE, 0, &hour, &min, &sec);
1133     jyear = (float) year
1134     + (month*31.+ (float) day)/365.
1135 mocchiut 1.44 + (hour*3600.+min*60.+(float)sec)/(24.*3600.*365.);
1136 mocchiut 1.15 //
1137 pam-mep 1.69 if ( debug ) printf(" %i compute magnetic dipole moment get dipole moment for year\n",procev);
1138 mocchiut 1.70 if ( debug ) printf(" %i jyear %f dimo %f \n",procev,jyear,dimo);
1139 mocchiut 1.9 feldcof_(&jyear, &dimo); // get dipole moment for year
1140 pam-mep 1.69 if ( debug ) printf(" %i compute magnetic dipole moment end\n",procev);
1141    
1142     GM_TimeAdjustCoefs(year, jyear, G0, G1, H1, &Model);
1143     GM_PoleLocation(Model, &Pole);
1144    
1145 mocchiut 1.9 } else {
1146     code = -56;
1147     goto closeandexit;
1148     };
1149 mocchiut 1.7 }
1150 mocchiut 1.15 coo = getCoo(atime, gltle->GetFromTime(), gltle->GetTle());
1151     //
1152     cOrbit orbits(*gltle->GetTle());
1153     //
1154     // synchronize with quaternions data
1155     //
1156 mocchiut 1.18 if ( isf && neventsm>0 ){
1157 mocchiut 1.15 //
1158     // First event
1159     //
1160     isf = false;
1161 mocchiut 1.64 // upperqtime = atime;
1162 mocchiut 1.15 lowerqtime = runinfo->RUNHEADER_TIME;
1163 mocchiut 1.44 for ( ik = 0; ik < neventsm; ik++){ //number of macrocommad packets
1164 mocchiut 1.43 if ( ch->GetEntry(ik) <= 0 ) throw -36;
1165 mocchiut 1.15 tmpSize = mcmdev->Records->GetEntries();
1166 mocchiut 1.64 // numrec = tmpSize;
1167 mocchiut 1.70 if ( debug ) cout << "packet number " << ik <<"\tnumber of subpackets is " << tmpSize << endl;
1168 mocchiut 1.44 for (Int_t j3 = 0;j3<tmpSize;j3++){ //number of subpackets
1169 mocchiut 1.15 mcmdrc = (pamela::McmdRecord*)mcmdev->Records->At(j3);
1170 mocchiut 1.38 if ( mcmdrc ){ // missing inclination bug [8RED 090116]
1171 mocchiut 1.70 if ( debug ) printf(" pluto \n");
1172 mocchiut 1.44 if ((int)mcmdrc->ID1 == 226 && mcmdrc->Mcmd_Block_crc_ok == 1){ //Check that it is Inclination Packet
1173 pam-mep 1.69 L_QQ_Q_l_upper->fill(mcmdrc->McmdData);
1174 mocchiut 1.38 for (UInt_t ui = 0; ui < 6; ui++){
1175     if (ui>0){
1176     if (L_QQ_Q_l_upper->time[ui]>L_QQ_Q_l_upper->time[0]){
1177 mocchiut 1.65 if ( debug ) printf(" here1 %i \n",ui);
1178 mocchiut 1.44 Double_t u_time = dbtime->DBabsTime((UInt_t)(L_QQ_Q_l_upper->time[ui]*1000-DeltaOBT*1000));
1179     Int_t recSize = recqtime.size();
1180 pam-mep 1.54 if(lowerqtime > recqtime[recSize-1]){
1181 pam-mep 1.69 // to avoid interpolation between bad quaternions arrays
1182     if(sqrt(pow(L_QQ_Q_l_upper->quat[ui][0],2)+pow(L_QQ_Q_l_upper->quat[ui][1],2)+pow(L_QQ_Q_l_upper->quat[ui][2],2)+pow(L_QQ_Q_l_upper->quat[ui][3],2))>0.99999){
1183     Int_t sizeqmcmd = qtime.size();
1184     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1185     qtime[sizeqmcmd]=u_time;
1186     q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][0];
1187     q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][1];
1188     q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][2];
1189     q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][3];
1190     qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
1191     lowerqtime = u_time;
1192     orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
1193     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,L_QQ_Q_l_upper->quat[ui][0],L_QQ_Q_l_upper->quat[ui][1],L_QQ_Q_l_upper->quat[ui][2],L_QQ_Q_l_upper->quat[ui][3]);
1194     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1195     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1196     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1197     }
1198 pam-mep 1.54 }
1199 mocchiut 1.44 for(Int_t mu = nt;mu<recSize;mu++){
1200     if(recqtime[mu]>lowerqtime && recqtime[mu]<u_time){
1201 pam-mep 1.69 if(sqrt(pow(recq0[mu],2)+pow(recq1[mu],2)+pow(recq2[mu],2)+pow(recq3[mu],2))>0.99999){
1202     nt=mu;
1203     Int_t sizeqmcmd = qtime.size();
1204     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1205     qtime[sizeqmcmd]=recqtime[mu];
1206     q0[sizeqmcmd]=recq0[mu];
1207     q1[sizeqmcmd]=recq1[mu];
1208     q2[sizeqmcmd]=recq2[mu];
1209     q3[sizeqmcmd]=recq3[mu];
1210     qmode[sizeqmcmd]=-10;
1211     orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
1212     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,recq0[mu],recq1[mu],recq2[mu],recq3[mu]);
1213     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1214     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1215     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1216     }
1217 mocchiut 1.44 }
1218     if(recqtime[mu]>=u_time){
1219 pam-mep 1.69 if(sqrt(pow(L_QQ_Q_l_upper->quat[ui][0],2)+pow(L_QQ_Q_l_upper->quat[ui][1],2)+pow(L_QQ_Q_l_upper->quat[ui][2],2)+pow(L_QQ_Q_l_upper->quat[ui][3],2))>0.99999){
1220     Int_t sizeqmcmd = qtime.size();
1221     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1222     qtime[sizeqmcmd]=u_time;
1223     q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][0];
1224     q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][1];
1225     q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][2];
1226     q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][3];
1227     qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
1228     lowerqtime = u_time;
1229     orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
1230     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,L_QQ_Q_l_upper->quat[ui][0],L_QQ_Q_l_upper->quat[ui][1],L_QQ_Q_l_upper->quat[ui][2],L_QQ_Q_l_upper->quat[ui][3]);
1231     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1232     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1233     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1234     break;
1235     }
1236 mocchiut 1.38 }
1237     }
1238     }
1239     }else{
1240 mocchiut 1.65 if ( debug ) printf(" here2 %i \n",ui);
1241 mocchiut 1.44 Double_t u_time = dbtime->DBabsTime((UInt_t)(L_QQ_Q_l_upper->time[0]*1000-DeltaOBT*1000));
1242     if(lowerqtime>u_time)nt=0;
1243     Int_t recSize = recqtime.size();
1244 pam-mep 1.54 if(lowerqtime > recqtime[recSize-1]){
1245 pam-mep 1.69 if(sqrt(pow(L_QQ_Q_l_upper->quat[ui][0],2)+pow(L_QQ_Q_l_upper->quat[ui][1],2)+pow(L_QQ_Q_l_upper->quat[ui][2],2)+pow(L_QQ_Q_l_upper->quat[ui][3],2))>0.99999){
1246     Int_t sizeqmcmd = qtime.size();
1247     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1248     qtime[sizeqmcmd]=u_time;
1249     q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][0];
1250     q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][1];
1251     q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][2];
1252     q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][3];
1253     qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
1254     lowerqtime = u_time;
1255     orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
1256     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,L_QQ_Q_l_upper->quat[0][0],L_QQ_Q_l_upper->quat[0][1],L_QQ_Q_l_upper->quat[0][2],L_QQ_Q_l_upper->quat[0][3]);
1257     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1258     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1259     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1260     }
1261 pam-mep 1.54 }
1262 mocchiut 1.44 for(Int_t mu = nt;mu<recSize;mu++){
1263     if(recqtime[mu]>lowerqtime && recqtime[mu]<u_time){
1264 pam-mep 1.69 if(sqrt(pow(recq0[mu],2)+pow(recq1[mu],2)+pow(recq2[mu],2)+pow(recq3[mu],2))>0.99999){
1265     nt=mu;
1266     Int_t sizeqmcmd = qtime.size();
1267     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1268     qtime[sizeqmcmd]=recqtime[mu];
1269     q0[sizeqmcmd]=recq0[mu];
1270     q1[sizeqmcmd]=recq1[mu];
1271     q2[sizeqmcmd]=recq2[mu];
1272     q3[sizeqmcmd]=recq3[mu];
1273     qmode[sizeqmcmd]=-10;
1274     orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
1275     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,recq0[mu],recq1[mu],recq2[mu],recq3[mu]);
1276     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1277     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1278     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1279     }
1280 mocchiut 1.44 }
1281     if(recqtime[mu]>=u_time){
1282 pam-mep 1.69 if(sqrt(pow(L_QQ_Q_l_upper->quat[0][0],2)+pow(L_QQ_Q_l_upper->quat[0][1],2)+pow(L_QQ_Q_l_upper->quat[0][2],2)+pow(L_QQ_Q_l_upper->quat[0][3],2))>0.99999){
1283     Int_t sizeqmcmd = qtime.size();
1284     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1285     qtime[sizeqmcmd]=u_time;
1286     q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][0];
1287     q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][1];
1288     q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][2];
1289     q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][3];
1290     qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
1291     lowerqtime = u_time;
1292     orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
1293     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,L_QQ_Q_l_upper->quat[0][0],L_QQ_Q_l_upper->quat[0][1],L_QQ_Q_l_upper->quat[0][2],L_QQ_Q_l_upper->quat[0][3]);
1294     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1295     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1296     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1297     CopyQ(L_QQ_Q_l_lower,L_QQ_Q_l_upper);
1298     break;
1299     }
1300 mocchiut 1.15 }
1301     }
1302 mocchiut 1.44 }
1303 mocchiut 1.38 }
1304 mocchiut 1.15 }
1305 mocchiut 1.44 }
1306 pam-mep 1.69 //if ( debug ) cout << "subpacket " << j3 << "\t qtime = " << qtime[qtime.size()-1] << endl;
1307 mocchiut 1.15 }
1308 mocchiut 1.44 }
1309 mocchiut 1.47
1310 pam-mep 1.69 if(qtime.size()==0){ // in case if no orientation information in data
1311 mocchiut 1.70 if ( debug ) cout << "qtime.size() = 0" << endl;
1312 mocchiut 1.65 for(UInt_t my=0;my<recqtime.size();my++){
1313 mocchiut 1.70 if(sqrt(pow(recq0[my],2)+pow(recq1[my],2)+pow(recq2[my],2)+pow(recq3[my],2))>0.99999){
1314     Int_t sizeqmcmd = qtime.size();
1315     inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1316     qtime[sizeqmcmd]=recqtime[my];
1317     q0[sizeqmcmd]=recq0[my];
1318     q1[sizeqmcmd]=recq1[my];
1319     q2[sizeqmcmd]=recq2[my];
1320     q3[sizeqmcmd]=recq3[my];
1321     qmode[sizeqmcmd]=-10;
1322     orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
1323     RYPang_upper->TransAngle(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,recq0[my],recq1[my],recq2[my],recq3[my]);
1324     qRoll[sizeqmcmd]=RYPang_upper->Kren;
1325     qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1326     qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1327     }
1328 mocchiut 1.65 }
1329 mocchiut 1.47 }
1330    
1331 mocchiut 1.49
1332 mocchiut 1.46 if ( debug ) printf(" puffi \n");
1333 mocchiut 1.44 Double_t tmin = 9999999999.;
1334     Double_t tmax = 0.;
1335     for(UInt_t tre = 0;tre<qtime.size();tre++){
1336     if(qtime[tre]>tmax)tmax = qtime[tre];
1337     if(qtime[tre]<tmin)tmin = qtime[tre];
1338     }
1339 pam-mep 1.69 // sorting quaternions by time
1340     Bool_t t = true;
1341     while(t){
1342     t=false;
1343     for(UInt_t i=0;i<qtime.size()-1;i++){
1344     if(qtime[i]>qtime[i+1]){
1345     Double_t tmpr = qtime[i];
1346     qtime[i]=qtime[i+1];
1347     qtime[i+1] = tmpr;
1348     tmpr = q0[i];
1349     q0[i]=q0[i+1];
1350     q0[i+1] = tmpr;
1351     tmpr = q1[i];
1352     q1[i]=q1[i+1];
1353     q1[i+1] = tmpr;
1354     tmpr = q2[i];
1355     q2[i]=q2[i+1];
1356     q2[i+1] = tmpr;
1357     tmpr = q3[i];
1358     q3[i]=q3[i+1];
1359     q3[i+1] = tmpr;
1360     tmpr = qRoll[i];
1361     qRoll[i]=qRoll[i+1];
1362     qRoll[i+1] = tmpr;
1363     tmpr = qYaw[i];
1364     qYaw[i]=qYaw[i+1];
1365     qYaw[i+1] = tmpr;
1366     tmpr = qPitch[i];
1367     qPitch[i]=qPitch[i+1];
1368     qPitch[i+1] = tmpr;
1369     t=true;
1370     }
1371     }
1372     }
1373    
1374     if ( debug ){
1375     cout << "we have loaded quaternions: size of quaternions set is "<< qtime.size() << endl;
1376     for(UInt_t i=0;i<qtime.size();i++) cout << qtime[i] << "\t";
1377     cout << endl << endl;
1378     Int_t lopu;
1379     cin >> lopu;
1380     }
1381 mocchiut 1.46
1382 mocchiut 1.44 } // if we processed first event
1383 pam-mep 1.60
1384 mocchiut 1.44
1385     //Filling Inclination information
1386     Double_t incli = 0;
1387 mocchiut 1.46 if ( qtime.size() > 1 ){
1388 pam-mep 1.69 if ( debug ) cout << "ok quaternions is exist and mu = " << must << endl;
1389     if ( debug ) cout << "qtimes[ " << qtime[0] << " , " << qtime[qtime.size()-1] << " ]\tatime = "<<atime<<endl;
1390 mocchiut 1.65 for(UInt_t mu = must;mu<qtime.size()-1;mu++){
1391     if ( debug ) printf(" ??grfuffi %i sixe %i must %i \n",mu,qtime.size()-1,must);
1392     if(qtime[mu+1]>qtime[mu]){
1393 mocchiut 1.70 if ( debug ) cout << "qtime[" << mu << "] = " << qtime[mu] << "\tqtime[" << mu+1 << "] = " << qtime[mu+1] << "\tatime = " << atime << endl;
1394 mocchiut 1.65 if(atime<=qtime[mu+1] && atime>=qtime[mu]){
1395 mocchiut 1.70 if ( debug ) cout << "here we have found proper quaternions for interpolation: mu = "<<mu<<endl;
1396 mocchiut 1.65 must = mu;
1397     incli = (qPitch[mu+1]-qPitch[mu])/(qtime[mu+1]-qtime[mu]);
1398     orbitalinfo->theta = incli*atime+qPitch[mu+1]-incli*qtime[mu+1];
1399     incli = (qRoll[mu+1]-qRoll[mu])/(qtime[mu+1]-qtime[mu]);
1400     orbitalinfo->etha = incli*atime+qRoll[mu+1]-incli*qtime[mu+1];
1401     incli = (qYaw[mu+1]-qYaw[mu])/(qtime[mu+1]-qtime[mu]);
1402     orbitalinfo->phi = incli*atime+qYaw[mu+1]-incli*qtime[mu+1];
1403    
1404     incli = (q0[mu+1]-q0[mu])/(qtime[mu+1]-qtime[mu]);
1405     orbitalinfo->q0 = incli*atime+q0[mu+1]-incli*qtime[mu+1];
1406     incli = (q1[mu+1]-q1[mu])/(qtime[mu+1]-qtime[mu]);
1407     orbitalinfo->q1 = incli*atime+q1[mu+1]-incli*qtime[mu+1];
1408     incli = (q2[mu+1]-q2[mu])/(qtime[mu+1]-qtime[mu]);
1409     orbitalinfo->q2 = incli*atime+q2[mu+1]-incli*qtime[mu+1];
1410     incli = (q3[mu+1]-q3[mu])/(qtime[mu+1]-qtime[mu]);
1411     orbitalinfo->q3 = incli*atime+q3[mu+1]-incli*qtime[mu+1];
1412 mocchiut 1.70 Float_t tg = (qtime[mu+1]-qtime[mu])/1000.;
1413     if(tg>=1) tg=0.00;
1414     orbitalinfo->TimeGap = TMath::Min(TMath::Abs(qtime[mu+1])-atime,TMath::Abs(atime-qtime[mu]))+tg;//qtime[mu+1]-qtime[mu];
1415     orbitalinfo->mode = qmode[mu+1];
1416 mocchiut 1.65
1417 mocchiut 1.70 //if(atime==qtime[mu] || atime==qtime[mu+1]) orbitalinfo->qkind = 0; else orbitalinfo->qkind=1;
1418 mocchiut 1.65 //if(qmode[mu+1]==-10) orbitalinfo->R10r = true;else orbitalinfo->R10r = false;
1419     if ( debug ) printf(" grfuffi4 %i \n",mu);
1420 mocchiut 1.70
1421 mocchiut 1.65 break;
1422     }
1423     }
1424     }
1425 mocchiut 1.46 }
1426 mocchiut 1.65 if ( debug ) printf(" grfuffi5 \n");
1427 mocchiut 1.30 //
1428     // ops no inclination information
1429     //
1430 mocchiut 1.65
1431 mocchiut 1.30 if ( orbitalinfo->q0< -999 || orbitalinfo->q1 < -999 || orbitalinfo->q2 < -999 || orbitalinfo->q3 < -999 || orbitalinfo->q0 != orbitalinfo->q0 || orbitalinfo->q1 != orbitalinfo->q1 || orbitalinfo->q2 != orbitalinfo->q2 || orbitalinfo->q3 != orbitalinfo->q3 ){
1432 pam-mep 1.69 if ( debug ) cout << "ops no iclination information" << endl;
1433 mocchiut 1.30 orbitalinfo->mode = 10;
1434     orbitalinfo->q0 = -1000.;
1435     orbitalinfo->q1 = -1000.;
1436     orbitalinfo->q2 = -1000.;
1437     orbitalinfo->q3 = -1000.;
1438     orbitalinfo->etha = -1000.;
1439     orbitalinfo->phi = -1000.;
1440     orbitalinfo->theta = -1000.;
1441 pam-mep 1.69 orbitalinfo->TimeGap = -1000.;
1442     //orbitalinfo->qkind = -1000;
1443 mocchiut 1.70
1444     // if ( debug ){
1445     // Int_t lopu;
1446     // cin >> lopu;
1447     // }
1448 mocchiut 1.65 if ( debug ) printf(" grfuffi6 \n");
1449     }
1450 mocchiut 1.30 //
1451 mocchiut 1.65 if ( debug ) printf(" filling \n");
1452 mocchiut 1.30 // #########################################################################################################################
1453 mocchiut 1.15 //
1454     // fill orbital positions
1455     //
1456 mocchiut 1.7 // Build coordinates in the right range. We want to convert,
1457     // longitude from (0, 2*pi) to (-180deg, 180deg). Altitude is
1458     // in meters.
1459     lon = (coo.m_Lon > M_PI) ? rad2deg(coo.m_Lon - 2*M_PI) : rad2deg(coo.m_Lon);
1460     lat = rad2deg(coo.m_Lat);
1461     alt = coo.m_Alt;
1462 pam-mep 1.69
1463     cOrbit orbits2(*gltle->GetTle());
1464     orbits2.getPosition((double) (atime - gltle->GetFromTime())/60., &eCi);
1465 mocchiut 1.70 // Float_t x=eCi.getPos().m_x;
1466     // Float_t y=eCi.getPos().m_y;
1467     // Float_t z=eCi.getPos().m_z;
1468    
1469 pam-mep 1.69 TVector3 V(eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z);
1470     TVector3 Pos(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z);
1471 mocchiut 1.70
1472 pam-mep 1.69 Float_t dlon=Pos.Phi()*TMath::RadToDeg()-lon;
1473 mocchiut 1.70
1474 pam-mep 1.69 Pos.RotateZ(-dlon*TMath::DegToRad());
1475     V.RotateZ(-dlon*TMath::DegToRad());
1476     Float_t diro;
1477     if(V.Z()>0) diro=1; else diro=-1;
1478 mocchiut 1.70
1479     // 10REDNEW
1480 pam-mep 1.69 Int_t errq=0;
1481     Int_t azim=0;;
1482 mocchiut 1.70 for(UInt_t mu = must;mu<RTtime2.size()-1;mu++){
1483     if(atime<=RTtime2[mu] && atime>=RTtime1[mu]){
1484     errq=RTerrq[mu];
1485     azim=RTazim[mu];
1486     }
1487     }
1488     orbitalinfo->errq = errq;
1489     orbitalinfo->azim = azim;
1490     orbitalinfo->qkind = 0;
1491    
1492 mocchiut 1.65 if ( debug ) printf(" coord done \n");
1493 mocchiut 1.7 if( lon<180 && lon>-180 && lat<90 && lat>-90 && alt>0 ){
1494 mocchiut 1.15 //
1495 mocchiut 1.7 orbitalinfo->lon = lon;
1496     orbitalinfo->lat = lat;
1497 pam-mep 1.69 orbitalinfo->alt = alt;
1498     orbitalinfo->V = V;
1499    
1500 mocchiut 1.70 // GMtype_CoordGeodetic location;
1501 pam-mep 1.69 location.lambda = lon;
1502     location.phi = lat;
1503     location.HeightAboveEllipsoid = alt;
1504    
1505     GM_GeodeticToSpherical(Ellip, location, &CoordSpherical);
1506     GM_SphericalToCartesian(CoordSpherical, &CoordCartesian);
1507     GM_EarthCartToDipoleCartCD(Pole, CoordCartesian, &DipoleCartesian);
1508     GM_CartesianToSpherical(DipoleCartesian, &DipoleSpherical);
1509     orbitalinfo->londip = DipoleSpherical.lambda;
1510     orbitalinfo->latdip = DipoleSpherical.phig;
1511    
1512     if(debug)cout<<"geodetic:\t"<<lon<<"\t"<<lat<<"\tgeomagnetic:\t"<<orbitalinfo->londip<<"\t"<<orbitalinfo->latdip<<endl;
1513    
1514 mocchiut 1.15 //
1515 mocchiut 1.7 // compute mag field components and L shell.
1516 mocchiut 1.15 //
1517 mocchiut 1.65 if ( debug ) printf(" call igrf feldg \n");
1518 mocchiut 1.7 feldg_(&lat, &lon, &alt, &bnorth, &beast, &bdown, &babs);
1519 mocchiut 1.65 if ( debug ) printf(" call igrf shellg \n");
1520 mocchiut 1.7 shellg_(&lat, &lon, &alt, &dimo, &xl, &icode, &bab1);
1521 mocchiut 1.65 if ( debug ) printf(" call igrf findb \n");
1522 mocchiut 1.7 findb0_(&stps, &bdel, &value, &bequ, &rr0);
1523 mocchiut 1.15 //
1524 mocchiut 1.65 if ( debug ) printf(" done igrf \n");
1525 mocchiut 1.7 orbitalinfo->Bnorth = bnorth;
1526     orbitalinfo->Beast = beast;
1527     orbitalinfo->Bdown = bdown;
1528     orbitalinfo->Babs = babs;
1529 pam-mep 1.60 orbitalinfo->M = dimo;
1530 mocchiut 1.7 orbitalinfo->BB0 = babs/bequ;
1531 mocchiut 1.15 orbitalinfo->L = xl;
1532 mocchiut 1.7 // Set Stormer vertical cutoff using L shell.
1533 mocchiut 1.57 orbitalinfo->cutoffsvl = 14.295 / (xl*xl); //
1534 pam-mep 1.69 if(debug)cout << "L = " << xl << "\tM = " << dimo << "\tvertical cutoff: "<< orbitalinfo->cutoffsvl << endl;
1535    
1536 mocchiut 1.57 /*
1537 mocchiut 1.65 ---------- Forwarded message ----------
1538     Date: Wed, 09 May 2012 12:16:47 +0200
1539     From: Alessandro Bruno <alessandro.bruno@ba.infn.it>
1540     To: Mirko Boezio <mirko.boezio@ts.infn.it>
1541     Cc: Francesco S. Cafagna <Francesco.Cafagna@ba.infn.it>
1542     Subject: Störmer vertical cutoff
1543    
1544     Ciao Mirko,
1545     volevo segnalarti che il valore dello Störmer vertical cutoff nel Level2 è
1546     sovrastimato di circa il 4%.
1547     Dopo un'approfondita analisi con l'IGRF-05 abbiamo ricavano un valore pari
1548     a: 14.295 / L^2 anzichè 14.9 / L^2, valore obsoleto in quanto riferito agli
1549     anni '50.
1550     */
1551 mocchiut 1.57 //14.9/(xl*xl);
1552 mocchiut 1.56 orbitalinfo->igrf_icode = icode;
1553 mocchiut 1.15 //
1554 mocchiut 1.65 }
1555 pamelaprod 1.11 //
1556 mocchiut 1.40 if ( debug ) printf(" pitch angle \n");
1557     //
1558 mocchiut 1.30 // pitch angles
1559     //
1560 pam-mep 1.69 if( orbitalinfo->TimeGap>0){
1561 mocchiut 1.30 //
1562 mocchiut 1.65 if ( debug ) printf(" timegap %f \n",orbitalinfo->TimeGap);
1563 mocchiut 1.49 Float_t Bx = -orbitalinfo->Bdown;
1564     Float_t By = orbitalinfo->Beast;
1565     Float_t Bz = orbitalinfo->Bnorth;
1566 pam-mep 1.69
1567     TMatrixD Qiji(3,3);
1568     TMatrixD Qij = PO->QuatoECI(orbitalinfo->q0,orbitalinfo->q1,orbitalinfo->q2,orbitalinfo->q3);
1569     TMatrixD Dij = PO->ECItoGEO(Qij,orbitalinfo->absTime,orbitalinfo->lat,orbitalinfo->lon);
1570    
1571     //10REDNEW
1572     /* If initial orientation data have reason to be inaccurate */
1573     Double_t tg = 0;
1574 mocchiut 1.72 if ( debug ) cout<<modf(orbitalinfo->TimeGap,&tg)<<endl;
1575 pam-mep 1.69 // if(orbitalinfo->TimeGap>0 && errq==0 && azim==0){ // 10RED CHECK (comparison between three metod of recovering orientation)
1576     if(((orbitalinfo->TimeGap>60.0 && TMath::Abs(orbitalinfo->etha)>0.5) || errq!=0 || modf(orbitalinfo->TimeGap,&tg)*1000>700 || modf(orbitalinfo->TimeGap,&tg)*1000==0.0 ) && azim==0){ //Standard condition to use this; One of these two cases should be commented
1577     /* found in Rotation Table this data for this time interval*/
1578     if(atime<RTtime1[0])
1579     orbitalinfo->azim = 5; //means that RotationTable no started yet
1580     else{
1581     for(UInt_t mu = must;mu<RTtime2.size()-1;mu++){
1582     if(atime<=RTtime2[mu] && atime>=RTtime1[mu]){
1583     // search for angle betwean velosity and direction to north in tangential to Earth surfase plane in satellite position
1584     Double_t tlat=orbitalinfo->lat;
1585 mocchiut 1.73 Double_t kar=(RTbank2[mu]-RTbank1[mu])/(RTtime2[mu]-RTtime1[mu]);
1586 pam-mep 1.69 Double_t bak=RTbank1[mu]-kar*RTtime1[mu];
1587     Double_t bank=kar*atime+bak;
1588     Float_t spitch = 0.00001; // temprary not zero to avoid problem with tranzition from Euler angles to orientation matrix
1589    
1590     //Estimations of pitch angle of satellite
1591     if(TMath::Abs(bank)>0.7){
1592     Float_t spitch1=TMath::DegToRad()*0.7*RTdir1[mu];
1593     Float_t spitch2=TMath::DegToRad()*0.7*RTdir2[mu];
1594     Float_t kva=(spitch2-spitch1)/(RTtime2[mu]-RTtime1[mu]);
1595     Float_t bva=spitch1-kva*RTtime1[mu];
1596     spitch=kva*atime+bva;
1597     }
1598 mocchiut 1.73
1599 pam-mep 1.69 //Calculate Yaw angle accordingly with fit, see picture FitYaw.jpg
1600     Double_t yaw=0.00001; // temprary not zero to avoid problem with tranzition from Euler angles to orientation matrix
1601     if(TMath::Abs(tlat)<70)
1602     yaw = -3.7e-8*tlat*tlat*tlat*tlat + 1.4e-7*tlat*tlat*tlat - 0.0005*tlat*tlat - 0.00025*tlat + 3.6;
1603     yaw = diro*yaw; //because should be different sign for ascending and descending orbits!
1604    
1605     if(TMath::Abs(bank)>0.5 && TMath::Abs(yaw-orbitalinfo->phi)<3.0) yaw=orbitalinfo->phi;
1606    
1607     // Qiji = PO->EulertoEci(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,bank,yaw,spitch); // 10RED CHECK
1608     Qij = PO->EulertoEci(eCi.getPos().m_x,eCi.getPos().m_y,eCi.getPos().m_z,eCi.getVel().m_x,eCi.getVel().m_y,eCi.getVel().m_z,bank,yaw,spitch); // STANDARD
1609     orbitalinfo->qkind = 1;
1610    
1611     break;
1612     }
1613     } // enf of loop for(UInt_t mu = must;mu<RTtime2.size()-1;mu++){
1614    
1615     //Qij = PO->GEOtoECI(Dij,orbitalinfo->absTime,orbitalinfo->lat,orbitalinfo->lon); // to convert from Dij to Qij
1616    
1617     } // end of if(atime<RTtime1[0]
1618     } // end of f(((orbitalinfo->TimeGap>60.0 && TMath...
1619    
1620     TMatrixD qij = PO->ColPermutation(Qij);
1621     TMatrixD Fij = PO->ECItoGreenwich(Qij,orbitalinfo->absTime);
1622 mocchiut 1.49 TMatrixD Gij = PO->ColPermutation(Fij);
1623 pam-mep 1.69 Dij = PO->ECItoGEO(Qij,orbitalinfo->absTime,orbitalinfo->lat,orbitalinfo->lon);
1624 mocchiut 1.30 TMatrixD Iij = PO->ColPermutation(Dij);
1625 pam-mep 1.69 TVector3 SP = PO->GetSunPosition(orbitalinfo->absTime);
1626     // go to Pamela reference frame from Resurs reference frame
1627     Float_t tmpy = SP.Y();
1628     SP.SetY(SP.Z());
1629     SP.SetZ(-tmpy);
1630     TVector3 SunZenith;
1631     SunZenith.SetMagThetaPhi(1,23.439281*TMath::DegToRad(),TMath::Pi()/2.);
1632     TVector3 SunMag = -SP;
1633     SunMag.Rotate(-45*TMath::DegToRad(),SunZenith);
1634     tmpy=SunMag.Y();
1635     SunMag.SetY(SunMag.Z());
1636     SunMag.SetZ(-tmpy);
1637    
1638 mocchiut 1.33 orbitalinfo->Iij.ResizeTo(Iij);
1639     orbitalinfo->Iij = Iij;
1640     //
1641 mocchiut 1.64 // A1 = Iij(0,2);
1642     // A2 = Iij(1,2);
1643     // A3 = Iij(2,2);
1644 mocchiut 1.49 //
1645 mocchiut 1.33 // orbitalinfo->pamzenitangle = (Float_t)PO->GetPitchAngle(1,0,0,A1,A2,A3); // Angle between zenit and Pamela's main axiz
1646     // orbitalinfo->pamBangle = (Float_t)PO->GetPitchAngle(A1,A2,A3,Bx,By,Bz); // Angle between Pamela's main axiz and B
1647 mocchiut 1.30 //
1648 mocchiut 1.65 if ( debug ) printf(" matrixes done \n");
1649 mocchiut 1.73 if ( !standalone ){
1650 mocchiut 1.65 if ( debug ) printf(" !standalone \n");
1651 mocchiut 1.30 //
1652 mocchiut 1.73 // Standard tracking algorithm
1653     //
1654 mocchiut 1.30 Int_t nn = 0;
1655 mocchiut 1.73 if ( verbose ) printf(" standard tracking \n");
1656 mocchiut 1.30 for(Int_t nt=0; nt < tof->ntrk(); nt++){
1657     //
1658     ToFTrkVar *ptt = tof->GetToFTrkVar(nt);
1659 pam-mep 1.69 if (debug) cout<<"tof->ntrk() = "<<tof->ntrk()<<"\tptt->trkseqno = "<<ptt->trkseqno<<"\ttrke->ntrk() = "<<trke->ntrk()<<endl;
1660 mocchiut 1.30 Double_t E11x = ptt->xtr_tof[0]; // tr->x[0];
1661     Double_t E11y = ptt->ytr_tof[0]; //tr->y[0];
1662     Double_t E11z = zin[0];
1663     Double_t E22x = ptt->xtr_tof[3];//tr->x[3];
1664     Double_t E22y = ptt->ytr_tof[3];//tr->y[3];
1665     Double_t E22z = zin[3];
1666 mocchiut 1.32 if ( (E11x < 100. && E11y < 100. && E22x < 100. && E22y < 100.) || ptt->trkseqno != -1 ){
1667 pam-mep 1.69 TrkTrack *mytrack = trke->GetStoredTrack(ptt->trkseqno);
1668     Float_t rig=1/mytrack->GetDeflection();
1669 mocchiut 1.30 Double_t norm = sqrt(pow(E22x-E11x,2)+pow(E22y-E11y,2)+pow(E22z-E11z,2));
1670 mocchiut 1.73 //
1671 mocchiut 1.30 Px = (E22x-E11x)/norm;
1672     Py = (E22y-E11y)/norm;
1673     Pz = (E22z-E11z)/norm;
1674     //
1675 mocchiut 1.33 t_orb->trkseqno = ptt->trkseqno;
1676     //
1677     TMatrixD Eij = PO->PamelatoGEO(Iij,Px,Py,Pz);
1678     t_orb->Eij.ResizeTo(Eij);
1679     t_orb->Eij = Eij;
1680     //
1681 pam-mep 1.51 TMatrixD Sij = PO->PamelatoGEO(Gij,Px,Py,Pz);
1682 mocchiut 1.33 t_orb->Sij.ResizeTo(Sij);
1683     t_orb->Sij = Sij;
1684 mocchiut 1.30 //
1685     t_orb->pitch = (Float_t)PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),Bx,By,Bz);
1686 mocchiut 1.33 //
1687 mocchiut 1.73 // SunPosition in instrumental reference frame
1688 pam-mep 1.69 TMatrixD Kij = PO->PamelatoGEO(qij,Px,Py,Pz);
1689     TMatrixD Lij = PO->PamelatoGEO(qij,0,0,1);
1690     t_orb->sunangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),-SP.X(),-SP.Y(),-SP.Z());
1691     t_orb->sunmagangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),SunMag.X(),SunMag.Y(),SunMag.Z());
1692 mocchiut 1.33 //
1693     //
1694 pam-mep 1.69 Double_t omega = PO->GetPitchAngle(-Eij(0,0),-Eij(1,0),-Eij(2,0),1,0,0) * TMath::DegToRad();
1695     TVector3 Bxy(0,By,Bz);
1696     TVector3 Exy(0,-Eij(1,0),-Eij(2,0));
1697     Double_t dzeta=Bxy.Angle(Exy);
1698     if (-Eij(1,0) < 0) dzeta=2.0*TMath::Pi() - dzeta;
1699 mocchiut 1.73
1700     if(debug) cout << "omega = "<<omega*TMath::RadToDeg()<<"\tdzeta = "<<dzeta*TMath::RadToDeg()<<endl;
1701 pam-mep 1.69
1702 mocchiut 1.73 // Formula from D.F. Smart *, M.A. Shea [2005]; A review of geomagnetic cutoff rigidities for earth-orbiting spacecraft
1703 pam-mep 1.69 if(rig>=0) t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1704 mocchiut 1.73 else t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(TMath::Pi()+dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1705 pam-mep 1.69 if (debug) cout << "R = " << rig << "\tcutoff = " << t_orb->cutoff << endl;
1706    
1707 mocchiut 1.33 //
1708 mocchiut 1.36 if ( t_orb->pitch != t_orb->pitch ) t_orb->pitch = -1000.;
1709     if ( t_orb->cutoff != t_orb->cutoff ) t_orb->cutoff = -1000.;
1710 pam-mep 1.69 if ( t_orb->sunangle != t_orb->sunangle ) t_orb->sunangle = -1000.;
1711     if ( t_orb->sunmagangle != t_orb->sunmagangle ) t_orb->sunmagangle = -1000.;
1712 mocchiut 1.36 //
1713 mocchiut 1.33 if ( debug ) printf(" orbitalinfo->cutoffsvl %f vitaly %f \n",orbitalinfo->cutoffsvl,t_orb->cutoff);
1714 mocchiut 1.30 //
1715     new(tor[nn]) OrbitalInfoTrkVar(*t_orb);
1716     nn++;
1717     //
1718     t_orb->Clear();
1719     //
1720 mocchiut 1.73 }
1721 mocchiut 1.30 //
1722 mocchiut 1.73 } // end standard tracking algorithm
1723    
1724     //
1725     // Code for extended tracking algorithm:
1726     //
1727     if ( hasNucleiTrk ){
1728     Int_t ttentry = 0;
1729     if ( verbose ) printf(" hasNucleiTrk \n");
1730     for(Int_t nt=0; nt < tcNucleiTof->GetEntries() ; nt++){
1731     //
1732     if ( verbose ) printf(" got1\n");
1733     ToFTrkVar *ptt = (ToFTrkVar*)(tcNucleiTof->At(nt));
1734     if (verbose) cout<<" tcNucleiTof->GetEntries() = "<<tcNucleiTof->GetEntries()<<"\tptt->trkseqno = "<<ptt->trkseqno<<endl;
1735     Double_t E11x = ptt->xtr_tof[0]; // tr->x[0];
1736     Double_t E11y = ptt->ytr_tof[0]; //tr->y[0];
1737     Double_t E11z = zin[0];
1738     Double_t E22x = ptt->xtr_tof[3];//tr->x[3];
1739     Double_t E22y = ptt->ytr_tof[3];//tr->y[3];
1740     Double_t E22z = zin[3];
1741     if ( (E11x < 100. && E11y < 100. && E22x < 100. && E22y < 100.) || ptt->trkseqno != -1 ){
1742     TrkTrack *mytrack = (TrkTrack*)(tcNucleiTrk->At(ptt->trkseqno));
1743     if ( verbose ) printf(" got tcNucleiTrk \n");
1744     Float_t rig=1/mytrack->GetDeflection();
1745     Double_t norm = sqrt(pow(E22x-E11x,2)+pow(E22y-E11y,2)+pow(E22z-E11z,2));
1746     //
1747     Px = (E22x-E11x)/norm;
1748     Py = (E22y-E11y)/norm;
1749     Pz = (E22z-E11z)/norm;
1750     //
1751     t_orb->trkseqno = ptt->trkseqno;
1752     //
1753     TMatrixD Eij = PO->PamelatoGEO(Iij,Px,Py,Pz);
1754     t_orb->Eij.ResizeTo(Eij);
1755     t_orb->Eij = Eij;
1756     //
1757     TMatrixD Sij = PO->PamelatoGEO(Gij,Px,Py,Pz);
1758     t_orb->Sij.ResizeTo(Sij);
1759     t_orb->Sij = Sij;
1760     //
1761     t_orb->pitch = (Float_t)PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),Bx,By,Bz);
1762     //
1763     // SunPosition in instrumental reference frame
1764     TMatrixD Kij = PO->PamelatoGEO(qij,Px,Py,Pz);
1765     TMatrixD Lij = PO->PamelatoGEO(qij,0,0,1);
1766     t_orb->sunangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),-SP.X(),-SP.Y(),-SP.Z());
1767     t_orb->sunmagangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),SunMag.X(),SunMag.Y(),SunMag.Z());
1768     //
1769     //
1770     Double_t omega = PO->GetPitchAngle(-Eij(0,0),-Eij(1,0),-Eij(2,0),1,0,0) * TMath::DegToRad();
1771     TVector3 Bxy(0,By,Bz);
1772     TVector3 Exy(0,-Eij(1,0),-Eij(2,0));
1773     Double_t dzeta=Bxy.Angle(Exy);
1774     if (-Eij(1,0) < 0) dzeta=2.0*TMath::Pi() - dzeta;
1775    
1776     if(debug) cout << "omega = "<<omega*TMath::RadToDeg()<<"\tdzeta = "<<dzeta*TMath::RadToDeg()<<endl;
1777    
1778     // Formula from D.F. Smart *, M.A. Shea [2005]; A review of geomagnetic cutoff rigidities for earth-orbiting spacecraft
1779     if(rig>=0) t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1780     else t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(TMath::Pi()+dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1781     if (debug) cout << "R = " << rig << "\tcutoff = " << t_orb->cutoff << endl;
1782    
1783     //
1784     if ( t_orb->pitch != t_orb->pitch ) t_orb->pitch = -1000.;
1785     if ( t_orb->cutoff != t_orb->cutoff ) t_orb->cutoff = -1000.;
1786     if ( t_orb->sunangle != t_orb->sunangle ) t_orb->sunangle = -1000.;
1787     if ( t_orb->sunmagangle != t_orb->sunmagangle ) t_orb->sunmagangle = -1000.;
1788     //
1789     if ( debug ) printf(" orbitalinfo->cutoffsvl %f vitaly %f \n",orbitalinfo->cutoffsvl,t_orb->cutoff);
1790     //
1791     TClonesArray &tt1 = *torbNucleiTrk;
1792     new(tt1[ttentry]) OrbitalInfoTrkVar(*t_orb);
1793     ttentry++;
1794     //
1795     t_orb->Clear();
1796     //
1797     }
1798     //
1799     }
1800     } // end standard tracking algorithm: nuclei
1801     if ( hasExtNucleiTrk ){
1802     Int_t ttentry = 0;
1803     if ( verbose ) printf(" hasExtNucleiTrk \n");
1804     for(Int_t nt=0; nt < tcExtNucleiTof->GetEntries() ; nt++){
1805     //
1806     if ( verbose ) printf(" got2\n");
1807     ToFTrkVar *ptt = (ToFTrkVar*)(tcExtNucleiTof->At(nt));
1808     if (verbose) cout<<" tcExtNucleiTof->GetEntries() = "<<tcExtNucleiTof->GetEntries()<<"\tptt->trkseqno = "<<ptt->trkseqno<<endl;
1809     Double_t E11x = ptt->xtr_tof[0]; // tr->x[0];
1810     Double_t E11y = ptt->ytr_tof[0]; //tr->y[0];
1811     Double_t E11z = zin[0];
1812     Double_t E22x = ptt->xtr_tof[3];//tr->x[3];
1813     Double_t E22y = ptt->ytr_tof[3];//tr->y[3];
1814     Double_t E22z = zin[3];
1815     if ( (E11x < 100. && E11y < 100. && E22x < 100. && E22y < 100.) || ptt->trkseqno != -1 ){
1816     ExtTrack *mytrack = (ExtTrack*)(tcExtNucleiTrk->At(ptt->trkseqno));
1817     if ( verbose ) printf(" got tcExtNucleiTrk \n");
1818     Float_t rig=1/mytrack->GetDeflection();
1819     Double_t norm = sqrt(pow(E22x-E11x,2)+pow(E22y-E11y,2)+pow(E22z-E11z,2));
1820     //
1821     Px = (E22x-E11x)/norm;
1822     Py = (E22y-E11y)/norm;
1823     Pz = (E22z-E11z)/norm;
1824     //
1825     t_orb->trkseqno = ptt->trkseqno;
1826     //
1827     TMatrixD Eij = PO->PamelatoGEO(Iij,Px,Py,Pz);
1828     t_orb->Eij.ResizeTo(Eij);
1829     t_orb->Eij = Eij;
1830     //
1831     TMatrixD Sij = PO->PamelatoGEO(Gij,Px,Py,Pz);
1832     t_orb->Sij.ResizeTo(Sij);
1833     t_orb->Sij = Sij;
1834     //
1835     t_orb->pitch = (Float_t)PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),Bx,By,Bz);
1836     //
1837     // SunPosition in instrumental reference frame
1838     TMatrixD Kij = PO->PamelatoGEO(qij,Px,Py,Pz);
1839     TMatrixD Lij = PO->PamelatoGEO(qij,0,0,1);
1840     t_orb->sunangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),-SP.X(),-SP.Y(),-SP.Z());
1841     t_orb->sunmagangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),SunMag.X(),SunMag.Y(),SunMag.Z());
1842     //
1843     //
1844     Double_t omega = PO->GetPitchAngle(-Eij(0,0),-Eij(1,0),-Eij(2,0),1,0,0) * TMath::DegToRad();
1845     TVector3 Bxy(0,By,Bz);
1846     TVector3 Exy(0,-Eij(1,0),-Eij(2,0));
1847     Double_t dzeta=Bxy.Angle(Exy);
1848     if (-Eij(1,0) < 0) dzeta=2.0*TMath::Pi() - dzeta;
1849    
1850     if(debug) cout << "omega = "<<omega*TMath::RadToDeg()<<"\tdzeta = "<<dzeta*TMath::RadToDeg()<<endl;
1851    
1852     // Formula from D.F. Smart *, M.A. Shea [2005]; A review of geomagnetic cutoff rigidities for earth-orbiting spacecraft
1853     if(rig>=0) t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1854     else t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(TMath::Pi()+dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1855     if (debug) cout << "R = " << rig << "\tcutoff = " << t_orb->cutoff << endl;
1856    
1857     //
1858     if ( t_orb->pitch != t_orb->pitch ) t_orb->pitch = -1000.;
1859     if ( t_orb->cutoff != t_orb->cutoff ) t_orb->cutoff = -1000.;
1860     if ( t_orb->sunangle != t_orb->sunangle ) t_orb->sunangle = -1000.;
1861     if ( t_orb->sunmagangle != t_orb->sunmagangle ) t_orb->sunmagangle = -1000.;
1862     //
1863     if ( debug ) printf(" orbitalinfo->cutoffsvl %f vitaly %f \n",orbitalinfo->cutoffsvl,t_orb->cutoff);
1864     //
1865     TClonesArray &tt2 = *torbExtNucleiTrk;
1866     new(tt2[ttentry]) OrbitalInfoTrkVar(*t_orb);
1867     ttentry++;
1868     //
1869     t_orb->Clear();
1870     //
1871     }
1872     //
1873     }
1874     } // end standard tracking algorithm: nuclei extended
1875     if ( hasExtTrk ){
1876     Int_t ttentry = 0;
1877     if ( verbose ) printf(" hasExtTrk \n");
1878     for(Int_t nt=0; nt < tcExtTof->GetEntries() ; nt++){
1879     //
1880     if ( verbose ) printf(" got3\n");
1881     ToFTrkVar *ptt = (ToFTrkVar*)(tcExtTof->At(nt));
1882     if (verbose) cout<<" tcExtTof->GetEntries() = "<<tcExtTof->GetEntries()<<"\tptt->trkseqno = "<<ptt->trkseqno<<endl;
1883     Double_t E11x = ptt->xtr_tof[0]; // tr->x[0];
1884     Double_t E11y = ptt->ytr_tof[0]; //tr->y[0];
1885     Double_t E11z = zin[0];
1886     Double_t E22x = ptt->xtr_tof[3];//tr->x[3];
1887     Double_t E22y = ptt->ytr_tof[3];//tr->y[3];
1888     Double_t E22z = zin[3];
1889     if ( (E11x < 100. && E11y < 100. && E22x < 100. && E22y < 100.) || ptt->trkseqno != -1 ){
1890     ExtTrack *mytrack = (ExtTrack*)(tcExtTrk->At(ptt->trkseqno));
1891     if ( verbose ) printf(" got tcExtTrk \n");
1892     Float_t rig=1/mytrack->GetDeflection();
1893     Double_t norm = sqrt(pow(E22x-E11x,2)+pow(E22y-E11y,2)+pow(E22z-E11z,2));
1894     //
1895     Px = (E22x-E11x)/norm;
1896     Py = (E22y-E11y)/norm;
1897     Pz = (E22z-E11z)/norm;
1898     //
1899     t_orb->trkseqno = ptt->trkseqno;
1900     //
1901     TMatrixD Eij = PO->PamelatoGEO(Iij,Px,Py,Pz);
1902     t_orb->Eij.ResizeTo(Eij);
1903     t_orb->Eij = Eij;
1904     //
1905     TMatrixD Sij = PO->PamelatoGEO(Gij,Px,Py,Pz);
1906     t_orb->Sij.ResizeTo(Sij);
1907     t_orb->Sij = Sij;
1908     //
1909     t_orb->pitch = (Float_t)PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),Bx,By,Bz);
1910     //
1911     // SunPosition in instrumental reference frame
1912     TMatrixD Kij = PO->PamelatoGEO(qij,Px,Py,Pz);
1913     TMatrixD Lij = PO->PamelatoGEO(qij,0,0,1);
1914     t_orb->sunangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),-SP.X(),-SP.Y(),-SP.Z());
1915     t_orb->sunmagangle=(Float_t)PO->GetPitchAngle(Kij(0,0),Kij(1,0),Kij(2,0),SunMag.X(),SunMag.Y(),SunMag.Z());
1916     //
1917     //
1918     Double_t omega = PO->GetPitchAngle(-Eij(0,0),-Eij(1,0),-Eij(2,0),1,0,0) * TMath::DegToRad();
1919     TVector3 Bxy(0,By,Bz);
1920     TVector3 Exy(0,-Eij(1,0),-Eij(2,0));
1921     Double_t dzeta=Bxy.Angle(Exy);
1922     if (-Eij(1,0) < 0) dzeta=2.0*TMath::Pi() - dzeta;
1923    
1924     if(debug) cout << "omega = "<<omega*TMath::RadToDeg()<<"\tdzeta = "<<dzeta*TMath::RadToDeg()<<endl;
1925    
1926     // Formula from D.F. Smart *, M.A. Shea [2005]; A review of geomagnetic cutoff rigidities for earth-orbiting spacecraft
1927     if(rig>=0) t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1928     else t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow(1+sqrt(1-sin(omega)*sin(TMath::Pi()+dzeta)*pow(cos(orbitalinfo->lat*TMath::DegToRad()),3)),2));
1929     if (debug) cout << "R = " << rig << "\tcutoff = " << t_orb->cutoff << endl;
1930    
1931     //
1932     if ( t_orb->pitch != t_orb->pitch ) t_orb->pitch = -1000.;
1933     if ( t_orb->cutoff != t_orb->cutoff ) t_orb->cutoff = -1000.;
1934     if ( t_orb->sunangle != t_orb->sunangle ) t_orb->sunangle = -1000.;
1935     if ( t_orb->sunmagangle != t_orb->sunmagangle ) t_orb->sunmagangle = -1000.;
1936     //
1937     if ( debug ) printf(" orbitalinfo->cutoffsvl %f vitaly %f \n",orbitalinfo->cutoffsvl,t_orb->cutoff);
1938     //
1939     TClonesArray &tt3 = *torbExtTrk;
1940     new(tt3[ttentry]) OrbitalInfoTrkVar(*t_orb);
1941     ttentry++;
1942     //
1943     t_orb->Clear();
1944     //
1945     }
1946     //
1947     }
1948     } // end standard tracking algorithm: extended
1949    
1950 mocchiut 1.30 } else {
1951 mocchiut 1.65 if ( debug ) printf(" mmm... mode %u standalone \n",orbitalinfo->mode);
1952     }
1953 mocchiut 1.30 //
1954 mocchiut 1.73 } else { // HERE IT MISS ALL CODE FOR EXTENDED TRACKING!
1955     if ( !standalone ){
1956 mocchiut 1.35 //
1957 mocchiut 1.73 if ( verbose ) printf(" no orb info for tracks \n");
1958 mocchiut 1.35 Int_t nn = 0;
1959     for(Int_t nt=0; nt < tof->ntrk(); nt++){
1960     //
1961     ToFTrkVar *ptt = tof->GetToFTrkVar(nt);
1962     if ( ptt->trkseqno != -1 ){
1963     //
1964     t_orb->trkseqno = ptt->trkseqno;
1965     //
1966     t_orb->Eij = 0;
1967     //
1968     t_orb->Sij = 0;
1969     //
1970     t_orb->pitch = -1000.;
1971     //
1972 pam-mep 1.69 t_orb->sunangle = -1000.;
1973     //
1974     t_orb->sunmagangle = -1000;
1975     //
1976 mocchiut 1.35 t_orb->cutoff = -1000.;
1977     //
1978     new(tor[nn]) OrbitalInfoTrkVar(*t_orb);
1979     nn++;
1980     //
1981     t_orb->Clear();
1982     //
1983 mocchiut 1.73 }
1984 mocchiut 1.35 //
1985 mocchiut 1.73 }
1986     //
1987     // Code for extended tracking algorithm:
1988     //
1989     if ( hasNucleiTrk ){
1990     Int_t ttentry = 0;
1991     if ( verbose ) printf(" hasNucleiTrk \n");
1992     for(Int_t nt=0; nt < tcNucleiTof->GetEntries() ; nt++){
1993     //
1994     ToFTrkVar *ptt = (ToFTrkVar*)(tcNucleiTof->At(nt));
1995     if ( ptt->trkseqno != -1 ){
1996     //
1997     t_orb->trkseqno = ptt->trkseqno;
1998     //
1999     t_orb->Eij = 0;
2000     //
2001     t_orb->Sij = 0;
2002     //
2003     t_orb->pitch = -1000.;
2004     //
2005     t_orb->sunangle = -1000.;
2006     //
2007     t_orb->sunmagangle = -1000;
2008     //
2009     t_orb->cutoff = -1000.;
2010     //
2011     TClonesArray &tz1 = *torbNucleiTrk;
2012     new(tz1[ttentry]) OrbitalInfoTrkVar(*t_orb);
2013     ttentry++;
2014     //
2015     t_orb->Clear();
2016     //
2017     }
2018     //
2019     }
2020     }
2021     if ( hasExtNucleiTrk ){
2022     Int_t ttentry = 0;
2023     if ( verbose ) printf(" hasExtNucleiTrk \n");
2024     for(Int_t nt=0; nt < tcExtNucleiTof->GetEntries() ; nt++){
2025     //
2026     if ( verbose ) printf(" got2\n");
2027     ToFTrkVar *ptt = (ToFTrkVar*)(tcExtNucleiTof->At(nt));
2028     if ( ptt->trkseqno != -1 ){
2029     //
2030     t_orb->trkseqno = ptt->trkseqno;
2031     //
2032     t_orb->Eij = 0;
2033     //
2034     t_orb->Sij = 0;
2035     //
2036     t_orb->pitch = -1000.;
2037     //
2038     t_orb->sunangle = -1000.;
2039     //
2040     t_orb->sunmagangle = -1000;
2041     //
2042     t_orb->cutoff = -1000.;
2043     //
2044     TClonesArray &tz2 = *torbExtNucleiTrk;
2045     new(tz2[ttentry]) OrbitalInfoTrkVar(*t_orb);
2046     ttentry++;
2047     //
2048     t_orb->Clear();
2049     //
2050     }
2051     //
2052     }
2053     }
2054     if ( hasExtTrk ){
2055     Int_t ttentry = 0;
2056     if ( verbose ) printf(" hasExtTrk \n");
2057     for(Int_t nt=0; nt < tcExtTof->GetEntries() ; nt++){
2058     //
2059     if ( verbose ) printf(" got3\n");
2060     ToFTrkVar *ptt = (ToFTrkVar*)(tcExtTof->At(nt));
2061     if ( ptt->trkseqno != -1 ){
2062     //
2063     t_orb->trkseqno = ptt->trkseqno;
2064     //
2065     t_orb->Eij = 0;
2066     //
2067     t_orb->Sij = 0;
2068     //
2069     t_orb->pitch = -1000.;
2070     //
2071     t_orb->sunangle = -1000.;
2072     //
2073     t_orb->sunmagangle = -1000;
2074     //
2075     t_orb->cutoff = -1000.;
2076     //
2077     TClonesArray &tz3 = *torbExtNucleiTrk;
2078     new(tz3[ttentry]) OrbitalInfoTrkVar(*t_orb);
2079     ttentry++;
2080     //
2081     t_orb->Clear();
2082     //
2083     }
2084     //
2085     }
2086     }
2087     }
2088     } // if( orbitalinfo->TimeGap>0){
2089 mocchiut 1.30 //
2090 mocchiut 1.15 // Fill the class
2091 pamelaprod 1.11 //
2092 mocchiut 1.1 OrbitalInfotr->Fill();
2093 mocchiut 1.15 //
2094 mocchiut 1.30 delete t_orb;
2095     //
2096 mocchiut 1.15 }; // loop over the events in the run
2097 mocchiut 1.1 //
2098     // Here you may want to clear some variables before processing another run
2099     //
2100 mocchiut 1.44
2101 mocchiut 1.57 if ( verbose ) printf(" Clear before new run \n");
2102 mocchiut 1.5 delete dbtime;
2103 mocchiut 1.57
2104 mocchiut 1.62 if ( mcmdrc ) mcmdrc->Clear();
2105 mocchiut 1.57 mcmdrc = 0;
2106    
2107     if ( verbose ) printf(" Clear before new run1 \n");
2108     if ( L_QQ_Q_l_lower ) delete L_QQ_Q_l_lower;
2109     if ( verbose ) printf(" Clear before new run2 \n");
2110 mocchiut 1.18 if ( L_QQ_Q_l_upper ) delete L_QQ_Q_l_upper;
2111 mocchiut 1.57 if ( verbose ) printf(" Clear before new run3 \n");
2112 mocchiut 1.18 if ( RYPang_upper ) delete RYPang_upper;
2113 mocchiut 1.57 if ( verbose ) printf(" Clear before new run4 \n");
2114 mocchiut 1.18 if ( RYPang_lower ) delete RYPang_lower;
2115 mocchiut 1.57
2116     if ( l0tr ) l0tr->Delete();
2117    
2118     if ( verbose ) printf(" End run \n");
2119    
2120 mocchiut 1.1 }; // process all the runs
2121 mocchiut 1.15
2122 mocchiut 1.1 if (verbose) printf("\n Finished processing data \n");
2123     //
2124     closeandexit:
2125     //
2126     // 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.
2127     //
2128     if ( !reprocall && reproc && code >= 0 ){
2129     if ( totfileentries > noaftrun ){
2130     if (verbose){
2131     printf("\n Post-processing: copying events from the old tree after the processed run\n");
2132     printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
2133     printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
2134     }
2135     for (UInt_t j = noaftrun; j < totfileentries; j++ ){
2136     //
2137     // Get entry from old tree
2138     //
2139 mocchiut 1.43 if ( OrbitalInfotrclone->GetEntry(j) <= 0 ) throw -36;
2140 mocchiut 1.1 //
2141     // copy orbitalinfoclone to OrbitalInfo
2142     //
2143 mocchiut 1.3 orbitalinfo->Clear();
2144 mocchiut 1.5 //
2145 mocchiut 1.1 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
2146     //
2147     // Fill entry in the new tree
2148     //
2149     OrbitalInfotr->Fill();
2150     };
2151     if (verbose) printf(" Finished successful copying!\n");
2152     };
2153 mocchiut 1.57 //if ( OrbitalInfotrclone ) OrbitalInfotrclone->Clear();
2154     //if ( OrbitalInfotrclone ) OrbitalInfotrclone->Delete();
2155 mocchiut 1.1 };
2156     //
2157     // Close files, delete old tree(s), write and close level2 file
2158     //
2159 pam-mep 1.69
2160 mocchiut 1.1 if ( l0File ) l0File->Close();
2161 mocchiut 1.23 if ( myfold ) gSystem->Unlink(tempname.str().c_str());
2162 mocchiut 1.5 //
2163 mocchiut 1.1 if ( OrbitalInfotr ) OrbitalInfotr->SetName("OrbitalInfo");
2164 mocchiut 1.30 //
2165 mocchiut 1.1 if ( file ){
2166     file->cd();
2167 mocchiut 1.63 if ( OrbitalInfotr ) OrbitalInfotr->Write("OrbitalInfo", TObject::kOverwrite); // 10 RED bug fixed
2168 mocchiut 1.1 };
2169     //
2170 mocchiut 1.57 if (verbose) printf("\n Exiting...\n");
2171    
2172 mocchiut 1.23 if ( myfold ) gSystem->Unlink(OrbitalInfofolder.str().c_str());
2173 mocchiut 1.1 //
2174     // the end
2175     //
2176 mocchiut 1.25 if ( dbc ){
2177     dbc->Close();
2178     delete dbc;
2179     };
2180 mocchiut 1.57 //
2181 mocchiut 1.1 if (verbose) printf("\n Exiting...\n");
2182 mocchiut 1.57 if ( tempfile ) tempfile->Close();
2183    
2184 mocchiut 1.30 if ( PO ) delete PO;
2185 mocchiut 1.57 if ( gltle ) delete gltle;
2186     if ( glparam ) delete glparam;
2187     if ( glparam2 ) delete glparam2;
2188     if ( glparam3 ) delete glparam3;
2189     if (verbose) printf("\n Exiting3...\n");
2190 mocchiut 1.4 if ( glroot ) delete glroot;
2191 mocchiut 1.57 if (verbose) printf("\n Exiting4...\n");
2192     if ( runinfo ) runinfo->Close();
2193 mocchiut 1.4 if ( runinfo ) delete runinfo;
2194 mocchiut 1.58
2195 mocchiut 1.73 if ( tcNucleiTrk ){
2196     tcNucleiTrk->Delete();
2197     delete tcNucleiTrk;
2198     tcNucleiTrk = NULL;
2199     }
2200     if ( tcExtNucleiTrk ){
2201     tcExtNucleiTrk->Delete();
2202     delete tcExtNucleiTrk;
2203     tcExtNucleiTrk = NULL;
2204     }
2205     if ( tcExtTrk ){
2206     tcExtTrk->Delete();
2207     delete tcExtTrk;
2208     tcExtTrk = NULL;
2209     }
2210    
2211     if ( tcNucleiTof ){
2212     tcNucleiTof->Delete();
2213     delete tcNucleiTof;
2214     tcNucleiTof = NULL;
2215     }
2216     if ( tcExtNucleiTof ){
2217     tcExtNucleiTof->Delete();
2218     delete tcExtNucleiTof;
2219     tcExtNucleiTof = NULL;
2220     }
2221     if ( tcExtTof ){
2222     tcExtTof->Delete();
2223     delete tcExtTof;
2224     tcExtTrk = NULL;
2225     }
2226    
2227    
2228 pam-mep 1.69 if ( tof ) delete tof;
2229     if ( trke ) delete trke;
2230    
2231 mocchiut 1.58 if ( debug ){
2232 mocchiut 1.57 cout << "1 0x" << OrbitalInfotr << endl;
2233     cout << "2 0x" << OrbitalInfotrclone << endl;
2234     cout << "3 0x" << l0tr << endl;
2235     cout << "4 0x" << tempOrbitalInfo << endl;
2236     cout << "5 0x" << ttof << endl;
2237 mocchiut 1.58 }
2238 mocchiut 1.57 //
2239 mocchiut 1.58 if ( debug ) file->ls();
2240 mocchiut 1.4 //
2241 mocchiut 1.1 if(code < 0) throw code;
2242     return(code);
2243     }
2244    
2245 mocchiut 1.7
2246     //
2247     // Returns the cCoordGeo structure holding the geographical
2248     // coordinates for the event (see sgp4.h).
2249     //
2250     // atime is the abstime of the event in UTC unix time.
2251     // tletime is the time of the tle in UTC unix time.
2252     // tle is the previous and nearest tle (compared to atime).
2253     cCoordGeo getCoo(UInt_t atime, UInt_t tletime, cTle *tle)
2254     {
2255     cEci eci;
2256     cOrbit orbit(*tle);
2257     orbit.getPosition((double) (atime - tletime)/60., &eci);
2258    
2259     return eci.toGeo();
2260     }
2261 mocchiut 1.15
2262     // function of copyng of quatrnions classes
2263    
2264     void CopyQ(Quaternions *Q1, Quaternions *Q2){
2265     for(UInt_t i = 0; i < 6; i++){
2266     Q1->time[i]=Q2->time[i];
2267     for (UInt_t j = 0; j < 4; j++)Q1->quat[i][j]=Q2->quat[i][j];
2268     }
2269     return;
2270     }
2271    
2272     // functions of copyng InclinationInfo classes
2273    
2274     void CopyAng(InclinationInfo *A1, InclinationInfo *A2){
2275     A1->Tangazh = A2->Tangazh;
2276     A1->Ryskanie = A2->Ryskanie;
2277     A1->Kren = A2->Kren;
2278     return;
2279     }
2280    
2281     UInt_t holeq(Double_t lower,Double_t upper,Quaternions *Qlower, Quaternions *Qupper, UInt_t f){
2282    
2283     UInt_t hole = 10;
2284 mocchiut 1.40 Bool_t R10l = false; // Sign of R10 mode in lower quaternions array
2285     Bool_t R10u = false; // Sign of R10 mode in upper quaternions array
2286     Bool_t insm = false; // Sign that we inside quaternions array
2287 mocchiut 1.64 // Bool_t mxtml = false; // Sign of mixt mode in lower quaternions array
2288     // Bool_t mxtmu = false; // Sign of mixt mode in upper quaternions array
2289 mocchiut 1.40 Bool_t npasm = false; // Sign of normall pass between R10 and non R10 or between non R10 and R10
2290 mocchiut 1.15 UInt_t NCQl = 6; // Number of correct quaternions in lower array
2291 mocchiut 1.64 // UInt_t NCQu = 6; // Number of correct quaternions in upper array
2292 mocchiut 1.15 if (f>0){
2293     insm = true;
2294     if(Qupper->time[f]-Qupper->time[f-1]==30) R10u = false;
2295     if(Qupper->time[f]-Qupper->time[f-1]<1) R10u = true;
2296     }else{
2297     insm = false;
2298     if((Qlower->time[5]-Qlower->time[0]<2)&&(Qlower->time[1]-Qlower->time[0]<2)) R10l = true;
2299     if((Qupper->time[5]-Qupper->time[0]<2)&&(Qupper->time[1]-Qupper->time[0]<2)) R10u = true;
2300     if((Qlower->time[5]-Qlower->time[0]==150)&&(Qlower->time[1]-Qlower->time[0]==30)) R10l = false;
2301     if((Qupper->time[5]-Qupper->time[0]==150)&&(Qupper->time[1]-Qupper->time[0]==30)) R10u = false;
2302     if((Qlower->time[5]-Qlower->time[0]<2)&&(Qlower->time[1]-Qlower->time[0]==30)){
2303 mocchiut 1.64 // mxtml = true;
2304 mocchiut 1.15 for(UInt_t i = 1; i < 6; i++){
2305     if(Qlower->time[i]-Qlower->time[0]==30*i) NCQl=i;
2306     }
2307     }
2308 mocchiut 1.64 // if((Qupper->time[5]-Qupper->time[0]<2)&&(Qupper->time[1]-Qupper->time[0]==30)){
2309     // mxtmu = true;
2310     // for(UInt_t i = 1; i < 6; i++){
2311     // if(Qupper->time[i]-Qupper->time[0]==30*i) NCQu=i;
2312     // }
2313     // }
2314 mocchiut 1.15 }
2315    
2316     if(((upper-lower==1.5)||(upper-lower==3.)||(upper-lower==30.)||(upper-lower==31.5)||(upper-lower==33.)||(upper-lower==181.5)||(upper-lower==210.)||(upper-lower==211.5))&&!insm) npasm = true;
2317    
2318    
2319     if (R10u&&insm) hole=0; // best event R10
2320     if ((upper-lower<=5)&&(!insm)&&R10l&&R10u) hole = 1; // when first of 6 quaternions in array is correct
2321     if (((!R10u)&&insm)||((!insm)&&(!R10u)&&(!R10l)&&((upper-lower==210+(6-NCQl)*30)||(upper-lower==30)))) hole = 2; //non R10
2322     if (npasm&&(!insm)&&((R10l&&!R10u)||(R10u&&!R10l))) hole = 3; //normall pass from R10 to non R10 or from non R10 to R10
2323     if ((!npasm)&&(upper-lower<=300)&&(!insm)&&((R10l&&!R10u)||(R10u&&!R10l))) hole = 4; // eliminable hole between R10 and non R10 or between non R10 and R10
2324     if ((upper-lower>=300)&&(!insm)&&((R10l&&!R10u)||(R10u&&!R10l))) hole = 5; //uneliminable hole between R10 and non R10 or between non R10 and R10
2325     if ((upper-lower>5)&&(upper-lower<=300)&&R10u&&R10l) hole = 6; // eliminable hole inside R10
2326     if ((upper-lower>300)&&R10u&&R10l) hole = 7; //uneliminable hole inside R10
2327     if ((upper-lower>210)&&(upper-lower<=1200)&&(!R10u)&&(!R10l)) hole = 8; //eliminable hole inside non R10
2328     if ((upper-lower>1200)&&!R10u&&!R10l) hole = 9; // uneliminable hole inside non R10
2329     return hole;
2330     }
2331    
2332 mocchiut 1.44 void inclresize(vector<Double_t>& t,vector<Float_t>& q0,vector<Float_t>& q1,vector<Float_t>& q2,vector<Float_t>& q3,vector<Int_t>& mode,vector<Float_t>& Roll,vector<Float_t>& Pitch,vector<Float_t>& Yaw){
2333     Int_t sizee = t.size()+1;
2334     t.resize(sizee);
2335     q0.resize(sizee);
2336     q1.resize(sizee);
2337     q2.resize(sizee);
2338     q3.resize(sizee);
2339     mode.resize(sizee);
2340     Roll.resize(sizee);
2341     Pitch.resize(sizee);
2342     Yaw.resize(sizee);
2343     }
2344    
2345 pam-mep 1.69 // geomagnetic calculation staff
2346    
2347     //void GM_ScanIGRF(TString PATH, GMtype_Data *G0, GMtype_Data *G1, GMtype_Data *H1)
2348     void GM_ScanIGRF(TSQLServer *dbc, GMtype_Data *G0, GMtype_Data *G1, GMtype_Data *H1)
2349     {
2350     GL_PARAM *glp = new GL_PARAM();
2351     Int_t parerror=glp->Query_GL_PARAM(1,304,dbc); // parameters stored in DB in GL_PRAM table
2352     if ( parerror<0 ) {
2353     throw -902;
2354     }
2355     /*This function scans inputs G0, G1, and H1 of the IGRF table into 3 data arrays*/
2356     // TString SATH="/data03/Malakhov/pam9Malakhov/installed10/calib/orb-param/";
2357     int i;
2358     double temp;
2359     char buffer[200];
2360     FILE *IGRF;
2361     IGRF = fopen((glp->PATH+glp->NAME).Data(), "r");
2362     // IGRF = fopen(PATH+"IGRF.tab", "r");
2363     G0->size = 25;
2364     G1->size = 25;
2365     H1->size = 25;
2366     for( i = 0; i < 4; i++)
2367     {
2368     fgets(buffer, 200, IGRF);
2369 mocchiut 1.44 }
2370 pam-mep 1.69 fscanf(IGRF, "g 1 0 %lf ", &G0->element[0]);
2371     for(i = 1; i <= 22; i++)
2372     {
2373     fscanf(IGRF ,"%lf ", &G0->element[i]);
2374 mocchiut 1.44 }
2375 pam-mep 1.69 fscanf(IGRF ,"%lf\n", &temp);
2376     G0->element[23] = temp * 5 + G0->element[22];
2377     G0->element[24] = G0->element[23] + 5 * temp;
2378     fscanf(IGRF, "g 1 1 %lf ", &G1->element[0]);
2379     for(i = 1; i <= 22; i++)
2380     {
2381     fscanf( IGRF, "%lf ", &G1->element[i]);
2382 mocchiut 1.44 }
2383 pam-mep 1.69 fscanf(IGRF, "%lf\n", &temp);
2384     G1->element[23] = temp * 5 + G1->element[22];
2385     G1->element[24] = temp * 5 + G1->element[23];
2386     fscanf(IGRF, "h 1 1 %lf ", &H1->element[0]);
2387     for(i = 1; i <= 22; i++)
2388     {
2389     fscanf( IGRF, "%lf ", &H1->element[i]);
2390 mocchiut 1.44 }
2391 pam-mep 1.69 fscanf(IGRF, "%lf\n", &temp);
2392     H1->element[23] = temp * 5 + H1->element[22];
2393     H1->element[24] = temp * 5 + H1->element[23];
2394     if ( glp ) delete glp;
2395     } /*GM_ScanIGRF*/
2396    
2397     void GM_SetEllipsoid(GMtype_Ellipsoid *Ellip)
2398     {
2399     /*This function sets the WGS84 reference ellipsoid to its default values*/
2400     Ellip->a = 6378.137; /*semi-major axis of the ellipsoid in */
2401     Ellip->b = 6356.7523142;/*semi-minor axis of the ellipsoid in */
2402     Ellip->fla = 1/298.257223563;/* flattening */
2403     Ellip->eps = sqrt(1- ( Ellip->b * Ellip->b) / (Ellip->a * Ellip->a )); /*first eccentricity */
2404     Ellip->epssq = (Ellip->eps * Ellip->eps); /*first eccentricity squared */
2405     Ellip->re = 6371.2;/* Earth's radius */
2406     } /*GM_SetEllipsoid*/
2407    
2408    
2409     void GM_EarthCartToDipoleCartCD(GMtype_Pole Pole, GMtype_CoordCartesian EarthCoord, GMtype_CoordCartesian *DipoleCoords)
2410     {
2411     /*This function converts from Earth centered cartesian coordinates to dipole centered cartesian coordinates*/
2412     double X, Y, Z, CosPhi, SinPhi, CosLambda, SinLambda;
2413     CosPhi = cos(TMath::DegToRad()*Pole.phi);
2414     SinPhi = sin(TMath::DegToRad()*Pole.phi);
2415     CosLambda = cos(TMath::DegToRad()*Pole.lambda);
2416     SinLambda = sin(TMath::DegToRad()*Pole.lambda);
2417     X = EarthCoord.x;
2418     Y = EarthCoord.y;
2419     Z = EarthCoord.z;
2420    
2421     /*These equations are taken from a document by Wallace H. Campbell*/
2422     DipoleCoords->x = X * CosPhi * CosLambda + Y * CosPhi * SinLambda - Z * SinPhi;
2423     DipoleCoords->y = -X * SinLambda + Y * CosLambda;
2424     DipoleCoords->z = X * SinPhi * CosLambda + Y * SinPhi * SinLambda + Z * CosPhi;
2425     } /*GM_EarthCartToDipoleCartCD*/
2426    
2427     void GM_GeodeticToSpherical(GMtype_Ellipsoid Ellip, GMtype_CoordGeodetic CoordGeodetic, GMtype_CoordSpherical *CoordSpherical)
2428     {
2429     double CosLat, SinLat, rc, xp, zp; /*all local variables */
2430     /*
2431     ** Convert geodetic coordinates, (defined by the WGS-84
2432     ** reference ellipsoid), to Earth Centered Earth Fixed Cartesian
2433     ** coordinates, and then to spherical coordinates.
2434     */
2435    
2436     CosLat = cos(TMath::DegToRad()*CoordGeodetic.phi);
2437     SinLat = sin(TMath::DegToRad()*CoordGeodetic.phi);
2438    
2439     /* compute the local radius of curvature on the WGS-84 reference ellipsoid */
2440    
2441     rc = Ellip.a / sqrt(1.0 - Ellip.epssq * SinLat * SinLat);
2442    
2443     /* compute ECEF Cartesian coordinates of specified point (for longitude=0) */
2444    
2445     xp = (rc + CoordGeodetic.HeightAboveEllipsoid) * CosLat;
2446     zp = (rc*(1.0 - Ellip.epssq) + CoordGeodetic.HeightAboveEllipsoid) * SinLat;
2447    
2448     /* compute spherical radius and angle lambda and phi of specified point */
2449    
2450     CoordSpherical->r = sqrt(xp * xp + zp * zp);
2451     CoordSpherical->phig = TMath::RadToDeg()*asin(zp / CoordSpherical->r); /* geocentric latitude */
2452     CoordSpherical->lambda = CoordGeodetic.lambda; /* longitude */
2453     } /*GM_GeodeticToSpherical*/
2454    
2455     void GM_PoleLocation(GMtype_Model Model, GMtype_Pole *Pole)
2456     {
2457     /*This function finds the location of the north magnetic pole in spherical coordinates. The equations are
2458     **from Wallace H. Campbell's Introduction to Geomagnetic Fields*/
2459    
2460     Pole->phi = TMath::RadToDeg()*-atan(sqrt(Model.h1 * Model.h1 + Model.g1 * Model.g1)/Model.g0);
2461     Pole->lambda = TMath::RadToDeg()*atan(Model.h1/Model.g1);
2462     } /*GM_PoleLocation*/
2463    
2464     void GM_SphericalToCartesian(GMtype_CoordSpherical CoordSpherical, GMtype_CoordCartesian *CoordCartesian)
2465     {
2466     /*This function converts spherical coordinates into Cartesian coordinates*/
2467     double CosPhi = cos(TMath::DegToRad()*CoordSpherical.phig);
2468     double SinPhi = sin(TMath::DegToRad()*CoordSpherical.phig);
2469     double CosLambda = cos(TMath::DegToRad()*CoordSpherical.lambda);
2470     double SinLambda = sin(TMath::DegToRad()*CoordSpherical.lambda);
2471    
2472     CoordCartesian->x = CoordSpherical.r * CosPhi * CosLambda;
2473     CoordCartesian->y = CoordSpherical.r * CosPhi * SinLambda;
2474     CoordCartesian->z = CoordSpherical.r * SinPhi;
2475     } /*GM_SphericalToCartesian*/
2476    
2477     void GM_TimeAdjustCoefs(Float_t year, Float_t jyear, GMtype_Data g0d, GMtype_Data g1d, GMtype_Data h1d, GMtype_Model *Model)
2478     {
2479     /*This function calls GM_LinearInterpolation for the coefficients to estimate the value of the
2480     **coefficient for the given date*/
2481     int index;
2482     double x;
2483     index = (year - GM_STARTYEAR) / 5;
2484     x = (jyear - GM_STARTYEAR) / 5;
2485     Model->g0 = GM_LinearInterpolation(index, index+1, g0d.element[index], g0d.element[index+1], x);
2486     Model->g1 = GM_LinearInterpolation(index, index+1, g1d.element[index], g1d.element[index+1], x);
2487     Model->h1 = GM_LinearInterpolation(index, index+1, h1d.element[index], h1d.element[index+1], x);
2488     } /*GM_TimeAdjustCoefs*/
2489    
2490     double GM_LinearInterpolation(double x1, double x2, double y1, double y2, double x)
2491     {
2492     /*This function takes a linear interpolation between two given points for x*/
2493     double weight, y;
2494     weight = (x - x1) / (x2 - x1);
2495     y = y1 * (1 - weight) + y2 * weight;
2496     return y;
2497     }/*GM_LinearInterpolation*/
2498 mocchiut 1.44
2499 pam-mep 1.69 void GM_CartesianToSpherical(GMtype_CoordCartesian CoordCartesian, GMtype_CoordSpherical *CoordSpherical)
2500     {
2501     /*This function converts a point from Cartesian coordinates into spherical coordinates*/
2502     double X, Y, Z;
2503    
2504     X = CoordCartesian.x;
2505     Y = CoordCartesian.y;
2506     Z = CoordCartesian.z;
2507    
2508     CoordSpherical->r = sqrt(X * X + Y * Y + Z * Z);
2509     CoordSpherical->phig = TMath::RadToDeg()*asin(Z / (CoordSpherical->r));
2510     CoordSpherical->lambda = TMath::RadToDeg()*atan2(Y, X);
2511     } /*GM_CartesianToSpherical*/

  ViewVC Help
Powered by ViewVC 1.1.23