source: trunk/LMDZ.COMMON/libf/evolution/orbit_param_criterion_mod.F90 @ 3046

Last change on this file since 3046 was 3040, checked in by jbclement, 17 months ago

Mars PCM:
The variable 'timeperi' (defined in "planete_h.F90" and computed in "iniorbit.F") is renamed into 'lsperi' and thus slightly changed to be coherent to the solar longitude of perihelion in radian. It can now be used out of the box by other subroutines/programs like the PEM.
JBC

File size: 8.1 KB
Line 
1MODULE orbit_param_criterion_mod
2
3!=======================================================================
4!
5! Purpose: Compute the maximum number of iteration for PEM based on the
6! stopping criterion given by the orbital parameters
7!
8! Author: RV, JBC
9!=======================================================================
10
11IMPLICIT NONE
12
13CONTAINS
14
15SUBROUTINE orbit_param_criterion(i_myear,year_iter_max)
16
17#ifdef CPP_IOIPSL
18    use IOIPSL,          only: getin
19#else
20    ! if not using IOIPSL, we still need to use (a local version of) getin
21    use ioipsl_getincom, only: getin
22#endif
23#ifndef CPP_STD
24    use planete_h,   only: e_elips, obliquit, lsperi
25#else
26    use planete_mod, only: e_elips, obliquit, lsperi
27#endif
28  use time_evol_mod,  only: year_bp_ini, var_obl, var_ecc, var_lsp, convert_years
29  use comconst_mod,   only: pi
30  use lask_param_mod, only: yearlask, obllask, ecclask, lsplask, ini_lask_param_mod, last_ilask
31
32  IMPLICIT NONE
33
34!--------------------------------------------------------
35! Input Variables
36!--------------------------------------------------------
37  integer, intent(in) :: i_myear ! Martian year of the simulation
38
39!--------------------------------------------------------
40! Output Variables
41!--------------------------------------------------------
42  integer, intent(out) :: year_iter_max ! Maximum number of iteration of the PEM before orb changes too much
43
44!--------------------------------------------------------
45! Local variables
46!--------------------------------------------------------
47  integer :: Year                                           ! Year of the simulation
48  real    :: Lsp                                            ! Ls perihelion
49  real    :: max_change_obl, max_change_ecc, max_change_lsp ! Maximum admissible change
50  real    :: max_obl, max_ecc, max_lsp                      ! Maximum value of orbit param given the admissible change
51  real    :: min_obl, min_ecc, min_lsp                      ! Minimum value of orbit param given the admissible change
52  integer :: max_obl_iter, max_ecc_iter, max_lsp_iter       ! Maximum year iteration before reaching an inadmissible value of orbit param
53  integer :: ilask, iylask                                  ! Loop variables
54  real    :: xa, xb, ya, yb, yc_obl, yc_ecc, yc_lsp         ! Variables for interpolation
55
56! **********************************************************************
57! 0. Initializations
58! **********************************************************************
59Year = year_bp_ini + i_myear ! We compute the current year
60Lsp = lsperi*180./pi         ! We convert in degree
61
62call ini_lask_param_mod ! Allocation of variables
63
64write(*,*) 'orbit_param_criterion, Current year =', Year
65write(*,*) 'orbit_param_criterion, Current obl. =', obliquit
66write(*,*) 'orbit_param_criterion, Current ecc. =', e_elips
67write(*,*) 'orbit_param_criterion, Current Lsp  =', Lsp
68
69! We read the file
70open(73,file = 'obl_ecc_lsp.asc')
71do ilask = 1,size(yearlask,1)
72    read(73,*) yearlask(ilask), obllask(ilask), ecclask(ilask), lsplask(ilask)
73    yearlask(ilask) = yearlask(ilask)*1000./convert_years ! Conversion from Laskar's kilo Earth years to PEM Martian years
74    if (yearlask(ilask) > Year) last_ilask = ilask + 1
75enddo
76close(73)
77write(*,*) 'Corresponding line in the obl_ecc_lsp.asc file =', last_ilask
78
79! Constant max change case
80max_change_obl = 0.6
81max_change_ecc = 2.8e-3
82max_change_lsp = 18.
83
84call getin('max_change_obl', max_change_obl)
85max_obl = obliquit + max_change_obl
86min_obl = obliquit - max_change_obl
87write(*,*) 'Maximum obliquity accepted =', max_obl
88write(*,*) 'Minimum obliquity accepted =', min_obl
89
90call getin('max_change_ecc', max_change_ecc)
91max_ecc = min(e_elips + max_change_ecc,1. - 1.e-10) ! ex < 1.
92min_ecc = max(e_elips - max_change_ecc,0.) ! ex >= 0.
93write(*,*) 'Maximum eccentricity accepted =', max_ecc
94write(*,*) 'Minimum eccentricity accepted =', min_ecc
95
96call getin('max_change_lsp', max_change_lsp)
97max_lsp = Lsp + max_change_lsp
98min_lsp = Lsp - max_change_lsp
99write(*,*) 'Maximum Lsp accepted =', max_lsp
100write(*,*) 'Minimum Lsp accepted =', min_lsp
101! End Constant max change case
102
103! If we do not want some orb parameter to change, they should not be a stopping criterion,
104! So the number of iteration corresponding is set to maximum
105max_obl_iter = 1000000
106max_ecc_iter = 1000000
107max_lsp_iter = 1000000
108
109        ! Tendency of the orbital parameter for the considered year gives the limitation between min and max
110do ilask = last_ilask,2,-1
111    xa = yearlask(ilask)
112    xb = yearlask(ilask - 1)
113    if (xa <= Year .and. Year < xb) then
114        ! Obliquity
115        if (var_obl) then
116            ya = obllask(ilask)
117            yb = obllask(ilask - 1)
118            if (ya < yb) then ! Increasing -> max is the limitation
119                yc_obl = max_obl
120            else ! Decreasing -> min is the limitation
121                yc_obl = min_obl
122            endif
123        endif
124
125        ! Eccentricity
126        if (var_ecc) then
127            ya = ecclask(ilask)
128            yb = ecclask(ilask - 1)
129            if (ya < yb) then ! Increasing -> max is the limitation
130                yc_ecc = max_ecc
131            else ! Decreasing -> min is the limitation
132                yc_ecc = min_ecc
133            endif
134        endif
135
136        ! Lsp
137        if (var_lsp) then
138            ya = lsplask(ilask)
139            yb = lsplask(ilask - 1)
140             if (ya < yb) then ! Increasing -> max is the limitation
141                if (yb - ya > 300.) then ! If modulo is "crossed" then it is the other way around
142                    yc_lsp = min_lsp
143                else
144                    yc_lsp = max_lsp
145                endif
146            else ! Decreasing -> min is the limitation
147                if (ya - yb > 300.) then ! If modulo is "crossed" then it is the other way around
148                    yc_lsp = max_lsp
149                else
150                    yc_lsp = min_lsp
151                endif
152            endif
153        endif
154        iylask = ilask
155        exit ! The loop is left as soon as the right interval is found
156    endif
157enddo
158
159if (ilask == 1) then
160    write(*,*) 'The year does not match with Laskar data in the file obl_ecc_lsp.asc.'
161    stop
162endif
163
164! Linear interpolation gives the maximum reachable year according to the limitation
165! Obliquity
166if (var_obl) then
167    do ilask = iylask,2,-1
168        ya = obllask(ilask)
169        yb = obllask(ilask - 1)
170        if (min(ya,yb) <= yc_obl .and. yc_obl < max(ya,yb)) then
171            xa = yearlask(ilask)
172            xb = yearlask(ilask - 1)
173            max_obl_iter = floor((yc_obl - ya)*(xb - xa)/(yb - ya) + xa) - Year
174            exit ! The loop is left as soon as the right interval is found
175        endif
176    enddo
177endif
178! Eccentricity
179if (var_ecc) then
180    do ilask = iylask,2,-1
181        ya = ecclask(ilask)
182        yb = ecclask(ilask - 1)
183        if (min(ya,yb) <= yc_ecc .and. yc_ecc < max(ya,yb)) then
184            xa = yearlask(ilask)
185            xb = yearlask(ilask - 1)
186            max_ecc_iter = floor((yc_ecc - ya)*(xb - xa)/(yb - ya) + xa) - Year
187            exit ! The loop is left as soon as the right interval is found
188        endif
189    enddo
190endif
191! Lsp
192if (var_lsp) then
193    do ilask = iylask,2,-1
194        ya = lsplask(ilask)
195        yb = lsplask(ilask - 1)
196        if (min(ya,yb) <= yc_lsp .and. yc_lsp < max(ya,yb)) then
197            xa = yearlask(ilask)
198            xb = yearlask(ilask - 1)
199            max_lsp_iter = floor((yc_lsp - ya)*(xb - xa)/(yb - ya) + xa) - Year
200            exit ! The loop is left as soon as the right interval is found
201        endif
202    enddo
203endif
204
205write(*,*) 'Max. number of iterations for the obl. parameter =', max_obl_iter
206write(*,*) 'Max. number of iterations for the ecc. parameter =', max_ecc_iter
207write(*,*) 'Max. number of iterations for the Lsp parameter  =', max_lsp_iter
208
209year_iter_max = min(max_obl_iter,max_ecc_iter,max_lsp_iter)
210if (year_iter_max > 0) then
211    write(*,*) 'So the max. number of iterations for the orbital criteria =', year_iter_max
212else
213    write(*,*) 'The max. number of iterations is not compatible with Laskar data in the file obl_ecc_lsp.asc.'
214    stop
215endif
216
217END SUBROUTINE orbit_param_criterion
218
219!********************************************************************************   
220
221END MODULE orbit_param_criterion_mod
222
Note: See TracBrowser for help on using the repository browser.