]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RALICE/AliTimestamp.cxx
Removal of AliDebug in the innermost raw-data decoding loop. Removal of unneeded...
[u/mrichter/AliRoot.git] / RALICE / AliTimestamp.cxx
CommitLineData
3ea81e9c 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16// $Id$
17
18///////////////////////////////////////////////////////////////////////////
19// Class AliTimestamp
20// Handling of timestamps for (astro)particle physics reserach.
21//
22// This class is derived from TTimeStamp and provides additional
23// facilities (e.g. Julian date) which are commonly used in the
24// field of (astro)particle physics.
25//
26// The Julian Date (JD) indicates the number of days since noon (UT) on
27// 01 jan -4712 (i.e. noon 01 jan 4713 BC), being day 0 of the Julian calendar.
28//
29// The Modified Julian Date (MJD) indicates the number of days since midnight
30// (UT) on 17-nov-1858, which corresponds to 2400000.5 days after day 0 of the
31// Julian calendar.
32//
33// The Truncated Julian Date (TJD) corresponds to 2440000.5 days after day 0
34// of the Julian calendar and consequently TJD=MJD-40000.
35// This TJD date indication was used by the Vela and Batse missions in
36// view of Gamma Ray Burst investigations.
37//
a4f7a3a1 38// The Julian Epoch (JE) indicates the fractional elapsed Julian year count
39// since the start of the Gregorian year count.
40// A Julian year is defined to be 365.25 days and starts at 01-jan 12:00:00 UT.
41// As such, the integer part of JE corresponds to the usual Gregorian year count,
42// apart from 01-jan before 12:00:00 UT.
43// So, 01-jan-1965 12:00:00 UT corresponds to JE=1965.0
44//
45// The Besselian Epoch (BE) indicates the fractional elapsed Besselian year count
46// since the start of the Gregorian year count.
47// A Besselian (or tropical) year is defined to be 365.242198781 days.
2fa2e816 48// The date 31-dec-1949 22:09:46.862 UT corresponds to BE=1950.0
a4f7a3a1 49//
50// The Besselian and Julian epochs are used in astronomical catalogs
51// to denote values of time varying observables like e.g. right ascension.
3ea81e9c 52//
53// Because of the fact that the Julian date indicators are all w.r.t. UT
54// they provide an absolute timescale irrespective of timezone or daylight
55// saving time (DST).
56//
a4f7a3a1 57// In view of astronomical observations and positioning it is convenient
58// to have also a UT equivalent related to stellar meridian transitions.
59// This is achieved by the Greenwich Sidereal Time (GST).
60// The GST is defined as the right ascension of the objects passing
61// the Greenwich meridian at 00:00:00 UT.
62// Due to the rotation of the Earth around the Sun, a sidereal day
63// lasts 86164.09 seconds (23h 56m 04.09s) compared to the mean solar
64// day of 86400 seconds (24h).
65// Furthermore, precession of the earth's spin axis results in the fact
66// that the zero point of right ascension (vernal equinox) gradually
67// moves along the celestial equator.
68// In addition, tidal friction and ocean and atmospheric effects will
69// induce seasonal variations in the earth's spin rate and polar motion
70// of the earth's spin axis.
2fa2e816 71// Taking the above effects into account leads to what is called
72// the Greenwich Mean Sidereal Time (GMST).
73// In case also the nutation of the earth's spin axis is taken into
74// account we speak of the Greenwich Apparent Sidereal Time (GAST).
a4f7a3a1 75//
a7dc0627 76// This AliTimestamp facility allows for picosecond precision, in view
77// of time of flight analyses for particle physics experiments.
78// For normal date/time indication the standard nanosecond precision
79// will in general be sufficient.
80// Note that when the fractional JD, MJD and TJD counts are used instead
3ea81e9c 81// of the integer (days,sec,ns) specification, the nanosecond precision
82// may be lost due to computer accuracy w.r.t. floating point operations.
83//
84// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
85// which corresponds to JD=2440587.5 or the start of MJD=40587 or TJD=587.
86// Using the corresponding MJD of this EPOCH allows construction of
87// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input (M/T)JD and time.
88// Obviously this TTimeStamp implementation would prevent usage of values
89// smaller than JD=2440587.5 or MJD=40587 or TJD=587.
0cfe76b5 90// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
91// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 92// However, this AliTimestamp facility provides support for the full range
93// of (M/T)JD values, but the setting of the corresponding TTimeStamp parameters
94// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 95// For these earlier/later (M/T)JD values, the standard TTimeStamp parameters will
3ea81e9c 96// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 97// This implies that for these earlier/later (M/T)JD values the TTimeStamp parameters
3ea81e9c 98// do not match the Julian parameters of AliTimestamp.
99// As such the standard TTimeStamp parameters do not appear on the print output
0cfe76b5 100// when invoking the Date() memberfunction for these earlier/later (M/T)JD values.
3ea81e9c 101//
102// Examples :
103// ==========
104//
105// Note : All TTimeStamp functionality is available as well.
106//
107// AliTimestamp t;
108//
109// t.Date();
110//
111// // Retrieve Julian Date
112// Int_t jd,jsec,jns;
113// t.GetJD(jd,jsec,jns);
114//
115// // Retrieve fractional Truncated Julian Date
116// Double_t tjd=t.GetTJD();
117//
118// // Retrieve fractional Julian Epoch
119// Double_t je=t.GetJE();
120//
121// // Set to a specific Modified Julian Date
122// Int_t mjd=50537;
123// Int_t mjsec=1528;
124// Int_t mjns=185643;
125// t.SetMJD(mjd,mjsec,mjns);
126//
127// t.Date();
128//
95cfc777 129// // Time intervals for e.g. trigger or TOF analysis
130// AliEvent evt;
ee26083f 131// AliTrack* tx=evt.GetTrack(5);
132// AliTimestamp* timex=tx->GetTimestamp();
133// Double_t dt=evt.GetDifference(timex,"ps");
134// AliTimestamp trig((AliTimestamp)evt);
135// trig.Add(0,0,2,173);
136// AliSignal* sx=evt.GetHit(23);
137// AliTimestamp* timex=sx->GetTimestamp();
138// Double_t dt=trig.GetDifference(timex,"ps");
95cfc777 139// Int_t d,s,ns,ps;
ee26083f 140// trig.GetDifference(timex,d,s,ns,ps);
95cfc777 141//
3ea81e9c 142// // Some practical conversion facilities
143// // Note : They don't influence the actual date/time settings
144// // and as such can also be invoked as AliTimestamp::Convert(...) etc...
145// Int_t y=1921;
146// Int_t m=7;
147// Int_t d=21;
148// Int_t hh=15;
149// Int_t mm=23;
150// Int_t ss=47;
151// Int_t ns=811743;
152// Double_t jdate=t.GetJD(y,m,d,hh,mm,ss,ns);
153//
154// Int_t days,secs,nsecs;
155// Double_t date=421.1949327;
156// t.Convert(date,days,secs,nsecs);
157//
158// days=875;
159// secs=23;
160// nsecs=9118483;
161// date=t.Convert(days,secs,nsecs);
162//
163// Double_t mjdate=40563.823744;
164// Double_t epoch=t.GetJE(mjdate,"mjd");
165//
166//--- Author: Nick van Eijndhoven 28-jan-2005 Utrecht University.
167//- Modified: NvE $Date$ Utrecht University.
168///////////////////////////////////////////////////////////////////////////
169
170#include "AliTimestamp.h"
171#include "Riostream.h"
172
173ClassImp(AliTimestamp) // Class implementation to enable ROOT I/O
174
175AliTimestamp::AliTimestamp() : TTimeStamp()
176{
177// Default constructor
178// Creation of an AliTimestamp object and initialisation of parameters.
179// All attributes are initialised to the current date/time as specified
180// in the docs of TTimeStamp.
181
182 FillJulian();
a7dc0627 183 fJps=0;
3ea81e9c 184}
185///////////////////////////////////////////////////////////////////////////
186AliTimestamp::AliTimestamp(TTimeStamp& t) : TTimeStamp(t)
187{
188// Creation of an AliTimestamp object and initialisation of parameters.
189// All attributes are initialised to the values of the input TTimeStamp.
190
191 FillJulian();
a7dc0627 192 fJps=0;
3ea81e9c 193}
194///////////////////////////////////////////////////////////////////////////
195AliTimestamp::~AliTimestamp()
196{
197// Destructor to delete dynamically allocated memory.
198}
199///////////////////////////////////////////////////////////////////////////
200AliTimestamp::AliTimestamp(const AliTimestamp& t) : TTimeStamp(t)
201{
202// Copy constructor
203
204 fMJD=t.fMJD;
205 fJsec=t.fJsec;
206 fJns=t.fJns;
a7dc0627 207 fJps=t.fJps;
3ea81e9c 208 fCalcs=t.fCalcs;
209 fCalcns=t.fCalcns;
210}
211///////////////////////////////////////////////////////////////////////////
e47fe004 212void AliTimestamp::Date(Int_t mode,Double_t offset)
3ea81e9c 213{
214// Print date/time info.
215//
2fa2e816 216// mode = 1 ==> Only the UT yy-mm-dd hh:mm:ss.sss and GMST info is printed
3ea81e9c 217// 2 ==> Only the Julian parameter info is printed
2fa2e816 218// 3 ==> Both the UT, GMST and Julian parameter info is printed
299f01aa 219// -1 ==> Only the UT yy-mm-dd hh:mm:ss.sss and GAST info is printed
220// -3 ==> Both the UT, GAST and Julian parameter info is printed
3ea81e9c 221//
2fa2e816 222// offset : Local time offset from UT (and also GMST) in fractional hours.
e47fe004 223//
224// When an offset value is specified, the corresponding local times
299f01aa 225// LT and LMST (or LAST) are printed as well.
e47fe004 226//
227// The default values are mode=3 and offset=0.
3ea81e9c 228//
229// Note : In case the (M/T)JD falls outside the TTimeStamp range,
145c9890 230// the yy-mm-dd info will be omitted.
a4f7a3a1 231
232 Int_t mjd,mjsec,mjns,mjps;
233 GetMJD(mjd,mjsec,mjns);
234 mjps=GetPs();
3ea81e9c 235
e47fe004 236 TString month[12]={"Jan","Feb","Mar","Apr","May","Jun","Jul","Aug","Sep","Oct","Nov","Dec"};
237 TString day[7]={"Mon","Tue","Wed","Thu","Fri","Sat","Sun"};
238 UInt_t y,m,d,wd;
a4f7a3a1 239 Int_t hh,mm,ss,ns,ps;
299f01aa 240 Double_t gast;
0cfe76b5 241
299f01aa 242 if (abs(mode)==1 || abs(mode)==3)
0cfe76b5 243 {
a4f7a3a1 244 if (mjd>=40587 && (mjd<65442 || (mjd==65442 && mjsec<8047)))
245 {
145c9890 246 GetDate(kTRUE,0,&y,&m,&d);
e47fe004 247 wd=GetDayOfWeek(kTRUE,0);
145c9890 248 cout << " " << day[wd-1].Data() << ", " << setfill('0') << setw(2) << d << " "
249 << setfill(' ') << month[m-1].Data() << " " << y << " ";
a4f7a3a1 250 }
251 else
252 {
145c9890 253 cout << " Time ";
a4f7a3a1 254 }
145c9890 255 GetUT(hh,mm,ss,ns,ps);
256 cout << setfill('0') << setw(2) << hh << ":"
257 << setw(2) << mm << ":" << setw(2) << ss << "."
258 << setw(9) << ns << setw(3) << ps << " (UT) ";
299f01aa 259 if (mode>0)
260 {
261 GetGMST(hh,mm,ss,ns,ps);
262 }
263 else
264 {
265 gast=GetGAST();
266 Convert(gast,hh,mm,ss,ns,ps);
267 }
145c9890 268 cout << setfill('0') << setw(2) << hh << ":"
269 << setw(2) << mm << ":" << setw(2) << ss << "."
299f01aa 270 << setw(9) << ns << setw(3) << ps;
271 if (mode>0)
272 {
273 cout << " (GMST)" << endl;
274 }
275 else
276 {
277 cout << " (GAST)" << endl;
278 }
2fa2e816 279
280 // Local time information
e47fe004 281 if (offset)
282 {
283 // Determine the new date by including the offset
284 AliTimestamp t2(*this);
285 t2.Add(offset);
286 Int_t mjd2,mjsec2,mjns2;
287 t2.GetMJD(mjd2,mjsec2,mjns2);
288 if (mjd2>=40587 && (mjd2<65442 || (mjd2==65442 && mjsec2<8047)))
289 {
290 t2.GetDate(kTRUE,0,&y,&m,&d);
291 wd=t2.GetDayOfWeek(kTRUE,0);
292 cout << " " << day[wd-1].Data() << ", " << setfill('0') << setw(2) << d << " "
293 << setfill(' ') << month[m-1].Data() << " " << y << " ";
294 }
295 else
296 {
297 cout << " Time ";
298 }
299 // Determine the local time by including the offset w.r.t. the original timestamp
2fa2e816 300 Double_t hlt=GetLT(offset);
299f01aa 301 Double_t hlst=0;
302 if (mode>0)
303 {
304 hlst=GetLMST(offset);
305 }
306 else
307 {
308 hlst=GetLAST(offset);
309 }
310 PrintTime(hlt,12); cout << " (LT) "; PrintTime(hlst,12);
311 if (mode>0)
312 {
313 cout << " (LMST)" << endl;
314 }
315 else
316 {
317 cout << " (LAST)" << endl;
318 }
e47fe004 319 }
0cfe76b5 320 }
299f01aa 321 if (abs(mode)==2 || abs(mode)==3)
3ea81e9c 322 {
323 Int_t jd,jsec,jns;
324 GetJD(jd,jsec,jns);
325 Int_t tjd,tjsec,tjns;
326 GetTJD(tjd,tjsec,tjns);
a4f7a3a1 327 cout << " Julian Epoch : " << setprecision(25) << GetJE()
328 << " Besselian Epoch : " << setprecision(25) << GetBE() << endl;
95cfc777 329 cout << " JD : " << jd << " sec : " << jsec << " ns : " << jns << " ps : " << fJps
3ea81e9c 330 << " Fractional : " << setprecision(25) << GetJD() << endl;
95cfc777 331 cout << " MJD : " << mjd << " sec : " << mjsec << " ns : " << mjns << " ps : " << fJps
3ea81e9c 332 << " Fractional : " << setprecision(25) << GetMJD() << endl;
95cfc777 333 cout << " TJD : " << tjd << " sec : " << tjsec << " ns : " << tjns << " ps : " << fJps
3ea81e9c 334 << " Fractional : " << setprecision(25) << GetTJD() << endl;
335 }
336}
337///////////////////////////////////////////////////////////////////////////
338Double_t AliTimestamp::GetJD(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns) const
339{
340// Provide the (fractional) Julian Date (JD) corresponding to the UT date
341// and time in the Gregorian calendar as specified by the input arguments.
342//
343// The input arguments represent the following :
344// y : year in UT (e.g. 1952, 2003 etc...)
345// m : month in UT (1=jan 2=feb etc...)
346// d : day in UT (1-31)
347// hh : elapsed hours in UT (0-23)
348// mm : elapsed minutes in UT (0-59)
349// ss : elapsed seconds in UT (0-59)
350// ns : remaining fractional elapsed second of UT in nanosecond
351//
352// This algorithm is valid for all AD dates in the Gregorian calendar
353// following the recipe of R.W. Sinnott Sky & Telescope 82, (aug. 1991) 183.
354// See also http://scienceworld.wolfram.com/astronomy/JulianDate.html
355//
356// In case of invalid input, a value of -1 is returned.
357//
358// Note :
359// ------
360// This memberfunction only provides the JD corresponding to the
361// UT input arguments. It does NOT set the corresponding Julian parameters
362// for the current AliTimestamp instance.
363// As such the TTimeStamp limitations do NOT apply to this memberfunction.
364// To set the Julian parameters for the current AliTimestamp instance,
365// please use the corresponding SET() memberfunctions of either AliTimestamp
366// or TTimeStamp.
367
368 if (y<0 || m<1 || m>12 || d<1 || d>31) return -1;
369 if (hh<0 || hh>23 || mm<0 || mm>59 || ss<0 || ss>59 || ns<0 || ns>1e9) return -1;
370
371 // The UT daytime in fractional hours
372 Double_t ut=double(hh)+double(mm)/60.+(double(ss)+double(ns)*1.e-9)/3600.;
373
374 Double_t JD=0;
375
376 JD=367*y-int(7*(y+int((m+9)/12))/4)
377 -int(3*(int((y+(m-9)/7)/100)+1)/4)
378 +int(275*m/9)+d+1721028.5+ut/24.;
379
380 return JD;
381}
382///////////////////////////////////////////////////////////////////////////
383Double_t AliTimestamp::GetMJD(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns) const
384{
385// Provide the (fractional) Modified Julian Date corresponding to the UT
386// date and time in the Gregorian calendar as specified by the input arguments.
387//
388// The input arguments represent the following :
389// y : year in UT (e.g. 1952, 2003 etc...)
390// m : month in UT (1=jan 2=feb etc...)
391// d : day in UT (1-31)
392// hh : elapsed hours in UT (0-23)
393// mm : elapsed minutes in UT (0-59)
394// ss : elapsed seconds in UT (0-59)
395// ns : remaining fractional elapsed second of UT in nanosecond
396//
397// This algorithm is valid for all AD dates in the Gregorian calendar
398// following the recipe of R.W. Sinnott Sky & Telescope 82, (aug. 1991) 183.
399// See also http://scienceworld.wolfram.com/astronomy/JulianDate.html
400//
401// In case of invalid input, a value of -1 is returned.
402//
403// Note :
404// ------
405// This memberfunction only provides the MJD corresponding to the
406// UT input arguments. It does NOT set the corresponding Julian parameters
407// for the current AliTimestamp instance.
408// As such the TTimeStamp limitations do NOT apply to this memberfunction.
409// To set the Julian parameters for the current AliTimestamp instance,
410// please use the corresponding SET() memberfunctions of either AliTimestamp
411// or TTimeStamp.
412
413 Double_t JD=GetJD(y,m,d,hh,mm,ss,ns);
414
415 if (JD<0) return JD;
416
417 Double_t MJD=JD-2400000.5;
418
419 return MJD;
420}
421///////////////////////////////////////////////////////////////////////////
422Double_t AliTimestamp::GetTJD(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns) const
423{
424// Provide the (fractional) Truncated Julian Date corresponding to the UT
425// date and time in the Gregorian calendar as specified by the input arguments.
426//
427// The input arguments represent the following :
428// y : year in UT (e.g. 1952, 2003 etc...)
429// m : month in UT (1=jan 2=feb etc...)
430// d : day in UT (1-31)
431// hh : elapsed hours in UT (0-23)
432// mm : elapsed minutes in UT (0-59)
433// ss : elapsed seconds in UT (0-59)
434// ns : remaining fractional elapsed second of UT in nanosecond
435//
436// This algorithm is valid for all AD dates in the Gregorian calendar
437// following the recipe of R.W. Sinnott Sky & Telescope 82, (aug. 1991) 183.
438// See also http://scienceworld.wolfram.com/astronomy/JulianDate.html
439//
440// In case of invalid input, a value of -1 is returned.
441//
442// Note :
443// ------
444// This memberfunction only provides the TJD corresponding to the
445// UT input arguments. It does NOT set the corresponding Julian parameters
446// for the current AliTimestamp instance.
447// As such the TTimeStamp limitations do NOT apply to this memberfunction.
448// To set the Julian parameters for the current AliTimestamp instance,
449// please use the corresponding SET() memberfunctions of either AliTimestamp
450// or TTimeStamp.
451
452 Double_t JD=GetJD(y,m,d,hh,mm,ss,ns);
453
454 if (JD<0) return JD;
455
456 Double_t TJD=JD-2440000.5;
457
458 return TJD;
459}
460///////////////////////////////////////////////////////////////////////////
461Double_t AliTimestamp::GetJE(Double_t date,TString mode) const
462{
463// Provide the Julian Epoch (JE) corresponding to the specified date.
464// The argument "mode" indicates the type of the argument "date".
465//
466// Available modes are :
467// mode = "jd" ==> date represents the Julian Date
468// = "mjd" ==> date represents the Modified Julian Date
469// = "tjd" ==> date represents the Truncated Julian Date
470//
471// The default is mode="jd".
472//
473// In case of invalid input, a value of -99999 is returned.
474//
475// Note :
476// ------
477// This memberfunction only provides the JE corresponding to the
478// input arguments. It does NOT set the corresponding Julian parameters
479// for the current AliTimestamp instance.
480// As such the TTimeStamp limitations do NOT apply to this memberfunction.
481// To set the Julian parameters for the current AliTimestamp instance,
482// please use the corresponding SET() memberfunctions of either AliTimestamp
483// or TTimeStamp.
484
485 if ((mode != "jd") && (mode != "mjd") && (mode != "tjd")) return -99999;
486
487 Double_t jd=date;
488 if (mode=="mjd") jd=date+2400000.5;
489 if (mode=="tjd") jd=date+2440000.5;
490
491 Double_t je=2000.+(jd-2451545.)/365.25;
492
493 return je;
494}
495///////////////////////////////////////////////////////////////////////////
a4f7a3a1 496Double_t AliTimestamp::GetBE(Double_t date,TString mode) const
497{
498// Provide the Besselian Epoch (JE) corresponding to the specified date.
499// The argument "mode" indicates the type of the argument "date".
500//
501// Available modes are :
502// mode = "jd" ==> date represents the Julian Date
503// = "mjd" ==> date represents the Modified Julian Date
504// = "tjd" ==> date represents the Truncated Julian Date
505//
506// The default is mode="jd".
507//
508// In case of invalid input, a value of -99999 is returned.
509//
510// Note :
511// ------
512// This memberfunction only provides the BE corresponding to the
513// input arguments. It does NOT set the corresponding Julian parameters
514// for the current AliTimestamp instance.
515// As such the TTimeStamp limitations do NOT apply to this memberfunction.
516// To set the Julian parameters for the current AliTimestamp instance,
517// please use the corresponding SET() memberfunctions of either AliTimestamp
518// or TTimeStamp.
519
520 if ((mode != "jd") && (mode != "mjd") && (mode != "tjd")) return -99999;
521
522 Double_t jd=date;
523 if (mode=="mjd") jd=date+2400000.5;
524 if (mode=="tjd") jd=date+2440000.5;
525
526 Double_t be=1900.+(jd-2415020.31352)/365.242198781;
527
528 return be;
529}
530///////////////////////////////////////////////////////////////////////////
3ea81e9c 531void AliTimestamp::Convert(Double_t date,Int_t& days,Int_t& secs,Int_t& ns) const
532{
533// Convert date as fractional day count into integer days, secs and ns.
534//
535// Note : Due to computer accuracy the ns value may become inaccurate.
536//
537// The arguments represent the following :
538// date : The input date as fractional day count
539// days : Number of elapsed days
540// secs : Remaining number of elapsed seconds
541// ns : Remaining fractional elapsed second in nanoseconds
542//
543// Note :
544// ------
545// This memberfunction only converts the input date into the corresponding
546// integer parameters. It does NOT set the corresponding Julian parameters
547// for the current AliTimestamp instance.
548// As such the TTimeStamp limitations do NOT apply to this memberfunction.
549// To set the Julian parameters for the current AliTimestamp instance,
550// please use the corresponding SET() memberfunctions of either AliTimestamp
551// or TTimeStamp.
552
553 days=int(date);
554 date=date-double(days);
555 Int_t daysecs=24*3600;
556 date=date*double(daysecs);
557 secs=int(date);
558 date=date-double(secs);
559 ns=int(date*1.e9);
560}
561///////////////////////////////////////////////////////////////////////////
562Double_t AliTimestamp::Convert(Int_t days,Int_t secs,Int_t ns) const
563{
564// Convert date in integer days, secs and ns into fractional day count.
565//
566// Note : Due to computer accuracy the ns precision may be lost.
567//
568// The input arguments represent the following :
569// days : Number of elapsed days
570// secs : Remaining number of elapsed seconds
571// ns : Remaining fractional elapsed second in nanoseconds
572//
573// Note :
574// ------
575// This memberfunction only converts the input integer parameters into the
576// corresponding fractional day count. It does NOT set the corresponding
577// Julian parameters for the current AliTimestamp instance.
578// As such the TTimeStamp limitations do NOT apply to this memberfunction.
579// To set the Julian parameters for the current AliTimestamp instance,
580// please use the corresponding SET() memberfunctions of either AliTimestamp
581// or TTimeStamp.
582
583 Double_t frac=double(secs)+double(ns)*1.e-9;
584 Int_t daysecs=24*3600;
585 frac=frac/double(daysecs);
586 Double_t date=double(days)+frac;
587 return date;
588}
589///////////////////////////////////////////////////////////////////////////
a4f7a3a1 590void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps) const
591{
592// Convert fractional hour count h into hh:mm:ss:ns:ps.
2fa2e816 593// The sign of the input value will be neglected, so h<0 will result in
594// the same output values as h>0.
a4f7a3a1 595//
596// Note : Due to computer accuracy the ps value may become inaccurate.
597//
598// Note :
599// ------
600// This memberfunction only converts the input "h" into the corresponding
601// integer parameters. It does NOT set the corresponding Julian parameters
602// for the current AliTimestamp instance.
603// As such the TTimeStamp limitations do NOT apply to this memberfunction.
604// To set the Julian parameters for the current AliTimestamp instance,
605// please use the corresponding SET() memberfunctions of either AliTimestamp
606// or TTimeStamp.
607
2fa2e816 608 // Neglect sign of h
609 h=fabs(h);
610
a4f7a3a1 611 hh=int(h);
612 h=h-double(hh);
e47fe004 613 h=h*60.;
614 mm=int(h);
615 h=h-double(mm);
616 h=h*60.;
a4f7a3a1 617 ss=int(h);
618 h=h-double(ss);
619 h=h*1.e9;
620 ns=int(h);
621 h=h-double(ns);
622 h=h*1000.;
623 ps=int(h);
624}
625///////////////////////////////////////////////////////////////////////////
2fa2e816 626void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Double_t& ss) const
627{
628// Convert fractional hour count h into hh:mm:ss.s.
629// The sign of the input value will be neglected, so h<0 will result in
630// the same output values as h>0.
631//
632// Notes :
633// -------
634// 1) This memberfunction only converts the input "h" into the corresponding
635// hh:mm:ss.s values. It does NOT set the corresponding Julian parameters
636// for the current AliTimestamp instance.
637// As such the TTimeStamp limitations do NOT apply to this memberfunction.
638// To set the Julian parameters for the current AliTimestamp instance,
639// please use the corresponding SET() memberfunctions of either AliTimestamp
640// or TTimeStamp.
641// 2) This facility can also be used to convert degrees in arcminutes etc...
642
643 // Neglect sign of h
644 h=fabs(h);
645
646 hh=int(h);
647 h=h-double(hh);
648 h=h*60.;
649 mm=int(h);
650 h=h-double(mm);
651 ss=h*60.;
652}
653///////////////////////////////////////////////////////////////////////////
a4f7a3a1 654Double_t AliTimestamp::Convert(Int_t hh,Int_t mm,Int_t ss,Int_t ns,Int_t ps) const
655{
656// Convert hh:mm:ss:ns:ps into fractional hour count.
2fa2e816 657// The sign of the input values will be neglected, so the output value
658// will always correspond to a positive hh:mm:ss:ns:ps specification.
a4f7a3a1 659//
660// Note : Due to computer accuracy the ps precision may be lost.
661//
662// Note :
663// ------
664// This memberfunction only converts the input integer parameters into the
665// corresponding fractional hour count. It does NOT set the corresponding
666// Julian parameters for the current AliTimestamp instance.
667// As such the TTimeStamp limitations do NOT apply to this memberfunction.
668// To set the Julian parameters for the current AliTimestamp instance,
669// please use the corresponding SET() memberfunctions of either AliTimestamp
670// or TTimeStamp.
671
2fa2e816 672 // Neglect the sign of the input values
673 hh=abs(hh);
674 mm=abs(mm);
675 ss=abs(ss);
676 ns=abs(ns);
677 ps=abs(ps);
678
a4f7a3a1 679 Double_t h=hh;
680 h+=double(mm)/60.+(double(ss)+double(ns)*1.e-9+double(ps)*1.e-12)/3600.;
681
682 return h;
683}
684///////////////////////////////////////////////////////////////////////////
2fa2e816 685Double_t AliTimestamp::Convert(Int_t hh,Int_t mm,Double_t ss) const
686{
687// Convert hh:mm:ss.s into fractional hour count.
688// The sign of the input values will be neglected, so the output value
689// will always correspond to a positive hh:mm:ss.s specification.
690//
691// Notes :
692// -------
693// 1) This memberfunction only converts the input hh:mm:ss.s data into the
694// corresponding fractional hour count. It does NOT set the corresponding
695// Julian parameters for the current AliTimestamp instance.
696// As such the TTimeStamp limitations do NOT apply to this memberfunction.
697// To set the Julian parameters for the current AliTimestamp instance,
698// please use the corresponding SET() memberfunctions of either AliTimestamp
699// or TTimeStamp.
700// 2) This facility can also be used to convert ddd:mm:ss.s into fractional degrees.
701
702 // Neglect the sign of the input values
703 hh=abs(hh);
704 mm=abs(mm);
705 ss=fabs(ss);
706
707 Double_t h=hh;
708 h+=double(mm)/60.+ss/3600.;
709
710 return h;
711}
712///////////////////////////////////////////////////////////////////////////
713void AliTimestamp::PrintTime(Double_t h,Int_t ndig) const
714{
715// Print a fractional hour count in hh:mm:ss.ssss format.
716// The range of the printed hour value is : -24 < hh < 24.
717// The argument "ndig" specifies the number of digits for the fractional
718// seconds (e.g. ndig=6 corresponds to microsecond precision).
719// No rounding will be performed, so a second count of 3.473 with ndig=1
720// will appear as 03.4 on the output.
721// Due to computer accuracy, precision on the picosecond level may get lost.
722//
723// The default is ndig=1.
724//
725// Note : The time info is printed without additional spaces or "endline".
726// This allows the print to be included in various composite output formats.
727
728 Int_t hh,mm,ss;
729 ULong64_t sfrac;
730 Double_t s;
731
732 while (h<-24)
733 {
734 h+=24.;
735 }
736 while (h>24)
737 {
738 h-=24.;
739 }
740
741 Convert(h,hh,mm,s);
742 ss=Int_t(s);
743 s-=Double_t(ss);
744 s*=pow(10.,ndig);
745 sfrac=ULong64_t(s);
746
747 if (h<0) cout << "-";
748 cout << setfill('0')
749 << setw(2) << hh << ":" << setw(2) << mm << ":"
750 << setw(2) << ss << "." << setw(ndig) << sfrac;
751}
752///////////////////////////////////////////////////////////////////////////
3ea81e9c 753void AliTimestamp::FillJulian()
754{
755// Calculation and setting of the Julian date/time parameters corresponding
756// to the current TTimeStamp date/time parameters.
757
758 UInt_t y,m,d,hh,mm,ss;
759
760 GetDate(kTRUE,0,&y,&m,&d);
761 GetTime(kTRUE,0,&hh,&mm,&ss);
762 Int_t ns=GetNanoSec();
763
764 Double_t mjd=GetMJD(y,m,d,hh,mm,ss,ns);
765
766 fMJD=int(mjd);
767 fJsec=GetSec()%(24*3600); // Daytime in elapsed seconds
768 fJns=ns; // Remaining fractional elapsed second in nanoseconds
769
770 // Store the TTimeStamp seconds and nanoseconds values
771 // for which this Julian calculation was performed.
772 fCalcs=GetSec();
773 fCalcns=GetNanoSec();
774}
775///////////////////////////////////////////////////////////////////////////
0cfe76b5 776void AliTimestamp::GetMJD(Int_t& mjd,Int_t& sec,Int_t& ns)
3ea81e9c 777{
778// Provide the Modified Julian Date (MJD) and time corresponding to the
779// currently stored AliTimestamp date/time parameters.
780//
781// The returned arguments represent the following :
782// mjd : The modified Julian date.
783// sec : The number of seconds elapsed within the MJD.
784// ns : The remaining fractional number of seconds (in ns) elapsed within the MJD.
785
786 if (fCalcs != GetSec() || fCalcns != GetNanoSec()) FillJulian();
787
788 mjd=fMJD;
789 sec=fJsec;
790 ns=fJns;
791}
792///////////////////////////////////////////////////////////////////////////
793Double_t AliTimestamp::GetMJD()
794{
795// Provide the (fractional) Modified Julian Date (MJD) corresponding to the
796// currently stored AliTimestamp date/time parameters.
797//
798// Due to computer accuracy the ns precision may be lost.
799// It is advised to use the (mjd,sec,ns) getter instead.
800
801 Int_t mjd=0;
802 Int_t sec=0;
803 Int_t ns=0;
804 GetMJD(mjd,sec,ns);
805
806 Double_t date=Convert(mjd,sec,ns);
807
808 return date;
809}
810///////////////////////////////////////////////////////////////////////////
811void AliTimestamp::GetTJD(Int_t& tjd,Int_t& sec, Int_t& ns)
812{
813// Provide the Truncated Julian Date (TJD) and time corresponding to the
814// currently stored AliTimestamp date/time parameters.
815//
816// The returned arguments represent the following :
817// tjd : The modified Julian date.
818// sec : The number of seconds elapsed within the MJD.
819// ns : The remaining fractional number of seconds (in ns) elapsed within the MJD.
820
821 Int_t mjd=0;
822 GetMJD(mjd,sec,ns);
823
824 tjd=mjd-40000;
825}
826///////////////////////////////////////////////////////////////////////////
827Double_t AliTimestamp::GetTJD()
828{
829// Provide the (fractional) Truncated Julian Date (TJD) corresponding to the
830// currently stored AliTimestamp date/time parameters.
831//
832// Due to computer accuracy the ns precision may be lost.
833// It is advised to use the (mjd,sec,ns) getter instead.
834
835 Int_t tjd=0;
836 Int_t sec=0;
837 Int_t ns=0;
838 GetTJD(tjd,sec,ns);
839
840 Double_t date=Convert(tjd,sec,ns);
841
842 return date;
843}
844///////////////////////////////////////////////////////////////////////////
845void AliTimestamp::GetJD(Int_t& jd,Int_t& sec, Int_t& ns)
846{
847// Provide the Julian Date (JD) and time corresponding to the currently
848// stored AliTimestamp date/time parameters.
849//
850// The returned arguments represent the following :
851// jd : The Julian date.
852// sec : The number of seconds elapsed within the JD.
853// ns : The remaining fractional number of seconds (in ns) elapsed within the JD.
854
855 Int_t mjd=0;
856 GetMJD(mjd,sec,ns);
857
858 jd=mjd+2400000;
859 sec+=12*3600;
860 if (sec >= 24*3600)
861 {
862 sec-=24*3600;
863 jd+=1;
864 }
865}
866///////////////////////////////////////////////////////////////////////////
867Double_t AliTimestamp::GetJD()
868{
869// Provide the (fractional) Julian Date (JD) corresponding to the currently
870// stored AliTimestamp date/time parameters.
871//
872// Due to computer accuracy the ns precision may be lost.
873// It is advised to use the (jd,sec,ns) getter instead.
874
875 Int_t jd=0;
876 Int_t sec=0;
877 Int_t ns=0;
878 GetJD(jd,sec,ns);
879
880 Double_t date=Convert(jd,sec,ns);
881
882 return date;
883}
884///////////////////////////////////////////////////////////////////////////
885Double_t AliTimestamp::GetJE()
886{
887// Provide the Julian Epoch (JE) corresponding to the currently stored
888// AliTimestamp date/time parameters.
889
890 Double_t jd=GetJD();
891 Double_t je=GetJE(jd);
892 return je;
893}
894///////////////////////////////////////////////////////////////////////////
a4f7a3a1 895Double_t AliTimestamp::GetBE()
896{
897// Provide the Besselian Epoch (BE) corresponding to the currently stored
898// AliTimestamp date/time parameters.
899
900 Double_t jd=GetJD();
901 Double_t be=GetBE(jd);
902 return be;
903}
904///////////////////////////////////////////////////////////////////////////
a7dc0627 905void AliTimestamp::SetMJD(Int_t mjd,Int_t sec,Int_t ns,Int_t ps)
3ea81e9c 906{
907// Set the Modified Julian Date (MJD) and time and update the TTimeStamp
908// parameters accordingly (if possible).
909//
910// Note :
911// ------
912// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
913// which corresponds to the start of MJD=40587.
914// Using the corresponding MJD of this EPOCH allows construction of
915// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
916// Obviously this TTimeStamp implementation would prevent usage of MJD values
917// smaller than 40587.
0cfe76b5 918// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
919// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 920// However, this AliTimestamp facility provides support for the full range
921// of (M)JD values, but the setting of the corresponding TTimeStamp parameters
922// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 923// For these earlier/later MJD values, the standard TTimeStamp parameters will
3ea81e9c 924// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 925// This implies that for these earlier/later MJD values the TTimeStamp parameters
3ea81e9c 926// do not match the Julian parameters of AliTimestamp.
927//
928// The input arguments represent the following :
929// mjd : The modified Julian date.
930// sec : The number of seconds elapsed within the MJD.
931// ns : The remaining fractional number of seconds (in ns) elapsed within the MJD.
a7dc0627 932// ps : The remaining fractional number of nanoseconds (in ps) elapsed within the MJD.
933//
934// Note : ps=0 is the default value.
3ea81e9c 935
a7dc0627 936 if (sec<0 || sec>=24*3600 || ns<0 || ns>=1e9 || ps<0 || ps>=1000)
3ea81e9c 937 {
938 cout << " *AliTimestamp::SetMJD* Invalid input."
939 << " sec : " << sec << " ns : " << ns << endl;
940 return;
941 }
942
943 fMJD=mjd;
944 fJsec=sec;
945 fJns=ns;
a7dc0627 946 fJps=ps;
3ea81e9c 947
0cfe76b5 948 Int_t epoch=40587; // MJD of the start of the epoch
949 Int_t limit=65442; // MJD of the latest possible TTimeStamp date/time
3ea81e9c 950
0cfe76b5 951 Int_t date,time;
5b917489 952 if (mjd<epoch || mjd>limit || (mjd==limit && sec>=8047))
3ea81e9c 953 {
954 Set(0,kFALSE,0,kFALSE);
0cfe76b5 955 date=GetDate();
956 time=GetTime();
957 Set(date,time,0,kTRUE,0);
3ea81e9c 958 }
959 else
960 {
961 // The elapsed time since start of EPOCH
962 Int_t days=mjd-epoch;
963 UInt_t secs=days*24*3600;
964 secs+=sec;
965 Set(secs,kFALSE,0,kFALSE);
0cfe76b5 966 date=GetDate();
967 time=GetTime();
3ea81e9c 968 Set(date,time,ns,kTRUE,0);
969 }
970
971 // Denote that the Julian and TTimeStamp parameters are synchronised,
972 // even in the case the MJD falls outside the TTimeStamp validity range.
973 // The latter still allows retrieval of Julian parameters for these
974 // earlier times.
975 fCalcs=GetSec();
976 fCalcns=GetNanoSec();
977}
978///////////////////////////////////////////////////////////////////////////
979void AliTimestamp::SetMJD(Double_t mjd)
980{
981// Set the Modified Julian Date (MJD) and time and update the TTimeStamp
982// parameters accordingly (if possible).
983//
984// Note :
985// ------
986// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
987// which corresponds to the start of MJD=40587.
988// Using the corresponding MJD of this EPOCH allows construction of
989// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
990// Obviously this TTimeStamp implementation would prevent usage of MJD values
991// smaller than 40587.
0cfe76b5 992// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
993// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 994// However, this AliTimestamp facility provides support for the full range
995// of (M)JD values, but the setting of the corresponding TTimeStamp parameters
996// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 997// For these earlier/later MJD values, the standard TTimeStamp parameters will
3ea81e9c 998// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 999// This implies that for these earlier/later MJD values the TTimeStamp parameters
3ea81e9c 1000// do not match the Julian parameters of AliTimestamp.
1001//
1002// Due to computer accuracy the ns precision may be lost.
1003// It is advised to use the (mjd,sec,ns) setting instead.
1004//
1005// The input argument represents the following :
1006// mjd : The modified Julian date as fractional day count.
1007
1008 Int_t days=0;
1009 Int_t secs=0;
1010 Int_t ns=0;
1011 Convert(mjd,days,secs,ns);
1012 SetMJD(days,secs,ns);
1013}
1014///////////////////////////////////////////////////////////////////////////
a7dc0627 1015void AliTimestamp::SetJD(Int_t jd,Int_t sec,Int_t ns,Int_t ps)
3ea81e9c 1016{
1017// Set the Julian Date (JD) and time and update the TTimeStamp
1018// parameters accordingly (if possible).
1019//
1020// Note :
1021// ------
1022// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
1023// which corresponds to JD=2440587.5 or the start of MJD=40587.
1024// Using the corresponding MJD of this EPOCH allows construction of
1025// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
1026// Obviously this TTimeStamp implementation would prevent usage of values
1027// smaller than JD=2440587.5.
0cfe76b5 1028// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
1029// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 1030// However, this AliTimestamp facility provides support for the full range
1031// of (M)JD values, but the setting of the corresponding TTimeStamp parameters
1032// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 1033// For these earlier/later JD values, the standard TTimeStamp parameters will
3ea81e9c 1034// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 1035// This implies that for these earlier/later (M)JD values the TTimeStamp parameters
3ea81e9c 1036// do not match the Julian parameters of AliTimestamp.
1037//
1038// The input arguments represent the following :
1039// jd : The Julian date.
1040// sec : The number of seconds elapsed within the JD.
1041// ns : The remaining fractional number of seconds (in ns) elapsed within the JD.
a7dc0627 1042// ps : The remaining fractional number of nanoseconds (in ps) elapsed within the JD.
1043//
1044// Note : ps=0 is the default value.
3ea81e9c 1045
1046 Int_t mjd=jd-2400000;
1047 sec-=12*3600;
1048 if (sec<0)
1049 {
1050 sec+=24*3600;
1051 mjd-=1;
1052 }
1053
a7dc0627 1054 SetMJD(mjd,sec,ns,ps);
3ea81e9c 1055}
1056///////////////////////////////////////////////////////////////////////////
1057void AliTimestamp::SetJD(Double_t jd)
1058{
1059// Set the Julian Date (JD) and time and update the TTimeStamp
1060// parameters accordingly (if possible).
1061//
1062// Note :
1063// ------
1064// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
1065// which corresponds to JD=2440587.5 or the start of MJD=40587.
1066// Using the corresponding MJD of this EPOCH allows construction of
1067// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
1068// Obviously this TTimeStamp implementation would prevent usage of values
1069// smaller than JD=2440587.5.
0cfe76b5 1070// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
1071// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 1072// However, this AliTimestamp facility provides support for the full range
1073// of (M)JD values, but the setting of the corresponding TTimeStamp parameters
1074// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 1075// For these earlier/later JD values, the standard TTimeStamp parameters will
3ea81e9c 1076// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 1077// This implies that for these earlier/later (M)JD values the TTimeStamp parameters
3ea81e9c 1078// do not match the Julian parameters of AliTimestamp.
1079//
1080// Due to computer accuracy the ns precision may be lost.
1081// It is advised to use the (jd,sec,ns) setting instead.
1082//
1083// The input argument represents the following :
1084// jd : The Julian date as fractional day count.
1085
1086 Int_t days=0;
1087 Int_t secs=0;
1088 Int_t ns=0;
1089 Convert(jd,days,secs,ns);
1090
1091 SetJD(days,secs,ns);
1092}
1093///////////////////////////////////////////////////////////////////////////
a7dc0627 1094void AliTimestamp::SetTJD(Int_t tjd,Int_t sec,Int_t ns,Int_t ps)
3ea81e9c 1095{
1096// Set the Truncated Julian Date (TJD) and time and update the TTimeStamp
1097// parameters accordingly (if possible).
1098//
1099// Note :
1100// ------
1101// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
1102// which corresponds to JD=2440587.5 or the start of TJD=587.
1103// Using the corresponding MJD of this EPOCH allows construction of
1104// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
1105// Obviously this TTimeStamp implementation would prevent usage of values
1106// smaller than TJD=587.
0cfe76b5 1107// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
1108// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 1109// However, this AliTimestamp facility provides support for the full range
1110// of (T)JD values, but the setting of the corresponding TTimeStamp parameters
1111// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 1112// For these earlier/later JD values, the standard TTimeStamp parameters will
3ea81e9c 1113// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 1114// This implies that for these earlier/later (T)JD values the TTimeStamp parameters
3ea81e9c 1115// do not match the Julian parameters of AliTimestamp.
1116//
1117// The input arguments represent the following :
1118// tjd : The Truncated Julian date.
1119// sec : The number of seconds elapsed within the JD.
1120// ns : The remaining fractional number of seconds (in ns) elapsed within the JD.
a7dc0627 1121// ps : The remaining fractional number of nanoseconds (in ps) elapsed within the JD.
1122//
1123// Note : ps=0 is the default value.
3ea81e9c 1124
1125 Int_t mjd=tjd+40000;
1126
5481c137 1127 SetMJD(mjd,sec,ns,ps);
3ea81e9c 1128}
1129///////////////////////////////////////////////////////////////////////////
1130void AliTimestamp::SetTJD(Double_t tjd)
1131{
1132// Set the Truncated Julian Date (TJD) and time and update the TTimeStamp
1133// parameters accordingly (if possible).
1134//
1135// Note :
1136// ------
1137// The TTimeStamp EPOCH starts at 01-jan-1970 00:00:00 UTC
1138// which corresponds to JD=2440587.5 or the start of TJD=587.
1139// Using the corresponding MJD of this EPOCH allows construction of
1140// the yy-mm-dd hh:mm:ss:ns TTimeStamp from a given input MJD and time.
1141// Obviously this TTimeStamp implementation would prevent usage of values
1142// smaller than TJD=587.
0cfe76b5 1143// Furthermore, due to a limitation on the "seconds since the EPOCH start" count
1144// in TTimeStamp, the latest accessible date/time is 19-jan-2038 02:14:08 UTC.
3ea81e9c 1145// However, this AliTimestamp facility provides support for the full range
1146// of (T)JD values, but the setting of the corresponding TTimeStamp parameters
1147// is restricted to the values allowed by the TTimeStamp implementation.
0cfe76b5 1148// For these earlier/later JD values, the standard TTimeStamp parameters will
3ea81e9c 1149// be set corresponding to the start of the TTimeStamp EPOCH.
0cfe76b5 1150// This implies that for these earlier/later (T)JD values the TTimeStamp parameters
3ea81e9c 1151// do not match the Julian parameters of AliTimestamp.
1152//
1153// Due to computer accuracy the ns precision may be lost.
1154// It is advised to use the (jd,sec,ns) setting instead.
1155//
1156// The input argument represents the following :
1157// tjd : The Truncated Julian date as fractional day count.
1158
1159 Int_t days=0;
1160 Int_t secs=0;
1161 Int_t ns=0;
1162 Convert(tjd,days,secs,ns);
1163
1164 SetTJD(days,secs,ns);
1165}
1166///////////////////////////////////////////////////////////////////////////
95cfc777 1167void AliTimestamp::SetNs(Int_t ns)
1168{
1169// Set the remaining fractional number of seconds in nanosecond precision.
1170// Notes :
1171// -------
1172// 1) The allowed range for the argument "ns" is [0,99999999].
1173// Outside that range no action is performed.
1174// 2) The ns fraction can also be entered directly via SetMJD() etc...
1175// 3) For additional accuracy see SetPs().
1176
1177 if (ns>=0 && ns<=99999999) fJns=ns;
1178}
1179///////////////////////////////////////////////////////////////////////////
1180Int_t AliTimestamp::GetNs() const
1181{
1182// Provide the remaining fractional number of seconds in nanosecond precision.
1183// This function allows trigger/timing analysis for (astro)particle physics
1184// experiments.
1185// Note : For additional accuracy see also GetPs().
1186
1187 return fJns;
1188}
1189///////////////////////////////////////////////////////////////////////////
1190void AliTimestamp::SetPs(Int_t ps)
1191{
1192// Set the remaining fractional number of nanoseconds in picoseconds.
1193// Notes :
1194// -------
1195// 1) The allowed range for the argument "ps" is [0,999].
1196// Outside that range no action is performed.
1197// 2) The ps fraction can also be entered directly via SetMJD() etc...
1198
1199 if (ps>=0 && ps<=999) fJps=ps;
1200}
1201///////////////////////////////////////////////////////////////////////////
1202Int_t AliTimestamp::GetPs() const
a7dc0627 1203{
1204// Provide remaining fractional number of nanoseconds in picoseconds.
95cfc777 1205// This function allows time of flight analysis for particle physics
a7dc0627 1206// experiments.
1207
1208 return fJps;
1209}
1210///////////////////////////////////////////////////////////////////////////
95cfc777 1211void AliTimestamp::Add(Int_t d,Int_t s,Int_t ns,Int_t ps)
1212{
1213// Add (or subtract) a certain time difference to the current timestamp.
ee26083f 1214// Subtraction can be achieved by entering negative values as input arguments.
95cfc777 1215//
25eefd00 1216// The time difference is entered via the following input arguments :
1217//
95cfc777 1218// d : elapsed number of days
25eefd00 1219// s : (remaining) elapsed number of seconds
1220// ns : (remaining) elapsed number of nanoseconds
1221// ps : (remaining) elapsed number of picoseconds
1222//
1223// The specified d, s, ns and ps values will be used in an additive
1224// way to determine the time difference.
1225// So, specification of d=1, s=100, ns=0, ps=0 will result in the
1226// same time difference addition as d=0, s=24*3600+100, ns=0, ps=0.
1227// However, by making use of the latter the user should take care
1228// of possible integer overflow problems in the input arguments,
1229// which obviously will provide incorrect results.
95cfc777 1230//
1231// Note : ps=0 is the default value.
1232
ee26083f 1233 Int_t days=0;
1234 Int_t secs=0;
1235 Int_t nsec=0;
1236 // Use Get functions to ensure updated Julian parameters.
1237 GetMJD(days,secs,nsec);
1238 Int_t psec=GetPs();
95cfc777 1239
25eefd00 1240 psec+=ps%1000;
1241 nsec+=ps/1000;
1242 while (psec<0)
95cfc777 1243 {
1244 nsec-=1;
1245 psec+=1000;
1246 }
25eefd00 1247 while (psec>999)
95cfc777 1248 {
1249 nsec+=1;
1250 psec-=1000;
1251 }
1252
25eefd00 1253 nsec+=ns%1000000000;
1254 secs+=ns/1000000000;
1255 while (nsec<0)
95cfc777 1256 {
1257 secs-=1;
1258 nsec+=1000000000;
1259 }
25eefd00 1260 while (nsec>999999999)
95cfc777 1261 {
1262 secs+=1;
1263 nsec-=1000000000;
1264 }
1265
25eefd00 1266 secs+=s%(24*3600);
1267 days+=s/(24*3600);
1268 while (secs<0)
95cfc777 1269 {
1270 days-=1;
1271 secs+=24*3600;
1272 }
25eefd00 1273 while (secs>=24*3600)
95cfc777 1274 {
1275 days+=1;
1276 secs-=24*3600;
1277 }
1278
1279 days+=d;
1280
6a7b0c73 1281 SetMJD(days,secs,nsec,psec);
95cfc777 1282}
1283///////////////////////////////////////////////////////////////////////////
e47fe004 1284void AliTimestamp::Add(Double_t hours)
1285{
1286// Add (or subtract) a certain time difference to the current timestamp.
1287// The time difference is specified as a (fractional) number of hours.
1288// Subtraction can be achieved by entering a negative value as input argument.
1289
1290 Int_t d,s,ns,ps;
1291 Double_t h=fabs(hours);
1292 d=int(h/24.);
1293 h-=double(d)*24.;
1294 h*=3600.;
1295 s=int(h);
1296 h-=double(s);
1297 h*=1.e9;
1298 ns=int(h);
1299 h-=double(ns);
1300 ps=int(h*1000.);
1301 if (hours>0) Add(d,s,ns,ps);
1302 if (hours<0) Add(-d,-s,-ns,-ps);
1303}
1304///////////////////////////////////////////////////////////////////////////
ee26083f 1305Int_t AliTimestamp::GetDifference(AliTimestamp* t,Int_t& d,Int_t& s,Int_t& ns,Int_t& ps)
a7dc0627 1306{
1307// Provide the time difference w.r.t the AliTimestamp specified on the input.
1308// This memberfunction supports both very small (i.e. time of flight analysis
1309// for particle physics experiments) and very long (i.e. investigation of
1310// astrophysical phenomena) timescales.
1311//
1312// The time difference is returned via the following output arguments :
1313// d : elapsed number of days
1314// s : remaining elapsed number of seconds
1315// ns : remaining elapsed number of nanoseconds
1316// ps : remaining elapsed number of picoseconds
1317//
95cfc777 1318// Note :
1319// ------
1320// The calculated time difference is the absolute value of the time interval.
1321// This implies that the values of d, s, ns and ps are always positive or zero.
1322//
a7dc0627 1323// The integer return argument indicates whether the AliTimestamp specified
1324// on the input argument occurred earlier (-1), simultaneously (0) or later (1).
1325
ee26083f 1326 if (!t) return 0;
1327
1328 // Ensure updated Julian parameters for this AliTimestamp instance
1329 if (fCalcs != GetSec() || fCalcns != GetNanoSec()) FillJulian();
1330
1331 // Use Get functions to ensure updated Julian parameters.
1332 t->GetMJD(d,s,ns);
1333 ps=t->GetPs();
1334
1335 d-=fMJD;
1336 s-=fJsec;
1337 ns-=fJns;
1338 ps-=fJps;
a7dc0627 1339
1340 if (!d && !s && !ns && !ps) return 0;
1341
1342 Int_t sign=0;
1343
95cfc777 1344 if (d>0) sign=1;
1345 if (d<0) sign=-1;
1346
1347 if (!sign && s>0) sign=1;
1348 if (!sign && s<0) sign=-1;
a7dc0627 1349
1350 if (!sign && ns>0) sign=1;
1351 if (!sign && ns<0) sign=-1;
1352
1353 if (!sign && ps>0) sign=1;
1354 if (!sign && ps<0) sign=-1;
1355
1356 // In case the input stamp was earlier, take the reverse difference
1357 // to simplify the algebra.
1358 if (sign<0)
1359 {
95cfc777 1360 d=-d;
a7dc0627 1361 s=-s;
1362 ns=-ns;
1363 ps=-ps;
1364 }
1365
1366 // Here we always have a positive time difference
95cfc777 1367 // and can now unambiguously correct for other negative values.
a7dc0627 1368 if (ps<0)
1369 {
1370 ns-=1;
1371 ps+=1000;
1372 }
1373
1374 if (ns<0)
1375 {
1376 s-=1;
95cfc777 1377 ns+=1000000000;
a7dc0627 1378 }
1379
95cfc777 1380 if (s<0)
1381 {
1382 d-=1;
1383 s+=24*3600;
1384 }
a7dc0627 1385
1386 return sign;
1387}
1388///////////////////////////////////////////////////////////////////////////
ee26083f 1389Int_t AliTimestamp::GetDifference(AliTimestamp& t,Int_t& d,Int_t& s,Int_t& ns,Int_t& ps)
1390{
1391// Provide the time difference w.r.t the AliTimestamp specified on the input.
1392// This memberfunction supports both very small (i.e. time of flight analysis
1393// for particle physics experiments) and very long (i.e. investigation of
1394// astrophysical phenomena) timescales.
1395//
1396// The time difference is returned via the following output arguments :
1397// d : elapsed number of days
1398// s : remaining elapsed number of seconds
1399// ns : remaining elapsed number of nanoseconds
1400// ps : remaining elapsed number of picoseconds
1401//
1402// Note :
1403// ------
1404// The calculated time difference is the absolute value of the time interval.
1405// This implies that the values of d, s, ns and ps are always positive or zero.
1406//
1407// The integer return argument indicates whether the AliTimestamp specified
1408// on the input argument occurred earlier (-1), simultaneously (0) or later (1).
1409
1410 return GetDifference(&t,d,s,ns,ps);
1411}
1412///////////////////////////////////////////////////////////////////////////
1413Double_t AliTimestamp::GetDifference(AliTimestamp* t,TString u,Int_t mode)
95cfc777 1414{
1415// Provide the time difference w.r.t the AliTimestamp specified on the input
1416// argument in the units as specified by the TString argument.
1417// A positive return value means that the AliTimestamp specified on the input
1418// argument occurred later, whereas a negative return value indicates an
1419// earlier occurence.
1420//
1421// The units may be specified as :
1422// u = "d" ==> Time difference returned as (fractional) day count
1423// "s" ==> Time difference returned as (fractional) second count
1424// "ns" ==> Time difference returned as (fractional) nanosecond count
1425// "ps" ==> Time difference returned as picosecond count
1426//
1427// It may be clear that for a time difference of several days, the picosecond
1428// and even the nanosecond accuracy may be lost.
1429// To cope with this, the "mode" argument has been introduced to allow
1430// timestamp comparison on only the specified units.
1431//
1432// The following operation modes are supported :
1433// mode = 1 : Full time difference is returned in specified units
1434// 2 : Time difference is returned in specified units by
1435// neglecting the elapsed time for the larger units than the
1436// ones specified.
1437// 3 : Time difference is returned in specified units by only
1438// comparing the timestamps on the level of the specified units.
1439//
1440// Example :
1441// ---------
1442// AliTimestamp t1; // Corresponding to days=3, secs=501, ns=31, ps=7
1443// AliTimestamp t2; // Corresponding to days=5, secs=535, ns=12, ps=15
1444//
1445// The statement : Double_t val=t1.GetDifference(t2,....)
1446// would return the following values :
1447// val=(2*24*3600)+34-(19*1e-9)+(8*1e-12) for u="s" and mode=1
1448// val=34-(19*1e-9)+(8*1e-12) for u="s" and mode=2
1449// val=34 for u="s" and mode=3
1450// val=-19 for u="ns" and mode=3
1451//
1452// The default is mode=1.
1453
ee26083f 1454 if (!t || mode<1 || mode>3) return 0;
95cfc777 1455
1456 Double_t dt=0;
1457
ee26083f 1458 // Ensure updated Julian parameters for this AliTimestamp instance
1459 if (fCalcs != GetSec() || fCalcns != GetNanoSec()) FillJulian();
1460
1461 Int_t dd=0;
1462 Int_t ds=0;
1463 Int_t dns=0;
1464 Int_t dps=0;
1465
1466 // Use Get functions to ensure updated Julian parameters.
1467 t->GetMJD(dd,ds,dns);
1468 dps=t->GetPs();
1469
1470 dd-=fMJD;
1471 ds-=fJsec;
1472 dns-=fJns;
1473 dps-=fJps;
95cfc777 1474
1475 // Time difference for the specified units only
1476 if (mode==3)
1477 {
1478 if (u=="d") dt=dd;
1479 if (u=="s") dt=ds;
1480 if (u=="ns") dt=dns;
1481 if (u=="ps") dt=dps;
1482 return dt;
1483 }
1484
1485 // Suppress elapsed time for the larger units than specified
1486 if (mode==2)
1487 {
1488 if (u=="s") dd=0;
1489 if (u=="ns")
1490 {
1491 dd=0;
1492 ds=0;
1493 }
1494 if (u=="ps")
1495 {
1496 dd=0;
1497 ds=0;
1498 dns=0;
1499 }
1500 }
1501
1502 // Compute the time difference as requested
1503 if (u=="s" || u=="d")
1504 {
1505 // The time difference in (fractional) seconds
1506 dt=double(dd*24*3600+ds)+(double(dns)*1e-9)+(double(dps)*1e-12);
1507 if (u=="d") dt=dt/double(24*3600);
1508 }
1509 if (u=="ns") dt=(double(dd*24*3600+ds)*1e9)+double(dns)+(double(dps)*1e-3);
1510 if (u=="ps") dt=(double(dd*24*3600+ds)*1e12)+(double(dns)*1e3)+double(dps);
1511
1512 return dt;
1513}
1514///////////////////////////////////////////////////////////////////////////
ee26083f 1515Double_t AliTimestamp::GetDifference(AliTimestamp& t,TString u,Int_t mode)
1516{
1517// Provide the time difference w.r.t the AliTimestamp specified on the input
1518// argument in the units as specified by the TString argument.
1519// A positive return value means that the AliTimestamp specified on the input
1520// argument occurred later, whereas a negative return value indicates an
1521// earlier occurence.
1522//
1523// The units may be specified as :
1524// u = "d" ==> Time difference returned as (fractional) day count
1525// "s" ==> Time difference returned as (fractional) second count
1526// "ns" ==> Time difference returned as (fractional) nanosecond count
1527// "ps" ==> Time difference returned as picosecond count
1528//
1529// It may be clear that for a time difference of several days, the picosecond
1530// and even the nanosecond accuracy may be lost.
1531// To cope with this, the "mode" argument has been introduced to allow
1532// timestamp comparison on only the specified units.
1533//
1534// The following operation modes are supported :
1535// mode = 1 : Full time difference is returned in specified units
1536// 2 : Time difference is returned in specified units by
1537// neglecting the elapsed time for the larger units than the
1538// ones specified.
1539// 3 : Time difference is returned in specified units by only
1540// comparing the timestamps on the level of the specified units.
1541//
1542// Example :
1543// ---------
1544// AliTimestamp t1; // Corresponding to days=3, secs=501, ns=31, ps=7
1545// AliTimestamp t2; // Corresponding to days=5, secs=535, ns=12, ps=15
1546//
1547// The statement : Double_t val=t1.GetDifference(t2,....)
1548// would return the following values :
1549// val=(2*24*3600)+34-(19*1e-9)+(8*1e-12) for u="s" and mode=1
1550// val=34-(19*1e-9)+(8*1e-12) for u="s" and mode=2
1551// val=34 for u="s" and mode=3
1552// val=-19 for u="ns" and mode=3
1553//
1554// The default is mode=1.
1555
1556 return GetDifference(&t,u,mode);
1557}
1558///////////////////////////////////////////////////////////////////////////
0cfe76b5 1559void AliTimestamp::SetUT(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns,Int_t ps)
1560{
1561// Set the AliTimestamp parameters corresponding to the UT date and time
1562// in the Gregorian calendar as specified by the input arguments.
1563// This facility is exact upto picosecond precision and as such is
1564// for scientific observations preferable above the corresponding
1565// Set function(s) of TTimestamp.
1566// The latter has a random spread in the sub-second part, which
1567// might be of use in generating distinguishable timestamps while
1568// still keeping second precision.
1569//
1570// The input arguments represent the following :
1571// y : year in UT (e.g. 1952, 2003 etc...)
1572// m : month in UT (1=jan 2=feb etc...)
1573// d : day in UT (1-31)
1574// hh : elapsed hours in UT (0-23)
1575// mm : elapsed minutes in UT (0-59)
1576// ss : elapsed seconds in UT (0-59)
1577// ns : remaining fractional elapsed second of UT in nanosecond
1578// ps : remaining fractional elapsed nanosecond of UT in picosecond
1579//
1580// Note : ns=0 and ps=0 are the default values.
1581//
1582// This facility first determines the elapsed days, seconds etc...
2fa2e816 1583// since the beginning of the specified UT year on basis of the
0cfe76b5 1584// input arguments. Subsequently it invokes the SetUT memberfunction
1585// for the elapsed timespan.
1586// As such this facility is valid for all AD dates in the Gregorian
1587// calendar with picosecond precision.
1588
1589 Int_t day=GetDayOfYear(d,m,y);
1590 Int_t secs=hh*3600+mm*60+ss;
1591 SetUT(y,day-1,secs,ns,ps);
1592}
1593///////////////////////////////////////////////////////////////////////////
1594void AliTimestamp::SetUT(Int_t y,Int_t d,Int_t s,Int_t ns,Int_t ps)
1595{
1596// Set the AliTimestamp parameters corresponding to the specified elapsed
1597// timespan since the beginning of the new UT year.
1598// This facility is exact upto picosecond precision and as such is
1599// for scientific observations preferable above the corresponding
1600// Set function(s) of TTimestamp.
1601// The latter has a random spread in the sub-second part, which
1602// might be of use in generating distinguishable timestamps while
1603// still keeping second precision.
1604//
1605// The UT year and elapsed time span is entered via the following input arguments :
1606//
1607// y : year in UT (e.g. 1952, 2003 etc...)
1608// d : elapsed number of days
1609// s : (remaining) elapsed number of seconds
1610// ns : (remaining) elapsed number of nanoseconds
1611// ps : (remaining) elapsed number of picoseconds
1612//
1613// The specified d, s, ns and ps values will be used in an additive
1614// way to determine the elapsed timespan.
1615// So, specification of d=1, s=100, ns=0, ps=0 will result in the
1616// same elapsed time span as d=0, s=24*3600+100, ns=0, ps=0.
1617// However, by making use of the latter the user should take care
1618// of possible integer overflow problems in the input arguments,
1619// which obviously will provide incorrect results.
1620//
1621// Note : ns=0 and ps=0 are the default values.
1622//
1623// This facility first sets the (M)JD corresponding to the start (01-jan 00:00:00)
1624// of the specified UT year following the recipe of R.W. Sinnott
1625// Sky & Telescope 82, (aug. 1991) 183.
1626// Subsequently the day and (sub)second parts are added to the AliTimestamp.
1627// As such this facility is valid for all AD dates in the Gregorian calendar.
1628
1629 Double_t jd=GetJD(y,1,1,0,0,0,0);
1630 SetJD(jd);
1631
1632 Int_t mjd,sec,nsec;
1633 GetMJD(mjd,sec,nsec);
1634 SetMJD(mjd,0,0,0);
1635 Add(d,s,ns,ps);
1636}
1637///////////////////////////////////////////////////////////////////////////
a4f7a3a1 1638void AliTimestamp::GetUT(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
1639{
1640// Provide the corrresponding UT as hh:mm:ss:ns:ps.
1641// This facility is based on the MJD, so the TTimeStamp limitations
1642// do not apply here.
1643
1644 Int_t mjd,sec,nsec,psec;
1645
1646 GetMJD(mjd,sec,nsec);
1647 psec=GetPs();
1648
1649 hh=sec/3600;
1650 sec=sec%3600;
1651 mm=sec/60;
1652 ss=sec%60;
1653 ns=nsec;
1654 ps=psec;
1655}
1656///////////////////////////////////////////////////////////////////////////
1657Double_t AliTimestamp::GetUT()
1658{
1659// Provide the corrresponding UT in fractional hours.
1660// This facility is based on the MJD, so the TTimeStamp limitations
1661// do not apply here.
1662
1663 Int_t hh,mm,ss,ns,ps;
1664
1665 GetUT(hh,mm,ss,ns,ps);
1666
1667 Double_t ut=Convert(hh,mm,ss,ns,ps);
1668
1669 return ut;
1670}
1671///////////////////////////////////////////////////////////////////////////
2fa2e816 1672void AliTimestamp::GetGMST(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
a4f7a3a1 1673{
2fa2e816 1674// Provide the corrresponding Greenwich Mean Sideral Time (GMST).
a4f7a3a1 1675// The algorithm used is the one described at p. 83 of the book
1676// Astronomy Methods by Hale Bradt.
1677// This facility is based on the MJD, so the TTimeStamp limitations
1678// do not apply here.
1679
1680 Int_t mjd,sec,nsec,psec;
1681
1682 // The current UT based timestamp data
1683 GetMJD(mjd,sec,nsec);
1684 psec=fJps;
1685
1686 // The basis for the daily corrections in units of Julian centuries w.r.t. J2000.
1687 // Note : Epoch J2000 starts at 01-jan-2000 12:00:00 UT.
1688 Double_t tau=(GetJD()-2451545.)/36525.;
1689
1690 // Syncronise sidereal time with current timestamp
1691 AliTimestamp sid;
1692 sid.SetMJD(mjd,sec,nsec,psec);
1693
2fa2e816 1694 // Add offset for GMST start value defined as 06:41:50.54841 at 01-jan 00:00:00 UT
a4f7a3a1 1695 sec=6*3600+41*60+50;
1696 nsec=548410000;
1697 psec=0;
1698 sid.Add(0,sec,nsec,psec);
1699
1700 // Daily correction for precession and polar motion
1701 Double_t addsec=8640184.812866*tau+0.093104*pow(tau,2)-6.2e-6*pow(tau,3);
1702 sec=int(addsec);
1703 addsec-=double(sec);
1704 nsec=int(addsec*1.e9);
1705 addsec-=double(nsec)*1.e-9;
1706 psec=int(addsec*1.e12);
1707 sid.Add(0,sec,nsec,psec);
1708
1709 sid.GetMJD(mjd,sec,nsec);
1710 psec=sid.GetPs();
1711
1712 hh=sec/3600;
1713 sec=sec%3600;
1714 mm=sec/60;
1715 ss=sec%60;
1716 ns=nsec;
1717 ps=psec;
1718}
1719///////////////////////////////////////////////////////////////////////////
2fa2e816 1720Double_t AliTimestamp::GetGMST()
a4f7a3a1 1721{
2fa2e816 1722// Provide the corrresponding Greenwich Mean Sideral Time (GMST)
a4f7a3a1 1723// in fractional hours.
1724// This facility is based on the MJD, so the TTimeStamp limitations
1725// do not apply here.
1726
1727 Int_t hh,mm,ss,ns,ps;
1728
2fa2e816 1729 GetGMST(hh,mm,ss,ns,ps);
a4f7a3a1 1730
1731 Double_t gst=Convert(hh,mm,ss,ns,ps);
1732
1733 return gst;
1734}
1735///////////////////////////////////////////////////////////////////////////
2fa2e816 1736Double_t AliTimestamp::GetGAST()
1737{
1738// Provide the corrresponding Greenwich Apparent Sideral Time (GAST)
1739// in fractional hours.
1740// In case a hh:mm:ss.sss format is needed, please invoke the Convert()
1741// memberfunction for conversion of the provided fractional hour value.
1742//
1743// The GAST is the GMST corrected for the shift of the vernal equinox
1744// due to nutation. The right ascension component of the nutation correction
1745// of the vernal equinox is called the "equation of the equinoxes".
1746// So we have :
1747//
1748// GAST = GMST + (equation of the equinoxes)
1749//
1750// With the right ascension and declination values to be (0,0) for the
1751// vernal equinox, we can workout the shift in right ascension of the vernal
1752// equinox due to nutation.
1753// The nutation model used is the new one as documented in :
1754// "The IAU Resolutions on Astronomical Reference Systems,
1755// Time Scales and Earth Rotation Models".
1756// This document is freely available as Circular 179 (2005) of the
1757// United States Naval Observatory (USNO).
1758// (See : http://aa.usno.navy.mil/publications/docs).
1759// The new expression for the equation of the equinoxes is based on a series
1760// expansion and is the most accurate one known to date.
1761// The components are documented on p.17 of the USNO Circular 179.
5b917489 1762//
1763// The change (dpsi) in a star's ecliptic longitude is evaluated using
1764// the formulas from the book "Astronomical Algorithms" by Jean Meeus.
2fa2e816 1765//
1766// Since GMST is based on the MJD, the TTimeStamp limitations
1767// do not apply here.
1768
1769 Double_t pi=acos(-1.);
1770
1771 Double_t gmst=GetGMST();
1772
1773 Double_t days; // Time difference in fractional Julian days w.r.t. the start of J2000.
1774 Double_t t; // Time difference in fractional Julian centuries w.r.t. the start of J2000.
1775 Double_t epsilon; // Mean obliquity of the ecliptic
5b917489 1776 Double_t lsun; // Mean longitude of the Sun
1777 Double_t lmoon; // Mean longitude of the Moon
2fa2e816 1778 Double_t omega; // Mean longitude of the Moon's mean ascending mode
1779 Double_t f; // Mean argument of latitude of the moon
1780 Double_t d; // Mean elongation of the Moon from the Sun
1781
1782 days=GetJD()-2451545.0;
1783 t=days/36525.;
1784
5b917489 1785 // Values in degrees
1786 lsun=280.4665+36000.7698*t;
1787 lmoon=218.3165+481267.8813*t;
1788
2fa2e816 1789 // Values in arcseconds
1790 epsilon=84381.406-46.836769*t-0.0001831*pow(t,2)+0.00200340*pow(t,3)
1791 -0.000000576*pow(t,4)-0.0000000434*pow(t,5);
1792 omega=450160.398036-6962890.5431*t+7.4722*pow(t,2)+0.007702*pow(t,3)-0.00005939*pow(t,4);
1793 f=335779.526232+1739527262.8478*t-12.7512*pow(t,2)-0.001037*pow(t,3)+0.00000417*pow(t,4);
1794 d=1072260.70369+1602961601.2090*t-6.3706*pow(t,2)+0.006593*pow(t,3)-0.00003169*pow(t,4);
1795
1796 // Convert to radians
5b917489 1797 lsun=lsun*pi/180.;
1798 lmoon=lmoon*pi/180.;
2fa2e816 1799 epsilon=epsilon*pi/(180.*3600.);
1800 omega=omega*pi/(180.*3600.);
1801 f=f*pi/(180.*3600.);
1802 d=d*pi/(180.*3600.);
1803
5b917489 1804 // Change in ecliptic longitude (in arcseconds) due to nutation
1805 Double_t dpsi=-17.2*sin(omega)-1.32*sin(2.*lsun)-0.23*sin(2.*lmoon)+0.21*sin(2.*omega);
2fa2e816 1806
5b917489 1807 // Right ascension shift of the vernal equinox (in arcseconds) due to nutation
1808 Double_t da;
2fa2e816 1809 da=dpsi*cos(epsilon)+0.00264096*sin(omega)+0.00006352*sin(2.*omega)
1810 +0.00001175*sin(2.*f-2.*d+3.*omega)+0.00001121*sin(2.*f-2.*d+omega)
1811 -0.00000455*sin(2.*f-2.*d+2.*omega)+0.00000202*sin(2.*f+3.*omega)+0.00000198*sin(2.*f+omega)
1812 -0.00000172*sin(3.*omega)-0.00000087*t*sin(omega);
1813
1814 // Convert to fractional hours
1815 da=da/(3600.*15.);
1816
1817 Double_t gast=gmst+da;
1818
000b4f12 1819 while (gast<0)
2fa2e816 1820 {
1821 gast+=24.;
1822 }
1823 while (gast>24.)
1824 {
1825 gast-=24.;
1826 }
1827
1828 return gast;
1829}
1830///////////////////////////////////////////////////////////////////////////
1831Double_t AliTimestamp::GetLT(Double_t offset)
1832{
1833// Provide the corresponding local time in fractional hours.
1834// The "offset" denotes the time difference in (fractional) hours w.r.t. UT.
1835// A mean solar day lasts 24h (i.e. 86400s).
1836//
1837// In case a hh:mm:ss format is needed, please use the Convert() facility.
1838
1839 // Current UT time in fractional hours
1840 Double_t h=GetUT();
1841
1842 h+=offset;
1843
000b4f12 1844 while (h<0)
2fa2e816 1845 {
1846 h+=24.;
1847 }
1848 while (h>24)
1849 {
1850 h-=24.;
1851 }
1852
1853 return h;
1854}
1855///////////////////////////////////////////////////////////////////////////
1856Double_t AliTimestamp::GetLMST(Double_t offset)
1857{
1858// Provide the corresponding Local Mean Sidereal Time (LMST) in fractional hours.
1859// The "offset" denotes the time difference in (fractional) hours w.r.t. GMST.
1860// A sidereal day corresponds to 23h 56m 04.09s (i.e. 86164.09s) mean solar time.
1861// The definition of GMST is such that a sidereal clock corresponds with
1862// 24 sidereal hours per revolution of the Earth.
1863// As such, local time offsets w.r.t. UT and GMST can be treated similarly.
1864//
1865// In case a hh:mm:ss format is needed, please use the Convert() facility.
1866
1867 // Current GMST time in fractional hours
1868 Double_t h=GetGMST();
1869
1870 h+=offset;
1871
000b4f12 1872 while (h<0)
2fa2e816 1873 {
1874 h+=24.;
1875 }
1876 while (h>24)
1877 {
1878 h-=24.;
1879 }
1880
1881 return h;
1882}
1883///////////////////////////////////////////////////////////////////////////
1884Double_t AliTimestamp::GetLAST(Double_t offset)
1885{
1886// Provide the corresponding Local Apparent Sidereal Time (LAST) in fractional hours.
1887// The "offset" denotes the time difference in (fractional) hours w.r.t. GAST.
1888// A sidereal day corresponds to 23h 56m 04.09s (i.e. 86164.09s) mean solar time.
1889// The definition of GMST and GAST is such that a sidereal clock corresponds with
1890// 24 sidereal hours per revolution of the Earth.
1891// As such, local time offsets w.r.t. UT, GMST and GAST can be treated similarly.
1892//
1893// In case a hh:mm:ss.sss format is needed, please use the Convert() facility.
1894
1895 // Current GAST time in fractional hours
1896 Double_t h=GetGAST();
1897
1898 h+=offset;
1899
000b4f12 1900 while (h<0)
2fa2e816 1901 {
1902 h+=24.;
1903 }
1904 while (h>24)
1905 {
1906 h-=24.;
1907 }
1908
1909 return h;
1910}
1911///////////////////////////////////////////////////////////////////////////
1912void AliTimestamp::SetLT(Double_t dt,Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_t ns,Int_t ps)
1913{
1914// Set the AliTimestamp parameters corresponding to the LT date and time
1915// in the Gregorian calendar as specified by the input arguments.
1916// This facility is exact upto picosecond precision and as such is
1917// for scientific observations preferable above the corresponding
1918// Set function(s) of TTimestamp.
1919// The latter has a random spread in the sub-second part, which
1920// might be of use in generating distinguishable timestamps while
1921// still keeping second precision.
1922//
1923// The input arguments represent the following :
1924//
1925// dt : the local time offset in fractional hours w.r.t. UT.
1926// y : year in LT (e.g. 1952, 2003 etc...)
1927// m : month in LT (1=jan 2=feb etc...)
1928// d : day in LT (1-31)
1929// hh : elapsed hours in LT (0-23)
1930// mm : elapsed minutes in LT (0-59)
1931// ss : elapsed seconds in LT (0-59)
1932// ns : remaining fractional elapsed second of LT in nanosecond
1933// ps : remaining fractional elapsed nanosecond of LT in picosecond
1934//
1935// Note : ns=0 and ps=0 are the default values.
1936//
1937// This facility first sets the UT as specified by the input arguments
1938// and then corrects the UT by subtracting the local time offset w.r.t. UT.
1939// As such this facility is valid for all AD dates in the Gregorian
1940// calendar with picosecond precision.
1941
1942 SetUT(y,m,d,hh,mm,ss,ns,ps);
1943 Add(-dt);
1944}
1945///////////////////////////////////////////////////////////////////////////
1946void AliTimestamp::SetLT(Double_t dt,Int_t y,Int_t d,Int_t s,Int_t ns,Int_t ps)
1947{
1948// Set the AliTimestamp parameters corresponding to the specified elapsed
1949// timespan since the beginning of the new LT year.
1950// This facility is exact upto picosecond precision and as such is
1951// for scientific observations preferable above the corresponding
1952// Set function(s) of TTimestamp.
1953// The latter has a random spread in the sub-second part, which
1954// might be of use in generating distinguishable timestamps while
1955// still keeping second precision.
1956//
1957// The LT year and elapsed time span is entered via the following input arguments :
1958//
1959// dt : the local time offset in fractional hours w.r.t. UT.
1960// y : year in LT (e.g. 1952, 2003 etc...)
1961// d : elapsed number of days
1962// s : (remaining) elapsed number of seconds
1963// ns : (remaining) elapsed number of nanoseconds
1964// ps : (remaining) elapsed number of picoseconds
1965//
1966// The specified d, s, ns and ps values will be used in an additive
1967// way to determine the elapsed timespan.
1968// So, specification of d=1, s=100, ns=0, ps=0 will result in the
1969// same elapsed time span as d=0, s=24*3600+100, ns=0, ps=0.
1970// However, by making use of the latter the user should take care
1971// of possible integer overflow problems in the input arguments,
1972// which obviously will provide incorrect results.
1973//
1974// Note : ns=0 and ps=0 are the default values.
1975//
1976// This facility first sets the UT as specified by the input arguments
1977// and then corrects the UT by subtracting the local time offset w.r.t. UT.
1978// As such this facility is valid for all AD dates in the Gregorian calendar.
1979
1980 SetUT(y,d,s,ns,ps);
1981 Add(-dt);
1982}
1983///////////////////////////////////////////////////////////////////////////
a4f7a3a1 1984Double_t AliTimestamp::GetJD(Double_t e,TString mode) const
1985{
1986// Provide the fractional Julian Date from epoch e.
1987// The sort of epoch may be specified via the "mode" parameter.
1988//
1989// mode = "J" ==> Julian epoch
1990// "B" ==> Besselian epoch
1991//
1992// The default value is mode="J".
1993
1994 Double_t jd=0;
1995
1996 if (mode=="J" || mode=="j") jd=(e-2000.0)*365.25+2451545.0;
1997
1998 if (mode=="B" || mode=="b") jd=(e-1900.0)*365.242198781+2415020.31352;
1999
2000 return jd;
2001}
2002///////////////////////////////////////////////////////////////////////////
2003Double_t AliTimestamp::GetMJD(Double_t e,TString mode) const
2004{
2005// Provide the fractional Modified Julian Date from epoch e.
2006// The sort of epoch may be specified via the "mode" parameter.
2007//
2008// mode = "J" ==> Julian epoch
2009// "B" ==> Besselian epoch
2010//
2011// The default value is mode="J".
2012
2013 Double_t mjd=GetJD(e,mode)-2400000.5;
2014
2015 return mjd;
2016}
2017///////////////////////////////////////////////////////////////////////////
2018Double_t AliTimestamp::GetTJD(Double_t e,TString mode) const
2019{
2020// Provide the fractional Truncated Julian Date from epoch e.
2021// The sort of epoch may be specified via the "mode" parameter.
2022//
2023// mode = "J" ==> Julian epoch
2024// "B" ==> Besselian epoch
2025//
2026// The default value is mode="J".
2027
2028 Double_t tjd=GetJD(e,mode)-2440000.5;
2029
2030 return tjd;
2031}
2032///////////////////////////////////////////////////////////////////////////