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