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

Last change on this file since 3532 was 3498, checked in by jbclement, 7 weeks ago

PEM:

  • Correction of the variable name for the ice table depth in "pemetat0.F90". So it is now got as intended from the "startpem.nc" file;
  • Renaming of the tendencies in the PEM with the prefix 'd_' instead of 'tend_';
  • Modification of the PEM time step type from integer to real. As a consequence, all time variables are now of real type. This change adds the possibility to consider fractions of year as time step.

JBC

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