source: trunk/LMDZ.COMMON/libf/evolution/recomp_orb_param_mod.F90 @ 3046

Last change on this file since 3046 was 3040, checked in by jbclement, 17 months ago

Mars PCM:
The variable 'timeperi' (defined in "planete_h.F90" and computed in "iniorbit.F") is renamed into 'lsperi' and thus slightly changed to be coherent to the solar longitude of perihelion in radian. It can now be used out of the box by other subroutines/programs like the PEM.
JBC

File size: 3.8 KB
Line 
1MODULE recomp_orb_param_mod
2!=======================================================================
3!
4! Purpose: Recompute orbit parameters based on values by Laskar et al.,
5!
6! Author: RV, JBC
7!=======================================================================
8IMPLICIT NONE
9
10CONTAINS
11
12SUBROUTINE recomp_orb_param(i_myear,year_iter)
13
14  use time_evol_mod, only: year_bp_ini, var_obl, var_ecc, var_lsp, convert_years
15#ifndef CPP_STD     
16    use comconst_mod, only: pi
17    use planete_h,    only: e_elips, obliquit, lsperi, periheli, aphelie, p_elips, peri_day, year_day
18#else
19    use planete_mod,  only: e_elips, obliquit, lsperi, periastr, apoastr, p_elips, peri_day, year_day
20    use comcstfi_mod, only: pi
21#endif
22
23  use lask_param_mod, only: yearlask, obllask, ecclask, lsplask, end_lask_param_mod, last_ilask
24
25  IMPLICIT NONE
26
27!--------------------------------------------------------
28! Input Variables
29!--------------------------------------------------------
30  integer, intent(in) :: i_myear   ! Number of simulated Martian years
31  integer, intent(in) :: year_iter ! Number of iterations of the PEM
32
33!--------------------------------------------------------
34! Output Variables
35!--------------------------------------------------------
36
37!--------------------------------------------------------
38! Local variables
39!--------------------------------------------------------
40  integer :: Year           ! Year of the simulation
41  real    :: Lsp            ! Ls perihelion
42  integer :: ilask          ! Loop variable
43  real    :: halfaxe        ! Million of km
44  real    :: unitastr
45  real    :: xa, xb, ya, yb ! Variables for interpolation
46
47! **********************************************************************
48! 0. Initializations
49! **********************************************************************
50#ifdef CPP_STD
51    real aphelie
52    real periheli
53
54    aphelie=apoastr
55    periheli=periastr
56#endif
57
58Year = year_bp_ini + i_myear ! We compute the current year
59Lsp = lsperi*180./pi         ! We convert in degree
60
61write(*,*) 'recomp_orb_param, Old year =', Year - year_iter
62write(*,*) 'recomp_orb_param, Old obl. =', obliquit
63write(*,*) 'recomp_orb_param, Old ecc. =', e_elips
64write(*,*) 'recomp_orb_param, Old Lsp  =', Lsp
65
66do ilask = last_ilask,2,-1
67    xa = yearlask(ilask)
68    xb = yearlask(ilask - 1)
69    if (xa <= Year .and. Year < xb) then
70        ! Obliquity
71        if (var_obl) then
72            ya = obllask(ilask)
73            yb = obllask(ilask - 1)
74            obliquit = (Year - xa)*(yb - ya)/(xb - xa) + ya
75        endif
76
77        ! Eccentricity
78        if (var_ecc) then
79            ya = ecclask(ilask)
80            yb = ecclask(ilask - 1)
81            e_elips = (Year - xa)*(yb - ya)/(xb - xa) + ya
82        endif
83
84        ! Lsp
85        if (var_lsp) then
86            ya = lsplask(ilask)
87            yb = lsplask(ilask - 1)
88            if (yb - ya > 300.) then ! If modulo is "crossed"
89                if (ya < yb) then ! Increasing
90                    yb = yb - 360.
91                else ! Decreasing
92                    yb = yb + 360.
93                endif
94            endif
95            Lsp = (Year - xa)*(yb - ya)/(xb - xa) + ya
96            Lsp = modulo(Lsp,360.)
97        endif
98        exit ! The loop is left as soon as the right interval is found
99    endif
100enddo
101halfaxe = 227.94
102Lsp = Lsp*pi/180.
103periheli = halfaxe*(1 - e_elips)
104aphelie =  halfaxe*(1 + e_elips)
105unitastr = 149.597927
106p_elips = 0.5*(periheli + aphelie)*(1. - e_elips*e_elips)/unitastr
107
108#ifndef CPP_STD 
109    call call_dayperi(Lsp,e_elips,peri_day,year_day)
110#endif
111
112write(*,*) 'recomp_orb_param, New year =', Year
113write(*,*) 'recomp_orb_param, New obl. =', obliquit
114write(*,*) 'recomp_orb_param, New ecc. =', e_elips
115write(*,*) 'recomp_orb_param, New Lsp  =', Lsp
116
117END SUBROUTINE recomp_orb_param
118
119!********************************************************************************   
120     
121END MODULE recomp_orb_param_mod
122
Note: See TracBrowser for help on using the repository browser.