/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz, 11/2003 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MAstro // ------ // //////////////////////////////////////////////////////////////////////////// #include "MAstro.h" #include #include // TVector3 #include "MTime.h" // MTime::GetGmst #include "MAstroCatalog.h" // FIXME: replace by MVector3! using namespace std; ClassImp(MAstro); Double_t MAstro::Trunc(Double_t val) { // dint(A) - truncate to nearest whole number towards zero (double) return val<0 ? TMath::Ceil(val) : TMath::Floor(val); } Double_t MAstro::Round(Double_t val) { // dnint(A) - round to nearest whole number (double) return val<0 ? TMath::Ceil(val-0.5) : TMath::Floor(val+0.5); } Double_t MAstro::Hms2Sec(Int_t deg, UInt_t min, Double_t sec, Char_t sgn) { const Double_t rc = TMath::Sign((60.0 * (60.0 * (Double_t)TMath::Abs(deg) + (Double_t)min) + sec), (Double_t)deg); return sgn=='-' ? -rc : rc; } Double_t MAstro::Dms2Rad(Int_t deg, UInt_t min, Double_t sec, Char_t sgn) { // pi/(180*3600): arcseconds to radians //#define DAS2R 4.8481368110953599358991410235794797595635330237270e-6 return Hms2Sec(deg, min, sec, sgn)*TMath::Pi()/(180*3600)/**DAS2R*/; } Double_t MAstro::Hms2Rad(Int_t hor, UInt_t min, Double_t sec, Char_t sgn) { // pi/(12*3600): seconds of time to radians //#define DS2R 7.2722052166430399038487115353692196393452995355905e-5 return Hms2Sec(hor, min, sec, sgn)*TMath::Pi()/(12*3600)/**DS2R*/; } Double_t MAstro::Dms2Deg(Int_t deg, UInt_t min, Double_t sec, Char_t sgn) { return Hms2Sec(deg, min, sec, sgn)/3600.; } Double_t MAstro::Hms2Deg(Int_t hor, UInt_t min, Double_t sec, Char_t sgn) { return Hms2Sec(hor, min, sec, sgn)/240.; } Double_t MAstro::Dms2Hor(Int_t deg, UInt_t min, Double_t sec, Char_t sgn) { return Hms2Sec(deg, min, sec, sgn)/54000.; } Double_t MAstro::Hms2Hor(Int_t hor, UInt_t min, Double_t sec, Char_t sgn) { return Hms2Sec(hor, min, sec, sgn)/3600.; } void MAstro::Day2Hms(Double_t day, Char_t &sgn, UShort_t &hor, UShort_t &min, UShort_t &sec) { /* Handle sign */ sgn = day<0?'-':'+'; /* Round interval and express in smallest units required */ Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds /* Separate into fields */ const Double_t ah = Trunc(a/3600.); a -= ah * 3600.; const Double_t am = Trunc(a/60.); a -= am * 60.; const Double_t as = Trunc(a); /* Return results */ hor = (UShort_t)ah; min = (UShort_t)am; sec = (UShort_t)as; } void MAstro::Rad2Hms(Double_t rad, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Day2Hms(rad/(TMath::Pi()*2), sgn, deg, min, sec); } void MAstro::Rad2Dms(Double_t rad, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Rad2Hms(rad*15, sgn, deg, min, sec); } void MAstro::Deg2Dms(Double_t d, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Day2Hms(d/24, sgn, deg, min, sec); } void MAstro::Deg2Hms(Double_t d, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Rad2Hms(d/360, sgn, deg, min, sec); } void MAstro::Hor2Dms(Double_t h, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Day2Hms(h*15/24, sgn, deg, min, sec); } void MAstro::Hor2Hms(Double_t h, Char_t &sgn, UShort_t °, UShort_t &min, UShort_t &sec) { Day2Hms(h/24, sgn, deg, min, sec); } void MAstro::Day2Hm(Double_t day, Char_t &sgn, UShort_t &hor, Double_t &min) { /* Handle sign */ sgn = day<0?'-':'+'; /* Round interval and express in smallest units required */ Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds /* Separate into fields */ const Double_t ah = Trunc(a/3600.); a -= ah * 3600.; /* Return results */ hor = (UShort_t)ah; min = a/60.; } void MAstro::Rad2Hm(Double_t rad, Char_t &sgn, UShort_t °, Double_t &min) { Day2Hm(rad/(TMath::Pi()*2), sgn, deg, min); } void MAstro::Rad2Dm(Double_t rad, Char_t &sgn, UShort_t °, Double_t &min) { Rad2Hm(rad*15, sgn, deg, min); } void MAstro::Deg2Dm(Double_t d, Char_t &sgn, UShort_t °, Double_t &min) { Day2Hm(d/24, sgn, deg, min); } void MAstro::Deg2Hm(Double_t d, Char_t &sgn, UShort_t °, Double_t &min) { Rad2Hm(d/360, sgn, deg, min); } void MAstro::Hor2Dm(Double_t h, Char_t &sgn, UShort_t °, Double_t &min) { Day2Hm(h*15/24, sgn, deg, min); } void MAstro::Hor2Hm(Double_t h, Char_t &sgn, UShort_t °, Double_t &min) { Day2Hm(h/24, sgn, deg, min); } // -------------------------------------------------------------------------- // // Interpretes a string ' - 12 30 00.0' or '+ 12 30 00.0' // as floating point value -12.5 or 12.5. If interpretation is // successfull kTRUE is returned, otherwise kFALSE. ret is not // touched if interpretation was not successfull. The successfull // interpreted part is removed from the TString. // Bool_t MAstro::String2Angle(TString &str, Double_t &ret) { Char_t sgn; Int_t d, len; UInt_t m; Float_t s; // Skip whitespaces before %c and after %f int n=sscanf(str.Data(), " %c %d %d %f %n", &sgn, &d, &m, &s, &len); if (n!=4 || (sgn!='+' && sgn!='-')) return kFALSE; str.Remove(0, len); ret = Dms2Deg(d, m, s, sgn); return kTRUE; } // -------------------------------------------------------------------------- // // Interpretes a string '-12:30:00.0', '12:30:00.0' or '+12:30:00.0' // as floating point value -12.5, 12.5 or 12.5. If interpretation is // successfull kTRUE is returned, otherwise kFALSE. ret is not // touched if interpretation was not successfull. // Bool_t MAstro::Coordinate2Angle(const TString &str, Double_t &ret) { Char_t sgn = str[0]=='-' ? '-' : '+'; Int_t d; UInt_t m; Float_t s; const int n=sscanf(str[0]=='+'||str[0]=='-' ? str.Data()+1 : str.Data(), "%d:%d:%f", &d, &m, &s); if (n!=3) return kFALSE; ret = Dms2Deg(d, m, s, sgn); return kTRUE; } // -------------------------------------------------------------------------- // // Return year y, month m and day d corresponding to Mjd. // void MAstro::Mjd2Ymd(UInt_t mjd, UShort_t &y, Byte_t &m, Byte_t &d) { // Express day in Gregorian calendar const ULong_t jd = mjd + 2400001; const ULong_t n4 = 4*(jd+((6*((4*jd-17918)/146097))/4+1)/2-37); const ULong_t nd10 = 10*(((n4-237)%1461)/4)+5; y = n4/1461L-4712; m = ((nd10/306+2)%12)+1; d = (nd10%306)/10+1; } // -------------------------------------------------------------------------- // // Return Mjd corresponding to year y, month m and day d. // Int_t MAstro::Ymd2Mjd(UShort_t y, Byte_t m, Byte_t d) { // Month lengths in days static int months[12] = { 31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31 }; // Validate month if (m<1 || m>12) return -1; // Allow for leap year months[1] = (y%4==0 && (y%100!=0 || y%400==0)) ? 29 : 28; // Validate day if (d<1 || d>months[m-1]) return -1; // Precalculate some values const Byte_t lm = 12-m; const ULong_t lm10 = 4712 + y - lm/10; // Perform the conversion return 1461L*lm10/4 + (306*((m+9)%12)+5)/10 - (3*((lm10+188)/100))/4 + d - 2399904; } // -------------------------------------------------------------------------- // // theta0, phi0 [rad]: polar angle/zenith distance, azimuth of 1st object // theta1, phi1 [rad]: polar angle/zenith distance, azimuth of 2nd object // AngularDistance [rad]: Angular distance between two objects // Double_t MAstro::AngularDistance(Double_t theta0, Double_t phi0, Double_t theta1, Double_t phi1) { TVector3 v0(1); v0.Rotate(phi0, theta0); TVector3 v1(1); v1.Rotate(phi1, theta1); return v0.Angle(v1); } // -------------------------------------------------------------------------- // // Calls MTime::GetGmst() Better use MTime::GetGmst() directly // Double_t MAstro::UT2GMST(Double_t ut1) { return MTime(ut1).GetGmst(); } // -------------------------------------------------------------------------- // // RotationAngle // // calculates the angle for the rotation of the sky coordinate system // with respect to the local coordinate system. This is identical // to the rotation angle of the sky image in the camera. // // sinl [rad]: sine of observers latitude // cosl [rad]: cosine of observers latitude // theta [rad]: polar angle/zenith distance // phi [rad]: rotation angle/azimuth // // Return sin/cos component of angle // // The convention is such, that the rotation angle is -pi/pi if // right ascension and local rotation angle are counted in the // same direction, 0 if counted in the opposite direction. // // (In other words: The rotation angle is 0 when the source culminates) // // Using vectors it can be done like: // TVector3 v, p; // v.SetMagThetaPhi(1, theta, phi); // p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0); // v = v.Cross(l)); // v.RotateZ(-phi); // v.Rotate(-theta) // rho = TMath::ATan2(v(2), v(1)); // // For more information see TDAS 00-11, eqs. (18) and (20) // void MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi, Double_t &sin, Double_t &cos) { const Double_t sint = TMath::Sin(theta); const Double_t cost = TMath::Cos(theta); const Double_t snlt = sinl*sint; const Double_t cslt = cosl*cost; const Double_t sinp = TMath::Sin(phi); const Double_t cosp = TMath::Cos(phi); const Double_t v1 = sint*sinp; const Double_t v2 = cslt - snlt*cosp; const Double_t denom = TMath::Sqrt(v1*v1 + v2*v2); sin = cosl*sinp / denom; // y-component cos = (snlt-cslt*cosp) / denom; // x-component } // -------------------------------------------------------------------------- // // RotationAngle // // calculates the angle for the rotation of the sky coordinate system // with respect to the local coordinate system. This is identical // to the rotation angle of the sky image in the camera. // // sinl [rad]: sine of observers latitude // cosl [rad]: cosine of observers latitude // theta [rad]: polar angle/zenith distance // phi [rad]: rotation angle/azimuth // // Return angle [rad] in the range -pi, pi // // The convention is such, that the rotation angle is -pi/pi if // right ascension and local rotation angle are counted in the // same direction, 0 if counted in the opposite direction. // // (In other words: The rotation angle is 0 when the source culminates) // // Using vectors it can be done like: // TVector3 v, p; // v.SetMagThetaPhi(1, theta, phi); // p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0); // v = v.Cross(l)); // v.RotateZ(-phi); // v.Rotate(-theta) // rho = TMath::ATan2(v(2), v(1)); // // For more information see TDAS 00-11, eqs. (18) and (20) // Double_t MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi) { const Double_t snlt = sinl*TMath::Sin(theta); const Double_t cslt = cosl*TMath::Cos(theta); const Double_t sinp = TMath::Sin(phi); const Double_t cosp = TMath::Cos(phi); return TMath::ATan2(cosl*sinp, snlt-cslt*cosp); } // -------------------------------------------------------------------------- // // Kepler - solve the equation of Kepler // Double_t MAstro::Kepler(Double_t m, Double_t ecc) { m *= TMath::DegToRad(); Double_t delta = 0; Double_t e = m; do { delta = e - ecc * sin(e) - m; e -= delta / (1 - ecc * cos(e)); } while (fabs(delta) > 1e-6); return e; } // -------------------------------------------------------------------------- // // GetMoonPhase - calculate phase of moon as a fraction: // Returns -1 if calculation failed // Double_t MAstro::GetMoonPhase(Double_t mjd) { /****** Calculation of the Sun's position. ******/ // date within epoch const Double_t epoch = 44238; // 1980 January 0.0 const Double_t day = mjd - epoch; if (day<0) { cout << "MAstro::GetMoonPhase - Day before Jan 1980" << endl; return -1; } // mean anomaly of the Sun const Double_t n = fmod(day*360/365.2422, 360); const Double_t elonge = 278.833540; // ecliptic longitude of the Sun at epoch 1980.0 const Double_t elongp = 282.596403; // ecliptic longitude of the Sun at perigee // convert from perigee co-ordinates to epoch 1980.0 const Double_t m = fmod(n + elonge - elongp + 360, 360); // solve equation of Kepler const Double_t eccent = 0.016718; // eccentricity of Earth's orbit const Double_t k = Kepler(m, eccent); const Double_t ec0 = sqrt((1 + eccent) / (1 - eccent)) * tan(k / 2); // true anomaly const Double_t ec = 2 * atan(ec0) * TMath::RadToDeg(); // Sun's geocentric ecliptic longitude const Double_t lambdasun = fmod(ec + elongp + 720, 360); /****** Calculation of the Moon's position. ******/ // Moon's mean longitude. const Double_t mmlong = 64.975464; // moon's mean lonigitude at the epoch const Double_t ml = fmod(13.1763966*day + mmlong + 360, 360); // Moon's mean anomaly. const Double_t mmlongp = 349.383063; // mean longitude of the perigee at the epoch const Double_t mm = fmod(ml - 0.1114041*day - mmlongp + 720, 360); // Evection. const Double_t ev = 1.2739 * sin((2 * (ml - lambdasun) - mm)*TMath::DegToRad()); // Annual equation. const Double_t sinm = TMath::Sin(m*TMath::DegToRad()); const Double_t ae = 0.1858 * sinm; // Correction term. const Double_t a3 = 0.37 * sinm; // Corrected anomaly. const Double_t mmp = (mm + ev - ae - a3)*TMath::DegToRad(); // Correction for the equation of the centre. const Double_t mec = 6.2886 * sin(mmp); // Another correction term. const Double_t a4 = 0.214 * sin(2 * mmp); // Corrected longitude. const Double_t lp = ml + ev + mec - ae + a4; // Variation. const Double_t v = 0.6583 * sin(2 * (lp - lambdasun)*TMath::DegToRad()); // True longitude. const Double_t lpp = lp + v; // Age of the Moon in degrees. const Double_t age = (lpp - lambdasun)*TMath::DegToRad(); // Calculation of the phase of the Moon. return (1 - TMath::Cos(age)) / 2; } // -------------------------------------------------------------------------- // // Calculate the Period to which the time belongs to. The Period is defined // as the number of synodic months ellapsed since the first full moon // after Jan 1st 1980 (which was @ MJD=44240.37917) // Double_t MAstro::GetMoonPeriod(Double_t mjd) { const Double_t synmonth = 29.53058868; // synodic month (new Moon to new Moon) const Double_t epoch0 = 44240.37917; // First full moon after 1980/1/1 const Double_t et = mjd-epoch0; // Ellapsed time return et/synmonth; } // -------------------------------------------------------------------------- // // To get the moon period as defined for MAGIC observation we take the // nearest integer mjd, eg: // 53257.8 --> 53258 // 53258.3 --> 53258 // Which is the time between 13h and 12:59h of the following day. To // this day-period we assign the moon-period at midnight. To get // the MAGIC definition we now substract 284. // // For MAGIC observation period do eg: // GetMagicPeriod(53257.91042) // or // MTime t; // t.SetMjd(53257.91042); // GetMagicPeriod(t.GetMjd()); // or // MTime t; // t.Set(2004, 1, 1, 12, 32, 11); // GetMagicPeriod(t.GetMjd()); // Int_t MAstro::GetMagicPeriod(Double_t mjd) { const Double_t mmjd = (Double_t)TMath::Nint(mjd); const Double_t period = GetMoonPeriod(mmjd); return (Int_t)TMath::Floor(period)-284; } // -------------------------------------------------------------------------- // // Returns the distance in x,y between two polar-vectors (eg. Alt/Az, Ra/Dec) // projected on aplain in a distance dist. For Magic this this the distance // of the camera plain (1700mm) dist also determins the unit in which // the TVector2 is returned. // // v0 is the reference vector (eg. the vector to the center of the camera) // v1 is the vector to which we determin the distance on the plain // // (see also MStarCamTrans::Loc0LocToCam()) // TVector2 MAstro::GetDistOnPlain(const TVector3 &v0, TVector3 v1, Double_t dist) { v1.RotateZ(-v0.Phi()); v1.RotateY(-v0.Theta()); v1.RotateZ(-TMath::Pi()/2); // exchange x and y v1 *= dist/v1.Z(); return v1.XYvector(); //TVector2(v1.Y(), -v1.X());//v1.XYvector(); }