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