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

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