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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.56 - (show annotations) (download)
Thu Apr 12 12:27:02 2012 UTC (13 years, 8 months ago) by mocchiut
Branch: MAIN
Changes since 1.55: +1 -0 lines
New TrigLevel2, ToFLevel2, OrbitalInfo classes, new order in DV processing (trigger first), ToF bugs fixed, new tof calibration, new variables in ToF and OrbitalInfo

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 //#include <TCanvas.h>
12 //#include <TH2F.h> //for test only. Vitaly.
13 //#include <TF1.h>
14
15 #include <TTree.h>
16 #include <TClassEdit.h>
17 #include <TObject.h>
18 #include <TList.h>
19 #include <TArrayI.h>
20 #include <TSystem.h>
21 #include <TSystemDirectory.h>
22 #include <TString.h>
23 #include <TFile.h>
24 #include <TClass.h>
25 #include <TSQLServer.h>
26 #include <TSQLRow.h>
27 #include <TSQLResult.h>
28 //
29 // RunInfo header
30 //
31 #include <RunInfo.h>
32 #include <GLTables.h>
33 //
34 // YODA headers
35 //
36 #include <PamelaRun.h>
37 #include <PscuHeader.h>
38 #include <PscuEvent.h>
39 #include <EventHeader.h>
40 #include <mcmd/McmdEvent.h>
41 #include <mcmd/McmdRecord.h>
42 //
43 // This program headers
44 //
45 #include <OrbitalInfo.h>
46 #include <OrbitalInfoVerl2.h>
47 #include <OrbitalInfoCore.h>
48 #include <InclinationInfo.h>
49
50
51 using namespace std;
52
53 //
54 // CORE ROUTINE
55 //
56 //
57 int OrbitalInfoCore(UInt_t run, TFile *file, GL_TABLES *glt, Int_t OrbitalInfoargc, char *OrbitalInfoargv[]){
58 //
59 Int_t i = 0;
60 TString host = glt->CGetHost();
61 TString user = glt->CGetUser();
62 TString psw = glt->CGetPsw();
63 TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
64 //
65 stringstream myquery;
66 myquery.str("");
67 myquery << "SET time_zone='+0:00'";
68 dbc->Query(myquery.str().c_str());
69 //
70 TString processFolder = Form("OrbitalInfoFolder_%u",run);
71 //
72 // Set these to true to have a very verbose output.
73 //
74 Bool_t debug = false;
75 //
76 Bool_t verbose = false;
77 //
78 Bool_t standalone = false;
79 //
80 if ( OrbitalInfoargc > 0 ){
81 i = 0;
82 while ( i < OrbitalInfoargc ){
83 if ( !strcmp(OrbitalInfoargv[i],"-processFolder") ) {
84 if ( OrbitalInfoargc < i+1 ){
85 throw -3;
86 };
87 processFolder = (TString)OrbitalInfoargv[i+1];
88 i++;
89 };
90 if ( (!strcmp(OrbitalInfoargv[i],"--debug")) || (!strcmp(OrbitalInfoargv[i],"-g")) ) {
91 verbose = true;
92 debug = true;
93 };
94 if ( (!strcmp(OrbitalInfoargv[i],"--verbose")) || (!strcmp(OrbitalInfoargv[i],"-v")) ) {
95 verbose = true;
96 };
97 if ( (!strcmp(OrbitalInfoargv[i],"--standalone")) ) {
98 standalone = true;
99 };
100 if ( (!strcmp(OrbitalInfoargv[i],"--calculate-pitch")) ) {
101 standalone = false;
102 };
103 i++;
104 };
105 };
106 //
107 const char* outDir = gSystem->DirName(gSystem->DirName(file->GetPath()));
108 //
109 TTree *OrbitalInfotr = 0;
110 UInt_t nevents = 0;
111 UInt_t neventsm = 0;
112 //
113 // variables needed to reprocess data
114 //
115 Long64_t maxsize = 10000000000LL;
116 TTree::SetMaxTreeSize(maxsize);
117 //
118 TString OrbitalInfoversion;
119 ItoRunInfo *runinfo = 0;
120 TArrayI *runlist = 0;
121 TTree *OrbitalInfotrclone = 0;
122 Bool_t reproc = false;
123 Bool_t reprocall = false;
124 UInt_t nobefrun = 0;
125 UInt_t noaftrun = 0;
126 UInt_t numbofrun = 0;
127 stringstream ftmpname;
128 TString fname;
129 UInt_t totfileentries = 0;
130 UInt_t idRun = 0;
131 UInt_t anni5 = 60 * 60 * 24 * 365 * 5 ;//1576800
132 //
133 // My variables. Vitaly.
134 //
135 // UInt_t oi = 0;
136 Int_t tmpSize = 0;
137 //
138 // variables needed to handle error signals
139 //
140 Int_t code = 0;
141 Int_t sgnl;
142 //
143 // OrbitalInfo classes
144 //
145 OrbitalInfo *orbitalinfo = new OrbitalInfo();
146 OrbitalInfo *orbitalinfoclone = new OrbitalInfo();
147
148 //
149 // define variables for opening and reading level0 file
150 //
151 TFile *l0File = 0;
152 TTree *l0tr = 0;
153 // TTree *l0trm = 0;
154 TChain *ch = 0;
155 // EM: open also header branch
156 TBranch *l0head = 0;
157 pamela::EventHeader *eh = 0;
158 pamela::PscuHeader *ph = 0;
159 pamela::McmdEvent *mcmdev = 0;
160 pamela::McmdRecord *mcmdrc = 0;
161 // end EM
162
163 // pamela::RunHeaderEvent *reh = new pamela::RunHeaderEvent;
164 // pamela::EventHeader *eH = new pamela::EventHeader;
165
166 //
167 // Define other basic variables
168 //
169 UInt_t procev = 0;
170 stringstream file2;
171 stringstream file3;
172 stringstream qy;
173 Int_t totevent = 0;
174 UInt_t atime = 0;
175 UInt_t re = 0;
176 UInt_t ik = 0;
177
178 // Position
179 Float_t lon, lat, alt;
180
181 //
182 // IGRF stuff
183 //
184 Float_t dimo = 0.0; // dipole moment (computed from dat files)
185 Float_t bnorth, beast, bdown, babs;
186 Float_t xl; // L value
187 Float_t icode; // code value for L accuracy (see fortran code)
188 Float_t bab1; // What's the difference with babs?
189 Float_t stps = 0.005; // step size for field line tracing
190 Float_t bdel = 0.01; // required accuracy
191 Float_t bequ; // equatorial b value (also called b_0)
192 Bool_t value = 0; // false if bequ is not the minimum b value
193 Float_t rr0; // equatorial radius normalized to earth radius
194
195 //
196 // Working filename
197 //
198 TString outputfile;
199 stringstream name;
200 name.str("");
201 name << outDir << "/";
202 //
203 // temporary file and folder
204 //
205 TFile *tempfile = 0;
206 TTree *tempOrbitalInfo = 0;
207 stringstream tempname;
208 stringstream OrbitalInfofolder;
209 Bool_t myfold = false;
210 tempname.str("");
211 tempname << outDir;
212 tempname << "/" << processFolder.Data();
213 OrbitalInfofolder.str("");
214 OrbitalInfofolder << tempname.str().c_str();
215 tempname << "/OrbitalInfotree_run";
216 tempname << run << ".root";
217 UInt_t totnorun = 0;
218 //
219 // DB classes
220 //
221 GL_ROOT *glroot = new GL_ROOT();
222 GL_TIMESYNC *dbtime = 0;
223 GL_TLE *gltle = new GL_TLE();
224 //
225 //Quaternions classes
226 //
227 Quaternions *L_QQ_Q_l_lower = new Quaternions();
228 InclinationInfo *RYPang_lower = new InclinationInfo();
229 Quaternions *L_QQ_Q_l_upper = new Quaternions();
230 InclinationInfo *RYPang_upper = new InclinationInfo();
231
232 cEci eCi;
233
234 // Initialize fortran routines!!!
235 Int_t ltp1 = 0;
236 Int_t ltp2 = 0;
237 Int_t ltp3 = 0;
238 // Int_t uno = 1;
239 // const char *niente = " ";
240 GL_PARAM *glparam = new GL_PARAM();
241 GL_PARAM *glparam2 = new GL_PARAM();
242 GL_PARAM *glparam3 = new GL_PARAM();
243
244 //
245 // Orientation variables. Vitaly
246 //
247 UInt_t evfrom = 0;
248 UInt_t jumped = 0;
249 Int_t itr = -1;
250 Double_t A1;
251 Double_t A2;
252 Double_t A3;
253 Double_t Px = 0;
254 Double_t Py = 0;
255 Double_t Pz = 0;
256 TTree *ttof = 0;
257 ToFLevel2 *tof = new ToFLevel2();
258 OrientationInfo *PO = new OrientationInfo();
259 Int_t nz = 6;
260 Float_t zin[6];
261 Int_t nevtofl2 = 0;
262 if ( verbose ) cout<<"Reading quaternions external file"<<endl;
263 cout.setf(ios::fixed,ios::floatfield);
264 /******Reading recovered quaternions...*********/
265 vector<Double_t> recqtime;
266 vector<Float_t> recq0;
267 vector<Float_t> recq1;
268 vector<Float_t> recq2;
269 vector<Float_t> recq3;
270 Float_t Norm = 1;
271 Int_t parerror=glparam->Query_GL_PARAM(1,303,dbc); // parameters stored in DB in GL_PRAM table
272 ifstream in((glparam->PATH+glparam->NAME).Data(),ios::in);
273 if ( parerror<0 ) {
274 code = parerror;
275 goto closeandexit;
276 }
277 while(!in.eof()){
278 recqtime.resize(recqtime.size()+1);
279 Int_t sizee = recqtime.size();
280 recq0.resize(sizee);
281 recq1.resize(sizee);
282 recq2.resize(sizee);
283 recq3.resize(sizee);
284 in>>recqtime[sizee-1];
285 in>>recq0[sizee-1];
286 in>>recq1[sizee-1];
287 in>>recq2[sizee-1];
288 in>>recq3[sizee-1];
289 in>>Norm;
290 }
291 if ( verbose ) cout<<"We have read recovered data"<<endl;
292
293
294 // IGRF stuff moved inside run loop!
295
296 for (Int_t ip=0;ip<nz;ip++){
297 zin[ip] = tof->GetZTOF(tof->GetToFPlaneID(ip));
298 };
299 //
300 if ( !standalone ){
301 //
302 // Does it contain the Tracker tree?
303 //
304 ttof = (TTree*)file->Get("ToF");
305 if ( !ttof ) {
306 if ( verbose ) printf(" OrbitalInfo - ERROR: no tof tree\n");
307 code = -900;
308 goto closeandexit;
309 };
310 ttof->SetBranchAddress("ToFLevel2",&tof);
311 nevtofl2 = ttof->GetEntries();
312 };
313 //
314 // Let's start!
315 //
316 // 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
317 // if run != 0 we must process only that run but first we have to check if the tree MyDetector2 already exist in the file
318 // if it exists we are reprocessing data and we must delete that entries, if not we must create it.
319 //
320 if ( run == 0 ) reproc = true;
321 //
322 //
323 // Output file is "outputfile"
324 //
325 if ( !file->IsOpen() ){
326 //printf(" OrbitalInfo - ERROR: cannot open file for writing\n");
327 throw -901;
328 };
329 //
330 // Retrieve GL_RUN variables from the level2 file
331 //
332 OrbitalInfoversion = OrbitalInfoInfo(false); // we should decide how to handle versioning system
333 //
334 // create an interface to RunInfo called "runinfo"
335 //
336 runinfo = new ItoRunInfo(file);
337 //
338 // open "Run" tree in level2 file, if not existing return an error (sngl != 0)
339 //
340 sgnl = 0;
341 sgnl = runinfo->Update(run, "ORB", OrbitalInfoversion);
342 //sgnl = runinfo->Read(run);
343
344 if ( sgnl ){
345 //printf("OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
346 code = sgnl;
347 goto closeandexit;
348 } else {
349 sgnl = 0;
350 };
351 //
352 // number of events in the file BEFORE the first event of our run
353 //
354 nobefrun = runinfo->GetFirstEntry();
355 //
356 // total number of events in the file
357 //
358 totfileentries = runinfo->GetFileEntries();
359 //
360 // first file entry AFTER the last event of our run
361 //
362 noaftrun = runinfo->GetLastEntry() + 1;
363 //
364 // number of run to be processed
365 //
366 numbofrun = runinfo->GetNoRun();
367 totnorun = runinfo->GetRunEntries();
368 //
369 // Try to access the OrbitalInfo tree in the file, if it exists we are reprocessing data if not we are processing a new run
370 //
371 OrbitalInfotrclone = (TTree*)file->Get("OrbitalInfo");
372 //
373 if ( !OrbitalInfotrclone ){
374 //
375 // tree does not exist, we are not reprocessing
376 //
377 reproc = false;
378 if ( run == 0 ){
379 if (verbose) printf(" OrbitalInfo - WARNING: you are reprocessing data but OrbitalInfo tree does not exist!\n");
380 }
381 if ( runinfo->IsReprocessing() && run != 0 ) {
382 if (verbose) printf(" OrbitalInfo - WARNING: it seems you are not reprocessing data but OrbitalInfo\n versioning information already exists in RunInfo.\n");
383 }
384 } else {
385 //
386 // tree exists, we are reprocessing data. Are we reprocessing a single run or all the file?
387 //
388 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
389 reproc = true;
390 //
391 //
392 if (verbose) printf("\n Preparing the pre-processing...\n");
393 //
394 if ( run == 0 || totnorun == 1 ){
395 //
396 // we are reprocessing all the file
397 // 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
398 //
399 reprocall = true;
400 //
401 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing all runs\n");
402 //
403 } else {
404 //
405 // 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
406 //
407 reprocall = false;
408 //
409 if (verbose) printf("\n OrbitalInfo - WARNING: Reprocessing run number %u \n",run);
410 //
411 // copying old tree to a new file
412 //
413 gSystem->MakeDirectory(OrbitalInfofolder.str().c_str());
414 myfold = true;
415 tempfile = new TFile(tempname.str().c_str(),"RECREATE");
416 tempOrbitalInfo = OrbitalInfotrclone->CloneTree(-1,"fast");
417 tempOrbitalInfo->SetName("OrbitalInfo-old");
418 tempfile->Write();
419 tempfile->Close();
420 }
421 //
422 // Delete the old tree from old file and memory
423 //
424 OrbitalInfotrclone->Delete("all");
425 //
426 if (verbose) printf(" ...done!\n");
427 //
428 };
429 //
430 // create mydetector tree mydect
431 //
432 file->cd();
433 OrbitalInfotr = new TTree("OrbitalInfo-new","PAMELA OrbitalInfo data");
434 OrbitalInfotr->SetAutoSave(900000000000000LL);
435 orbitalinfo->Set();//ELENA **TEMPORANEO?**
436 OrbitalInfotr->Branch("OrbitalInfo","OrbitalInfo",&orbitalinfo);
437 //
438 if ( reproc && !reprocall ){
439 //
440 // open new file and retrieve also tree informations
441 //
442 tempfile = new TFile(tempname.str().c_str(),"READ");
443 OrbitalInfotrclone = (TTree*)tempfile->Get("OrbitalInfo-old");
444 OrbitalInfotrclone->SetAutoSave(900000000000000LL);
445 OrbitalInfotrclone->SetBranchAddress("OrbitalInfo",&orbitalinfoclone);
446 //
447 if ( nobefrun > 0 ){
448 if (verbose){
449 printf("\n Pre-processing: copying events from the old tree before the processed run\n");
450 printf(" Copying %u events in the file which are before the beginning of the run %u \n",nobefrun,run);
451 printf(" Start copying at event number 0, end copying at event number %u \n",nobefrun);
452 }
453 for (UInt_t j = 0; j < nobefrun; j++){
454 //
455 if ( OrbitalInfotrclone->GetEntry(j) <= 0 ) throw -36;
456 //
457 // copy orbitalinfoclone to mydec
458 //
459 orbitalinfo->Clear();
460 //
461 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
462 //
463 // Fill entry in the new tree
464 //
465 OrbitalInfotr->Fill();
466 //
467 };
468 if (verbose) printf(" Finished successful copying!\n");
469 };
470 };
471 //
472 //
473 // Get the list of run to be processed, if only one run has to be processed the list will contain one entry only.
474 //
475 runlist = runinfo->GetRunList();
476 //
477 // Loop over the run to be processed
478 //
479 for (UInt_t irun=0; irun < numbofrun; irun++){
480 //
481 // retrieve the first run ID to be processed using the RunInfo list
482 //
483
484 idRun = runlist->At(irun);
485 if (verbose){
486 printf("\n\n\n ####################################################################### \n");
487 printf(" PROCESSING RUN NUMBER %i \n",(int)idRun);
488 printf(" ####################################################################### \n\n\n");
489 }
490 //
491 runinfo->ID_ROOT_L0 = 0;
492 //
493 // store in the runinfo class the GL_RUN variables for our run
494 //
495 sgnl = 0;
496 sgnl = runinfo->GetRunInfo(idRun);
497 if ( sgnl ){
498 if ( debug ) printf("\n OrbitalInfo - ERROR: RunInfo exited with non-zero status\n");
499 code = sgnl;
500 goto closeandexit;
501 } else {
502 sgnl = 0;
503 };
504 //
505 // now you can access that variables using the RunInfo class this way runinfo->ID_REG_RUN
506 //
507 if ( runinfo->ID_ROOT_L0 == 0 ){
508 if ( debug ) printf("\n OrbitalInfo - ERROR: no run with ID_RUN = %u \n\n Exiting... \n\n",idRun);
509 code = -5;
510 goto closeandexit;
511 };
512 //
513 // prepare the timesync for the db
514 //
515 dbtime = new GL_TIMESYNC(runinfo->ID_ROOT_L0,"ID",dbc);
516
517 //
518 // Search in the DB the path and name of the LEVEL0 file to be processed.
519 //
520 glroot->Query_GL_ROOT(runinfo->ID_ROOT_L0,dbc);
521 //
522 ftmpname.str("");
523 ftmpname << glroot->PATH.Data() << "/";
524 ftmpname << glroot->NAME.Data();
525 fname = ftmpname.str().c_str();
526 ftmpname.str("");
527 //
528 // print nout informations
529 //
530 totevent = runinfo->NEVENTS;
531 evfrom = runinfo->EV_FROM;
532 //cout<<"totevents = "<<totevent<<"\n";
533 if (verbose){
534 printf("\n LEVEL0 data file: %s \n",fname.Data());
535 printf(" RUN HEADER absolute time is: %u \n",runinfo->RUNHEADER_TIME);
536 printf(" RUN TRAILER absolute time is: %u \n",runinfo->RUNTRAILER_TIME);
537 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);
538 }//
539 //
540 // if ( !totevent ) goto closeandexit;
541 // Open Level0 file
542 l0File = new TFile(fname.Data());
543 if ( !l0File ) {
544 if ( debug ) printf(" OrbitalInfo - ERROR: problems opening Level0 file\n");
545 code = -6;
546 goto closeandexit;
547 };
548 l0tr = (TTree*)l0File->Get("Physics");
549 if ( !l0tr ) {
550 if ( debug ) printf(" OrbitalInfo - ERROR: no Physics tree in Level0 file\n");
551 l0File->Close();
552 code = -7;
553 goto closeandexit;
554 };
555 // EM: open header branch as well
556 l0head = l0tr->GetBranch("Header");
557 if ( !l0head ) {
558 if ( debug ) printf(" OrbitalInfo - ERROR: no Header branch in Level0 tree\n");
559 l0File->Close();
560 code = -8;
561 goto closeandexit;
562 };
563 l0tr->SetBranchAddress("Header", &eh);
564 // end EM
565 nevents = l0head->GetEntries();
566 //
567 if ( nevents < 1 && totevent ) {
568 if ( debug ) printf(" OrbitalInfo - ERROR: Level0 file is empty\n\n");
569 l0File->Close();
570 code = -11;
571 goto closeandexit;
572 };
573 //
574 if ( runinfo->EV_TO > nevents-1 && totevent ) {
575 if ( debug ) printf(" OrbitalInfo - ERROR: too few entries in the registry tree\n");
576 l0File->Close();
577 code = -12;
578 goto closeandexit;
579 };
580
581 //
582 // open IGRF files and do it only once if we are processing a full level2 file
583 //
584 if ( irun == 0 ){
585 if ( l0head->GetEntry(runinfo->EV_FROM) <= 0 ) throw -36;
586 //
587 // absolute time of first event of the run (it should not matter a lot)
588 //
589 ph = eh->GetPscuHeader();
590 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
591
592 parerror=glparam->Query_GL_PARAM(atime-anni5,301,dbc); // parameters stored in DB in GL_PRAM table
593 if ( parerror<0 ) {
594 code = parerror;
595 goto closeandexit;
596 };
597 ltp1 = (Int_t)(glparam->PATH+glparam->NAME).Length();
598 if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
599 //
600 parerror=glparam2->Query_GL_PARAM(atime,301,dbc); // parameters stored in DB in GL_PRAM table
601 if ( parerror<0 ) {
602 code = parerror;
603 goto closeandexit;
604 };
605 ltp2 = (Int_t)(glparam2->PATH+glparam->NAME).Length();
606 if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam2->PATH+glparam2->NAME).Data());
607 //
608 parerror=glparam3->Query_GL_PARAM(atime,302,dbc); // parameters stored in DB in GL_PRAM table
609 if ( parerror<0 ) {
610 code = parerror;
611 goto closeandexit;
612 };
613 ltp3 = (Int_t)(glparam3->PATH+glparam2->NAME).Length();
614 if ( verbose ) printf(" Reading Earth's Magnetic Field parameter file: %s \n",(glparam3->PATH+glparam3->NAME).Data());
615 //
616 initize_((char *)(glparam->PATH+glparam->NAME).Data(),&ltp1,(char *)(glparam2->PATH+glparam2->NAME).Data(),&ltp2,(char *)(glparam3->PATH+glparam3->NAME).Data(),&ltp3);
617 //
618 }
619 //
620 // End IGRF stuff//
621 //
622
623 //
624 // TTree *tp = (TTree*)l0File->Get("RunHeader");
625 // tp->SetBranchAddress("Header", &eH);
626 // tp->SetBranchAddress("RunHeader", &reh);
627 // tp->GetEntry(0);
628 // ph = eH->GetPscuHeader();
629 // ULong_t TimeSync = reh->LAST_TIME_SYNC_INFO;
630 // ULong_t ObtSync = reh->OBT_TIME_SYNC;
631 // if ( debug ) printf(" 1 TimeSync %lu ObtSync %lu DeltaOBT %lu\n",TimeSync,ObtSync,TimeSync-ObtSync);
632 //
633 ULong_t TimeSync = (ULong_t)dbtime->GetTimesync();
634 ULong_t ObtSync = (ULong_t)(dbtime->GetObt0()/1000);
635 ULong_t DeltaOBT = TimeSync - ObtSync;
636
637 if ( debug ) printf(" 2 TimeSync %lu ObtSync %lu DeltaOBT %lu\n",(ULong_t)(dbtime->GetTimesync()/1000),(ULong_t)dbtime->GetObt0(),TimeSync-ObtSync);
638 //
639 // Read MCMDs from up to 11 files, 5 before and 5 after the present one in order to have some kind of inclination information
640 //
641 ch = new TChain("Mcmd","Mcmd");
642 //
643 // look in the DB to find the closest files to this run
644 //
645 TSQLResult *pResult = 0;
646 TSQLRow *Row = 0;
647 stringstream myquery;
648 UInt_t l0fid[10];
649 Int_t i = 0;
650 memset(l0fid,0,10*sizeof(Int_t));
651 //
652 myquery.str("");
653 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;";
654 //
655 pResult = dbc->Query(myquery.str().c_str());
656 //
657 i = 9;
658 if( pResult ){
659 //
660 Row = pResult->Next();
661 //
662 while ( Row ){
663 //
664 // store infos and exit
665 //
666 l0fid[i] = (UInt_t)atoll(Row->GetField(0));
667 i--;
668 Row = pResult->Next();
669 //
670 };
671 pResult->Delete();
672 };
673 //
674 myquery.str("");
675 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;";
676 //
677 pResult = dbc->Query(myquery.str().c_str());
678 //
679 i = 0;
680 if( pResult ){
681 //
682 Row = pResult->Next();
683 //
684 while ( Row ){
685 //
686 // store infos and exit
687 //
688 l0fid[i] = (UInt_t)atoll(Row->GetField(0));
689 i++;
690 Row = pResult->Next();
691 //
692 };
693 pResult->Delete();
694 };
695 //
696 i = 0;
697 UInt_t previd = 0;
698 while ( i < 10 ){
699 if ( l0fid[i] && previd != l0fid[i] ){
700 previd = l0fid[i];
701 myquery.str("");
702 myquery << "select PATH,NAME from GL_ROOT where ID=" << l0fid[i] << " ;";
703 //
704 pResult = dbc->Query(myquery.str().c_str());
705 //
706 if( pResult ){
707 //
708 Row = pResult->Next();
709 //
710 if ( debug ) printf(" Using inclination informations from file: %s \n",(((TString)gSystem->ExpandPathName(Row->GetField(0)))+"/"+(TString)Row->GetField(1)).Data());
711 ch->Add(((TString)gSystem->ExpandPathName(Row->GetField(0)))+"/"+(TString)Row->GetField(1));
712 //
713 pResult->Delete();
714 };
715 };
716 i++;
717 };
718 //
719 // l0trm = (TTree*)l0File->Get("Mcmd");
720 // ch->ls();
721 ch->SetBranchAddress("Mcmd",&mcmdev);
722 // printf(" entries %llu \n", ch->GetEntries());
723 // l0trm = ch->GetTree();
724 // neventsm = l0trm->GetEntries();
725 neventsm = ch->GetEntries();
726 if ( debug ) printf(" entries %u \n", neventsm);
727 // neventsm = 0;
728 //
729 if (neventsm == 0){
730 if ( debug ) printf("InclinationInfo - WARNING: No quaternions in this File");
731 // l0File->Close();
732 code = 900;
733 // goto closeandexit;
734 }
735 //
736
737 // l0trm->SetBranchAddress("Mcmd", &mcmdev);
738 // l0trm->SetBranchAddress("Header", &eh);
739 //
740 //
741 //
742
743 // UInt_t mctren = 0;
744 // UInt_t mcreen = 0;
745 UInt_t numrec = 0;
746 //
747 Double_t upperqtime = 0;
748 Double_t lowerqtime = 0;
749
750 // Double_t incli = 0;
751 // oi = 0;
752 // UInt_t ooi = 0;
753 //
754 // init quaternions information from mcmd-packets
755 //
756 Bool_t isf = true;
757 // Int_t fgh = 0;
758
759 vector<Float_t> q0;
760 vector<Float_t> q1;
761 vector<Float_t> q2;
762 vector<Float_t> q3;
763 vector<Double_t> qtime;
764 vector<Float_t> qPitch;
765 vector<Float_t> qRoll;
766 vector<Float_t> qYaw;
767 vector<Int_t> qmode;
768
769 Int_t nt = 0;
770
771 //init sine-function interpolation
772
773 //cout<<"Sine coeficient initialisation..."<<endl;
774 vector<Sine> q0sine;
775 vector<Sine> q1sine;
776 vector<Sine> q2sine;
777 vector<Sine> q3sine;
778 vector<Sine> Yawsine;
779
780 /*TH2F* q0testing = new TH2F();
781 TH2F* q1testing = new TH2F();
782 TH2F* q2testing = new TH2F();
783 TH2F* q3testing = new TH2F();
784 TH2F* Pitchtesting = new TH2F();
785 */
786 UInt_t must = 0;
787
788 //
789 // run over all the events of the run
790 //
791 if (verbose) printf("\n Ready to start! \n\n Processed events: \n\n");
792 //
793 //
794 for ( re = runinfo->EV_FROM; re < (runinfo->EV_FROM+runinfo->NEVENTS); re++){
795 //
796 if ( procev%1000 == 0 && procev > 0 && verbose ) printf(" %iK \n",procev/1000);
797 if ( debug ) printf(" %i \n",procev);
798 //
799 if ( l0head->GetEntry(re) <= 0 ) throw -36;
800 //
801 // absolute time of this event
802 //
803 ph = eh->GetPscuHeader();
804 atime = dbtime->DBabsTime(ph->GetOrbitalTime());
805 if ( debug ) printf(" %i absolute time \n",procev);
806 //
807 // paranoid check
808 //
809 if ( (atime > (runinfo->RUNTRAILER_TIME+1)) || (atime < (runinfo->RUNHEADER_TIME-1)) ) {
810 if (verbose) printf(" OrbitalInfo - WARNING: event at time outside the run time window, skipping it\n");
811 jumped++;
812 // debug = true;
813 continue;
814 }
815
816 //
817 // retrieve tof informations
818 //
819 if ( !reprocall ){
820 itr = nobefrun + (re - evfrom - jumped);
821 //itr = re-(46438+200241);
822 } else {
823 itr = runinfo->GetFirstEntry() + (re - evfrom - jumped);
824 };
825 //
826 if ( !standalone ){
827 if ( itr > nevtofl2 ){
828 if ( verbose ) printf(" OrbitalInfo - ERROR: no tof events with entry = %i in Level2 file\n",itr);
829 if ( debug ) printf(" nobefrun %u re %u evfrom %u jumped %u reprocall %i \n",nobefrun,re,evfrom,jumped,reprocall);
830 l0File->Close();
831 code = -901;
832 goto closeandexit;
833 };
834 //
835 tof->Clear();
836 //
837 if ( ttof->GetEntry(itr) <= 0 ) throw -36;
838 //
839 };
840 //
841 procev++;
842 //
843 // start processing
844 //
845 if ( debug ) printf(" %i start processing \n",procev);
846 orbitalinfo->Clear();
847 //
848 OrbitalInfoTrkVar *t_orb = new OrbitalInfoTrkVar();
849 if( !(orbitalinfo->OrbitalInfoTrk) ) orbitalinfo->OrbitalInfoTrk = new TClonesArray("OrbitalInfoTrkVar",2);
850 TClonesArray &tor = *orbitalinfo->OrbitalInfoTrk;
851 //
852 // Fill OBT, pkt_num and absTime
853 //
854 orbitalinfo->pkt_num = ph->GetCounter();
855 orbitalinfo->OBT = ph->GetOrbitalTime();
856 orbitalinfo->absTime = atime;
857 if ( debug ) printf(" %i pktnum obt abstime \n",procev);
858 //
859 // Propagate the orbit from the tle time to atime, using SGP(D)4.
860 //
861 if ( debug ) printf(" %i sgp4 \n",procev);
862 cCoordGeo coo;
863 Float_t jyear=0.;
864 //
865 if(atime >= gltle->GetToTime()) {
866 if ( !gltle->Query(atime, dbc) ){
867 //
868 // Compute the magnetic dipole moment.
869 //
870 if ( debug ) printf(" %i compute magnetic dipole moment \n",procev);
871 UInt_t year, month, day, hour, min, sec;
872 //
873 TTimeStamp t = TTimeStamp(atime, kTRUE);
874 t.GetDate(kTRUE, 0, &year, &month, &day);
875 t.GetTime(kTRUE, 0, &hour, &min, &sec);
876 jyear = (float) year
877 + (month*31.+ (float) day)/365.
878 + (hour*3600.+min*60.+(float)sec)/(24.*3600.*365.);
879 //
880 if ( debug ) printf(" %i compute magnetic dipole moment get dipole moment for year\n",procev);
881 feldcof_(&jyear, &dimo); // get dipole moment for year
882 if ( debug ) printf(" %i compute magnetic dipole moment end\n",procev);
883 } else {
884 code = -56;
885 goto closeandexit;
886 };
887 }
888 coo = getCoo(atime, gltle->GetFromTime(), gltle->GetTle());
889 //
890 cOrbit orbits(*gltle->GetTle());
891 //
892 if ( debug ) printf(" I am Here \n");
893 //
894 // synchronize with quaternions data
895 //
896 if ( isf && neventsm>0 ){
897 //
898 // First event
899 //
900 isf = false;
901 upperqtime = atime;
902 lowerqtime = runinfo->RUNHEADER_TIME;
903 for ( ik = 0; ik < neventsm; ik++){ //number of macrocommad packets
904 if ( ch->GetEntry(ik) <= 0 ) throw -36;
905 tmpSize = mcmdev->Records->GetEntries();
906 numrec = tmpSize;
907 for (Int_t j3 = 0;j3<tmpSize;j3++){ //number of subpackets
908 if ( debug ) printf(" ik %i j3 %i eh eh eh \n",ik,j3);
909 mcmdrc = (pamela::McmdRecord*)mcmdev->Records->At(j3);
910 if ( mcmdrc ){ // missing inclination bug [8RED 090116]
911 if ( debug ) printf(" pluto \n");
912 if ((int)mcmdrc->ID1 == 226 && mcmdrc->Mcmd_Block_crc_ok == 1){ //Check that it is Inclination Packet
913 L_QQ_Q_l_upper->fill(mcmdrc->McmdData);
914 for (UInt_t ui = 0; ui < 6; ui++){
915 if (ui>0){
916 if (L_QQ_Q_l_upper->time[ui]>L_QQ_Q_l_upper->time[0]){
917 if ( debug ) printf(" here1 %i \n",ui);
918 Double_t u_time = dbtime->DBabsTime((UInt_t)(L_QQ_Q_l_upper->time[ui]*1000-DeltaOBT*1000));
919 Int_t recSize = recqtime.size();
920 if(lowerqtime > recqtime[recSize-1]){
921 Int_t sizeqmcmd = qtime.size();
922 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
923 qtime[sizeqmcmd]=u_time;
924 q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][0];
925 q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][1];
926 q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][2];
927 q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][3];
928 qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
929 lowerqtime = u_time;
930 orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
931 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]);
932 qRoll[sizeqmcmd]=RYPang_upper->Kren;
933 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
934 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
935 }
936 for(Int_t mu = nt;mu<recSize;mu++){
937 if(recqtime[mu]>lowerqtime && recqtime[mu]<u_time){
938 nt=mu;
939 Int_t sizeqmcmd = qtime.size();
940 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
941 qtime[sizeqmcmd]=recqtime[mu];
942 q0[sizeqmcmd]=recq0[mu];
943 q1[sizeqmcmd]=recq1[mu];
944 q2[sizeqmcmd]=recq2[mu];
945 q3[sizeqmcmd]=recq3[mu];
946 qmode[sizeqmcmd]=-10;
947 orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
948 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]);
949 qRoll[sizeqmcmd]=RYPang_upper->Kren;
950 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
951 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
952 }
953 if(recqtime[mu]>=u_time){
954 Int_t sizeqmcmd = qtime.size();
955 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
956 qtime[sizeqmcmd]=u_time;
957 q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][0];
958 q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][1];
959 q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][2];
960 q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[ui][3];
961 qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
962 lowerqtime = u_time;
963 orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
964 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]);
965 qRoll[sizeqmcmd]=RYPang_upper->Kren;
966 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
967 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
968 break;
969 }
970 }
971 }
972 }else{
973 if ( debug ) printf(" here2 %i \n",ui);
974 Double_t u_time = dbtime->DBabsTime((UInt_t)(L_QQ_Q_l_upper->time[0]*1000-DeltaOBT*1000));
975 if(lowerqtime>u_time)nt=0;
976 Int_t recSize = recqtime.size();
977 if(lowerqtime > recqtime[recSize-1]){
978 Int_t sizeqmcmd = qtime.size();
979 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
980 qtime[sizeqmcmd]=u_time;
981 q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][0];
982 q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][1];
983 q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][2];
984 q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][3];
985 qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
986 lowerqtime = u_time;
987 orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
988 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]);
989 qRoll[sizeqmcmd]=RYPang_upper->Kren;
990 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
991 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
992 }
993 for(Int_t mu = nt;mu<recSize;mu++){
994 if(recqtime[mu]>lowerqtime && recqtime[mu]<u_time){
995 nt=mu;
996 Int_t sizeqmcmd = qtime.size();
997 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
998 qtime[sizeqmcmd]=recqtime[mu];
999 q0[sizeqmcmd]=recq0[mu];
1000 q1[sizeqmcmd]=recq1[mu];
1001 q2[sizeqmcmd]=recq2[mu];
1002 q3[sizeqmcmd]=recq3[mu];
1003 qmode[sizeqmcmd]=-10;
1004 orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
1005 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]);
1006 qRoll[sizeqmcmd]=RYPang_upper->Kren;
1007 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1008 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1009 }
1010 if(recqtime[mu]>=u_time){
1011 Int_t sizeqmcmd = qtime.size();
1012 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1013 qtime[sizeqmcmd]=u_time;
1014 q0[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][0];
1015 q1[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][1];
1016 q2[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][2];
1017 q3[sizeqmcmd]=L_QQ_Q_l_upper->quat[0][3];
1018 qmode[sizeqmcmd]=holeq(lowerqtime,qtime[sizeqmcmd],L_QQ_Q_l_lower,L_QQ_Q_l_upper,ui);
1019 lowerqtime = u_time;
1020 orbits.getPosition((double) (u_time - gltle->GetFromTime())/60., &eCi);
1021 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]);
1022 qRoll[sizeqmcmd]=RYPang_upper->Kren;
1023 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1024 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1025 CopyQ(L_QQ_Q_l_lower,L_QQ_Q_l_upper);
1026 break;
1027 }
1028 }
1029 }
1030 }
1031 }
1032 }
1033 if ( debug ) printf(" ciccio \n");
1034 }
1035 }
1036
1037 if(qtime.size()==0){
1038 for(UInt_t my=0;my<recqtime.size();my++){
1039 Int_t sizeqmcmd = qtime.size();
1040 inclresize(qtime,q0,q1,q2,q3,qmode,qRoll,qPitch,qYaw);
1041 qtime[sizeqmcmd]=recqtime[my];
1042 q0[sizeqmcmd]=recq0[my];
1043 q1[sizeqmcmd]=recq1[my];
1044 q2[sizeqmcmd]=recq2[my];
1045 q3[sizeqmcmd]=recq3[my];
1046 qmode[sizeqmcmd]=-10;
1047 orbits.getPosition((double) (qtime[sizeqmcmd] - gltle->GetFromTime())/60., &eCi);
1048 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]);
1049 qRoll[sizeqmcmd]=RYPang_upper->Kren;
1050 qYaw[sizeqmcmd]=RYPang_upper->Ryskanie;
1051 qPitch[sizeqmcmd]=RYPang_upper->Tangazh;
1052 }
1053 }
1054
1055 if ( debug ) printf(" fuffi \n");
1056
1057 //sineparam(q0sine,qtime,q0,qRoll,qPitch,0.60);
1058 //sineparam(q1sine,qtime,q1,qRoll,qPitch,0.82);
1059 //sineparam(q2sine,qtime,q2,qRoll,qPitch,0.82);
1060 //sineparam(q3sine,qtime,q3,qRoll,qPitch,0.60);
1061 //sineparam(Yawsine,qtime,qYaw,qRoll,qPitch,4);
1062
1063 if ( debug ) printf(" puffi \n");
1064 Double_t tmin = 9999999999.;
1065 Double_t tmax = 0.;
1066 for(UInt_t tre = 0;tre<qtime.size();tre++){
1067 if(qtime[tre]>tmax)tmax = qtime[tre];
1068 if(qtime[tre]<tmin)tmin = qtime[tre];
1069 }
1070 if ( debug ) printf(" gnfuffi \n");
1071
1072 //q0testing->SetName("q0testing");
1073 //q1testing->SetName("q1testing");
1074 //q2testing->SetName("q2testing");
1075 //q3testing->SetName("q3testing");
1076
1077 // Int_t ss=10.*(tmax-tmin);
1078 //q0testing->SetBins(ss,tmin,tmax,1000,-1.,1.);
1079 //Pitchtesting->SetBins(ss,tmin,tmax,1000,-40.,40.);
1080
1081 // for(Int_t tre = 0;tre<qtime.size();tre++){
1082 //cout<<"q0["<<tre<<" = "<<q0[tre]<<endl;
1083 //q0testing->Fill(qtime[tre],q0[tre]);
1084 //q1testing->Fill(qtime[tre],q1[tre]);
1085 //Pitchtesting->Fill(qtime[tre],qPitch[tre],100);
1086 //if(qmode[tre] == -10)Pitchtesting->Fill(qtime[tre],10,100);
1087 //q2testing->Fill(qtime[tre],q2[tre],100);
1088 //q3testing->Fill(qtime[tre],q3[tre],100);
1089 // }
1090
1091 //for(Int_t tre=0;tre<q0sine.size();tre++)cout<<q1sine[tre].A<<"*sin("<<q1sine[tre].b<<"x+"<<q1sine[tre].c<<")\t time start: "<<q1sine[tre].startPoint<<"\ttime end: "<<q1sine[tre].finishPoint<<endl;
1092 //for(Int_t tre=0;tre<q0sine.size();tre++)cout<<q1sine[tre].A<<"*sin("<<q1sine[tre].b<<"x+"<<q1sine[tre].c<<")\t time start: "<<q0sine[tre].startPoint<<"\ttime end: "<<q0sine[tre].finishPoint<<endl;
1093 } // if we processed first event
1094
1095 //Filling Inclination information
1096 Double_t incli = 0;
1097 if ( qtime.size() > 1 ){
1098 for(UInt_t mu = must;mu<qtime.size()-1;mu++){
1099 if ( debug ) printf(" ??grfuffi %i sixe %i must %i \n",mu,qtime.size()-1,must);
1100 if(qtime[mu+1]>qtime[mu]){
1101 if ( debug ) printf(" grfuffi2 %i \n",mu);
1102 if(atime<=qtime[mu+1] && atime>=qtime[mu]){
1103 must = mu;
1104 incli = (qPitch[mu+1]-qPitch[mu])/(qtime[mu+1]-qtime[mu]);
1105 orbitalinfo->theta = incli*atime+qPitch[mu+1]-incli*qtime[mu+1];
1106 incli = (qRoll[mu+1]-qRoll[mu])/(qtime[mu+1]-qtime[mu]);
1107 orbitalinfo->etha = incli*atime+qRoll[mu+1]-incli*qtime[mu+1];
1108 incli = (qYaw[mu+1]-qYaw[mu])/(qtime[mu+1]-qtime[mu]);
1109 orbitalinfo->phi = incli*atime+qYaw[mu+1]-incli*qtime[mu+1];
1110
1111 incli = (q0[mu+1]-q0[mu])/(qtime[mu+1]-qtime[mu]);
1112 orbitalinfo->q0 = incli*atime+q0[mu+1]-incli*qtime[mu+1];
1113 incli = (q1[mu+1]-q1[mu])/(qtime[mu+1]-qtime[mu]);
1114 orbitalinfo->q1 = incli*atime+q1[mu+1]-incli*qtime[mu+1];
1115 incli = (q2[mu+1]-q2[mu])/(qtime[mu+1]-qtime[mu]);
1116 orbitalinfo->q2 = incli*atime+q2[mu+1]-incli*qtime[mu+1];
1117 incli = (q3[mu+1]-q3[mu])/(qtime[mu+1]-qtime[mu]);
1118 orbitalinfo->q3 = incli*atime+q3[mu+1]-incli*qtime[mu+1];
1119
1120 orbitalinfo->TimeGap = qtime[mu+1]-qtime[mu];
1121 orbitalinfo->mode = qmode[mu+1];
1122 //if(qmode[mu+1]==-10) orbitalinfo->R10r = true;else orbitalinfo->R10r = false;
1123 //reserved for next versions Vitaly.
1124 /*if(qmode[mu+1]==-10 || qmode[mu+1]==0 || qmode[mu+1]==1 || qmode[mu+1]==3 || qmode[mu+1]==4 || qmode[mu+1]==6){
1125 //linear interpolation
1126 incli = (q0[mu+1]-q0[mu])/(qtime[mu+1]-qtime[mu]);
1127 orbitalinfo->q0 = incli*atime+q0[mu+1]-incli*qtime[mu+1];
1128 incli = (q1[mu+1]-q1[mu])/(qtime[mu+1]-qtime[mu]);
1129 orbitalinfo->q1 = incli*atime+q1[mu+1]-incli*qtime[mu+1];
1130 incli = (q2[mu+1]-q2[mu])/(qtime[mu+1]-qtime[mu]);
1131 orbitalinfo->q2 = incli*atime+q2[mu+1]-incli*qtime[mu+1];
1132 incli = (q3[mu+1]-q3[mu])/(qtime[mu+1]-qtime[mu]);
1133 orbitalinfo->q3 = incli*atime+q3[mu+1]-incli*qtime[mu+1];
1134 }else{
1135 //sine interpolation
1136 for(UInt_t mt=0;mt<q0sine.size();mt++){
1137 if(atime<=q0sine[mt].finishPoint && atime>=q0sine[mt].startPoint){
1138 if(!q0sine[mt].NeedFit)orbitalinfo->q0=q0sine[mt].A*sin(q0sine[mt].b*atime+q0sine[mt].c);else{
1139 incli = (q0[mu+1]-q0[mu])/(qtime[mu+1]-qtime[mu]);
1140 orbitalinfo->q0 = incli*atime+q0[mu+1]-incli*qtime[mu+1];
1141 }
1142 }
1143 if(atime<=q1sine[mt].finishPoint && atime>=q1sine[mt].startPoint){
1144 if(!q1sine[mt].NeedFit)orbitalinfo->q1=q1sine[mt].A*sin(q1sine[mt].b*atime+q1sine[mt].c);else{
1145 incli = (q1[mu+1]-q1[mu])/(qtime[mu+1]-qtime[mu]);
1146 orbitalinfo->q1 = incli*atime+q1[mu+1]-incli*qtime[mu+1];
1147 }
1148 }
1149 if(atime<=q2sine[mt].finishPoint && atime>=q2sine[mt].startPoint){
1150 if(!q2sine[mt].NeedFit)orbitalinfo->q2=q0sine[mt].A*sin(q2sine[mt].b*atime+q2sine[mt].c);else{
1151 incli = (q2[mu+1]-q2[mu])/(qtime[mu+1]-qtime[mu]);
1152 orbitalinfo->q2 = incli*atime+q2[mu+1]-incli*qtime[mu+1];
1153 }
1154 }
1155 if(atime<=q3sine[mt].finishPoint && atime>=q3sine[mt].startPoint){
1156 if(!q3sine[mt].NeedFit)orbitalinfo->q3=q0sine[mt].A*sin(q3sine[mt].b*atime+q3sine[mt].c);else{
1157 incli = (q3[mu+1]-q3[mu])/(qtime[mu+1]-qtime[mu]);
1158 orbitalinfo->q3 = incli*atime+q3[mu+1]-incli*qtime[mu+1];
1159 }
1160 }
1161 if(atime<=Yawsine[mt].finishPoint && atime>=Yawsine[mt].startPoint){
1162 if(!Yawsine[mt].NeedFit)orbitalinfo->phi=Yawsine[mt].A*sin(Yawsine[mt].b*atime+Yawsine[mt].c);else{
1163 incli = (qYaw[mu+1]-qYaw[mu])/(qtime[mu+1]-qtime[mu]);
1164 orbitalinfo->phi = incli*atime+qYaw[mu+1]-incli*qtime[mu+1];
1165 }
1166 }
1167 }
1168 }*/
1169 //q0testing->Fill(atime,orbitalinfo->q0,100);
1170 //q1testing->Fill(atime,orbitalinfo->q1,100);
1171 //Pitchtesting->Fill(atime,orbitalinfo->etha);
1172 //q2testing->Fill(atime,orbitalinfo->q2);
1173 //q3testing->Fill(atime,orbitalinfo->q3);
1174 break;
1175 }
1176 }
1177 }
1178 }
1179 //
1180 // ops no inclination information
1181 //
1182
1183 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 ){
1184 orbitalinfo->mode = 10;
1185 orbitalinfo->q0 = -1000.;
1186 orbitalinfo->q1 = -1000.;
1187 orbitalinfo->q2 = -1000.;
1188 orbitalinfo->q3 = -1000.;
1189 orbitalinfo->etha = -1000.;
1190 orbitalinfo->phi = -1000.;
1191 orbitalinfo->theta = -1000.;
1192 };
1193 //
1194 // #########################################################################################################################
1195 //
1196 // fill orbital positions
1197 //
1198 // Build coordinates in the right range. We want to convert,
1199 // longitude from (0, 2*pi) to (-180deg, 180deg). Altitude is
1200 // in meters.
1201 lon = (coo.m_Lon > M_PI) ? rad2deg(coo.m_Lon - 2*M_PI) : rad2deg(coo.m_Lon);
1202 lat = rad2deg(coo.m_Lat);
1203 alt = coo.m_Alt;
1204 //
1205 if( lon<180 && lon>-180 && lat<90 && lat>-90 && alt>0 ){
1206 //
1207 orbitalinfo->lon = lon;
1208 orbitalinfo->lat = lat;
1209 orbitalinfo->alt = alt ;
1210 //
1211 // compute mag field components and L shell.
1212 //
1213 feldg_(&lat, &lon, &alt, &bnorth, &beast, &bdown, &babs);
1214 shellg_(&lat, &lon, &alt, &dimo, &xl, &icode, &bab1);
1215 findb0_(&stps, &bdel, &value, &bequ, &rr0);
1216 //
1217 orbitalinfo->Bnorth = bnorth;
1218 orbitalinfo->Beast = beast;
1219 orbitalinfo->Bdown = bdown;
1220 orbitalinfo->Babs = babs;
1221 orbitalinfo->BB0 = babs/bequ;
1222 orbitalinfo->L = xl;
1223 // Set Stormer vertical cutoff using L shell.
1224 orbitalinfo->cutoffsvl = 14.9/(xl*xl);
1225 orbitalinfo->igrf_icode = icode;
1226 //
1227 };
1228 //
1229 if ( debug ) printf(" pitch angle \n");
1230 //
1231 // pitch angles
1232 //
1233 //if ( orbitalinfo->mode != 10 && orbitalinfo->mode != 5 && orbitalinfo->mode !=7 && orbitalinfo->mode != 9 ){
1234 if( orbitalinfo->TimeGap>0 && orbitalinfo->TimeGap<2000000){
1235 //
1236 Float_t Bx = -orbitalinfo->Bdown;
1237 Float_t By = orbitalinfo->Beast;
1238 Float_t Bz = orbitalinfo->Bnorth;
1239 //
1240 TMatrixD Fij = PO->ECItoGreenwich(PO->QuatoECI(orbitalinfo->q0,orbitalinfo->q1,orbitalinfo->q2,orbitalinfo->q3),orbitalinfo->absTime);
1241 TMatrixD Gij = PO->ColPermutation(Fij);
1242 TMatrixD Dij = PO->GreenwichtoGEO(orbitalinfo->lat,orbitalinfo->lon,Fij);
1243 TMatrixD Iij = PO->ColPermutation(Dij);
1244 //
1245 orbitalinfo->Iij.ResizeTo(Iij);
1246 orbitalinfo->Iij = Iij;
1247 //
1248 A1 = Iij(0,2);
1249 A2 = Iij(1,2);
1250 A3 = Iij(2,2);
1251 //
1252 // orbitalinfo->pamzenitangle = (Float_t)PO->GetPitchAngle(1,0,0,A1,A2,A3); // Angle between zenit and Pamela's main axiz
1253 // orbitalinfo->pamBangle = (Float_t)PO->GetPitchAngle(A1,A2,A3,Bx,By,Bz); // Angle between Pamela's main axiz and B
1254 //
1255 if ( !standalone && tof->ntrk() > 0 ){
1256 //
1257 Int_t nn = 0;
1258 for(Int_t nt=0; nt < tof->ntrk(); nt++){
1259 //
1260 ToFTrkVar *ptt = tof->GetToFTrkVar(nt);
1261 Double_t E11x = ptt->xtr_tof[0]; // tr->x[0];
1262 Double_t E11y = ptt->ytr_tof[0]; //tr->y[0];
1263 Double_t E11z = zin[0];
1264 Double_t E22x = ptt->xtr_tof[3];//tr->x[3];
1265 Double_t E22y = ptt->ytr_tof[3];//tr->y[3];
1266 Double_t E22z = zin[3];
1267 if ( (E11x < 100. && E11y < 100. && E22x < 100. && E22y < 100.) || ptt->trkseqno != -1 ){
1268 Double_t norm = sqrt(pow(E22x-E11x,2)+pow(E22y-E11y,2)+pow(E22z-E11z,2));
1269 // Double_t MyAzim = TMath::RadToDeg()*atan(TMath::Abs(E22y-E11y)/TMath::Abs(E22x-E11x));
1270 // if(E22x-E11x>=0 && E22y-E11y <0) MyAzim = 360. - MyAzim;
1271 // if(E22x-E11x>=0 && E22y-E11y >=0) MyAzim = MyAzim;
1272 // if(E22x-E11x<0 && E22y-E11y >0) MyAzim = 180. - MyAzim;
1273 // if(E22x-E11x<0 && E22y-E11y <0) MyAzim = 180. + MyAzim;
1274 Px = (E22x-E11x)/norm;
1275 Py = (E22y-E11y)/norm;
1276 Pz = (E22z-E11z)/norm;
1277 //
1278 t_orb->trkseqno = ptt->trkseqno;
1279 //
1280 TMatrixD Eij = PO->PamelatoGEO(Iij,Px,Py,Pz);
1281 t_orb->Eij.ResizeTo(Eij);
1282 t_orb->Eij = Eij;
1283 //
1284 TMatrixD Sij = PO->PamelatoGEO(Gij,Px,Py,Pz);
1285 t_orb->Sij.ResizeTo(Sij);
1286 t_orb->Sij = Sij;
1287 //
1288 t_orb->pitch = (Float_t)PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),Bx,By,Bz);
1289 //
1290 //
1291 Double_t omega = PO->GetPitchAngle(Eij(0,0),Eij(1,0),Eij(2,0),cos(orbitalinfo->lon+TMath::Pi()/2)-sin(orbitalinfo->lon+TMath::Pi()/2),cos(orbitalinfo->lon+TMath::Pi()/2)+sin(orbitalinfo->lon+TMath::Pi()/2),1);
1292 //
1293 t_orb->cutoff = 59.3/(pow(orbitalinfo->L,2)*pow((1+sqrt(1-pow(orbitalinfo->L,-3/2)*cos(omega))),2));
1294 //
1295 if ( t_orb->pitch != t_orb->pitch ) t_orb->pitch = -1000.;
1296 if ( t_orb->cutoff != t_orb->cutoff ) t_orb->cutoff = -1000.;
1297 //
1298 if ( debug ) printf(" orbitalinfo->cutoffsvl %f vitaly %f \n",orbitalinfo->cutoffsvl,t_orb->cutoff);
1299 //
1300 new(tor[nn]) OrbitalInfoTrkVar(*t_orb);
1301 nn++;
1302 //
1303 t_orb->Clear();
1304 //
1305 };
1306 //
1307 };
1308 } else {
1309 if ( debug ) printf(" mmm... mode %u standalone %i ntrk %i \n",orbitalinfo->mode,standalone,tof->ntrk());
1310 };
1311 //
1312 } else {
1313 if ( !standalone && tof->ntrk() > 0 ){
1314 //
1315 Int_t nn = 0;
1316 for(Int_t nt=0; nt < tof->ntrk(); nt++){
1317 //
1318 ToFTrkVar *ptt = tof->GetToFTrkVar(nt);
1319 if ( ptt->trkseqno != -1 ){
1320 //
1321 t_orb->trkseqno = ptt->trkseqno;
1322 //
1323 t_orb->Eij = 0;
1324 //
1325 t_orb->Sij = 0;
1326 //
1327 t_orb->pitch = -1000.;
1328 //
1329 t_orb->cutoff = -1000.;
1330 //
1331 new(tor[nn]) OrbitalInfoTrkVar(*t_orb);
1332 nn++;
1333 //
1334 t_orb->Clear();
1335 //
1336 };
1337 //
1338 };
1339 };
1340 };
1341 //
1342 // Fill the class
1343 //
1344 OrbitalInfotr->Fill();
1345 //
1346 delete t_orb;
1347 //
1348 }; // loop over the events in the run
1349 //
1350 // Here you may want to clear some variables before processing another run
1351 //
1352
1353 //gStyle->SetOptStat(000000);
1354 //gStyle->SetPalette(1);
1355
1356 /*TCanvas* c1 = new TCanvas("c1","",1200,800);
1357 //c1->Divide(1,4);
1358 c1->cd(1);
1359 //q0testing->Draw("colz");
1360 //c1->cd(2);
1361 //q1testing->Draw("colz");
1362 //c1->cd(3);
1363 Pitchtesting->Draw("colz");
1364 //c1->cd(4);
1365 //q3testing->Draw("colz");
1366 c1->SaveAs("9.Rollhyst.png");
1367 delete c1;*/
1368
1369 delete dbtime;
1370 if ( L_QQ_Q_l_upper ) delete L_QQ_Q_l_upper;
1371 if ( L_QQ_Q_l_lower ) delete L_QQ_Q_l_lower;
1372 if ( RYPang_upper ) delete RYPang_upper;
1373 if ( RYPang_lower ) delete RYPang_lower;
1374 }; // process all the runs
1375
1376 if (verbose) printf("\n Finished processing data \n");
1377 //
1378 closeandexit:
1379 //
1380 // 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.
1381 //
1382 if ( !reprocall && reproc && code >= 0 ){
1383 if ( totfileentries > noaftrun ){
1384 if (verbose){
1385 printf("\n Post-processing: copying events from the old tree after the processed run\n");
1386 printf(" Copying %i events in the file which are after the end of the run %i \n",(int)(totfileentries-noaftrun),(int)run);
1387 printf(" Start copying at event number %i end copying at event number %i \n",(int)noaftrun,(int)totfileentries);
1388 }
1389 for (UInt_t j = noaftrun; j < totfileentries; j++ ){
1390 //
1391 // Get entry from old tree
1392 //
1393 if ( OrbitalInfotrclone->GetEntry(j) <= 0 ) throw -36;
1394 //
1395 // copy orbitalinfoclone to OrbitalInfo
1396 //
1397 orbitalinfo->Clear();
1398 //
1399 memcpy(&orbitalinfo,&orbitalinfoclone,sizeof(orbitalinfoclone));
1400 //
1401 // Fill entry in the new tree
1402 //
1403 OrbitalInfotr->Fill();
1404 };
1405 if (verbose) printf(" Finished successful copying!\n");
1406 };
1407 };
1408 //
1409 // Close files, delete old tree(s), write and close level2 file
1410 //
1411 if ( l0File ) l0File->Close();
1412 if ( tempfile ) tempfile->Close();
1413 if ( myfold ) gSystem->Unlink(tempname.str().c_str());
1414 //
1415 if ( runinfo ) runinfo->Close();
1416 if ( OrbitalInfotr ) OrbitalInfotr->SetName("OrbitalInfo");
1417 if ( tof ) tof->Delete();
1418 if ( ttof ) ttof->Delete();
1419 //
1420 if ( file ){
1421 file->cd();
1422 file->Write();
1423 };
1424 //
1425 if ( myfold ) gSystem->Unlink(OrbitalInfofolder.str().c_str());
1426 //
1427 // the end
1428 //
1429 if ( dbc ){
1430 dbc->Close();
1431 delete dbc;
1432 };
1433 if (verbose) printf("\n Exiting...\n");
1434 if(OrbitalInfotr)OrbitalInfotr->Delete();
1435 //
1436 if ( PO ) delete PO;
1437 if ( orbitalinfo ) delete orbitalinfo;
1438 if ( orbitalinfoclone ) delete orbitalinfoclone;
1439 if ( glroot ) delete glroot;
1440 if ( runinfo ) delete runinfo;
1441 //
1442 if(code < 0) throw code;
1443 return(code);
1444 }
1445
1446
1447 //
1448 // Returns the cCoordGeo structure holding the geographical
1449 // coordinates for the event (see sgp4.h).
1450 //
1451 // atime is the abstime of the event in UTC unix time.
1452 // tletime is the time of the tle in UTC unix time.
1453 // tle is the previous and nearest tle (compared to atime).
1454 cCoordGeo getCoo(UInt_t atime, UInt_t tletime, cTle *tle)
1455 {
1456 cEci eci;
1457 cOrbit orbit(*tle);
1458 orbit.getPosition((double) (atime - tletime)/60., &eci);
1459
1460 return eci.toGeo();
1461 }
1462
1463 // function of copyng of quatrnions classes
1464
1465 void CopyQ(Quaternions *Q1, Quaternions *Q2){
1466 for(UInt_t i = 0; i < 6; i++){
1467 Q1->time[i]=Q2->time[i];
1468 for (UInt_t j = 0; j < 4; j++)Q1->quat[i][j]=Q2->quat[i][j];
1469 }
1470 return;
1471 }
1472
1473 // functions of copyng InclinationInfo classes
1474
1475 void CopyAng(InclinationInfo *A1, InclinationInfo *A2){
1476 A1->Tangazh = A2->Tangazh;
1477 A1->Ryskanie = A2->Ryskanie;
1478 A1->Kren = A2->Kren;
1479 return;
1480 }
1481
1482 UInt_t holeq(Double_t lower,Double_t upper,Quaternions *Qlower, Quaternions *Qupper, UInt_t f){
1483
1484 UInt_t hole = 10;
1485 Bool_t R10l = false; // Sign of R10 mode in lower quaternions array
1486 Bool_t R10u = false; // Sign of R10 mode in upper quaternions array
1487 Bool_t insm = false; // Sign that we inside quaternions array
1488 Bool_t mxtml = false; // Sign of mixt mode in lower quaternions array
1489 Bool_t mxtmu = false; // Sign of mixt mode in upper quaternions array
1490 Bool_t npasm = false; // Sign of normall pass between R10 and non R10 or between non R10 and R10
1491 UInt_t NCQl = 6; // Number of correct quaternions in lower array
1492 UInt_t NCQu = 6; // Number of correct quaternions in upper array
1493 if (f>0){
1494 insm = true;
1495 if(Qupper->time[f]-Qupper->time[f-1]==30) R10u = false;
1496 if(Qupper->time[f]-Qupper->time[f-1]<1) R10u = true;
1497 }else{
1498 insm = false;
1499 if((Qlower->time[5]-Qlower->time[0]<2)&&(Qlower->time[1]-Qlower->time[0]<2)) R10l = true;
1500 if((Qupper->time[5]-Qupper->time[0]<2)&&(Qupper->time[1]-Qupper->time[0]<2)) R10u = true;
1501 if((Qlower->time[5]-Qlower->time[0]==150)&&(Qlower->time[1]-Qlower->time[0]==30)) R10l = false;
1502 if((Qupper->time[5]-Qupper->time[0]==150)&&(Qupper->time[1]-Qupper->time[0]==30)) R10u = false;
1503 if((Qlower->time[5]-Qlower->time[0]<2)&&(Qlower->time[1]-Qlower->time[0]==30)){
1504 mxtml = true;
1505 for(UInt_t i = 1; i < 6; i++){
1506 if(Qlower->time[i]-Qlower->time[0]==30*i) NCQl=i;
1507 }
1508 }
1509 if((Qupper->time[5]-Qupper->time[0]<2)&&(Qupper->time[1]-Qupper->time[0]==30)){
1510 mxtmu = true;
1511 for(UInt_t i = 1; i < 6; i++){
1512 if(Qupper->time[i]-Qupper->time[0]==30*i) NCQu=i;
1513 }
1514 }
1515 }
1516
1517 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;
1518
1519
1520 if (R10u&&insm) hole=0; // best event R10
1521 if ((upper-lower<=5)&&(!insm)&&R10l&&R10u) hole = 1; // when first of 6 quaternions in array is correct
1522 if (((!R10u)&&insm)||((!insm)&&(!R10u)&&(!R10l)&&((upper-lower==210+(6-NCQl)*30)||(upper-lower==30)))) hole = 2; //non R10
1523 if (npasm&&(!insm)&&((R10l&&!R10u)||(R10u&&!R10l))) hole = 3; //normall pass from R10 to non R10 or from non R10 to R10
1524 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
1525 if ((upper-lower>=300)&&(!insm)&&((R10l&&!R10u)||(R10u&&!R10l))) hole = 5; //uneliminable hole between R10 and non R10 or between non R10 and R10
1526 if ((upper-lower>5)&&(upper-lower<=300)&&R10u&&R10l) hole = 6; // eliminable hole inside R10
1527 if ((upper-lower>300)&&R10u&&R10l) hole = 7; //uneliminable hole inside R10
1528 if ((upper-lower>210)&&(upper-lower<=1200)&&(!R10u)&&(!R10l)) hole = 8; //eliminable hole inside non R10
1529 if ((upper-lower>1200)&&!R10u&&!R10l) hole = 9; // uneliminable hole inside non R10
1530 return hole;
1531 }
1532
1533 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){
1534 Int_t sizee = t.size()+1;
1535 t.resize(sizee);
1536 q0.resize(sizee);
1537 q1.resize(sizee);
1538 q2.resize(sizee);
1539 q3.resize(sizee);
1540 mode.resize(sizee);
1541 Roll.resize(sizee);
1542 Pitch.resize(sizee);
1543 Yaw.resize(sizee);
1544 }
1545
1546 //Find fitting sine functions for q0,q1,q2,q3 and Yaw-angle;
1547 void sineparam(vector<Sine>& qsine, vector<Double_t>& qtime, vector<Float_t>& q, vector<Float_t>& Roll, vector<Float_t>& Pitch, Float_t limsin){
1548 UInt_t mulast = 0;
1549 UInt_t munow = 0;
1550 UInt_t munext = 0;
1551 Bool_t increase = false;
1552 Bool_t decrease = false;
1553 Bool_t Max_is_defined = false;
1554 Bool_t Start_point_is_defined = false;
1555 Bool_t Period_is_defined = false;
1556 Bool_t Large_gap = false;
1557 Bool_t normal_way = true;
1558 Bool_t small_gap_on_ridge = false;
1559 Double_t t1 = 0;
1560 Double_t t1A = 0;
1561 Int_t sinesize = 0;
1562 Int_t nfi = 0;
1563 for(UInt_t mu = 0;mu<qtime.size();mu++){
1564 //cout<<"Roll["<<mu<<"] = "<<Roll[mu]<<endl;
1565 if(TMath::Abs(Roll[mu])<1. && TMath::Abs(Pitch[mu])<1. && TMath::Abs(q[mu])<limsin){
1566 //cout<<"q["<<mu<<endl<<"] = "<<q[mu]<<endl;
1567 if(mulast!=0 && munow!=0 && munext!=0){mulast=munow;munow=munext;munext=mu;}
1568 if(munext==0 && munow!=0)munext=mu;
1569 if(munow==0 && mulast!=0)munow=mu;
1570 if(mulast==0)mulast=mu;
1571
1572 //cout<<"mulast = "<<mulast<<"\tmunow = "<<munow<<"\tmunext = "<<munext<<endl;
1573 //Int_t ref;
1574 //cin>>ref;
1575 if(TMath::Abs(q[munow])>TMath::Abs(q[mulast]) && TMath::Abs(q[munow])>TMath::Abs(q[munext]) && q[mulast]*q[munext]>0 && qtime[munext]-qtime[mulast]>400)small_gap_on_ridge = true;
1576 if(munext>mulast && (qtime[munext]-qtime[mulast]>=2000 || qtime[munext]-qtime[mulast]<0)){if(Large_gap){normal_way = false;Large_gap = false;}else{Large_gap = true;normal_way = false;}}else normal_way = true;
1577 //if(normal_way)cout<<"Normal_Way"<<endl;
1578 if(Large_gap || small_gap_on_ridge){
1579 //cout<<"Large gap..."<<endl;
1580 //if(small_gap_on_ridge)cout<<"small gap..."<<endl;
1581 //cout<<"q["<<mulast<<"] = "<<q[mulast]<<"\tq["<<munow<<"] = "<<q[munow]<<"\tq["<<munext<<"] = "<<q[munext]<<endl;
1582 //cout<<"qtime["<<mulast<<"] = "<<qtime[mulast]<<"\tqtime["<<munow<<"] = "<<qtime[munow]<<"\tqtime["<<munext<<"] = "<<qtime[munext]<<endl;
1583 increase = false;
1584 decrease = false;
1585 if(nfi>0){
1586 qsine.resize(qsine.size()-1);
1587 sinesize = qsine.size();
1588 //cout<<"nfi was larger then zero"<<endl;
1589 }else{
1590 //cout<<"nfi was not larger then zero :( nfi = "<<nfi<<endl;
1591 //cout<<"qsine.size = "<<qsine.size()<<endl;
1592 if(!Period_is_defined){
1593 //cout<<"Period was defined"<<endl;
1594 if(qsine.size()>1){
1595 qsine[sinesize-1].b = qsine[sinesize-2].b;
1596 qsine[sinesize-1].c = qsine[sinesize-2].c;
1597 }else{
1598 qsine[sinesize-1].b = TMath::Pi()/1591.54;
1599 qsine[sinesize-1].c = qsine[sinesize-1].startPoint;
1600 }
1601 }
1602 if(!Max_is_defined){
1603 //cout<<"Max was already defined"<<endl;
1604 if(qsine.size()>1)qsine[sinesize-1].A = qsine[sinesize-2].A;else qsine[sinesize-1].A = limsin;
1605 }
1606 qsine[sinesize-1].NeedFit = true;
1607 }
1608 qsine[sinesize-1].finishPoint = qtime[munow];
1609 //cout<<"finish point before large gap = "<<qtime[munow]<<endl;
1610 nfi = 0;
1611 Max_is_defined = false;
1612 Start_point_is_defined = false;
1613 Period_is_defined = false;
1614 small_gap_on_ridge = false;
1615 }
1616 //cout<<"Slope "<<increase<<"\t"<<decrease<<endl;
1617 //cout<<"mulast = "<<mulast<<"\tmunow = "<<munow<<"\tmunext = "<<munext<<endl;
1618 if((munext>munow) && (munow>mulast) && normal_way){
1619 if(!increase && !decrease){
1620 //cout<<"Normal way have started"<<endl;
1621 qsine.resize(qsine.size()+1);
1622 sinesize = qsine.size();
1623 qsine[sinesize-1].startPoint=qtime[mulast];
1624 if(q[munext]>q[munow] && q[munow]>q[mulast]) increase = true;
1625 if(q[munext]<q[munow] && q[munow]<q[mulast]) decrease = true;
1626 }
1627 //if(TMath::Abs(q[munow])>TMath::Abs(q[mulast]) && TMath::Abs(q[munow])>TMath::Abs(q[munext]) && TMath::Abs(q[munow])>limsin && qtime[munow]-qtime[mulast]>=400 || qtime[munext]-qtime[munow]>=400){small_gap_on_ridge = true;mu--;continue;}
1628 if(TMath::Abs(q[munow])>TMath::Abs(q[mulast]) && TMath::Abs(q[munow])>TMath::Abs(q[munext]) && TMath::Abs(q[munow])>0.9*limsin && qtime[munow]-qtime[mulast]<400 && qtime[munext]-qtime[munow]<400){
1629 //cout<<"Max point is qtime = "<<qtime[munow]<<"\tq = "<<q[munow]<<endl;
1630 if(q[munow]>q[mulast]){
1631 increase = false;
1632 decrease = true;
1633 }
1634 if(q[munow]<q[mulast]){
1635 increase = true;
1636 decrease = false;
1637 }
1638 if(Max_is_defined && !Start_point_is_defined){
1639 Double_t qPer = qtime[munow]-t1A;
1640 if(qPer>1000){
1641 //cout<<"qsine["<<sinesize-1<<"] = "<<qPer<<" = "<<qtime[munow]<<" - "<<t1A<<"\tlim = "<<limsin<<endl;
1642 qsine[sinesize-1].b=TMath::Pi()/qPer;
1643 if(decrease)qsine[sinesize-1].c=-qsine[sinesize-1].b*t1A;
1644 if(increase)qsine[sinesize-1].c=-qsine[sinesize-1].b*(t1A-qPer);
1645 Period_is_defined = true;
1646 }
1647 }
1648 Max_is_defined = true;
1649 qsine[sinesize-1].A = TMath::Abs(q[munow]);
1650 if(Start_point_is_defined && Period_is_defined){
1651 qsine[sinesize-1].finishPoint = qtime[munow];
1652 nfi++;
1653 qsine[sinesize-1].NeedFit = false;
1654 Max_is_defined = false;
1655 Start_point_is_defined = false;
1656 Period_is_defined = false;
1657 qsine.resize(qsine.size()+1);
1658 sinesize = qsine.size();
1659 qsine[sinesize-1].startPoint = qtime[munow];
1660 }
1661 if(!Start_point_is_defined) t1A=qtime[munow];
1662 }
1663 //if((q[munow]>=0 && q[mulast]<=0) || (q[munow]<=0 && q[mulast]>=0))cout<<"cross zero point diference = "<<qtime[munext] - qtime[mulast]<<"\tqlast = "<<qtime[mulast]<<"\tqnow = "<<qtime[munow]<<"\tqnext = "<<qtime[munext]<<endl;
1664 if(((q[munow]>=0 && q[mulast]<=0) || (q[munow]<=0 && q[mulast]>=0)) && qtime[munow]-qtime[mulast]<2000 && qtime[munext]-qtime[munow]<2000){
1665 Double_t tcrosszero = 0;
1666 //cout<<"cross zero point...qtime = "<<qtime[munow]<<endl;
1667 if(q[munow]==0.) tcrosszero = qtime[munow];else
1668 if(q[mulast]==0.)tcrosszero = qtime[mulast];else{
1669 Double_t k_ = (q[munow]-q[mulast])/(qtime[munow]-qtime[mulast]);
1670 Double_t b_ = q[munow]-k_*qtime[munow];
1671 tcrosszero = -b_/k_;
1672 }
1673 if(Start_point_is_defined){
1674 //cout<<"Start Point allready defined"<<endl;
1675 Double_t qPer = tcrosszero - t1;
1676 qsine[sinesize-1].b = TMath::Pi()/qPer;
1677 //cout<<"qsine["<<sinesize-1<<"].b = "<<TMath::Pi()/qPer<<endl;
1678 Period_is_defined = true;
1679 Float_t x0 = 0;
1680 if(decrease)x0 = t1;
1681 if(increase)x0 = tcrosszero;
1682 qsine[sinesize-1].c = -qsine[sinesize-1].b*x0;
1683 if(Max_is_defined){
1684 //cout<<"Max was previous defined"<<endl;
1685 qsine[sinesize-1].finishPoint = qtime[munow];
1686 nfi++;
1687 qsine[sinesize-1].NeedFit = false;
1688 Max_is_defined = false;
1689 t1 = tcrosszero;
1690 Start_point_is_defined = true;
1691 Period_is_defined = false;
1692 qsine.resize(qsine.size()+1);
1693 sinesize = qsine.size();
1694 qsine[sinesize-1].startPoint = qtime[munow];
1695 }
1696 }else{
1697 t1 = tcrosszero;
1698 Start_point_is_defined = true;
1699 }
1700 }
1701 }
1702 }
1703 }
1704
1705 //cout<<"FINISH SINE INTERPOLATION FUNCTION..."<<endl<<endl;
1706 }

  ViewVC Help
Powered by ViewVC 1.1.23