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

Last change on this file was 3742, checked in by jbclement, 2 months ago

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

File size: 4.2 KB
RevLine 
[3039]1MODULE recomp_orb_param_mod
[2855]2!=======================================================================
3!
4! Purpose: Recompute orbit parameters based on values by Laskar et al.,
5!
[3022]6! Author: RV, JBC
[2855]7!=======================================================================
[2835]8
[3076]9implicit none
[2835]10
[3076]11!=======================================================================
12contains
13!=======================================================================
14
[3498]15SUBROUTINE recomp_orb_param(i_myear,i_myear_leg)
[2835]16
[3742]17use time_evol_mod,  only: year_bp_ini, var_obl, var_ecc, var_lsp
[3076]18use lask_param_mod, only: yearlask, obllask, ecclask, lsplask, end_lask_param_mod, last_ilask
19#ifndef CPP_STD
[3742]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
[2842]23#else
[3742]24    use comcstfi_mod,     only: pi
25    use planete_mod,      only: e_elips, obliquit, lsperi, periastr, apoastr, p_elips, peri_day, year_day
[2842]26#endif
[2849]27
[3076]28implicit none
[2835]29
30!--------------------------------------------------------
31! Input Variables
32!--------------------------------------------------------
[3498]33real, intent(in) :: i_myear   ! Number of simulated Martian years
34real, intent(in) :: i_myear_leg ! Number of iterations of the PEM
[2835]35
36!--------------------------------------------------------
37! Output Variables
38!--------------------------------------------------------
39
40!--------------------------------------------------------
[3076]41! Local variables
[2835]42!--------------------------------------------------------
[3498]43real    :: Year           ! Year of the simulation
[3076]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
[2835]50
[3039]51! **********************************************************************
52! 0. Initializations
53! **********************************************************************
[3076]54#ifdef CPP_STD
55    real :: aphelie, periheli
[2835]56
[3076]57    aphelie = apoastr
58    periheli = periastr
[2985]59#endif
[2835]60
[3040]61Year = year_bp_ini + i_myear ! We compute the current year
62Lsp = lsperi*180./pi         ! We convert in degree
[2835]63
[3498]64write(*,*) 'recomp_orb_param, Old year =', Year - i_myear_leg
[3039]65write(*,*) 'recomp_orb_param, Old obl. =', obliquit
66write(*,*) 'recomp_orb_param, Old ecc. =', e_elips
67write(*,*) 'recomp_orb_param, Old Lsp  =', Lsp
[3553]68write(*,*) 'recomp_orb_param, New year =', Year
[2835]69
[3047]70found_year = .false.
[3039]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
[2980]81
[3039]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
[2980]88
[3039]89        ! Lsp
90        if (var_lsp) then
91            ya = lsplask(ilask)
92            yb = lsplask(ilask - 1)
[3047]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
[3039]97                    yb = yb + 360.
[3022]98                endif
99            endif
[3047]100            Lsp = modulo((Year - xa)*(yb - ya)/(xb - xa) + ya,360.)
[3039]101        endif
[3047]102        found_year = .true.
[3039]103        exit ! The loop is left as soon as the right interval is found
104    endif
105enddo
[3047]106
[3069]107if (.not. found_year) error stop 'The new year could not be found in the file "obl_ecc_lsp.asc".'
[3047]108
[3039]109halfaxe = 227.94
[3040]110Lsp = Lsp*pi/180.
[3039]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
[3076]116#ifndef CPP_STD
[3040]117    call call_dayperi(Lsp,e_elips,peri_day,year_day)
[2985]118#endif
[2835]119
[3039]120write(*,*) 'recomp_orb_param, New obl. =', obliquit
121write(*,*) 'recomp_orb_param, New ecc. =', e_elips
122write(*,*) 'recomp_orb_param, New Lsp  =', Lsp
[2835]123
[3039]124END SUBROUTINE recomp_orb_param
[2835]125
[3076]126!********************************************************************************
127
[3039]128END MODULE recomp_orb_param_mod
129
Note: See TracBrowser for help on using the repository browser.