source: tags/Mars-V1.1/mastro/MAstro.cc

Last change on this file was 8324, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 31.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 <TArrayD.h> // TArrayD
36#include <TVector3.h> // TVector3
37
38#include "MTime.h" // MTime::GetGmst
39#include "MString.h"
40
41#include "MAstroCatalog.h" // FIXME: replace by MVector3!
42
43using namespace std;
44
45ClassImp(MAstro);
46
47Double_t MAstro::Trunc(Double_t val)
48{
49 // dint(A) - truncate to nearest whole number towards zero (double)
50 return val<0 ? TMath::Ceil(val) : TMath::Floor(val);
51}
52
53Double_t MAstro::Round(Double_t val)
54{
55 // dnint(A) - round to nearest whole number (double)
56 return val<0 ? TMath::Ceil(val-0.5) : TMath::Floor(val+0.5);
57}
58
59Double_t MAstro::Hms2Sec(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
60{
61 const Double_t rc = TMath::Sign((60.0 * (60.0 * (Double_t)TMath::Abs(deg) + (Double_t)min) + sec), (Double_t)deg);
62 return sgn=='-' ? -rc : rc;
63}
64
65Double_t MAstro::Dms2Rad(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
66{
67 // pi/(180*3600): arcseconds to radians
68 //#define DAS2R 4.8481368110953599358991410235794797595635330237270e-6
69 return Hms2Sec(deg, min, sec, sgn)*TMath::Pi()/(180*3600)/**DAS2R*/;
70}
71
72Double_t MAstro::Hms2Rad(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
73{
74 // pi/(12*3600): seconds of time to radians
75//#define DS2R 7.2722052166430399038487115353692196393452995355905e-5
76 return Hms2Sec(hor, min, sec, sgn)*TMath::Pi()/(12*3600)/**DS2R*/;
77}
78
79Double_t MAstro::Dms2Deg(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
80{
81 return Hms2Sec(deg, min, sec, sgn)/3600.;
82}
83
84Double_t MAstro::Hms2Deg(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
85{
86 return Hms2Sec(hor, min, sec, sgn)/240.;
87}
88
89Double_t MAstro::Dms2Hor(Int_t deg, UInt_t min, Double_t sec, Char_t sgn)
90{
91 return Hms2Sec(deg, min, sec, sgn)/54000.;
92}
93
94Double_t MAstro::Hms2Hor(Int_t hor, UInt_t min, Double_t sec, Char_t sgn)
95{
96 return Hms2Sec(hor, min, sec, sgn)/3600.;
97}
98
99void MAstro::Day2Hms(Double_t day, Char_t &sgn, UShort_t &hor, UShort_t &min, UShort_t &sec)
100{
101 /* Handle sign */
102 sgn = day<0?'-':'+';
103
104 /* Round interval and express in smallest units required */
105 Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds
106
107 /* Separate into fields */
108 const Double_t ah = Trunc(a/3600.);
109 a -= ah * 3600.;
110 const Double_t am = Trunc(a/60.);
111 a -= am * 60.;
112 const Double_t as = Trunc(a);
113
114 /* Return results */
115 hor = (UShort_t)ah;
116 min = (UShort_t)am;
117 sec = (UShort_t)as;
118}
119
120void MAstro::Rad2Hms(Double_t rad, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
121{
122 Day2Hms(rad/(TMath::Pi()*2), sgn, deg, min, sec);
123}
124
125void MAstro::Rad2Dms(Double_t rad, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
126{
127 Rad2Hms(rad*15, sgn, deg, min, sec);
128}
129
130void MAstro::Deg2Dms(Double_t d, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
131{
132 Day2Hms(d/24, sgn, deg, min, sec);
133}
134
135void MAstro::Deg2Hms(Double_t d, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
136{
137 Day2Hms(d/360, sgn, deg, min, sec);
138}
139
140void MAstro::Hor2Dms(Double_t h, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
141{
142 Day2Hms(h*15/24, sgn, deg, min, sec);
143}
144
145void MAstro::Hor2Hms(Double_t h, Char_t &sgn, UShort_t &deg, UShort_t &min, UShort_t &sec)
146{
147 Day2Hms(h/24, sgn, deg, min, sec);
148}
149
150void MAstro::Day2Hm(Double_t day, Char_t &sgn, UShort_t &hor, Double_t &min)
151{
152 /* Handle sign */
153 sgn = day<0?'-':'+';
154
155 /* Round interval and express in smallest units required */
156 Double_t a = Round(86400. * TMath::Abs(day)); // Days to seconds
157
158 /* Separate into fields */
159 const Double_t ah = Trunc(a/3600.);
160 a -= ah * 3600.;
161
162 /* Return results */
163 hor = (UShort_t)ah;
164 min = a/60.;
165}
166
167void MAstro::Rad2Hm(Double_t rad, Char_t &sgn, UShort_t &deg, Double_t &min)
168{
169 Day2Hm(rad/(TMath::Pi()*2), sgn, deg, min);
170}
171
172void MAstro::Rad2Dm(Double_t rad, Char_t &sgn, UShort_t &deg, Double_t &min)
173{
174 Rad2Hm(rad*15, sgn, deg, min);
175}
176
177void MAstro::Deg2Dm(Double_t d, Char_t &sgn, UShort_t &deg, Double_t &min)
178{
179 Day2Hm(d/24, sgn, deg, min);
180}
181
182void MAstro::Deg2Hm(Double_t d, Char_t &sgn, UShort_t &deg, Double_t &min)
183{
184 Rad2Hm(d/360, sgn, deg, min);
185}
186
187void MAstro::Hor2Dm(Double_t h, Char_t &sgn, UShort_t &deg, Double_t &min)
188{
189 Day2Hm(h*15/24, sgn, deg, min);
190}
191
192void MAstro::Hor2Hm(Double_t h, Char_t &sgn, UShort_t &deg, Double_t &min)
193{
194 Day2Hm(h/24, sgn, deg, min);
195}
196
197TString MAstro::GetStringDeg(Double_t deg, const char *fmt)
198{
199 Char_t sgn;
200 UShort_t d, m, s;
201 Deg2Dms(deg, sgn, d, m, s);
202
203 MString str;
204 str.Print(fmt, sgn, d, m ,s);
205 return str;
206}
207
208TString MAstro::GetStringHor(Double_t deg, const char *fmt)
209{
210 Char_t sgn;
211 UShort_t h, m, s;
212 Hor2Hms(deg, sgn, h, m, s);
213
214 MString str;
215 str.Print(fmt, sgn, h, m ,s);
216 return str;
217}
218
219// --------------------------------------------------------------------------
220//
221// Interpretes a string ' - 12 30 00.0' or '+ 12 30 00.0'
222// as floating point value -12.5 or 12.5. If interpretation is
223// successfull kTRUE is returned, otherwise kFALSE. ret is not
224// touched if interpretation was not successfull. The successfull
225// interpreted part is removed from the TString.
226//
227Bool_t MAstro::String2Angle(TString &str, Double_t &ret)
228{
229 Char_t sgn;
230 Int_t d, len;
231 UInt_t m;
232 Float_t s;
233
234 // Skip whitespaces before %c and after %f
235 int n=sscanf(str.Data(), " %c %d %d %f %n", &sgn, &d, &m, &s, &len);
236
237 if (n!=4 || (sgn!='+' && sgn!='-'))
238 return kFALSE;
239
240 str.Remove(0, len);
241
242 ret = Dms2Deg(d, m, s, sgn);
243 return kTRUE;
244}
245
246// --------------------------------------------------------------------------
247//
248// Interpretes a string '-12:30:00.0', '12:30:00.0' or '+12:30:00.0'
249// as floating point value -12.5, 12.5 or 12.5. If interpretation is
250// successfull kTRUE is returned, otherwise kFALSE. ret is not
251// touched if interpretation was not successfull.
252//
253Bool_t MAstro::Coordinate2Angle(const TString &str, Double_t &ret)
254{
255 Char_t sgn = str[0]=='-' ? '-' : '+';
256 Int_t d;
257 UInt_t m;
258 Float_t s;
259
260 const int n=sscanf(str[0]=='+'||str[0]=='-' ? str.Data()+1 : str.Data(), "%d:%d:%f", &d, &m, &s);
261
262 if (n!=3)
263 return kFALSE;
264
265 ret = Dms2Deg(d, m, s, sgn);
266 return kTRUE;
267}
268
269// --------------------------------------------------------------------------
270//
271// Returns val=-12.5 as string '-12:30:00'
272//
273TString MAstro::Angle2Coordinate(Double_t val)
274{
275 Char_t sgn;
276 UShort_t d,m,s;
277
278 Deg2Dms(val, sgn, d, m, s);
279
280 return Form("%c%02d:%02d:%02d", sgn, d, m, s);
281}
282
283// --------------------------------------------------------------------------
284//
285// Return year y, month m and day d corresponding to Mjd.
286//
287void MAstro::Mjd2Ymd(UInt_t mjd, UShort_t &y, Byte_t &m, Byte_t &d)
288{
289 // Express day in Gregorian calendar
290 const ULong_t jd = mjd + 2400001;
291 const ULong_t n4 = 4*(jd+((6*((4*jd-17918)/146097))/4+1)/2-37);
292 const ULong_t nd10 = 10*(((n4-237)%1461)/4)+5;
293
294 y = n4/1461L-4712;
295 m = ((nd10/306+2)%12)+1;
296 d = (nd10%306)/10+1;
297}
298
299// --------------------------------------------------------------------------
300//
301// Return Mjd corresponding to year y, month m and day d.
302//
303Int_t MAstro::Ymd2Mjd(UShort_t y, Byte_t m, Byte_t d)
304{
305 // Month lengths in days
306 static int months[12] = { 31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31 };
307
308 // Validate month
309 if (m<1 || m>12)
310 return -1;
311
312 // Allow for leap year
313 months[1] = (y%4==0 && (y%100!=0 || y%400==0)) ? 29 : 28;
314
315 // Validate day
316 if (d<1 || d>months[m-1])
317 return -1;
318
319 // Precalculate some values
320 const Byte_t lm = 12-m;
321 const ULong_t lm10 = 4712 + y - lm/10;
322
323 // Perform the conversion
324 return 1461L*lm10/4 + (306*((m+9)%12)+5)/10 - (3*((lm10+188)/100))/4 + d - 2399904;
325}
326
327// --------------------------------------------------------------------------
328//
329// theta0, phi0 [rad]: polar angle/zenith distance, azimuth of 1st object
330// theta1, phi1 [rad]: polar angle/zenith distance, azimuth of 2nd object
331// AngularDistance [rad]: Angular distance between two objects
332//
333Double_t MAstro::AngularDistance(Double_t theta0, Double_t phi0, Double_t theta1, Double_t phi1)
334{
335 TVector3 v0(1);
336 v0.Rotate(phi0, theta0);
337
338 TVector3 v1(1);
339 v1.Rotate(phi1, theta1);
340
341 return v0.Angle(v1);
342}
343
344// --------------------------------------------------------------------------
345//
346// Calls MTime::GetGmst() Better use MTime::GetGmst() directly
347//
348Double_t MAstro::UT2GMST(Double_t ut1)
349{
350 return MTime(ut1).GetGmst();
351}
352
353// --------------------------------------------------------------------------
354//
355// RotationAngle
356//
357// calculates the angle for the rotation of the sky coordinate system
358// with respect to the local coordinate system. This is identical
359// to the rotation angle of the sky image in the camera.
360//
361// sinl [rad]: sine of observers latitude
362// cosl [rad]: cosine of observers latitude
363// theta [rad]: polar angle/zenith distance
364// phi [rad]: rotation angle/azimuth
365//
366// Return sin/cos component of angle
367//
368// The convention is such, that the rotation angle is -pi/pi if
369// right ascension and local rotation angle are counted in the
370// same direction, 0 if counted in the opposite direction.
371//
372// (In other words: The rotation angle is 0 when the source culminates)
373//
374// Using vectors it can be done like:
375// TVector3 v, p;
376// v.SetMagThetaPhi(1, theta, phi);
377// p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0);
378// v = v.Cross(l));
379// v.RotateZ(-phi);
380// v.Rotate(-theta)
381// rho = TMath::ATan2(v(2), v(1));
382//
383// For more information see TDAS 00-11, eqs. (18) and (20)
384//
385void MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi, Double_t &sin, Double_t &cos)
386{
387 const Double_t sint = TMath::Sin(theta);
388 const Double_t cost = TMath::Cos(theta);
389
390 const Double_t snlt = sinl*sint;
391 const Double_t cslt = cosl*cost;
392
393 const Double_t sinp = TMath::Sin(phi);
394 const Double_t cosp = TMath::Cos(phi);
395
396 const Double_t v1 = sint*sinp;
397 const Double_t v2 = cslt - snlt*cosp;
398
399 const Double_t denom = TMath::Sqrt(v1*v1 + v2*v2);
400
401 sin = cosl*sinp / denom; // y-component
402 cos = (snlt-cslt*cosp) / denom; // x-component
403}
404
405// --------------------------------------------------------------------------
406//
407// RotationAngle
408//
409// calculates the angle for the rotation of the sky coordinate system
410// with respect to the local coordinate system. This is identical
411// to the rotation angle of the sky image in the camera.
412//
413// sinl [rad]: sine of observers latitude
414// cosl [rad]: cosine of observers latitude
415// theta [rad]: polar angle/zenith distance
416// phi [rad]: rotation angle/azimuth
417//
418// Return angle [rad] in the range -pi, pi
419//
420// The convention is such, that the rotation angle is -pi/pi if
421// right ascension and local rotation angle are counted in the
422// same direction, 0 if counted in the opposite direction.
423//
424// (In other words: The rotation angle is 0 when the source culminates)
425//
426// Using vectors it can be done like:
427// TVector3 v, p;
428// v.SetMagThetaPhi(1, theta, phi);
429// p.SetMagThetaPhi(1, TMath::Pi()/2-latitude, 0);
430// v = v.Cross(l));
431// v.RotateZ(-phi);
432// v.Rotate(-theta)
433// rho = TMath::ATan2(v(2), v(1));
434//
435// For more information see TDAS 00-11, eqs. (18) and (20)
436//
437Double_t MAstro::RotationAngle(Double_t sinl, Double_t cosl, Double_t theta, Double_t phi)
438{
439 const Double_t snlt = sinl*TMath::Sin(theta);
440 const Double_t cslt = cosl*TMath::Cos(theta);
441
442 const Double_t sinp = TMath::Sin(phi);
443 const Double_t cosp = TMath::Cos(phi);
444
445 return TMath::ATan2(cosl*sinp, snlt-cslt*cosp);
446}
447
448// --------------------------------------------------------------------------
449//
450// Estimates the time at which a source culminates.
451//
452// ra: right ascension [rad]
453// elong: observers longitude [rad]
454// mjd: modified julian date (utc)
455//
456// return time in [-12;12]
457//
458Double_t MAstro::EstimateCulminationTime(Double_t mjd, Double_t elong, Double_t ra)
459{
460 // startime at 1.1.2000 for greenwich 0h
461 const Double_t gmt0 = 6.664520;
462
463 // difference of startime for greenwich for two calendar days [h]
464 const Double_t d0 = 0.06570982224;
465
466 // mjd of greenwich 1.1.2000 0h
467 const Double_t mjd0 = 51544;
468
469 // mjd today
470 const Double_t mjd1 = TMath::Floor(mjd);
471
472 // scale between star-time and sun-time
473 const Double_t scale = 1;//1.00273790926;
474
475 const Double_t UT = (ra-elong)*RadToHor() - (gmt0 + d0 * (mjd1-mjd0))/scale;
476
477 return fmod(2412 + UT, 24) - 12;
478}
479
480// --------------------------------------------------------------------------
481//
482// Kepler - solve the equation of Kepler
483//
484Double_t MAstro::Kepler(Double_t m, Double_t ecc)
485{
486 m *= TMath::DegToRad();
487
488 Double_t delta = 0;
489 Double_t e = m;
490 do {
491 delta = e - ecc * sin(e) - m;
492 e -= delta / (1 - ecc * cos(e));
493 } while (fabs(delta) > 1e-6);
494
495 return e;
496}
497
498// --------------------------------------------------------------------------
499//
500// GetMoonPhase - calculate phase of moon as a fraction:
501// Returns -1 if calculation failed
502//
503Double_t MAstro::GetMoonPhase(Double_t mjd)
504{
505 /****** Calculation of the Sun's position. ******/
506
507 // date within epoch
508 const Double_t epoch = 44238; // 1980 January 0.0
509 const Double_t day = mjd - epoch;
510 if (day<0)
511 {
512 cout << "MAstro::GetMoonPhase - Day before Jan 1980" << endl;
513 return -1;
514 }
515
516 // mean anomaly of the Sun
517 const Double_t n = fmod(day*360/365.2422, 360);
518
519 const Double_t elonge = 278.833540; // ecliptic longitude of the Sun at epoch 1980.0
520 const Double_t elongp = 282.596403; // ecliptic longitude of the Sun at perigee
521
522 // convert from perigee co-ordinates to epoch 1980.0
523 const Double_t m = fmod(n + elonge - elongp + 360, 360);
524
525 // solve equation of Kepler
526 const Double_t eccent = 0.016718; // eccentricity of Earth's orbit
527 const Double_t k = Kepler(m, eccent);
528 const Double_t ec0 = sqrt((1 + eccent) / (1 - eccent)) * tan(k / 2);
529 // true anomaly
530 const Double_t ec = 2 * atan(ec0) * TMath::RadToDeg();
531
532 // Sun's geocentric ecliptic longitude
533 const Double_t lambdasun = fmod(ec + elongp + 720, 360);
534
535
536 /****** Calculation of the Moon's position. ******/
537
538 // Moon's mean longitude.
539 const Double_t mmlong = 64.975464; // moon's mean lonigitude at the epoch
540 const Double_t ml = fmod(13.1763966*day + mmlong + 360, 360);
541 // Moon's mean anomaly.
542 const Double_t mmlongp = 349.383063; // mean longitude of the perigee at the epoch
543 const Double_t mm = fmod(ml - 0.1114041*day - mmlongp + 720, 360);
544 // Evection.
545 const Double_t ev = 1.2739 * sin((2 * (ml - lambdasun) - mm)*TMath::DegToRad());
546 // Annual equation.
547 const Double_t sinm = TMath::Sin(m*TMath::DegToRad());
548 const Double_t ae = 0.1858 * sinm;
549 // Correction term.
550 const Double_t a3 = 0.37 * sinm;
551 // Corrected anomaly.
552 const Double_t mmp = (mm + ev - ae - a3)*TMath::DegToRad();
553 // Correction for the equation of the centre.
554 const Double_t mec = 6.2886 * sin(mmp);
555 // Another correction term.
556 const Double_t a4 = 0.214 * sin(2 * mmp);
557 // Corrected longitude.
558 const Double_t lp = ml + ev + mec - ae + a4;
559 // Variation.
560 const Double_t v = 0.6583 * sin(2 * (lp - lambdasun)*TMath::DegToRad());
561 // True longitude.
562 const Double_t lpp = lp + v;
563 // Age of the Moon in degrees.
564 const Double_t age = (lpp - lambdasun)*TMath::DegToRad();
565
566 // Calculation of the phase of the Moon.
567 return (1 - TMath::Cos(age)) / 2;
568}
569
570// --------------------------------------------------------------------------
571//
572// Calculate the Period to which the time belongs to. The Period is defined
573// as the number of synodic months ellapsed since the first full moon
574// after Jan 1st 1980 (which was @ MJD=44240.37917)
575//
576Double_t MAstro::GetMoonPeriod(Double_t mjd)
577{
578 const Double_t synmonth = 29.53058868; // synodic month (new Moon to new Moon)
579 const Double_t epoch0 = 44240.37917; // First full moon after 1980/1/1
580
581 const Double_t et = mjd-epoch0; // Ellapsed time
582 return et/synmonth;
583}
584
585// --------------------------------------------------------------------------
586//
587// To get the moon period as defined for MAGIC observation we take the
588// nearest integer mjd, eg:
589// 53257.8 --> 53258
590// 53258.3 --> 53258
591// Which is the time between 13h and 12:59h of the following day. To
592// this day-period we assign the moon-period at midnight. To get
593// the MAGIC definition we now substract 284.
594//
595// For MAGIC observation period do eg:
596// GetMagicPeriod(53257.91042)
597// or
598// MTime t;
599// t.SetMjd(53257.91042);
600// GetMagicPeriod(t.GetMjd());
601// or
602// MTime t;
603// t.Set(2004, 1, 1, 12, 32, 11);
604// GetMagicPeriod(t.GetMjd());
605//
606// To get a floating point magic period use
607// GetMoonPeriod(mjd)-284
608//
609Int_t MAstro::GetMagicPeriod(Double_t mjd)
610{
611 const Double_t mmjd = (Double_t)TMath::Nint(mjd);
612 const Double_t period = GetMoonPeriod(mmjd);
613
614 return (Int_t)TMath::Floor(period)-284;
615}
616
617// --------------------------------------------------------------------------
618//
619// Returns right ascension and declination [rad] of the sun at the
620// given mjd (ra, dec).
621//
622// returns the mean longitude [rad].
623//
624// from http://xoomer.alice.it/vtomezzo/sunriset/formulas/index.html
625//
626Double_t MAstro::GetSunRaDec(Double_t mjd, Double_t &ra, Double_t &dec)
627{
628 const Double_t T = (mjd-51544.5)/36525;// + (h-12)/24.0;
629
630 const Double_t T2 = T<0 ? -T*T : T*T;
631 const Double_t T3 = T*T*T;
632
633 // Find the ecliptic longitude of the Sun
634
635 // Geometric mean longitude of the Sun
636 const Double_t L = 280.46646 + 36000.76983*T + 0.0003032*T2;
637
638 // mean anomaly of the Sun
639 Double_t g = 357.52911 + 35999.05029*T - 0.0001537*T2;
640 g *= TMath::DegToRad();
641
642 // Longitude of the moon's ascending node
643 Double_t omega = 125.04452 - 1934.136261*T + 0.0020708*T2 + T3/450000;
644 omega *= TMath::DegToRad();
645
646 const Double_t coso = cos(omega);
647 const Double_t sino = sin(omega);
648
649 // Equation of the center
650 const Double_t C = (1.914602 - 0.004817*T - 0.000014*T2)*sin(g) +
651 (0.019993 - 0.000101*T)*sin(2*g) + 0.000289*sin(3*g);
652
653 // True longitude of the sun
654 const Double_t tlong = L + C;
655
656 // Apperent longitude of the Sun (ecliptic)
657 Double_t lambda = tlong - 0.00569 - 0.00478*sino;
658 lambda *= TMath::DegToRad();
659
660 // Obliquity of the ecliptic
661 Double_t obliq = 23.4392911 - 0.01300416667*T - 0.00000016389*T2 + 0.00000050361*T3 + 0.00255625*coso;
662 obliq *= TMath::DegToRad();
663
664 // Find the RA and DEC of the Sun
665 const Double_t sinl = sin(lambda);
666
667 ra = atan2(cos(obliq) * sinl, cos(lambda));
668 dec = asin(sin(obliq) * sinl);
669
670 return L*TMath::DegToRad();
671}
672
673// --------------------------------------------------------------------------
674//
675// Returns right ascension and declination [rad] of the moon at the
676// given mjd (ra, dec).
677//
678void MAstro::GetMoonRaDec(Double_t mjd, Double_t &ra, Double_t &dec)
679{
680 // Mean Moon orbit elements as of 1990.0
681 const Double_t l0 = 318.351648 * TMath::DegToRad();
682 const Double_t P0 = 36.340410 * TMath::DegToRad();
683 const Double_t N0 = 318.510107 * TMath::DegToRad();
684 const Double_t i = 5.145396 * TMath::DegToRad();
685
686 Double_t sunra, sundec, g;
687 {
688 const Double_t T = (mjd-51544.5)/36525;// + (h-12)/24.0;
689 const Double_t T2 = T<0 ? -T*T : T*T;
690
691 GetSunRaDec(mjd, sunra, sundec);
692
693 // mean anomaly of the Sun
694 g = 357.52911 + 35999.05029*T - 0.0001537*T2;
695 g *= TMath::DegToRad();
696 }
697
698 const Double_t sing = sin(g)*TMath::DegToRad();
699
700 const Double_t D = (mjd-47891) * TMath::DegToRad();
701 const Double_t l = 13.1763966*D + l0;
702 const Double_t MMoon = l -0.1114041*D - P0; // Moon's mean anomaly M
703 const Double_t N = N0 -0.0529539*D; // Moon's mean ascending node longitude
704
705 const Double_t C = l-sunra;
706 const Double_t Ev = 1.2739 * sin(2*C-MMoon) * TMath::DegToRad();
707 const Double_t Ae = 0.1858 * sing;
708 const Double_t A3 = 0.37 * sing;
709 const Double_t MMoon2 = MMoon+Ev-Ae-A3; // corrected Moon anomaly
710
711 const Double_t Ec = 6.2886 * sin(MMoon2) * TMath::DegToRad(); // equation of centre
712 const Double_t A4 = 0.214 * sin(2*MMoon2)* TMath::DegToRad();
713 const Double_t l2 = l+Ev+Ec-Ae+A4; // corrected Moon's longitude
714
715 const Double_t V = 0.6583 * sin(2*(l2-sunra)) * TMath::DegToRad();
716 const Double_t l3 = l2+V; // true orbital longitude;
717
718 const Double_t N2 = N -0.16*sing;
719
720 ra = fmod( N2 + atan2( sin(l3-N2)*cos(i), cos(l3-N2) ), TMath::TwoPi() );
721 dec = asin(sin(l3-N2)*sin(i) );
722}
723
724// --------------------------------------------------------------------------
725//
726// Return Euqation of time in hours for given mjd
727//
728Double_t MAstro::GetEquationOfTime(Double_t mjd)
729{
730 Double_t ra, dec;
731 const Double_t L = fmod(GetSunRaDec(mjd, ra, dec), TMath::TwoPi());
732
733 if (L-ra>TMath::Pi())
734 ra += TMath::TwoPi();
735
736 return 24*(L - ra)/TMath::TwoPi();
737}
738
739// --------------------------------------------------------------------------
740//
741// Returns noon time (the time of the highest altitude of the sun)
742// at the given mjd and at the given observers longitude [deg]
743//
744// The maximum altitude reached at noon time is
745// altmax = 90.0 + dec - latit;
746// if (dec > latit)
747// altmax = 90.0 + latit - dec;
748// dec=Declination of the sun
749//
750Double_t MAstro::GetNoonTime(Double_t mjd, Double_t longit)
751{
752 const Double_t equation = GetEquationOfTime(TMath::Floor(mjd));
753 return 12. + equation - longit/15;
754}
755
756// --------------------------------------------------------------------------
757//
758// Returns the time (in hours) between noon (the sun culmination)
759// and the sun being at height alt[deg] (90=zenith, 0=horizont)
760//
761// civil twilight: 0deg to -6deg
762// nautical twilight: -6deg to -12deg
763// astronom twilight: -12deg to -18deg
764//
765// latit is the observers latitude in rad
766//
767// returns -1 in case the sun doesn't reach this altitude.
768// (eg. alt=0: Polarnight or -day)
769//
770// To get the sun rise/set:
771// double timediff = MAstro::GetTimeFromNoonToAlt(mjd, latit*TMath::DegToRad(), par[0]);
772// double noon = MAstro::GetNoonTime(mjd, longit);
773// double N = TMath::Floor(mjd)+noon/24.;
774// double risetime = N-timediff/24.;
775// double settime = N+timediff/24.;
776//
777Double_t MAstro::GetTimeFromNoonToAlt(Double_t mjd, Double_t latit, Double_t alt)
778{
779 Double_t ra, dec;
780 GetSunRaDec(mjd, ra, dec);
781
782 const Double_t h = alt*TMath::DegToRad();
783
784 const Double_t arg = (sin(h) - sin(latit)*sin(dec))/(cos(latit)*cos(dec));
785
786 return TMath::Abs(arg)>1 ? -1 : 12*acos(arg)/TMath::Pi();
787}
788
789// --------------------------------------------------------------------------
790//
791// Returns the time of the sunrise/set calculated before and after
792// the noon of floor(mjd) (TO BE IMPROVED)
793//
794// Being longit and latit the longitude and latitude of the observer
795// in deg and alt the hight above or below the horizont in deg.
796//
797// civil twilight: 0deg to -6deg
798// nautical twilight: -6deg to -12deg
799// astronom twilight: -12deg to -18deg
800//
801// A TArrayD(2) is returned with the the mjd of the sunrise in
802// TArray[0] and the mjd of the sunset in TArrayD[1].
803//
804TArrayD MAstro::GetSunRiseSet(Double_t mjd, Double_t longit, Double_t latit, Double_t alt)
805{
806 const Double_t timediff = MAstro::GetTimeFromNoonToAlt(mjd, latit*TMath::DegToRad(), alt);
807 const Double_t noon = MAstro::GetNoonTime(mjd, longit);
808
809 const Double_t N = TMath::Floor(mjd)+noon/24.;
810
811 const Double_t rise = timediff<0 ? N-0.5 : N-timediff/24.;
812 const Double_t set = timediff<0 ? N+0.5 : N+timediff/24.;
813
814 TArrayD rc(2);
815 rc[0] = rise;
816 rc[1] = set;
817 return rc;
818}
819
820// --------------------------------------------------------------------------
821//
822// Returns the distance in x,y between two polar-vectors (eg. Alt/Az, Ra/Dec)
823// projected on aplain in a distance dist. For Magic this this the distance
824// of the camera plain (1700mm) dist also determins the unit in which
825// the TVector2 is returned.
826//
827// v0 is the reference vector (eg. the vector to the center of the camera)
828// v1 is the vector to which we determin the distance on the plain
829//
830// (see also MStarCamTrans::Loc0LocToCam())
831//
832TVector2 MAstro::GetDistOnPlain(const TVector3 &v0, TVector3 v1, Double_t dist)
833{
834 v1.RotateZ(-v0.Phi());
835 v1.RotateY(-v0.Theta());
836 v1.RotateZ(-TMath::Pi()/2); // exchange x and y
837 v1 *= dist/v1.Z();
838
839 return v1.XYvector(); //TVector2(v1.Y(), -v1.X());//v1.XYvector();
840}
841
842// --------------------------------------------------------------------------
843//
844// Calculate the absolute misspointing from the nominal zenith angle nomzd
845// and the deviations in zd (devzd) and az (devaz).
846// All values given in deg, the return value, too.
847//
848Double_t MAstro::GetDevAbs(Double_t nomzd, Double_t devzd, Double_t devaz)
849{
850 const Double_t pzd = nomzd * TMath::DegToRad();
851 const Double_t azd = devzd * TMath::DegToRad();
852 const Double_t aaz = devaz * TMath::DegToRad();
853
854 const double el = TMath::Pi()/2-pzd;
855
856 const double dphi2 = aaz/2.;
857 const double cos2 = TMath::Cos(dphi2)*TMath::Cos(dphi2);
858 const double sin2 = TMath::Sin(dphi2)*TMath::Sin(dphi2);
859 const double d = TMath::Cos(azd)*cos2 - TMath::Cos(2*el)*sin2;
860
861 return TMath::ACos(d)*TMath::RadToDeg();
862}
863
864// --------------------------------------------------------------------------
865//
866// Returned is the offset (number of days) which must be added to
867// March 1st of the given year, eg:
868//
869// Int_t offset = GetDayOfEaster(2004);
870//
871// MTime t;
872// t.Set(year, 3, 1);
873// t.SetMjd(t.GetMjd()+offset);
874//
875// cout << t << endl;
876//
877// If the date coudn't be calculated -1 is returned.
878//
879// The minimum value returned is 21 corresponding to March 22.
880// The maximum value returned is 55 corresponding to April 25.
881//
882// --------------------------------------------------------------------------
883//
884// Gauss'sche Formel zur Berechnung des Osterdatums
885// Wann wird Ostern gefeiert? Wie erfährt man das Osterdatum für ein
886// bestimmtes Jahr, ohne in einen Kalender zu schauen?
887//
888// Ostern ist ein "bewegliches" Fest. Es wird am ersten Sonntag nach dem
889// ersten Frühlingsvollmond gefeiert. Damit ist der 22. März der früheste
890// Termin, der 25. April der letzte, auf den Ostern fallen kann. Von
891// diesem Termin hängen auch die Feste Christi Himmelfahrt, das 40 Tage
892// nach Ostern, und Pfingsten, das 50 Tage nach Ostern gefeiert wird, ab.
893//
894// Von Carl Friedrich Gauß (Mathematiker, Astronom und Physiker;
895// 1777-1855) stammt ein Algorithmus, der es erlaubt ohne Kenntnis des
896// Mondkalenders die Daten der Osterfeste für die Jahre 1700 bis 2199 zu
897// bestimmen.
898//
899// Gib eine Jahreszahl zwischen 1700 und 2199 ein:
900//
901// Und so funktioniert der Algorithmus:
902//
903// Es sei:
904//
905// J die Jahreszahl
906// a der Divisionsrest von J/19
907// b der Divisionsrest von J/4
908// c der Divisionsrest von J/7
909// d der Divisionsrest von (19*a + M)/30
910// e der Divisionsrest von (2*b + 4*c + 6*d + N)/7
911//
912// wobei M und N folgende Werte annehmen:
913//
914// für die Jahre M N
915// 1583-1599 22 2
916// 1600-1699 22 2
917// 1700-1799 23 3
918// 1800-1899 23 4
919// 1900-1999 24 5
920// 2000-2099 24 5
921// 2100-2199 24 6
922// 2200-2299 25 0
923// 2300-2399 26 1
924// 2400-2499 25 1
925//
926// Dann fällt Ostern auf den
927// (22 + d + e)ten März
928//
929// oder den
930// (d + e - 9)ten April
931//
932// Beachte:
933// Anstelle des 26. Aprils ist immer der 19. April zu setzen,
934// anstelle des 25. Aprils immer dann der 18. April, wenn d=28 und a>10.
935//
936// Literatur:
937// Schüler-Rechenduden
938// Bibliographisches Institut
939// Mannheim, 1966
940//
941// --------------------------------------------------------------------------
942//
943// Der Ostersonntag ist ein sog. unregelmäßiger Feiertag. Alle anderen
944// unregelmäßigen Feiertage eines Jahres leiten sich von diesem Tag ab:
945//
946// * Aschermittwoch ist 46 Tage vor Ostern.
947// * Pfingsten ist 49 Tage nach Ostern.
948// * Christi Himmelfahrt ist 10 Tage vor Pfingsten.
949// * Fronleichnam ist 11 Tage nach Pfingsten.
950//
951// Man muß also nur den Ostersonntag ermitteln, um alle anderen
952// unregelmäßigen Feiertage zu berechnen. Doch wie geht das?
953//
954// Dazu etwas Geschichte:
955//
956// Das 1. Kirchenkonzil im Jahre 325 hat festgelegt:
957//
958// * Ostern ist stets am ersten Sonntag nach dem ersten Vollmond des
959// Frühlings.
960// * Stichtag ist der 21. März, die "Frühlings-Tagundnachtgleiche".
961//
962// Am 15.10.1582 wurde von Papst Gregor XIII. der bis dahin gültige
963// Julianische Kalender reformiert. Der noch heute gültige Gregorianische
964// Kalender legt dabei folgendes fest:
965//
966// Ein Jahr hat 365 Tage und ein Schaltjahr wird eingefügt, wenn das Jahr
967// durch 4 oder durch 400, aber nicht durch 100 teilbar ist. Hieraus
968// ergeben sich die zwei notwendigen Konstanten, um den Ostersonntag zu
969// berechnen:
970//
971// 1. Die Jahreslänge von und bis zum Zeitpunkt der
972// Frühlings-Tagundnachtgleiche: 365,2422 mittlere Sonnentage
973// 2. Ein Mondmonat: 29,5306 mittlere Sonnentage
974//
975// Mit der "Osterformel", von Carl Friedrich Gauß (1777-1855) im Jahre 1800
976// entwickelt, läßt sich der Ostersonntag für jedes Jahr von 1583 bis 8202
977// berechnen.
978//
979// Der früheste mögliche Ostertermin ist der 22. März. (Wenn der Vollmond
980// auf den 21. März fällt und der 22. März ein Sonntag ist.)
981//
982// Der späteste mögliche Ostertermin ist der 25. April. (Wenn der Vollmond
983// auf den 21. März fällt und der 21. März ein Sonntag ist.)
984//
985Int_t MAstro::GetEasterOffset(UShort_t year)
986{
987 if (year<1583 || year>2499)
988 {
989 cout << "MAstro::GetDayOfEaster - Year " << year << " not between 1700 and 2199" << endl;
990 return -1;
991 }
992
993 Int_t M=0;
994 Int_t N=0;
995 switch (year/100)
996 {
997 case 15:
998 case 16: M=22; N=2; break;
999 case 17: M=23; N=3; break;
1000 case 18: M=23; N=4; break;
1001 case 19:
1002 case 20: M=24; N=5; break;
1003 case 21: M=24; N=6; break;
1004 case 22: M=25; N=0; break;
1005 case 23: M=26; N=1; break;
1006 case 24: M=25; N=1; break;
1007 }
1008
1009 const Int_t a = year%19;
1010 const Int_t b = year%4;
1011 const Int_t c = year%7;
1012 const Int_t d = (19*a + M)%30;
1013 const Int_t e = (2*b + 4*c + 6*d + N)%7;
1014
1015 if (e==6 && d==28 && a>10)
1016 return 48;
1017
1018 if (d+e==35)
1019 return 49;
1020
1021 return d + e + 21;
1022}
Note: See TracBrowser for help on using the repository browser.