]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - RALICE/AliTimestamp.cxx
Call VEventHandler::Init(Option_t*) in SlaveBegin() and VEventHandler::Init(TTree...
[u/mrichter/AliRoot.git] / RALICE / AliTimestamp.cxx
index bba3c90db1b0f6847d5640aece2b368249aa72c0..df6aa74617d23374d9b11c7c26527f8ecfb675a0 100644 (file)
@@ -45,6 +45,7 @@
 // The Besselian Epoch (BE) indicates the fractional elapsed Besselian year count
 // since the start of the Gregorian year count.
 // A Besselian (or tropical) year is defined to be 365.242198781 days.
+// The date 31-dec-1949 22:09:46.862 UT corresponds to BE=1950.0
 //
 // The Besselian and Julian epochs are used in astronomical catalogs
 // to denote values of time varying observables like e.g. right ascension.
 // In addition, tidal friction and ocean and atmospheric effects will
 // induce seasonal variations in the earth's spin rate and polar motion
 // of the earth's spin axis.
-// To obtain a sidereal time measure, the above efects are taken
-// into account via corrections in the UT to GST conversion.
+// Taking the above effects into account leads to what is called
+// the Greenwich Mean Sidereal Time (GMST).
+// In case also the nutation of the earth's spin axis is taken into
+// account we speak of the Greenwich Apparent Sidereal Time (GAST).
 //
 // This AliTimestamp facility allows for picosecond precision, in view
 // of time of flight analyses for particle physics experiments.
 
 #include "AliTimestamp.h"
 #include "Riostream.h"
-#include <cmath>
 
 ClassImp(AliTimestamp) // Class implementation to enable ROOT I/O
  
@@ -207,44 +209,116 @@ AliTimestamp::AliTimestamp(const AliTimestamp& t) : TTimeStamp(t)
  fCalcns=t.fCalcns;
 }
 ///////////////////////////////////////////////////////////////////////////
