source: trunk/MESOSCALE/LMD_MM_MARS/SRC/PYTHON/time.F @ 261

Last change on this file since 261 was 261, checked in by aslmd, 13 years ago

MESOSCALE: user manual. a decent complete version. needs to be read again. also moved python graphical scripts to /home/aymeric/Science/MODELESMESOSCALELMD_MM_MARS/SRC/PYTHON so that it is available to users (and not only developers)

File size: 3.3 KB
Line 
1cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
2      real function sol2ls(sol)
3
4c  convert a given martian day number (sol)
5c  into corresponding solar longitude, Ls (in degr.),
6c  where sol=0=Ls=0 is the
7c  northern hemisphere spring equinox.
8
9      implicit none
10
11c     input
12      real    sol
13c     output
14      real    ls
15
16c     local variables
17      double precision    year_day,peri_day,timeperi,e_elips
18      double precision    pi,radtodeg
19c number of martian days (sols) in a martian year
20      parameter (year_day=668.6d0)
21c perihelion date (in sols)
22      parameter (peri_day=485.35d0)
23c orbital eccentricity
24      parameter (e_elips=0.09340d0)
25      parameter (pi=3.14159265358979d0)
26c  radtodeg: 180/pi
27      parameter (radtodeg=57.2957795130823d0)
28c  timeperi: 2*pi*( 1 - Ls(perihelion)/ 360 ); Ls(perihelion)=250.99
29      parameter (timeperi=1.90258341759902d0)
30
31      double precision    zanom,xref,zx0,zdx,zteta,zz
32c  xref: mean anomaly, zx0: eccentric anomaly, zteta: true anomaly     
33      integer iter
34
35      zz=(sol-peri_day)/year_day
36      zanom=2.*pi*(zz-nint(zz))
37      xref=dabs(zanom)
38
39c  The equation zx0 - e * sin (zx0) = xref, solved by Newton
40      zx0=xref+e_elips*dsin(xref)
41      do 110 iter=1,10
42         zdx=-(zx0-e_elips*dsin(zx0)-xref)/(1.-e_elips*dcos(zx0))
43         if(dabs(zdx).le.(1.d-7)) then ! typically, 2 or 3 iterations are enough
44           goto 120
45         endif
46         zx0=zx0+zdx
47  110 continue
48  120 continue
49      zx0=zx0+zdx
50      if(zanom.lt.0.) zx0=-zx0
51
52c compute true anomaly zteta, now that eccentric anomaly zx0 is known
53      zteta=2.*datan(dsqrt((1.+e_elips)/(1.-e_elips))*dtan(zx0/2.))
54
55c compute Ls
56      ls=real(zteta-timeperi)
57      if(ls.lt.0.) ls=ls+2.*real(pi)
58      if(ls.gt.2.*pi) ls=ls-2.*real(pi)
59c convert Ls in deg.
60      ls=real(radtodeg)*ls
61
62      sol2ls = ls
63      return
64      end
65
66cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
67      real function ls2sol(ls)
68
69c  Returns solar longitude, Ls (in deg.), from day number (in sol),
70c  where sol=0=Ls=0 at the northern hemisphere spring equinox
71
72      implicit none
73
74c  Arguments:
75      real ls
76
77c  Local:
78      double precision xref,zx0,zteta,zz
79c       xref: mean anomaly, zteta: true anomaly, zx0: eccentric anomaly
80      double precision year_day
81      double precision peri_day,timeperi,e_elips
82      double precision pi,degrad
83      parameter (year_day=668.6d0) ! number of sols in a amartian year
84c      data peri_day /485.0/
85      parameter (peri_day=485.35d0) ! date (in sols) of perihelion
86c  timeperi: 2*pi*( 1 - Ls(perihelion)/ 360 ); Ls(perihelion)=250.99
87      parameter (timeperi=1.90258341759902d0)
88      parameter (e_elips=0.0934d0)  ! eccentricity of orbit
89      parameter (pi=3.14159265358979d0)
90      parameter (degrad=57.2957795130823d0)
91
92      if (abs(ls).lt.1.0e-5) then
93         if (ls.ge.0.0) then
94            ls2sol = 0.0
95         else
96            ls2sol = real(year_day)
97         end if
98         return
99      end if
100
101      zteta = ls/degrad + timeperi
102      zx0 = 2.0*datan(dtan(0.5*zteta)/dsqrt((1.+e_elips)/(1.-e_elips)))
103      xref = zx0-e_elips*dsin(zx0)
104      zz = xref/(2.*pi)
105      ls2sol = real(zz*year_day + peri_day)
106      if (ls2sol.lt.0.0) ls2sol = ls2sol + real(year_day)
107      if (ls2sol.ge.year_day) ls2sol = ls2sol - real(year_day)
108
109      return
110      end
Note: See TracBrowser for help on using the repository browser.