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

Last change on this file since 3032 was 3023, checked in by jbclement, 2 years ago

Mars Deftank:
Complete rewriting (much simpler) of script "run_pem1d_1" in deftank/pem/ to launch chained simulations of 1D PCM/PEM. To be adapted to 3D in the future.
JBC

File size: 4.2 KB
Line 
1      MODULE 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      IMPLICIT NONE
9
10      CONTAINS
11
12      SUBROUTINE recomp_orb_param(final_iter)
13
14      USE temps_mod_evol, ONLY: year_bp_ini, year_PEM, var_obl, var_ex, var_lsp
15#ifndef CPP_STD     
16      USE comconst_mod, ONLY: pi
17      USE planete_h, ONLY: e_elips, obliquit, timeperi, periheli,aphelie,p_elips,peri_day,year_day
18#else
19      use planete_mod, only: e_elips, obliquit, timeperi,periastr,apoastr,p_elips,peri_day,year_day
20      USE comcstfi_mod, only: pi
21
22#endif
23
24      USE lask_param_mod, only: yearlask,oblask,exlask,lsplask, &
25                                end_lask_param_mod,last_ilask
26
27      IMPLICIT NONE
28
29!--------------------------------------------------------
30! Input Variables
31!--------------------------------------------------------
32
33      integer,intent(in) :: final_iter ! Number of year iteration of the PEM
34
35!--------------------------------------------------------
36! Output Variables
37!--------------------------------------------------------
38
39!--------------------------------------------------------
40! Local variables
41!--------------------------------------------------------
42
43      real :: Year     ! Year of the simulation
44      real :: Lsp      ! time of peri in ls
45      integer ilask    ! Loop variable
46      real :: halfaxe  ! Million of km
47      real :: unitastr
48
49      real xa,xb,ya,yb
50
51
52      ! **********************************************************************
53      ! 0. Initializations
54      ! **********************************************************************
55#ifdef CPP_STD
56      real aphelie
57      real periheli
58
59      aphelie=apoastr
60      periheli=periastr
61#endif
62
63          Year = year_bp_ini + year_PEM + final_iter ! We compute the current year
64          Lsp = 360. - timeperi*360./(2.*pi)         ! We convert in degree
65
66!          Year=-953200
67          print*, "recomp_orb_param, Old year=", year_bp_ini+year_PEM
68          print*, "recomp_orb_param, New year=", Year
69          print*, "recomp_orb_param, Old obl=", obliquit
70          print*, "recomp_orb_param, Old ex=", e_elips
71          print*, "recomp_orb_param, Old lsp=", Lsp
72          print*, "recomp_orb_param, Old timeperi=", timeperi
73
74        do ilask = last_ilask,2,-1
75          xa = yearlask(ilask)
76          xb = yearlask(ilask - 1)
77          if (xa <= Year .and. Year < xb) then
78            ! Obliquity
79            if (var_obl) then
80              ya = oblask(ilask)
81              yb = oblask(ilask - 1)
82              obliquit = (Year - xa)*(yb - ya)/(xb - xa) + ya
83            endif
84
85            ! Excentricity
86            if (var_ex) then
87              ya = exlask(ilask)
88              yb = exlask(ilask - 1)
89              e_elips = (Year - xa)*(yb - ya)/(xb - xa) + ya
90            endif
91
92            ! Lsp
93            if (var_lsp) then
94              ya = lsplask(ilask)
95              yb = lsplask(ilask - 1)
96              if (yb - ya > 300.) then ! If modulo is "crossed"
97                if (ya < yb) then ! Increasing
98                  yb = yb - 360.
99                else ! Decreasing
100                  yb = yb + 360.
101                endif
102              endif
103              Lsp = (Year - xa)*(yb - ya)/(xb - xa) + ya
104              Lsp = modulo(Lsp,360.)
105            endif
106            exit ! The loop is left as soon as the right interval is found
107          endif
108        enddo
109        halfaxe = 227.94
110        timeperi = Lsp*2.*pi/360.
111        periheli = halfaxe*(1 - e_elips)
112        aphelie =  halfaxe*(1 + e_elips)
113        unitastr = 149.597927
114        p_elips = 0.5*(periheli + aphelie)*(1. - e_elips*e_elips)/unitastr
115#ifndef CPP_STD 
116        call call_dayperi(timeperi,e_elips,peri_day,year_day)
117#endif
118
119       print*, "recomp_orb_param, Final year of the PEM run:", Year
120       print*, "recomp_orb_param, New obl=", obliquit
121       print*, "recomp_orb_param, New ex=", e_elips
122       print*, "recomp_orb_param, New timeperi=", timeperi
123
124      END SUBROUTINE recomp_orb_param
125
126!********************************************************************************   
127     
128      END MODULE recomp_orb_param_mod
Note: See TracBrowser for help on using the repository browser.