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

Last change on this file since 3470 was 3399, checked in by jbclement, 4 months ago

PEM:
Small corection for the 1D related to r3386 and r3369 + Making the computation of maximum number of iterations due to orbital variations more robust.
JBC

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