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

Last change on this file since 3790 was 3742, checked in by jbclement, 8 weeks ago

PEM:
Resulting modifications due to changes in the Mars PCM with r3739 and r3741.
JBC

File size: 4.2 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!=======================================================================
8
9implicit none
10
11!=======================================================================
12contains
13!=======================================================================
14
15SUBROUTINE recomp_orb_param(i_myear,i_myear_leg)
16
17use time_evol_mod,  only: year_bp_ini, var_obl, var_ecc, var_lsp
18use lask_param_mod, only: yearlask, obllask, ecclask, lsplask, end_lask_param_mod, last_ilask
19#ifndef CPP_STD
20    use comcstfi_h,       only: pi
21    use planete_h,        only: e_elips, obliquit, lsperi, periheli, aphelie, p_elips, peri_day, year_day
22    use call_dayperi_mod, only: call_dayperi
23#else
24    use comcstfi_mod,     only: pi
25    use planete_mod,      only: e_elips, obliquit, lsperi, periastr, apoastr, p_elips, peri_day, year_day
26#endif
27
28implicit none
29
30!--------------------------------------------------------
31! Input Variables
32!--------------------------------------------------------
33real, intent(in) :: i_myear   ! Number of simulated Martian years
34real, intent(in) :: i_myear_leg ! Number of iterations of the PEM
35
36!--------------------------------------------------------
37! Output Variables
38!--------------------------------------------------------
39
40!--------------------------------------------------------
41! Local variables
42!--------------------------------------------------------
43real    :: Year           ! Year of the simulation
44real    :: Lsp            ! Ls perihelion
45integer :: ilask          ! Loop variable
46real    :: halfaxe        ! Million of km
47real    :: unitastr
48real    :: xa, xb, ya, yb ! Variables for interpolation
49logical :: found_year     ! Flag variable
50
51! **********************************************************************
52! 0. Initializations
53! **********************************************************************
54#ifdef CPP_STD
55    real :: aphelie, periheli
56
57    aphelie = apoastr
58    periheli = periastr
59#endif
60
61Year = year_bp_ini + i_myear ! We compute the current year
62Lsp = lsperi*180./pi         ! We convert in degree
63
64write(*,*) 'recomp_orb_param, Old year =', Year - i_myear_leg
65write(*,*) 'recomp_orb_param, Old obl. =', obliquit
66write(*,*) 'recomp_orb_param, Old ecc. =', e_elips
67write(*,*) 'recomp_orb_param, Old Lsp  =', Lsp
68write(*,*) 'recomp_orb_param, New year =', Year
69
70found_year = .false.
71do ilask = last_ilask,2,-1
72    xa = yearlask(ilask)
73    xb = yearlask(ilask - 1)
74    if (xa <= Year .and. Year < xb) then
75        ! Obliquity
76        if (var_obl) then
77            ya = obllask(ilask)
78            yb = obllask(ilask - 1)
79            obliquit = (Year - xa)*(yb - ya)/(xb - xa) + ya
80        endif
81
82        ! Eccentricity
83        if (var_ecc) then
84            ya = ecclask(ilask)
85            yb = ecclask(ilask - 1)
86            e_elips = (Year - xa)*(yb - ya)/(xb - xa) + ya
87        endif
88
89        ! Lsp
90        if (var_lsp) then
91            ya = lsplask(ilask)
92            yb = lsplask(ilask - 1)
93            if (abs(yb - ya) > 300.) then ! If modulo is "crossed" through the interval
94                if (ya < yb) then ! Lsp should be decreasing
95                    ya = ya + 360.
96                else ! Lsp should be increasing
97                    yb = yb + 360.
98                endif
99            endif
100            Lsp = modulo((Year - xa)*(yb - ya)/(xb - xa) + ya,360.)
101        endif
102        found_year = .true.
103        exit ! The loop is left as soon as the right interval is found
104    endif
105enddo
106
107if (.not. found_year) error stop 'The new year could not be found in the file "obl_ecc_lsp.asc".'
108
109halfaxe = 227.94
110Lsp = Lsp*pi/180.
111periheli = halfaxe*(1 - e_elips)
112aphelie =  halfaxe*(1 + e_elips)
113unitastr = 149.597927
114p_elips = 0.5*(periheli + aphelie)*(1. - e_elips*e_elips)/unitastr
115
116#ifndef CPP_STD
117    call call_dayperi(Lsp,e_elips,peri_day,year_day)
118#endif
119
120write(*,*) 'recomp_orb_param, New obl. =', obliquit
121write(*,*) 'recomp_orb_param, New ecc. =', e_elips
122write(*,*) 'recomp_orb_param, New Lsp  =', Lsp
123
124END SUBROUTINE recomp_orb_param
125
126!********************************************************************************
127
128END MODULE recomp_orb_param_mod
129
Note: See TracBrowser for help on using the repository browser.