-void AliTimestamp::Date(Int_t mode)
+void AliTimestamp::Date(Int_t mode,Double_t offset)
 {
 // Print date/time info.
 //
-// mode = 1 ==> Only the TTimeStamp yy-mm-dd hh:mm:ss:ns and GMST info is printed
+// mode = 1 ==> Only the UT yy-mm-dd hh:mm:ss.sss and GMST info is printed
 //        2 ==> Only the Julian parameter info is printed
-//        3 ==> Both the TTimeStamp, GMST and Julian parameter info is printed
+//        3 ==> Both the UT, GMST and Julian parameter info is printed
+//       -1 ==> Only the UT yy-mm-dd hh:mm:ss.sss and GAST info is printed
+//       -3 ==> Both the UT, GAST and Julian parameter info is printed
 //
-// The default is mode=3.
+// offset : Local time offset from UT (and also GMST) in fractional hours.
+//
+// When an offset value is specified, the corresponding local times
+// LT and LMST (or LAST) are printed as well.
+//
+// The default values are mode=3 and offset=0.
 //
 // Note : In case the (M/T)JD falls outside the TTimeStamp range,
-//        the TTimeStamp info will be replaced by UT hh:mm:ss:ns:ps info.
+//        the yy-mm-dd info will be omitted.
 
  Int_t mjd,mjsec,mjns,mjps;
  GetMJD(mjd,mjsec,mjns);
  mjps=GetPs();
 
+ TString month[12]={"Jan","Feb","Mar","Apr","May","Jun","Jul","Aug","Sep","Oct","Nov","Dec"};
+ TString day[7]={"Mon","Tue","Wed","Thu","Fri","Sat","Sun"};
+ UInt_t y,m,d,wd;
  Int_t hh,mm,ss,ns,ps;
+ Double_t gast;
  
- if (mode==1 || mode==3)
+ if (abs(mode)==1 || abs(mode)==3)
  {
   if (mjd>=40587 && (mjd<65442 || (mjd==65442 && mjsec<8047)))
   {
-   cout << " " << AsString() << endl;
+   GetDate(kTRUE,0,&y,&m,&d);
+   wd=GetDayOfWeek(kTRUE,0);
+   cout << " " << day[wd-1].Data() << ", " << setfill('0') << setw(2) << d << " "
+        << setfill(' ') << month[m-1].Data() << " " << y << " ";
+  }
+  else
+  {
+   cout << " Time ";
+  }
+  GetUT(hh,mm,ss,ns,ps);
+  cout << setfill('0') << setw(2) << hh << ":"
+       << setw(2) << mm << ":" << setw(2) << ss << "."
+       << setw(9) << ns << setw(3) << ps << " (UT)  ";
+  if (mode>0)
+  {
+   GetGMST(hh,mm,ss,ns,ps);
+  }
+  else
+  {
+   gast=GetGAST();
+   Convert(gast,hh,mm,ss,ns,ps);
+  }
+  cout << setfill('0') << setw(2) << hh << ":"
+       << setw(2) << mm << ":" << setw(2) << ss << "."
+       << setw(9) << ns << setw(3) << ps;
+  if (mode>0)
+  {
+   cout << " (GMST)" << endl;
   }
   else
   {
-   GetUT(hh,mm,ss,ns,ps);
-   cout << " UT : " << setfill('0') << setw(2) << hh << ":"
-                    << setw(2) << mm << ":" << setw(2) << ss
-                    << " ns : " << ns << " ps : " << ps << " ";
+   cout << " (GAST)" << endl;
+  }
+
+  // Local time information
+  if (offset)
+  {
+   // Determine the new date by including the offset
+   AliTimestamp t2(*this);
+   t2.Add(offset);
+   Int_t mjd2,mjsec2,mjns2;
+   t2.GetMJD(mjd2,mjsec2,mjns2);
+   if (mjd2>=40587 && (mjd2<65442 || (mjd2==65442 && mjsec2<8047)))
+   {
+    t2.GetDate(kTRUE,0,&y,&m,&d);
+    wd=t2.GetDayOfWeek(kTRUE,0);
+    cout << " " << day[wd-1].Data() << ", " << setfill('0') << setw(2) << d << " "
+         << setfill(' ') << month[m-1].Data() << " " << y << " ";
+   }
+   else
+   {
+    cout << " Time ";
+   }
+   // Determine the local time by including the offset w.r.t. the original timestamp
+   Double_t hlt=GetLT(offset);
+   Double_t hlst=0;
+   if (mode>0)
+   {
+    hlst=GetLMST(offset);
+   }
+   else
+   {
+    hlst=GetLAST(offset);
+   }
+   PrintTime(hlt,12); cout << " (LT)  "; PrintTime(hlst,12);
+   if (mode>0)
+   {
+    cout << " (LMST)" << endl;
+   }
+   else
+   {
+    cout << " (LAST)" << endl;
+   }
   }
-  GetGST(hh,mm,ss,ns,ps);
-  cout << " GST : " << setfill('0') << setw(2) << hh << ":"
-                    << setw(2) << mm << ":" << setw(2) << ss
-                    << " ns : " << ns << " ps : " << ps << endl;
  }
- if (mode==2 || mode==3)
+ if (abs(mode)==2 || abs(mode)==3)
  {
   Int_t jd,jsec,jns;
   GetJD(jd,jsec,jns);
@@ -516,6 +590,8 @@ Double_t AliTimestamp::Convert(Int_t days,Int_t secs,Int_t ns) const
 void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps) const
 {
 // Convert fractional hour count h into hh:mm:ss:ns:ps.
+// The sign of the input value will be neglected, so h<0 will result in
+// the same output values as h>0.
 //
 // Note : Due to computer accuracy the ps value may become inaccurate.
 //
@@ -529,9 +605,15 @@ void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,In
 // please use the corresponding SET() memberfunctions of either AliTimestamp
 // or TTimeStamp.
  
+ // Neglect sign of h
+ h=fabs(h);
+
  hh=int(h);
  h=h-double(hh);
- h=h*3600.;
+ h=h*60.;
+ mm=int(h);
+ h=h-double(mm);
+ h=h*60.;
  ss=int(h);
  h=h-double(ss);
  h=h*1.e9;
@@ -541,9 +623,39 @@ void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,In
  ps=int(h);
 }
 ///////////////////////////////////////////////////////////////////////////
+void AliTimestamp::Convert(Double_t h,Int_t& hh,Int_t& mm,Double_t& ss) const
+{
+// Convert fractional hour count h into hh:mm:ss.s.
+// The sign of the input value will be neglected, so h<0 will result in
+// the same output values as h>0.
+//
+// Notes :
+// -------
+// 1) This memberfunction only converts the input "h" into the corresponding
+//    hh:mm:ss.s values. It does NOT set the corresponding Julian parameters
+//    for the current AliTimestamp instance.
+//    As such the TTimeStamp limitations do NOT apply to this memberfunction.
+//    To set the Julian parameters for the current AliTimestamp instance,
+//    please use the corresponding SET() memberfunctions of either AliTimestamp
+//    or TTimeStamp.
+// 2) This facility can also be used to convert degrees in arcminutes etc...
+ // Neglect sign of h
+ h=fabs(h);
+ hh=int(h);
+ h=h-double(hh);
+ h=h*60.;
+ mm=int(h);
+ h=h-double(mm);
+ ss=h*60.;
+}
+///////////////////////////////////////////////////////////////////////////
 Double_t AliTimestamp::Convert(Int_t hh,Int_t mm,Int_t ss,Int_t ns,Int_t ps) const
 {
 // Convert hh:mm:ss:ns:ps into fractional hour count. 
+// The sign of the input values will be neglected, so the output value
+// will always correspond to a positive hh:mm:ss:ns:ps specification.
 //
 // Note : Due to computer accuracy the ps precision may be lost.
 //
@@ -557,12 +669,87 @@ Double_t AliTimestamp::Convert(Int_t hh,Int_t mm,Int_t ss,Int_t ns,Int_t ps) con
 // please use the corresponding SET() memberfunctions of either AliTimestamp
 // or TTimeStamp.
 
+ // Neglect the sign of the input values
+ hh=abs(hh);
+ mm=abs(mm);
+ ss=abs(ss);
+ ns=abs(ns);
+ ps=abs(ps);
+
  Double_t h=hh;
  h+=double(mm)/60.+(double(ss)+double(ns)*1.e-9+double(ps)*1.e-12)/3600.;
 
  return h;
 }
 ///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::Convert(Int_t hh,Int_t mm,Double_t ss) const
+{
+// Convert hh:mm:ss.s into fractional hour count. 
+// The sign of the input values will be neglected, so the output value
+// will always correspond to a positive hh:mm:ss.s specification.
+//
+// Notes :
+// -------
+// 1) This memberfunction only converts the input hh:mm:ss.s data into the
+//    corresponding fractional hour count. It does NOT set the corresponding
+//    Julian parameters for the current AliTimestamp instance.
+//    As such the TTimeStamp limitations do NOT apply to this memberfunction.
+//    To set the Julian parameters for the current AliTimestamp instance,
+//    please use the corresponding SET() memberfunctions of either AliTimestamp
+//    or TTimeStamp.
+// 2) This facility can also be used to convert ddd:mm:ss.s into fractional degrees.
+
+ // Neglect the sign of the input values
+ hh=abs(hh);
+ mm=abs(mm);
+ ss=fabs(ss);
+
+ Double_t h=hh;
+ h+=double(mm)/60.+ss/3600.;
+
+ return h;
+}
+///////////////////////////////////////////////////////////////////////////
+void AliTimestamp::PrintTime(Double_t h,Int_t ndig) const
+{
+// Print a fractional hour count in hh:mm:ss.ssss format.
+// The range of the printed hour value is : -24 < hh < 24.
+// The argument "ndig" specifies the number of digits for the fractional
+// seconds (e.g. ndig=6 corresponds to microsecond precision).
+// No rounding will be performed, so a second count of 3.473 with ndig=1
+// will appear as 03.4 on the output.
+// Due to computer accuracy, precision on the picosecond level may get lost.
+//
+// The default is ndig=1.
+//
+// Note : The time info is printed without additional spaces or "endline".
+//        This allows the print to be included in various composite output formats.
+
+ Int_t hh,mm,ss;
+ ULong64_t sfrac;
+ Double_t s;
+
+ while (h<-24)
+ {
+  h+=24.;
+ }
+ while (h>24)
+ {
+  h-=24.;
+ }
+   
+ Convert(h,hh,mm,s);
+ ss=Int_t(s);
+ s-=Double_t(ss);
+ s*=pow(10.,ndig);
+ sfrac=ULong64_t(s);
+
+ if (h<0) cout << "-";
+ cout << setfill('0')
+      << setw(2) << hh << ":" << setw(2) << mm << ":"
+      << setw(2) << ss << "." << setw(ndig) << sfrac;
+}
+///////////////////////////////////////////////////////////////////////////
 void AliTimestamp::FillJulian()
 {
 // Calculation and setting of the Julian date/time parameters corresponding
@@ -762,7 +949,7 @@ void AliTimestamp::SetMJD(Int_t mjd,Int_t sec,Int_t ns,Int_t ps)
  Int_t limit=65442; // MJD of the latest possible TTimeStamp date/time
  
  Int_t date,time;
- if (mjd<epoch || (mjd>=limit && sec>=8047))
+ if (mjd<epoch || mjd>limit || (mjd==limit && sec>=8047))
  {
   Set(0,kFALSE,0,kFALSE);
   date=GetDate();
@@ -1094,6 +1281,27 @@ void AliTimestamp::Add(Int_t d,Int_t s,Int_t ns,Int_t ps)
  SetMJD(days,secs,nsec,psec);
 }
 ///////////////////////////////////////////////////////////////////////////
+void AliTimestamp::Add(Double_t hours)
+{
+// Add (or subtract) a certain time difference to the current timestamp.
+// The time difference is specified as a (fractional) number of hours.
+// Subtraction can be achieved by entering a negative value as input argument.
+
+ Int_t d,s,ns,ps;
+ Double_t h=fabs(hours);
+ d=int(h/24.);
+ h-=double(d)*24.;
+ h*=3600.;
+ s=int(h);
+ h-=double(s);
+ h*=1.e9;
+ ns=int(h);
+ h-=double(ns);
+ ps=int(h*1000.);
+ if (hours>0) Add(d,s,ns,ps);
+ if (hours<0) Add(-d,-s,-ns,-ps);
+}
+///////////////////////////////////////////////////////////////////////////
 Int_t AliTimestamp::GetDifference(AliTimestamp* t,Int_t& d,Int_t& s,Int_t& ns,Int_t& ps)
 {
 // Provide the time difference w.r.t the AliTimestamp specified on the input.
@@ -1372,7 +1580,7 @@ void AliTimestamp::SetUT(Int_t y,Int_t m,Int_t d,Int_t hh,Int_t mm,Int_t ss,Int_
 // Note : ns=0 and ps=0 are the default values.
 //
 // This facility first determines the elapsed days, seconds etc...
-// since the beginning of the specified UT year on bais of the
+// since the beginning of the specified UT year on basis of the
 // input arguments. Subsequently it invokes the SetUT memberfunction
 // for the elapsed timespan.
 // As such this facility is valid for all AD dates in the Gregorian
@@ -1461,9 +1669,9 @@ Double_t AliTimestamp::GetUT()
  return ut;
 }
 ///////////////////////////////////////////////////////////////////////////
-void AliTimestamp::GetGST(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
+void AliTimestamp::GetGMST(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
 {
-// Provide the corrresponding Greenwich Sideral Time (GST).
+// Provide the corrresponding Greenwich Mean Sideral Time (GMST).
 // The algorithm used is the one described at p. 83 of the book
 // Astronomy Methods by Hale Bradt.
 // This facility is based on the MJD, so the TTimeStamp limitations
@@ -1483,7 +1691,7 @@ void AliTimestamp::GetGST(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
  AliTimestamp sid;
  sid.SetMJD(mjd,sec,nsec,psec);
 
- // Add offset for GST start value defined as 06:41:50.54841 at 01-jan 00:00:00 UT
+ // Add offset for GMST start value defined as 06:41:50.54841 at 01-jan 00:00:00 UT
  sec=6*3600+41*60+50;
  nsec=548410000;
  psec=0;
@@ -1509,22 +1717,212 @@ void AliTimestamp::GetGST(Int_t& hh,Int_t& mm,Int_t& ss,Int_t& ns,Int_t& ps)
  ps=psec;
 }
 ///////////////////////////////////////////////////////////////////////////
-Double_t AliTimestamp::GetGST()
+Double_t AliTimestamp::GetGMST()
 {
-// Provide the corrresponding Greenwich Sideral Time (GMST)
+// Provide the corrresponding Greenwich Mean Sideral Time (GMST)
 // in fractional hours.
 // This facility is based on the MJD, so the TTimeStamp limitations
 // do not apply here.
 
  Int_t hh,mm,ss,ns,ps;
 
- GetGST(hh,mm,ss,ns,ps);
+ GetGMST(hh,mm,ss,ns,ps);
 
  Double_t gst=Convert(hh,mm,ss,ns,ps);
 
  return gst;
 }
 ///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::GetGAST()
+{
+// Provide the corrresponding Greenwich Apparent Sideral Time (GAST)
+// in fractional hours.
+// In case a hh:mm:ss.sss format is needed, please invoke the Convert()
+// memberfunction for conversion of the provided fractional hour value.
+//
+// The GAST is the GMST corrected for the shift of the vernal equinox
+// due to nutation. The right ascension component of the nutation correction
+// of the vernal equinox is called the "equation of the equinoxes".
+// So we have :
+//
+//      GAST = GMST + (equation of the equinoxes)
+//
+// The equation of the equinoxes is determined via the Almanac() memberfunction.
+// 
+// Since GMST is based on the MJD, the TTimeStamp limitations do not apply here.
+
+ Double_t da=Almanac();
+
+ // Convert to fractional hours
+ da/=3600.;
+
+ Double_t gast=GetGMST()+da;
+
+ while (gast<0)
+ {
+  gast+=24.;
+ }
+ while (gast>24.)
+ {
+  gast-=24.;
+ }
+
+ return gast;
+}
+///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::GetLT(Double_t offset)
+{
+// Provide the corresponding local time in fractional hours.
+// The "offset" denotes the time difference in (fractional) hours w.r.t. UT.
+// A mean solar day lasts 24h (i.e. 86400s).
+//
+// In case a hh:mm:ss format is needed, please use the Convert() facility. 
+
+ // Current UT time in fractional hours
+ Double_t h=GetUT();
+ h+=offset;
+
+ while (h<0)
+ {
+  h+=24.;
+ }
+ while (h>24)
+ {
+  h-=24.;
+ }
+
+ return h;
+}
+///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::GetLMST(Double_t offset)
+{
+// Provide the corresponding Local Mean Sidereal Time (LMST) in fractional hours.
+// The "offset" denotes the time difference in (fractional) hours w.r.t. GMST.
+// A sidereal day corresponds to 23h 56m 04.09s (i.e. 86164.09s) mean solar time.
+// The definition of GMST is such that a sidereal clock corresponds with
+// 24 sidereal hours per revolution of the Earth.
+// As such, local time offsets w.r.t. UT and GMST can be treated similarly. 
+//
+// In case a hh:mm:ss format is needed, please use the Convert() facility. 
+
+ // Current GMST time in fractional hours
+ Double_t h=GetGMST();
+
+ h+=offset;
+
+ while (h<0)
+ {
+  h+=24.;
+ }
+ while (h>24)
+ {
+  h-=24.;
+ }
+
+ return h;
+}
+///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::GetLAST(Double_t offset)
+{
+// Provide the corresponding Local Apparent Sidereal Time (LAST) in fractional hours.
+// The "offset" denotes the time difference in (fractional) hours w.r.t. GAST.
+// A sidereal day corresponds to 23h 56m 04.09s (i.e. 86164.09s) mean solar time.
+// The definition of GMST and GAST is such that a sidereal clock corresponds with
+// 24 sidereal hours per revolution of the Earth.
+// As such, local time offsets w.r.t. UT, GMST and GAST can be treated similarly. 
+//
+// In case a hh:mm:ss.sss format is needed, please use the Convert() facility. 
+
+ // Current GAST time in fractional hours
+ Double_t h=GetGAST();
+
+ h+=offset;
+
+ while (h<0)
+ {
+  h+=24.;
+ }
+ while (h>24)
+ {
+  h-=24.;
+ }
+
+ return h;
+}
+///////////////////////////////////////////////////////////////////////////
+void 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)
+{
+// Set the AliTimestamp parameters corresponding to the LT date and time
+// in the Gregorian calendar as specified by the input arguments.
+// This facility is exact upto picosecond precision and as such is
+// for scientific observations preferable above the corresponding
+// Set function(s) of TTimestamp.
+// The latter has a random spread in the sub-second part, which
+// might be of use in generating distinguishable timestamps while
+// still keeping second precision.
+//
+// The input arguments represent the following :
+//
+// dt : the local time offset in fractional hours w.r.t. UT.
+// y  : year in LT (e.g. 1952, 2003 etc...)
+// m  : month in LT (1=jan  2=feb etc...)
+// d  : day in LT (1-31)
+// hh : elapsed hours in LT (0-23) 
+// mm : elapsed minutes in LT (0-59)
+// ss : elapsed seconds in LT (0-59)
+// ns : remaining fractional elapsed second of LT in nanosecond
+// ps : remaining fractional elapsed nanosecond of LT in picosecond
+//
+// Note : ns=0 and ps=0 are the default values.
+//
+// This facility first sets the UT as specified by the input arguments
+// and then corrects the UT by subtracting the local time offset w.r.t. UT.
+// As such this facility is valid for all AD dates in the Gregorian
+// calendar with picosecond precision.
+
+ SetUT(y,m,d,hh,mm,ss,ns,ps);
+ Add(-dt);
+}
+///////////////////////////////////////////////////////////////////////////
+void AliTimestamp::SetLT(Double_t dt,Int_t y,Int_t d,Int_t s,Int_t ns,Int_t ps)
+{
+// Set the AliTimestamp parameters corresponding to the specified elapsed
+// timespan since the beginning of the new LT year.
+// This facility is exact upto picosecond precision and as such is
+// for scientific observations preferable above the corresponding
+// Set function(s) of TTimestamp.
+// The latter has a random spread in the sub-second part, which
+// might be of use in generating distinguishable timestamps while
+// still keeping second precision.
+//
+// The LT year and elapsed time span is entered via the following input arguments :
+//
+// dt : the local time offset in fractional hours w.r.t. UT.
+// y  : year in LT (e.g. 1952, 2003 etc...)
+// d  : elapsed number of days 
+// s  : (remaining) elapsed number of seconds
+// ns : (remaining) elapsed number of nanoseconds
+// ps : (remaining) elapsed number of picoseconds
+//
+// The specified d, s, ns and ps values will be used in an additive
+// way to determine the elapsed timespan.
+// So, specification of d=1, s=100, ns=0, ps=0 will result in the
+// same elapsed time span as d=0, s=24*3600+100, ns=0, ps=0.
+// However, by making use of the latter the user should take care
+// of possible integer overflow problems in the input arguments,
+// which obviously will provide incorrect results. 
+//
+// Note : ns=0 and ps=0 are the default values.
+//
+// This facility first sets the UT as specified by the input arguments
+// and then corrects the UT by subtracting the local time offset w.r.t. UT.
+// As such this facility is valid for all AD dates in the Gregorian calendar.
+
+ SetUT(y,d,s,ns,ps);
+ Add(-dt);
+}
+///////////////////////////////////////////////////////////////////////////
 Double_t AliTimestamp::GetJD(Double_t e,TString mode) const
 {
 // Provide the fractional Julian Date from epoch e.
@@ -1574,3 +1972,157 @@ Double_t AliTimestamp::GetTJD(Double_t e,TString mode) const
  return tjd;
 }
 ///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::Almanac(Double_t* dpsi,Double_t* deps,Double_t* eps)
+{
+// Determination of some astronomical observables which may be needed
+// for further calculations like e.g. precession of coordinates.
+//
+// The standard returned value is the "equation of the equinoxes" 
+// (i.e. the nutational shift of the RA of the vernal equinox) in seconds.
+// The memberfunction arguments provide the possibility of retrieving
+// optional returned values. The corresponding observables are :
+//
+// dpsi : Nutational shift in ecliptic longitude in arcseconds
+// deps : Nutational shift in ecliptic obliquity in arcseconds
+// eps  : Mean obliquity of the ecliptic in arcseconds
+//
+// All shifts are determined for the current timestamp with
+// J2000.0 (i.e. 01-jan-2000 12:00:00 UT) as the reference epoch.
+//
+// Invokation example :
+// --------------------
+// AliTimestamp t;
+// Double_t da,dpsi,deps,eps;
+// da=t.Almanac(&dpsi,&deps,&eps); 
+//
+// The nutation model used is the new one as documented in :
+//   "The IAU Resolutions on Astronomical Reference Systems,
+//    Time Scales and Earth Rotation Models".
+// This document is freely available as Circular 179 (2005) of the
+// United States Naval Observatory (USNO).
+// (See : http://aa.usno.navy.mil/publications/docs).
+//
+// The change in ecliptic longitude (dpsi) and ecliptic obliquity (deps)
+// are evaluated using the IAU 2000A nutation series expansion
+// as provided in the USNO Circular 179.
+// The new expression for the equation of the equinoxes is based on a series
+// expansion and is the most accurate one known to date.
+// The components are documented on p.17 of the USNO Circular 179.
+//
+// In the current implementation only the first 28 terms of the nutation series
+// are used. This provides an accuracy of about 0.01 arcsec corresponding to 0.001 sec.
+// In case a better accuracy is required, the series can be extended.
+// The total series expansion consists of 1365 terms.
+// 
+// Since all calculations are based on the JD, the TTimeStamp limitations
+// do not apply here.
+
+ Double_t pi=acos(-1.);
+
+ Double_t t;       // Time difference in fractional Julian centuries w.r.t. the start of J2000.
+ Double_t epsilon; // Mean obliquity of the ecliptic
+ Double_t l;       // Mean anomaly of the Moon
+ Double_t lp;      // Mean anomaly of the Sun
+ Double_t f;       // Mean argument of latitude of the moon
+ Double_t d;       // Mean elongation of the Moon from the Sun
+ Double_t om;      // Mean longitude of the Moon's mean ascending mode
+
+ t=(GetJD()-2451545.0)/36525.;
+
+ // Values of epsilon and the fundamental luni-solar arguments in arcseconds
+ epsilon=84381.406-46.836769*t-0.0001831*pow(t,2)+0.00200340*pow(t,3)
+         -0.000000576*pow(t,4)-0.0000000434*pow(t,5);
+ l=485868.249036+1717915923.2178*t+31.8792*pow(t,2)+0.051635*pow(t,3)-0.00024470*pow(t,4);
+ lp=1287104.79305+129596581.0481*t-0.5532*pow(t,2)+0.000136*pow(t,3)-0.00001149*pow(t,4);
+ f=335779.526232+1739527262.8478*t-12.7512*pow(t,2)-0.001037*pow(t,3)+0.00000417*pow(t,4);
+ d=1072260.70369+1602961601.2090*t-6.3706*pow(t,2)+0.006593*pow(t,3)-0.00003169*pow(t,4);
+ om=450160.398036-6962890.5431*t+7.4722*pow(t,2)+0.007702*pow(t,3)-0.00005939*pow(t,4);
+
+ if (eps) *eps=epsilon;
+
+ // Convert to radians
+ epsilon=epsilon*pi/(180.*3600.);
+ f=f*pi/(180.*3600.);
+ d=d*pi/(180.*3600.);
+ l=l*pi/(180.*3600.);
+ lp=lp*pi/(180.*3600.);
+ om=om*pi/(180.*3600.);
+
+ //The IAU 2000A nutation series expansion.
+ Double_t phi[28]={om,2.*(f-d+om),2.*(f+om),2.*om,lp,lp+2.*(f-d+om),l,
+                   2.*f+om,l+2.*(f+om),2.*(f-d+om)-lp,2.*(f-d)+om,2.*(f+om)-l,2.*d-l,l+om,
+                   om-l,2.*(f+d+om)-l,l+2.*f+om,2.*(f-l)+om,2.*d,2.*(f+d+om),2.*(f-d+om-lp),
+                   2.*(d-l),2.*(l+d+om),l+2.*(f-d+om),2.*f+om-l,2.*l,2.*f,lp+om};
+ Double_t s[28]={-17.2064161,-1.3170907,-0.2276413, 0.2074554, 0.1475877,-0.0516821, 0.0711159,
+                  -0.0387298,-0.0301461, 0.0215829, 0.0128227, 0.0123457, 0.0156994, 0.0063110,
+                  -0.0057976,-0.0059641,-0.0051613, 0.0045893, 0.0063384,-0.0038571, 0.0032481,
+                  -0.0047722,-0.0031046, 0.0028593, 0.0020441, 0.0029243, 0.0025887,-0.0014053};
+ Double_t sd[28]={-0.0174666,-0.0001675,-0.0000234, 0.0000207,-0.0003633, 0.0001226, 0.0000073,
+                  -0.0000367,-0.0000036,-0.0000494, 0.0000137, 0.0000011, 0.0000010, 0.0000063,
+                  -0.0000063,-0.0000011,-0.0000042, 0.0000050, 0.0000011,-0.0000001, 0.0000000,
+                   0.0000000,-0.0000001, 0.0000000, 0.0000021, 0.0000000, 0.0000000,-0.0000025};
+ Double_t cp[28]={ 0.0033386,-0.0013696, 0.0002796,-0.0000698, 0.0011817,-0.0000524,-0.0000872,
+                   0.0000380, 0.0000816, 0.0000111, 0.0000181, 0.0000019,-0.0000168, 0.0000027,
+                  -0.0000189, 0.0000149, 0.0000129, 0.0000031,-0.0000150, 0.0000158, 0.0000000,
+                  -0.0000018, 0.0000131,-0.0000001, 0.0000010,-0.0000074,-0.0000066, 0.0000079};
+ Double_t c[28]= { 9.2052331, 0.5730336, 0.0978459,-0.0897492, 0.0073871, 0.0224386,-0.0006750,
+                   0.0200728, 0.0129025,-0.0095929,-0.0068982,-0.0053311,-0.0001235,-0.0033228,
+                   0.0031429, 0.0025543, 0.0026366,-0.0024236,-0.0001220, 0.0016452,-0.0013870,
+                   0.0000477, 0.0013238,-0.0012338,-0.0010758,-0.0000609,-0.0000550, 0.0008551};
+ Double_t cd[28]={ 0.0009086,-0.0003015,-0.0000485, 0.0000470,-0.0000184,-0.0000677, 0.0000000,
+                   0.0000018,-0.0000063, 0.0000299,-0.0000009, 0.0000032, 0.0000000, 0.0000000,
+                   0.0000000,-0.0000011, 0.0000000,-0.0000010, 0.0000000,-0.0000011, 0.0000000,
+                   0.0000000,-0.0000011, 0.0000010, 0.0000000, 0.0000000, 0.0000000,-0.0000002};
+ Double_t sp[28]={ 0.0015377,-0.0004587, 0.0001374,-0.0000291,-0.0001924,-0.0000174, 0.0000358,
+                   0.0000318, 0.0000367, 0.0000132, 0.0000039,-0.0000004, 0.0000082,-0.0000009,
+                  -0.0000075, 0.0000066, 0.0000078, 0.0000020, 0.0000029, 0.0000068, 0.0000000,
+                  -0.0000025, 0.0000059,-0.0000003,-0.0000003, 0.0000013, 0.0000011,-0.0000045};
+ Double_t dp=0,de=0,da=0;
+ for (Int_t i=0; i<28; i++)
+ {
+  dp+=(s[i]+sd[i]*t)*sin(phi[i])+cp[i]*cos(phi[i]);
+  de+=(c[i]+cd[i]*t)*cos(phi[i])+sp[i]*sin(phi[i]);
+ }
+
+ da=dp*cos(epsilon)+0.00264096*sin(om)+0.00006352*sin(2.*om)
+     +0.00001175*sin(2.*f-2.*d+3.*om)+0.00001121*sin(2.*f-2.*d+om)
+     -0.00000455*sin(2.*f-2.*d+2.*om)+0.00000202*sin(2.*f+3.*om)+0.00000198*sin(2.*f+om)
+     -0.00000172*sin(3.*om)-0.00000087*t*sin(om);
+
+ if (dpsi) *dpsi=dp;
+ if (deps) *deps=de;
+
+ // Convert to seconds
+ da/=15.;
+
+ return da;
+}
+///////////////////////////////////////////////////////////////////////////
+void AliTimestamp::SetEpoch(Double_t e,TString mode)
+{
+// Set the timestamp parameters according to the epoch as specified by
+// the input argument "e".
+// Via the input argument "mode" the user can specify the type of epoch
+//
+// mode = "B" ==> Besselian epoch
+//        "J" ==> Julian epoch
+
+ Double_t jd=GetJD(e,mode);
+ SetJD(jd);
+}
+///////////////////////////////////////////////////////////////////////////
+Double_t AliTimestamp::GetEpoch(TString mode)
+{
+// Provide the corresponding epoch value.
+// Via the input argument "mode" the user can specify the type of epoch
+//
+// mode = "B" ==> Besselian epoch
+//        "J" ==> Julian epoch
+
+ Double_t e=0;
+ if (mode=="B" || mode=="b") e=GetBE();
+ if (mode=="J" || mode=="j") e=GetJE();
+ return e;
+}
+///////////////////////////////////////////////////////////////////////////