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

Last change on this file since 2835 was 2835, checked in by romain.vande, 2 years ago

Mars PEM:
Introduction of the possibility to follow an orbital forcing.
Introduction of new control parameters.
Cleaning of the PEM (removing unused files, add comments and new files)

A file named run_PEM.def can be added to the run.def. It contains the following variables:

_ evol_orbit_pem: Boolean. Do you want to follow an orbital forcing predefined (read in ob_ex_lsp.asc for example)? (default=false)
_ year_bp_ini: Integer. Number of year before present to start the pem run if evol_orbit_pem=.true. , default=0
_ Max_iter_pem: Integer. Maximal number of iteration if none of the stopping criterion is reached and if evol_orbit_pem=.false., default=99999999
_ dt_pem: Integer. Time step of the PEM in year, default=1
_ alpha_criterion: Real. Acceptance rate of sublimating ice surface change, default=0.2
_ soil_pem: Boolean. Do you want to run with subsurface physical processes in the PEM? default=.true.

RV

File size: 2.6 KB
Line 
1      MODULE recomp_orb_param_mod
2
3      IMPLICIT NONE
4
5      CONTAINS
6
7      SUBROUTINE recomp_orb_param(final_iter)
8
9      USE temps_mod_evol, ONLY: year_bp_ini, year_PEM
10      USE planete_h, ONLY: e_elips, obliquit, timeperi
11      USE comcstfi_h, only: pi
12      USE lask_param_mod, only: yearlask,oblask,exlask,lsplask, &
13                                end_lask_param_mod,last_ilask
14
15      IMPLICIT NONE
16
17!--------------------------------------------------------
18! Input Variables
19!--------------------------------------------------------
20
21      integer,intent(in) :: final_iter ! Number of year iteration of the PEM
22
23!--------------------------------------------------------
24! Output Variables
25!--------------------------------------------------------
26
27!--------------------------------------------------------
28! Local variables
29!--------------------------------------------------------
30
31      real :: Year                      ! Year of the simulation
32      real :: timeperi_ls               ! time of peri in ls
33      integer ilask                     ! Loop variable
34
35
36      ! **********************************************************************
37      ! 0. Initializations
38      ! **********************************************************************
39
40
41          Year=year_bp_ini+year_PEM+final_iter
42
43          timeperi_ls=timeperi*360/2/pi
44
45!          Year=-953200
46          print *, "recomp_orb_param, Old year=", year_bp_ini+year_PEM
47          print *, "recomp_orb_param, New year=", Year
48          print *, "recomp_orb_param, Old obl=", obliquit
49          print *, "recomp_orb_param, Old ex=", e_elips
50          print *, "recomp_orb_param, Old lsp=", timeperi_ls
51          print *, "recomp_orb_param, Old timeperi=", timeperi
52
53        do ilask=last_ilask,1,-1
54           if(yearlask(ilask) .GT.Year) then
55              obliquit=oblask(ilask+1)+(oblask(ilask)-oblask(ilask+1))*(Year-yearlask(ilask+1))/(yearlask(ilask)-yearlask(ilask+1))
56              e_elips=exlask(ilask+1)+(exlask(ilask)-exlask(ilask+1))*(Year-yearlask(ilask+1))/(yearlask(ilask)-yearlask(ilask+1))
57              timeperi_ls=lsplask(ilask+1)+(lsplask(ilask)-lsplask(ilask+1))*(Year-yearlask(ilask+1))/(yearlask(ilask)-yearlask(ilask+1))
58              exit
59           endif
60        enddo
61        timeperi=timeperi_ls*2*pi/360
62
63       print *, "recomp_orb_param, Final year of the PEM run:", Year
64       print *, "recomp_orb_param, New obl=", obliquit
65       print *, "recomp_orb_param, New ex=", e_elips
66       print *, "recomp_orb_param, New timeperi=", timeperi
67
68      END SUBROUTINE recomp_orb_param
69
70!********************************************************************************   
71     
72      END MODULE recomp_orb_param_mod
Note: See TracBrowser for help on using the repository browser.