source: trunk/MagicSoft/Mars/mastro/MAstro.cc@ 4795

Last change on this file since 4795 was 4566, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 16.2 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 11/2003 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MAstro
28// ------
29//
30////////////////////////////////////////////////////////////////////////////
31#include "MAstro.h"
32
33#include <iostream>
34
35#include <TVector3.h> // TVector3
36
37#include "MTime.h" // MTime::GetGmst
38
39using namespace std;
40
41ClassImp(MAstro);
42
43Double_t MAstro::Trunc(Double_t val)
44{
45 // dint(A) - truncate to nearest whole number towards zero (double)
46 return val<0 ? TMath::Ceil(val) : TMath::Floor(val);
47}
48
49Double_t MAstro::Round(Double_t val)
50{
51 // dnint(A) - round to nearest whole number (double)
52 return val<0 ? TMath::Ceil(val-0.5) : TMath::Floor(val+0.5);
53}
54
55Double_t MAstro::Hms2Sec(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
56{
57 const Double_t rc = TMath::Sign((60.0 * (60.0 * (Double_t)TMath::Abs(deg) + (Double_t)min) + sec), (Double_t)deg);
58 return sgn=='-' ? -rc : rc;
59}
60
61Double_t MAstro::Dms2Rad(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
62{
63 // pi/(180*3600): arcseconds to radians
64 //#define DAS2R 4.8481368110953599358991410235794797595635330237270e-6
65 return Hms2Sec(deg, min, sec, sgn)*TMath::Pi()/(180*3600)/**DAS2R*/;
66}
67
68Double_t MAstro::Hms2Rad(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
69{
70 // pi/(12*3600): seconds of time to radians
71//#define DS2R 7.2722052166430399038487115353692196393452995355905e-5
72 return Hms2Sec(hor, min, sec, sgn)*TMath::Pi()/(12*3600)/**DS2R*/;
73}
74
75Double_t MAstro::Dms2Deg(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
76{
77 return Hms2Sec(deg, min, sec, sgn)/3600.;
78}
79
80Double_t MAstro::Hms2Deg(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
81{
82 return Hms2Sec(hor, min, sec, sgn)/240.;
83}
84
85Double_t MAstro::Dms2Hor(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
86{
87 return Hms2Sec(deg, min, sec, sgn)/54000.;
88}
89
90Double_t MAstro::Hms2Hor(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
91{
92 return Hms2Sec(hor, min, sec, sgn)/3600.;
93}
94
95void MAstro::Day2Hms(Double_t day, Char_t &sgn, UShort_t &hor, UShort_t &min, UShort_t &sec)
96{
97 /* Handle sign */
98 sgn = day<0?'-':'+';
99
100 /* Round interval and express in smallest units required */
101 Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds
102
103 /* Separate into fields */
104 const Double_t ah = Trunc(a/3600.);
105 a -= ah * 3600.;
106 const Double_t am = Trunc(a/60.);
107 a -= am * 60.;
108 const Double_t as = Trunc(a);
109
110 /* Return results */
111 hor = (UShort_t)ah;
112 min = (UShort_t)am;
113 sec = (UShort_t)as;
114}
115
116void MAstro::Rad2Hms(Double_t rad, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
117{
118 Day2Hms(rad/(TMath::Pi()*2), sgn, deg, min, sec);
119}
120
121void MAstro::Rad2Dms(Double_t rad, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
122{
123 Rad2Hms(rad*15, sgn, deg, min, sec);
124}
125
126void MAstro::Deg2Dms(Double_t d, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
127{
128 Day2Hms(d/24, sgn, deg, min, sec);
129}
130
131void MAstro::Deg2Hms(Double_t d, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
132{
133 Rad2Hms(d/360, sgn, deg, min, sec);
134}
135
136void MAstro::Hor2Dms(Double_t h, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
137{
138 Day2Hms(h*15/24, sgn, deg, min, sec);
139}
140
141void MAstro::Hor2Hms(Double_t h, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
142{
143 Day2Hms(h/24, sgn, deg, min, sec);
144}
145
146void MAstro::Day2Hm(Double_t day, Char_t &sgn, UShort_t &hor, Double_t &min)
147{
148 /* Handle sign */
149 sgn = day<0?'-':'+';
150
151 /* Round interval and express in smallest units required */
152 Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds
153
154 /* Separate into fields */
155 const Double_t ah = Trunc(a/3600.);
156 a -= ah * 3600.;
157
158 /* Return results */
159 hor = (UShort_t)ah;
160 min = a/60.;
161}
162
163void MAstro::Rad2Hm(Double_t rad, Char_t &sgn, UShort_t &deg, Double_t &min)
164{
165 Day2Hm(rad/(TMath::Pi()*2), sgn, deg, min);
166}
167
168void MAstro::Rad2Dm(Double_t rad, Char_t &sgn, UShort_t &deg, Double_t &min)
169{
170 Rad2Hm(rad*15, sgn, deg, min);
171}
172
173void MAstro::Deg2Dm(Double_t d, Char_t &sgn, UShort_t &deg, Double_t &min)
174{
175 Day2Hm(d/24, sgn, deg, min);
176}
177
178void MAstro::Deg2Hm(Double_t d, Char_t &sgn, UShort_t &deg, Double_t &min)
179{
180 Rad2Hm(d/360, sgn, deg, min);
181}
182
183void MAstro::Hor2Dm(Double_t h, Char_t &sgn, UShort_t &deg, Double_t &min)
184{
185 Day2Hm(h*15/24, sgn, deg, min);
186}
187
188void MAstro::Hor2Hm(Double_t h, Char_t &sgn, UShort_t &deg, Double_t &min)
189{
190 Day2Hm(h/24, sgn, deg, min);
191}
192
193// --------------------------------------------------------------------------
194//
195// Interpretes a string ' - 12 30 00.0' or '+ 12 30 00.0'
196// as floating point value -12.5 or 12.5. If interpretation is
197// successfull kTRUE is returned, otherwise kFALSE. ret is not
198// touched if interpretation was not successfull. The successfull
199// interpreted part is removed from the TString.
200//
201Bool_t MAstro::String2Angle(TString &str, Double_t &ret)
202{
203 Char_t sgn;
204 Int_t d, len;
205 UInt_t m;
206 Float_t s;
207
208 // Skip whitespaces before %c and after %f
209 int n=sscanf(str.Data(), " %c %d %d %f %n", &sgn, &d, &m, &s, &len);
210
211 if (n!=4 || (sgn!='+' && sgn!='-'))
212 return kFALSE;
213
214 str.Remove(0, len);
215
216 ret = Dms2Deg(d, m, s, sgn);
217 return kTRUE;
218}
219
220// --------------------------------------------------------------------------
221//
222// Interpretes a string '-12:30:00.0', '12:30:00.0' or '+12:30:00.0'
223// as floating point value -12.5, 12.5 or 12.5. If interpretation is
224// successfull kTRUE is returned, otherwise kFALSE. ret is not
225// touched if interpretation was not successfull.
226//
227Bool_t MAstro::Coordinate2Angle(const TString &str, Double_t &ret)
228{
229 Char_t sgn = str[0]=='-' ? '-' : '+';
230 Int_t d;
231 UInt_t m;
232 Float_t s;
233
234 const int n=sscanf(str[0]=='+'||str[0]=='-' ? str.Data()+1 : str.Data(), "%d:%d:%f", &d, &m, &s);
235
236 if (n!=3)
237 return kFALSE;
238
239 ret = Dms2Deg(d, m, s, sgn);
240 return kTRUE;
241}
242
243// --------------------------------------------------------------------------
244//
245// Return year y, month m and day d corresponding to Mjd.
246//
247void MAstro::Mjd2Ymd(UInt_t mjd, UShort_t &y, Byte_t &m, Byte_t &d)
248{
249 // Express day in Gregorian calendar
250 const ULong_t jd = mjd + 2400001;
251 const ULong_t n4 = 4*(jd+((6*((4*jd-17918)/146097))/4+1)/2-37);
252 const ULong_t nd10 = 10*(((n4-237)%1461)/4)+5;
253
254 y = n4/1461L-4712;
255 m = ((nd10/306+2)%12)+1;
256 d = (nd10%306)/10+1;
257}
258
259// --------------------------------------------------------------------------
260//
261// Return Mjd corresponding to year y, month m and day d.
262//
263Int_t MAstro::Ymd2Mjd(UShort_t y, Byte_t m, Byte_t d)
264{
265 // Month lengths in days
266 static int months[12] = { 31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31 };
267
268 // Validate month
269 if (m<1 || m>12)
270 return -1;
271
272 // Allow for leap year
273 months[1] = (y%4==0 && (y%100!=0 || y%400==0)) ? 29 : 28;
274
275 // Validate day
276 if (d<1 || d>months[m-1])
277 return -1;
278
279 // Precalculate some values
280 const Byte_t lm = 12-m;
281 const ULong_t lm10 = 4712 + y - lm/10;
282
283 // Perform the conversion
284 return 1461L*lm10/4 + (306*((m+9)%12)+5)/10 - (3*((lm10+188)/100))/4 + d - 2399904;
285}
286
287// --------------------------------------------------------------------------
288//
289// theta0, phi0 [rad]: polar angle/zenith distance, azimuth of 1st object
290// theta1, phi1 [rad]: polar angle/zenith distance, azimuth of 2nd object
291// AngularDistance [rad]: Angular distance between two objects
292//
293Double_t MAstro::AngularDistance(Double_t theta0, Double_t phi0, Double_t theta1, Double_t phi1)
294{
295 TVector3 v0(1);
296 v0.Rotate(phi0, theta0);
297
298 TVector3 v1(1);
299 v1.Rotate(phi1, theta1);
300
301 return v0.Angle(v1);
302}
303
304// --------------------------------------------------------------------------
305//
306// Calls MTime::GetGmst() Better use MTime::GetGmst() directly
307//
308Double_t MAstro::UT2GMST(Double_t ut1)
309{
310 return MTime(ut1).GetGmst();
311}
312
313// --------------------------------------------------------------------------
314//
315// RotationAngle
316//
317// calculates the angle for the rotation of the sky coordinate system
318// with respect to the local coordinate system. This is identical
319// to the rotation angle of the sky image in the camera.
320//
321// sinl [rad]: sine of observers latitude
322// cosl [rad]: cosine of observers latitude
323// theta [rad]: polar angle/zenith distance
324// phi [rad]: rotation angle/azimuth
325//
326// Return sin/cos component of angle
327//
328// The convention is such, that the rotation angle is -pi/pi if
329// right ascension and local rotation angle are counted in the
330// same direction, 0 if counted in the opposite direction.
331//
332// (In other words: The rotation angle is 0 when the source culminates)
333//
334// Using vectors it can be done like:
335// TVector3 v, p;
336// v.SetMagThetaPhi(1, theta, phi);
337// p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0);
338// v = v.Cross(l));
339// v.RotateZ(-phi);
340// v.Rotate(-theta)
341// rho = TMath::ATan2(v(2), v(1));
342//
343// For more information see TDAS 00-11, eqs. (18) and (20)
344//
345void MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi, Double_t &sin, Double_t &cos)
346{
347 const Double_t sint = TMath::Sin(theta);
348 const Double_t cost = TMath::Cos(theta);
349
350 const Double_t snlt = sinl*sint;
351 const Double_t cslt = cosl*cost;
352
353 const Double_t sinp = TMath::Sin(phi);
354 const Double_t cosp = TMath::Cos(phi);
355
356 const Double_t v1 = sint*sinp;
357 const Double_t v2 = cslt - snlt*cosp;
358
359 const Double_t denom = TMath::Sqrt(v1*v1 + v2*v2);
360
361 sin = cosl*sinp / denom; // y-component
362 cos = (snlt-cslt*cosp) / denom; // x-component
363}
364
365// --------------------------------------------------------------------------
366//
367// RotationAngle
368//
369// calculates the angle for the rotation of the sky coordinate system
370// with respect to the local coordinate system. This is identical
371// to the rotation angle of the sky image in the camera.
372//
373// sinl [rad]: sine of observers latitude
374// cosl [rad]: cosine of observers latitude
375// theta [rad]: polar angle/zenith distance
376// phi [rad]: rotation angle/azimuth
377//
378// Return angle [rad] in the range -pi, pi
379//
380// The convention is such, that the rotation angle is -pi/pi if
381// right ascension and local rotation angle are counted in the
382// same direction, 0 if counted in the opposite direction.
383//
384// (In other words: The rotation angle is 0 when the source culminates)
385//
386// Using vectors it can be done like:
387// TVector3 v, p;
388// v.SetMagThetaPhi(1, theta, phi);
389// p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0);
390// v = v.Cross(l));
391// v.RotateZ(-phi);
392// v.Rotate(-theta)
393// rho = TMath::ATan2(v(2), v(1));
394//
395// For more information see TDAS 00-11, eqs. (18) and (20)
396//
397Double_t MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi)
398{
399 const Double_t snlt = sinl*TMath::Sin(theta);
400 const Double_t cslt = cosl*TMath::Cos(theta);
401
402 const Double_t sinp = TMath::Sin(phi);
403 const Double_t cosp = TMath::Cos(phi);
404
405 return TMath::ATan2(cosl*sinp, snlt-cslt*cosp);
406}
407
408
409// --------------------------------------------------------------------------
410//
411// Kepler - solve the equation of Kepler
412//
413Double_t MAstro::Kepler(Double_t m, Double_t ecc)
414{
415 m *= TMath::DegToRad();
416
417 Double_t delta = 0;
418 Double_t e = m;
419 do {
420 delta = e - ecc * sin(e) - m;
421 e -= delta / (1 - ecc * cos(e));
422 } while (fabs(delta) > 1e-6);
423
424 return e;
425}
426
427// --------------------------------------------------------------------------
428//
429// GetMoonPhase - calculate phase of moon as a fraction:
430// Returns -1 if calculation failed
431//
432Double_t MAstro::GetMoonPhase(Double_t mjd)
433{
434 /****** Calculation of the Sun's position. ******/
435
436 // date within epoch
437 const Double_t epoch = 44238; // 1980 January 0.0
438 const Double_t day = mjd - epoch;
439 if (day<0)
440 {
441 cout << "MAstro::GetMoonPhase - Day before Jan 1980" << endl;
442 return -1;
443 }
444
445 // mean anomaly of the Sun
446 const Double_t n = fmod(day*360/365.2422, 360);
447
448 const Double_t elonge = 278.833540; // ecliptic longitude of the Sun at epoch 1980.0
449 const Double_t elongp = 282.596403; // ecliptic longitude of the Sun at perigee
450
451 // convert from perigee co-ordinates to epoch 1980.0
452 const Double_t m = fmod(n + elonge - elongp + 360, 360);
453
454 // solve equation of Kepler
455 const Double_t eccent = 0.016718; // eccentricity of Earth's orbit
456 const Double_t k = Kepler(m, eccent);
457 const Double_t ec0 = sqrt((1 + eccent) / (1 - eccent)) * tan(k / 2);
458 // true anomaly
459 const Double_t ec = 2 * atan(ec0) * TMath::RadToDeg();
460
461 // Sun's geocentric ecliptic longitude
462 const Double_t lambdasun = fmod(ec + elongp + 720, 360);
463
464
465 /****** Calculation of the Moon's position. ******/
466
467 // Moon's mean longitude.
468 const Double_t mmlong = 64.975464; // moon's mean lonigitude at the epoch
469 const Double_t ml = fmod(13.1763966*day + mmlong + 360, 360);
470 // Moon's mean anomaly.
471 const Double_t mmlongp = 349.383063; // mean longitude of the perigee at the epoch
472 const Double_t mm = fmod(ml - 0.1114041*day - mmlongp + 720, 360);
473 // Evection.
474 const Double_t ev = 1.2739 * sin((2 * (ml - lambdasun) - mm)*TMath::DegToRad());
475 // Annual equation.
476 const Double_t sinm = TMath::Sin(m*TMath::DegToRad());
477 const Double_t ae = 0.1858 * sinm;
478 // Correction term.
479 const Double_t a3 = 0.37 * sinm;
480 // Corrected anomaly.
481 const Double_t mmp = (mm + ev - ae - a3)*TMath::DegToRad();
482 // Correction for the equation of the centre.
483 const Double_t mec = 6.2886 * sin(mmp);
484 // Another correction term.
485 const Double_t a4 = 0.214 * sin(2 * mmp);
486 // Corrected longitude.
487 const Double_t lp = ml + ev + mec - ae + a4;
488 // Variation.
489 const Double_t v = 0.6583 * sin(2 * (lp - lambdasun)*TMath::DegToRad());
490 // True longitude.
491 const Double_t lpp = lp + v;
492 // Age of the Moon in degrees.
493 const Double_t age = (lpp - lambdasun)*TMath::DegToRad();
494
495 // Calculation of the phase of the Moon.
496 return (1 - TMath::Cos(age)) / 2;
497}
498
499// --------------------------------------------------------------------------
500//
501// Calculate the Period to which the time belongs to. The Period is defined
502// as the number of synodic months ellapsed since the first full moon
503// after Jan 1st 1980 (which was @ MJD=44240.37917)
504//
505Double_t MAstro::GetMoonPeriod(Double_t mjd)
506{
507 const Double_t synmonth = 29.53058868; // synodic month (new Moon to new Moon)
508 const Double_t epoch0 = 44240.37917; // First full moon after 1980/1/1
509
510 const Double_t et = mjd-epoch0; // Ellapsed time
511 return et/synmonth;
512}
513
514// --------------------------------------------------------------------------
515//
516// To get the moon period as defined for MAGIC observation we take the
517// nearest integer mjd, eg:
518// 53257.8 --> 53258
519// 53258.3 --> 53258
520// Which is the time between 13h and 12:59h of the following day. To
521// this day-period we assign the moon-period at midnight. To get
522// the MAGIC definition we now substract 284.
523//
524// For MAGIC observation period do eg:
525// GetMagicPeriod(53257.91042)
526// or
527// MTime t;
528// t.SetMjd(53257.91042);
529// GetMagicPeriod(t.GetMjd());
530// or
531// MTime t;
532// t.Set(2004, 1, 1, 12, 32, 11);
533// GetMagicPeriod(t.GetMjd());
534//
535Int_t MAstro::GetMagicPeriod(Double_t mjd)
536{
537 const Double_t mmjd = (Double_t)TMath::Nint(mjd);
538 const Double_t period = GetMoonPeriod(mmjd);
539
540 return (Int_t)TMath::Floor(period)-284;
541}
Note: See TracBrowser for help on using the repository browser.