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

Last change on this file since 3070 was 3069, checked in by jbclement, 2 years ago

Mars PCM:
Related to commit r3066, correction of a bug to write/read a restart/start in 1D and more adaptations of the code.
JBC

File size: 8.6 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                                 ! Variables for interpolation
55  logical :: found_obl, found_ecc, found_lsp                ! Flag variables for orbital parameters
56  logical :: crossed_modulo_d, crossed_modulo_i             ! Flag variables for modulo of Lsp
57
58! **********************************************************************
59! 0. Initializations
60! **********************************************************************
61Year = year_bp_ini + i_myear ! We compute the current year
62Lsp = lsperi*180./pi         ! We convert in degree
63
64call ini_lask_param_mod ! Allocation of variables
65
66write(*,*) 'orbit_param_criterion, Current year =', Year
67write(*,*) 'orbit_param_criterion, Current obl. =', obliquit
68write(*,*) 'orbit_param_criterion, Current ecc. =', e_elips
69write(*,*) 'orbit_param_criterion, Current Lsp  =', Lsp
70
71! We read the file
72open(73,file = 'obl_ecc_lsp.asc')
73last_ilask = 0
74do ilask = 1,size(yearlask,1)
75    read(73,*) yearlask(ilask), obllask(ilask), ecclask(ilask), lsplask(ilask)
76    yearlask(ilask) = yearlask(ilask)*1000./convert_years ! Conversion from Laskar's kilo Earth years to PEM Martian years
77    if (yearlask(ilask) > Year) last_ilask = ilask + 1
78enddo
79close(73)
80if (last_ilask == 0 .or. last_ilask == size(yearlask,1) + 1) then
81    error stop 'The current year could not be found in the file "obl_ecc_lsp.asc".'
82else
83    write(*,*) 'The current year in the "obl_ecc_lsp.asc" file is at line:', last_ilask
84endif
85
86! Constant max change case
87max_change_obl = 0.6
88max_change_ecc = 2.8e-3
89max_change_lsp = 18.
90
91call getin('max_change_obl',max_change_obl)
92max_obl = obliquit + max_change_obl
93min_obl = obliquit - max_change_obl
94write(*,*) 'Maximum obliquity accepted =', max_obl
95write(*,*) 'Minimum obliquity accepted =', min_obl
96
97call getin('max_change_ecc',max_change_ecc)
98max_ecc = min(e_elips + max_change_ecc,1. - 1.e-10) ! ex < 1.
99min_ecc = max(e_elips - max_change_ecc,0.) ! ex >= 0.
100write(*,*) 'Maximum eccentricity accepted =', max_ecc
101write(*,*) 'Minimum eccentricity accepted =', min_ecc
102
103call getin('max_change_lsp',max_change_lsp)
104max_lsp = Lsp + max_change_lsp
105min_lsp = Lsp - max_change_lsp
106write(*,*) 'Maximum Lsp accepted =', max_lsp
107write(*,*) 'Minimum Lsp accepted =', min_lsp
108! End Constant max change case
109
110! If we do not want some orb parameter to change, they should not be a stopping criterion,
111! So the number of iteration corresponding is set to maximum
112max_obl_iter = 1000000
113max_ecc_iter = 1000000
114max_lsp_iter = 1000000
115
116! The maximum reachable year is found by interpolation according to the maximum admissible change of orbital parameters
117found_obl = .false.
118found_ecc = .false.
119found_lsp = .false.
120if (.not. var_obl) found_obl = .true.
121if (.not. var_ecc) found_ecc = .true.
122if (.not. var_lsp) found_lsp = .true.
123crossed_modulo_d = .false.
124crossed_modulo_i = .false.
125
126do ilask = last_ilask,2,-1
127    xa = yearlask(ilask)
128    xb = yearlask(ilask - 1)
129
130    ! Obliquity
131    if (var_obl .and. .not. found_obl) then
132        ya = obllask(ilask)
133        yb = obllask(ilask - 1)
134        if (yb < min_obl) then ! The minimum accepted is overtaken
135            max_obl_iter = floor((min_obl - ya)*(xb - xa)/(yb - ya) + xa) - Year
136            found_obl = .true.
137        else if (max_obl < yb) then ! The maximum accepted is overtaken
138            max_obl_iter = floor((max_obl - ya)*(xb - xa)/(yb - ya) + xa) - Year
139            found_obl = .true.
140        endif
141    endif
142
143    ! Eccentricity
144    if (var_ecc .and. .not. found_ecc) then
145        ya = ecclask(ilask)
146        yb = ecclask(ilask - 1)
147        if (yb < min_ecc) then ! The minimum accepted is overtaken
148            max_ecc_iter = floor((min_ecc - ya)*(xb - xa)/(yb - ya) + xa) - Year
149            found_ecc = .true.
150        else if (max_ecc < yb) then ! The maximum accepted is overtaken
151            max_ecc_iter = floor((max_ecc - ya)*(xb - xa)/(yb - ya) + xa) - Year
152            found_ecc = .true.
153        endif
154    endif
155
156    ! Lsp
157    if (var_lsp .and. .not. found_lsp) then
158        ya = lsplask(ilask)
159        yb = lsplask(ilask - 1)
160        if (abs(yb - ya) > 300.) then ! If modulo is "crossed" through the interval
161            if (ya < yb) then ! Lsp should be decreasing
162                if (.not. crossed_modulo_i) then
163                    yb = yb - 360.
164                    if (min_lsp < 0.) crossed_modulo_d = .true.
165                else
166                    crossed_modulo_i = .false.
167                endif
168            else ! Lsp should be increasing
169                if (.not. crossed_modulo_d) then
170                    yb = yb + 360.
171                    if (max_lsp > 360.) crossed_modulo_i = .true.
172                else
173                    crossed_modulo_d = .false.
174                endif
175            endif
176        else
177            if (crossed_modulo_d) then ! If decreasing Lsp "crossed" the modulo but min_lsp has not been met yet
178                ya = ya - 360.
179                yb = yb - 360.
180            else if (crossed_modulo_i) then ! If increasing Lsp "crossed" the modulo but max_lsp has not been met yet
181                ya = ya + 360.
182                yb = yb + 360.
183            endif
184        endif
185        if (yb < min_lsp) then ! The minimum accepted is overtaken
186            max_lsp_iter = floor((min_lsp - ya)*(xb - xa)/(yb - ya) + xa) - Year
187            found_lsp = .true.
188        else if (max_lsp < yb) then ! The maximum accepted is overtaken
189            max_lsp_iter = floor((max_lsp - ya)*(xb - xa)/(yb - ya) + xa) - Year
190            found_lsp = .true.
191        endif
192    endif
193
194    if (found_obl .and. found_ecc .and. found_lsp) exit ! The loop is left as soon as everything is found
195enddo
196
197if (.not. found_obl) write(*,*) 'The maximum reachable year for obl. could not be found in the file "obl_ecc_lsp.asc".'
198if (.not. found_ecc) write(*,*) 'The maximum reachable year for ecc. could not be found in the file "obl_ecc_lsp.asc".'
199if (.not. found_lsp) write(*,*) 'The maximum reachable year for Lsp could not be found in the file "obl_ecc_lsp.asc".'
200
201write(*,*) 'Max. number of iterations for the obl. parameter =', max_obl_iter
202write(*,*) 'Max. number of iterations for the ecc. parameter =', max_ecc_iter
203write(*,*) 'Max. number of iterations for the Lsp parameter  =', max_lsp_iter
204
205year_iter_max = min(max_obl_iter,max_ecc_iter,max_lsp_iter)
206write(*,*) 'So the max. number of iterations for the orbital criteria =', year_iter_max
207
208END SUBROUTINE orbit_param_criterion
209
210!********************************************************************************   
211
212END MODULE orbit_param_criterion_mod
213
Note: See TracBrowser for help on using the repository browser.