source: trunk/LMDZ.PLUTO/libf/phypluto/calcmilank.F90 @ 3504

Last change on this file since 3504 was 3247, checked in by afalco, 9 months ago

Added missing files from pluto.old

File size: 1.8 KB
Line 
1      subroutine calcmilank(tim,obliq,perid,eccen)
2         
3      use planete_mod, only: year_day
4      use comcstfi_mod, only: pi
5
6      implicit none
7
8!==================================================================
9!     Purpose
10!     -------
11!     calculate the milankovitch parametersdepeding on the
12!     date (tim, in Myr)
13!     Obliquity
14!     Day of perihelion from the vernal equinox
15!     Eccentricity
16
17!     Inputs
18!     ------
19!     tim                   Date (Myr)
20!     
21!     Outputs
22!     -------
23!     obliq                 new value for obliquity
24!     perid                 new value for peri_day
25!     eccen                 new value for eccentricity
26!     
27!     Authors
28!     -------
29!     Tanguy Bertrand (new equations in 2022)
30!     
31!==================================================================
32
33!------------------------------------------------------------------
34!     Arguments
35
36      REAL,intent(in) :: tim
37      REAL,intent(out) :: perid,obliq,eccen
38      real teta,xo,xref,psi,phi
39!------------------------------------------------------------------
40
41      write(*,*) 'TB22 tim  =',tim
42
43      psi=72.8+91.*tim
44      phi=19.5+130.2*tim
45
46      ! update new obliquity
47      obliq=115.5+11.8*sin(phi*pi/180.)
48
49      ! update new eccentricity
50      eccen=0.244+0.022*cos(psi*pi/180.)+0.005*cos(3*psi*pi/180.)
51
52      ! update new solar longitude of perihelion
53      teta=pi/180.*modulo(19.5+130.2*tim-24.*sin(psi*pi/180.),360.)
54      xo=2.*atan(((1-eccen)/(1+eccen))**0.5*tan(teta/2.))
55      xref=xo-eccen*sin(xo)
56      perid=year_day*(1.-xref/(2.*pi))
57
58      if (perid.ge.year_day) then
59         perid = perid - year_day
60      endif
61      write(*,*) 'new ecc, lsp, perid=',eccen,360.-teta*180./pi,perid
62      write(*,*) 'new obliquity=',obliq
63
64      end subroutine calcmilank
65
Note: See TracBrowser for help on using the repository browser.