source: trunk/LMDZ.COMMON/libf/evolution/orbit.F90 @ 4066

Last change on this file since 4066 was 4065, checked in by jbclement, 5 weeks ago

PEM:
Major refactor following the previous ones (r3989 and r3991) completing the large structural reorganization and cleanup of the PEM codebase. This revision introduces newly designed modules, standardizes interfaces with explicit ini/end APIs and adds native NetCDF I/O together with explicit PCM/PEM adapters. In detail:

  • Some PEM models were corrected or improved:
    • Frost/perennial ice semantics are clarified via renaming;
    • Soil temperature remapping clarified, notably by removing the rescaling of temperature deviation;
    • Geothermal flux for the PCM is computed based on the PEM state;
  • New explicit PEM/PCM adapters ("set_*"/"build4PCM_*") to decouple PEM internal representation from PCM file layouts and reconstruct consistent fields returned to the PCM;
  • New explicit build/teardown routines that centralize allocation and initialization ordering, reducing accidental use of uninitialized data and making the model lifecycle explicit;
  • Add native read/write helpers for NetCDF that centralize all low-level NetCDF interactions with major improvements (and more simplicity) compared to legacy PEM/PCM I/O (see the modules "io_netcdf" and "output"). They support reading, creation and writing of "diagevol.nc" (renamed from "diagpem.nc") and start/restart files;
  • Provide well-focused modules ("numerics"/"maths"/"utility"/"display") to host commonly-used primitives:
    • "numerics" defines numerical types and constants for reproducibility, portability across compilers and future transitions (e.g. quadruple precision experiments);
    • "display" provides a single controlled interface for runtime messages, status output and diagnostics, avoiding direct 'print'/'write' to enable silent mode, log redirection, and MPI-safe output in the future.
    • "utility" (new module) hosts generic helpers used throughout the code (e.g. "int2str" or "real2str");
  • Add modules "clim_state_init"/"clim_state_rec" which provide robust read/write logic for "start/startfi/startpem", including 1D fallbacks, mesh conversions and dimension checks. NetCDF file creation is centralized and explicit. Restart files are now self-consistent and future-proof, requiring changes only to affected variables;
  • Add module "atmosphere" which computes pressure fields, reconstructs potential temperature and air mass. It also holds the whole logic to define sigma or hybrid coordinates for altitudes;
  • Add module "geometry" to centrilize dimensions logic and grid conversions routines (including 2 new ones "dyngrd2vect"/"vect2dyngrd");
  • Add module "slopes" to isolate slopes handling;
  • Add module "surface" to isolate surface management. Notably, albedo and emissivity are now fully reconstructed following the PCM settings;
  • Add module "allocation" to check the dimension initialization and centrilize allocation/deallocation;
  • Finalize module decomposition and renaming to consolidate domain-based modules, purpose-based routines and physics/process-based variables;
  • The main program now drives a clearer sequence of conceptual steps (initialization / reading / evolution / update / build / writing) and fails explicitly instead of silently defaulting;
  • Ice table logic is made restart-safe;
  • 'Open'/'read' intrinsic logic is made safe and automatic;
  • Improve discoverability and standardize the data handling (private vs protected vs public);
  • Apply consistent documentation/header style already introduced;
  • Update deftank scripts to reflect new names and launch wrappers;

This revision is a structural milestone aiming to be behavior-preserving where possible. It has been tested via compilation and short integration runs. However, due to extensive renames, moves, and API changes, full validation is still ongoing.
Note: the revision includes one (possibly two) easter egg hidden in the code for future archaeologists of the PEM. No physics were harmed.
JBC

File size: 22.3 KB
RevLine 
[3989]1MODULE orbit
[3991]2!-----------------------------------------------------------------------
3! NAME
4!     orbit
5! DESCRIPTION
6!     Compute orbital-parameter-based stopping criteria and update
7!     planetary orbit parameters.
8! AUTHORS & DATE
9!     R. Vandemeulebrouck
10!     JB Clement, 2023-2025
11! NOTES
[2980]12!
[3991]13!-----------------------------------------------------------------------
[2980]14
[4065]15! DEPENDENCIES
16! ------------
17use numerics, only: dp, di, k4, largest_nb, minieps
18
[3991]19! DECLARATION
20! -----------
[3076]21implicit none
[2835]22
[4065]23! PARAMETERS
24! ----------
25! File-related:
26character(15),                       parameter :: orbitfile_name = 'obl_ecc_lsp.asc'
27integer(di),                         protected :: iyear_lask! Index of the line in the file year_obl_lask.asc corresponding to the closest lower year to the current year
28real(dp), dimension(:), allocatable, protected :: year_lask ! Current_year before present from Laskar et al.
29real(dp), dimension(:), allocatable, protected :: obl_lask  ! Obliquity in Laskar's data [deg]
30real(dp), dimension(:), allocatable, protected :: ecc_lask  ! Eccentricity in Laskar's data
31real(dp), dimension(:), allocatable, protected :: lsp_lask  ! Ls perihelie in Laskar's data [deg]
32! Planet-related:
33real(dp),                            protected :: obliquity    ! Planet obliquity [deg]
34real(dp),                            protected :: eccentricity ! Orbit eccentricity
35real(dp),                            protected :: ls_peri      ! Solar longitude of the perihelion [deg]
36real(dp),                            protected :: perihelion   ! Perihelion [Mkm]
37real(dp),                            protected :: semimaj_axis ! Semi-major axis [Mkm]
38real(dp),                            protected :: aphelion     ! Aphelion [Mkm]
39real(dp),                            protected :: date_peri    ! Date of perihelion [sol]
40real(dp),                            protected :: sol_len_s    ! Length of a sol [s]
41real(dp),                            protected :: yr_len_sol   ! Length of a year [sol]
42! Evolution-related:
43logical(k4), protected :: evol_orbit     ! Flag to follow the orbital parameters
44logical(k4), protected :: evol_obl       ! Flag the PEM to follow obliquity
45logical(k4), protected :: evol_ecc       ! Flag the PEM to follow eccenticity
46logical(k4), protected :: evol_lsp       ! Flag the PEM to follow Ls perihelie
47real(dp),    protected :: max_change_obl ! Maximum admissible change for obliquity
48real(dp),    protected :: max_change_ecc ! Maximum admissible change for eccentricity
49real(dp),    protected :: max_change_lsp ! Maximum admissible change for Lsp
[3989]50
[3076]51contains
[3991]52!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
53
[3076]54!=======================================================================
[4065]55SUBROUTINE set_orbit_config(evol_orbit_in,evol_obl_in,evol_ecc_in,evol_lsp_in,max_change_obl_in,max_change_ecc_in,max_change_lsp_in)
[3991]56!-----------------------------------------------------------------------
57! NAME
[4065]58!     set_orbit_config
[3991]59!
60! DESCRIPTION
[4065]61!     Setter for 'orbit' configuration parameters.
[3991]62!
63! AUTHORS & DATE
[4065]64!     JB Clement, 02/2026
[3991]65!
66! NOTES
[4065]67!
[3991]68!-----------------------------------------------------------------------
[3039]69
[3991]70! DEPENDENCIES
71! ------------
[4065]72use utility, only: real2str, bool2str
73use display, only: print_msg
[3989]74
[3991]75! DECLARATION
76! -----------
[3989]77implicit none
78
[3991]79! ARGUMENTS
80! ---------
[4065]81logical(k4), intent(in) :: evol_orbit_in, evol_obl_in, evol_ecc_in, evol_lsp_in
82real(dp),    intent(in) :: max_change_obl_in, max_change_ecc_in, max_change_lsp_in
[3989]83
[4065]84! CODE
85! ----
86evol_orbit = evol_orbit_in
87evol_obl = evol_obl_in
88evol_ecc = evol_ecc_in
89evol_lsp = evol_lsp_in
90max_change_obl = max_change_obl_in
91max_change_ecc = max_change_ecc_in
92max_change_lsp = max_change_lsp_in
93call print_msg('evol_orbit     = '//bool2str(evol_orbit))
94call print_msg('evol_obl       = '//bool2str(evol_obl))
95call print_msg('evol_ecc       = '//bool2str(evol_ecc))
96call print_msg('evol_lsp       = '//bool2str(evol_lsp))
97call print_msg('max_change_obl = '//real2str(max_change_obl))
98call print_msg('max_change_ecc = '//real2str(max_change_ecc))
99call print_msg('max_change_lsp = '//real2str(max_change_lsp))
100
101END SUBROUTINE set_orbit_config
102!=======================================================================
103
104!=======================================================================
105SUBROUTINE init_orbit(obliquity_in,perihelion_in,aphelion_in,date_peri_in,yr_len_sol_in,sol_len_s_in)
106!-----------------------------------------------------------------------
107! NAME
108!     init_orbit
109!
110! DESCRIPTION
111!     Initialize orbital characteristics of the planet.
112!
113! AUTHORS & DATE
114!     JB Clement, 12/2025
115!
116! NOTES
117!     Taken from 'iniorbit' in the Mars PCM.
118!-----------------------------------------------------------------------
119
120! DEPENDENCIES
121! ------------
122use maths,   only: pi
123use display, only: print_msg
124use utility, only: real2str
125
126! DECLARATION
127! -----------
128implicit none
129
130! ARGUMENTS
131! ---------
132real(dp), intent(in) :: obliquity_in, perihelion_in, aphelion_in, date_peri_in, yr_len_sol_in, sol_len_s_in
133
[3991]134! LOCAL VARIABLES
135! ---------------
[4065]136real(dp)    :: zxref, zanom, zz, zx0, zdx
137real(dp)    :: lsp ! Solar longitude of the perihelion [rad]
138integer(di) :: iter
[3989]139
[3991]140! CODE
141! ----
[4065]142! Setters
143obliquity = obliquity_in
144perihelion = perihelion_in
145aphelion = aphelion_in
146date_peri = date_peri_in
147sol_len_s = sol_len_s_in
148yr_len_sol = yr_len_sol_in
149
150! Compute eccentricity and semi-major axis
151eccentricity = (aphelion - perihelion)/(aphelion + perihelion)
152semimaj_axis = 0.5_dp*(aphelion + perihelion)
153
154! Compute mean anomaly zanom
155zz = (yr_len_sol - date_peri)/yr_len_sol
156zanom = 2._dp*pi*(zz - nint(zz))
157zxref = abs(zanom)
158
159! Solve equation zx0 - e*sin(zx0) = zxref for eccentric anomaly zx0 using Newton method
160zx0 = zxref + eccentricity*sin(zxref)
161do iter = 1,100
162    zdx = -(zx0 - eccentricity*sin(zx0) - zxref)/(1._dp - eccentricity*cos(zx0))
163    if (abs(zdx) <= 1.e-12_dp) exit
164    zx0 = zx0 + zdx
165end do
166
167zx0 = zx0 + zdx
168if (zanom < 0.) zx0 = -zx0
169
170! Compute solar longitude of the perihelion
171lsp = -2._dp*atan(sqrt((1._dp + eccentricity)/(1._dp - eccentricity))*tan(zx0/2._dp))
172lsp = modulo(lsp,2._dp*pi)
173ls_peri = lsp*180._dp/pi ! Conversion in degree
174
175call print_msg('Obliquity          [deg] = '//real2str(obliquity))
176call print_msg('Eccentricity       [-]   = '//real2str(eccentricity))
177call print_msg('Ls of perihelion   [deg] = '//real2str(ls_peri))
178call print_msg('Perihelion         [Mkm] = '//real2str(perihelion))
179call print_msg('Aphelion           [Mkm] = '//real2str(aphelion))
180call print_msg('Semi-major axis    [Mkm] = '//real2str(semimaj_axis))
181call print_msg('Date of perihelion [sol] = '//real2str(date_peri))
182call print_msg('Length of a sol    [s]   = '//real2str(sol_len_s))
183call print_msg('Length of a year   [sol] = '//real2str(yr_len_sol))
184
185END SUBROUTINE init_orbit
186!=======================================================================
187
188!=======================================================================
189SUBROUTINE lsp2datep(lsp,datep)
190!-----------------------------------------------------------------------
191! NAME
192!     lsp2datep
193!
194! DESCRIPTION
195!     Initialize orbital data of the planet.
196!
197! AUTHORS & DATE
198!     JB Clement, 12/2025
199!
200! NOTES
201!     Taken from 'lsp2solp' in the Mars PCM.
202!-----------------------------------------------------------------------
203
204! DEPENDENCIES
205! ------------
206use maths, only: pi
207
208! DECLARATION
209! -----------
210implicit none
211
212! ARGUMENTS
213! ---------
214real(dp), intent(in)  :: lsp   ! Ls at perihelion [deg]
215real(dp), intent(out) :: datep ! Date at perihelion [sol]
216
217! LOCAL VARIABLES
218! ---------------
219real(dp) :: zx0 ! Eccentric anomaly at Ls = 0
220
221! CODE
222! ----
223! Compute orbit geometrical parameters
224zx0 = -2._dp*atan(tan(0.5_dp*lsp*pi/180._dp)*sqrt((1._dp - eccentricity)/(1._dp + eccentricity)))
225if (zx0 <= 0._dp) zx0 = zx0 + 2._dp*pi
226
227! Compute sol at perihelion
228datep = yr_len_sol*(1._dp - (zx0 - eccentricity*sin(zx0))/(2._dp*pi))
229
230END SUBROUTINE lsp2datep
231!=======================================================================
232
233!=======================================================================
234SUBROUTINE ini_orbit()
235!-----------------------------------------------------------------------
236! NAME
237!     ini_orbit
238!
239! DESCRIPTION
240!     Initialize the parameters of module 'orbit'.
241!
242! AUTHORS & DATE
243!     R. Vandemeulebrouck
244!     JB Clement, 2023-2025
245!
246! NOTES
247!
248!-----------------------------------------------------------------------
249
250! DEPENDENCIES
251! ------------
252use stoppage, only: stop_clean
253
254! DECLARATION
255! -----------
256implicit none
257
258! LOCAL VARIABLES
259! ---------------
260integer(di) :: n ! Number of lines in the file
261integer(di) :: ierr, funit
262logical(k4) :: here
263
264! CODE
265! ----
266inquire(file = orbitfile_name,exist = here)
267if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//orbitfile_name//'"!',1)
268
269! Get the number of lines for Laskar's orbital file
270open(newunit = funit,file = orbitfile_name,status = 'old',action = 'read',iostat = ierr)
271if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "'//orbitfile_name//'"!',ierr)
[3989]272n = 0
273do
[4065]274    read(funit,*,iostat = ierr)
[3989]275    if (ierr /= 0) exit
276    n = n + 1
[4065]277end do
278close(funit)
[3989]279
[4065]280! Allocation of module arrays
281if (.not. allocated(year_lask)) allocate(year_lask(n))
282if (.not. allocated(obl_lask)) allocate(obl_lask(n))
283if (.not. allocated(ecc_lask)) allocate(ecc_lask(n))
284if (.not. allocated(lsp_lask)) allocate(lsp_lask(n))
[3989]285
[4065]286END SUBROUTINE ini_orbit
[3991]287!=======================================================================
[3989]288
289!=======================================================================
290SUBROUTINE end_orbit
[3991]291!-----------------------------------------------------------------------
292! NAME
293!     end_orbit
294!
295! DESCRIPTION
296!     Deallocate orbital data arrays.
297!
298! AUTHORS & DATE
299!     R. Vandemeulebrouck
300!     JB Clement, 2023-2025
301!
302! NOTES
303!
304!-----------------------------------------------------------------------
[3989]305
[3991]306! DECLARATION
307! -----------
[3989]308implicit none
309
[3991]310! CODE
311! ----
[3989]312if (allocated(year_lask)) deallocate(year_lask)
[4065]313if (allocated(obl_lask)) deallocate(obl_lask)
314if (allocated(ecc_lask)) deallocate(ecc_lask)
315if (allocated(lsp_lask)) deallocate(lsp_lask)
[3989]316
317END SUBROUTINE end_orbit
318!=======================================================================
319
[3991]320!=======================================================================
[4065]321SUBROUTINE read_orbitpm(n_yr_sim)
[3991]322!-----------------------------------------------------------------------
323! NAME
[4065]324!     read_orbitpm
[3991]325!
326! DESCRIPTION
[4065]327!     Read Laskar's orbital file and locate index for
328!     closest lower year relative to current simulation year.
[3991]329!
330! AUTHORS & DATE
331!     JB Clement, 2023-2025
332!
333! NOTES
334!
335!-----------------------------------------------------------------------
[2835]336
[4065]337! DEPENDENCIES
338! ------------
339use evolution, only: pem_ini_date, r_plnt2earth_yr
340use stoppage,  only: stop_clean
341use display,   only: print_msg
342use utility,   only: real2str, int2str
343
[3991]344! DECLARATION
345! -----------
[3076]346implicit none
[2835]347
[3991]348! ARGUMENTS
349! ---------
[4065]350real(dp), intent(in) :: n_yr_sim
[2835]351
[3991]352! LOCAL VARIABLES
353! ---------------
[4065]354integer(di) :: ierr, funit, i, n
355logical(k4) :: here
356real(dp)    :: current_year
[2835]357
[3991]358! CODE
359! ----
[4065]360call print_msg('> Reading "'//orbitfile_name//'"')
361! Compute the current year
362current_year = pem_ini_date + n_yr_sim
363call print_msg('Current year = '//real2str(current_year))
[2835]364
[4065]365inquire(file = orbitfile_name,exist = here)
366if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//orbitfile_name//'"!',1)
[2835]367
[4065]368! Read the Laskars's orbital file
369open(newunit = funit,file = orbitfile_name,status = 'old',action = 'read',iostat = ierr)
370if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "'//orbitfile_name//'"!',ierr)
371n = size(year_lask,1) ! Number of lines in the file
372iyear_lask = 0 ! Line number for closest lower year relative to current simulation year
373do i = 1,n
374    read(funit,*) year_lask(i), obl_lask(i), ecc_lask(i), lsp_lask(i)
375    year_lask(i) = year_lask(i)*1000._dp/r_plnt2earth_yr ! Conversion from Laskar's kilo Earth years to planetary years
376    if (year_lask(i) > current_year) iyear_lask = i + 1
377end do
378close(funit)
379if (iyear_lask == 0 .or. iyear_lask == n + 1) then
380    call stop_clean(__FILE__,__LINE__,'the current year could not be found in the file "'//orbitfile_name//'".',1)
381else
382    call print_msg('The current year in "'//orbitfile_name//'" is at line '//int2str(iyear_lask)//'.')
383end if
[2980]384
[4065]385END SUBROUTINE read_orbitpm
386!=======================================================================
387
388!=======================================================================
389SUBROUTINE compute_maxyr_orbit(n_yr_sim,nmax_yr_runorb)
390!-----------------------------------------------------------------------
391! NAME
392!     compute_maxyr_orbit
393!
394! DESCRIPTION
395!     Determine maximum number of years before orbital params change
396!     beyond admissible thresholds.
397!
398! AUTHORS & DATE
399!     R. Vandemeulebrouck
400!     JB Clement, 2023-2025
401!
402! NOTES
403!     Handles Lsp modulo crossings.
404!-----------------------------------------------------------------------
405
406! DEPENDENCIES
407! ------------
408use maths,     only: pi
409use evolution, only: pem_ini_date
410use display,   only: print_msg
411use utility,   only: real2str
412
413! DECLARATION
414! -----------
415implicit none
416
417! ARGUMENTS
418! ---------
419real(dp), intent(in)  :: n_yr_sim       ! Current year of the simulation
420real(dp), intent(out) :: nmax_yr_runorb ! Maximum number of years of the PEM before orbital parameters change too much
421
422! LOCAL VARIABLES
423! ---------------
424real(dp)                            :: current_year
425real(dp)                            :: max_obl, max_ecc, max_lsp
426real(dp)                            :: min_obl, min_ecc, min_lsp
427real(dp)                            :: nmax_yr_runobl, nmax_yr_runecc, nmax_yr_runlsp ! Maximum year iteration before reaching an inadmissible value of orbit param
428real(dp)                            :: xa, xb, ya, yb
429logical(k4)                         :: found_obl, found_ecc, found_lsp
430real(dp)                            :: lsp_corr ! Lsp corrected if the "modulo is crossed"
431real(dp)                            :: delta
432real(dp), dimension(:), allocatable :: lsplask_unwrap
433integer(di)                         :: i
434
435! CODE
436! ----
437if (.not. evol_orbit) then
438    nmax_yr_runorb = largest_nb ! Default huge value
439    return
440end if
441call print_msg('> Computing the accepted maximum number of years due to orbital parameters')
442current_year = pem_ini_date + n_yr_sim ! We compute the current year
443
[3791]444! Unwrap the Lsp changes in Laskar's file
[3989]445allocate(lsplask_unwrap(iyear_lask))
446lsplask_unwrap(iyear_lask) = lsp_lask(iyear_lask)
447do i = iyear_lask,2,-1
448    delta = lsp_lask(i - 1) - lsp_lask(i)
[4065]449    if (delta > 180._dp) then
450        lsp_corr = lsp_lask(i - 1) - 360._dp
451    else if (delta < -180._dp) then
452        lsp_corr = lsp_lask(i - 1) + 360._dp
[3791]453    else
[3989]454        lsp_corr = lsp_lask(i - 1)
[4065]455    end if
[3989]456    lsplask_unwrap(i - 1) = lsplask_unwrap(i) + lsp_corr - lsp_lask(i)
[4065]457end do
[3791]458
459! Correct the current Lsp according to the unwrapping
[4065]460delta = ls_peri - lsplask_unwrap(iyear_lask)
461if (delta > 180._dp) then
462    ls_peri = ls_peri - 360._dp
463else if (delta < -180._dp) then
464    ls_peri = ls_peri + 360._dp
465end if
[3791]466
[4065]467! Maximum change accepted for orbital parameters
468max_obl = obliquity + max_change_obl
469min_obl = obliquity - max_change_obl
470call print_msg('Obl. (current|accepted min|accepted max): '//real2str(obliquity)//' | '//real2str(min_obl)//real2str(max_obl))
[2835]471
[4065]472max_ecc = min(eccentricity + max_change_ecc,1. - minieps) ! ecc < 1.
473min_ecc = max(eccentricity - max_change_ecc,0._dp) ! ecc >= 0.
474call print_msg('Ecc. (current|accepted min|accepted max): '//real2str(eccentricity)//' | '//real2str(min_ecc)//' | '//real2str(max_ecc))
[2835]475
[4065]476max_lsp = ls_peri + max_change_lsp
477min_lsp = ls_peri - max_change_lsp
478call print_msg('Lsp  (current|accepted min|accepted max): '//real2str(ls_peri)//' | '//real2str(min_lsp)//' | '//real2str(max_lsp))
[2835]479
[4065]480! Initialization
481nmax_yr_runobl = largest_nb ! Default huge value
482nmax_yr_runecc = largest_nb ! Default huge value
483nmax_yr_runlsp =largest_nb ! Default huge value
[3047]484found_obl = .false.
485found_ecc = .false.
486found_lsp = .false.
[4065]487if (.not. evol_obl) found_obl = .true.
488if (.not. evol_ecc) found_ecc = .true.
489if (.not. evol_lsp) found_lsp = .true.
[3047]490
[4065]491! The maximum reachable year is found by interpolation according to the maximum admissible change of orbital parameters
[3989]492do i = iyear_lask,2,-1
493    xa = year_lask(i)
494    xb = year_lask(i - 1)
[3047]495
496    ! Obliquity
[4065]497    if (evol_obl .and. .not. found_obl) then
[3989]498        ya = obl_lask(i)
499        yb = obl_lask(i - 1)
[3047]500        if (yb < min_obl) then ! The minimum accepted is overtaken
[4065]501            nmax_yr_runobl = (min_obl - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]502            found_obl = .true.
503        else if (max_obl < yb) then ! The maximum accepted is overtaken
[4065]504            nmax_yr_runobl = (max_obl - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]505            found_obl = .true.
[4065]506        end if
507    end if
[2835]508
[3047]509    ! Eccentricity
[4065]510    if (evol_ecc .and. .not. found_ecc) then
[3989]511        ya = ecc_lask(i)
512        yb = ecc_lask(i - 1)
[3047]513        if (yb < min_ecc) then ! The minimum accepted is overtaken
[4065]514            nmax_yr_runecc = (min_ecc - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]515            found_ecc = .true.
516        else if (max_ecc < yb) then ! The maximum accepted is overtaken
[4065]517            nmax_yr_runecc = (max_ecc - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]518            found_ecc = .true.
[4065]519        end if
520    end if
[2980]521
[3047]522    ! Lsp
[4065]523    if (evol_lsp .and. .not. found_lsp) then
[3989]524        ya = lsplask_unwrap(i)
525        yb = lsplask_unwrap(i - 1)
[3047]526        if (yb < min_lsp) then ! The minimum accepted is overtaken
[4065]527            nmax_yr_runlsp = (min_lsp - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]528            found_lsp = .true.
529        else if (max_lsp < yb) then ! The maximum accepted is overtaken
[4065]530            nmax_yr_runlsp = (max_lsp - ya)*(xb - xa)/(yb - ya) + xa - current_year
[3047]531            found_lsp = .true.
[4065]532        end if
533    end if
[3047]534
535    if (found_obl .and. found_ecc .and. found_lsp) exit ! The loop is left as soon as everything is found
[4065]536end do
[3039]537
[3791]538deallocate(lsplask_unwrap)
539
[4065]540if (.not. found_obl) call print_msg('The maximum reachable year for obl. could not be found in "'//orbitfile_name//'".')
541if (.not. found_ecc) call print_msg('The maximum reachable year for ecc. could not be found in "'//orbitfile_name//'".')
542if (.not. found_lsp) call print_msg('The maximum reachable year for Lsp could not be found in "'//orbitfile_name//'".')
[3039]543
[4065]544call print_msg('Number of years accepted for obl. = '//real2str(nmax_yr_runobl))
545call print_msg('Number of years accepted for ecc. = '//real2str(nmax_yr_runecc))
546call print_msg('Number of years accepted for Lsp  = '//real2str(nmax_yr_runlsp))
[2835]547
[4065]548nmax_yr_runorb = min(nmax_yr_runobl,nmax_yr_runecc,nmax_yr_runlsp)
549if (nmax_yr_runorb < 1._dp) nmax_yr_runorb = 1._dp ! nmax_yr_runorb should be at least equal to 1
[2835]550
[4065]551call print_msg('So the max. number of years for the orbital criterion = '//real2str(nmax_yr_runorb))
552
553END SUBROUTINE compute_maxyr_orbit
[3991]554!=======================================================================
[2835]555
[3991]556!=======================================================================
[4065]557SUBROUTINE update_orbit(n_yr_sim,n_yr_run)
[3991]558!-----------------------------------------------------------------------
559! NAME
[4065]560!     update_orbit
[3991]561!
562! DESCRIPTION
[4065]563!     Update orbital parameters from Laskar values at a new year.
[3991]564!
565! AUTHORS & DATE
566!     R. Vandemeulebrouck
567!     JB Clement, 2023-2025
568!
569! NOTES
570!
571!-----------------------------------------------------------------------
[2835]572
[4065]573use evolution,        only: pem_ini_date
574use maths,            only: pi
[3989]575use call_dayperi_mod, only: call_dayperi
[4065]576use stoppage,         only: stop_clean
577use display,          only: print_msg
578use utility,          only: real2str
[3039]579
[3991]580! DECLARATION
581! -----------
[3989]582implicit none
583
[3991]584! ARGUMENTS
585! ---------
[4065]586real(dp), intent(in) :: n_yr_sim ! Number of simulated years
587real(dp), intent(in) :: n_yr_run ! Number of iterations of the PEM
[3989]588
[3991]589! LOCAL VARIABLES
590! ---------------
[4065]591integer(di) :: i
592real(dp)    :: current_year, old_year, obl_old, ecc_old, lsp_old
593real(dp)    :: xa, xb, ya, yb ! Variables for interpolation
594logical(k4) :: found_year
[3989]595
[3991]596! CODE
597! ----
[4065]598call print_msg('> Updating the orbital parameters for the new year')
599current_year = pem_ini_date + n_yr_sim ! We compute the current year
600old_year = current_year - n_yr_run
601obl_old = obliquity
602ecc_old = eccentricity
603lsp_old = ls_peri
[3989]604
605found_year = .false.
[4065]606if (current_year < 0._dp) then
[3989]607    do i = iyear_lask,2,-1
608        xa = year_lask(i)
609        xb = year_lask(i - 1)
[4065]610        if (xa <= current_year .and. current_year < xb) then
[3989]611            ! Obliquity
[4065]612            if (evol_obl) then
[3989]613                ya = obl_lask(i)
614                yb = obl_lask(i - 1)
[4065]615                obliquity = (current_year - xa)*(yb - ya)/(xb - xa) + ya
616            end if
[3989]617
618            ! Eccentricity
[4065]619            if (evol_ecc) then
[3989]620                ya = ecc_lask(i)
621                yb = ecc_lask(i - 1)
[4065]622                eccentricity = (current_year - xa)*(yb - ya)/(xb - xa) + ya
623            end if
[3989]624
625            ! Lsp
[4065]626            if (evol_lsp) then
[3989]627                ya = lsp_lask(i)
628                yb = lsp_lask(i - 1)
[4065]629                if (abs(yb - ya) > 300._dp) then ! If modulo is "crossed" through the interval
[3989]630                    if (ya < yb) then ! Lsp should be decreasing
[4065]631                        ya = ya + 360._dp
[3989]632                    else ! Lsp should be increasing
[4065]633                        yb = yb + 360._dp
634                    end if
635                end if
636                ls_peri = modulo((current_year - xa)*(yb - ya)/(xb - xa) + ya,360._dp)
637            end if
[3989]638            found_year = .true.
639            exit ! The loop is left as soon as the right interval is found
[4065]640        end if
641    end do
642else if (current_year >= 0._dp .and. current_year < 100._dp) then ! For present orbital characteristics
643    if (evol_obl) obliquity = obl_lask(1)
644    if (evol_ecc) eccentricity = ecc_lask(1)
645    if (evol_lsp) ls_peri = lsp_lask(1)
[3989]646    found_year = .true.
[4065]647end if
[3989]648
[4065]649if (.not. found_year) call stop_clean(__FILE__,__LINE__,'the new year could not be found in the file "'//orbitfile_name//'".',1)
[3989]650
[4065]651call lsp2datep(ls_peri,date_peri)
652perihelion = semimaj_axis*(1._dp - eccentricity)
653aphelion = semimaj_axis*(1._dp + eccentricity)
[3989]654
[4065]655call print_msg('Year (ini|now) = '//real2str(old_year)//' | '//real2str(current_year))
656call print_msg('Obl. (ini|now) = '//real2str(obl_old)//' | '//real2str(obliquity))
657call print_msg('Ecc. (ini|now) = '//real2str(ecc_old)//' | '//real2str(eccentricity))
658call print_msg('Lsp  (ini|now) = '//real2str(lsp_old)//' | '//real2str(ls_peri))
[3989]659
[4065]660END SUBROUTINE update_orbit
[3991]661!=======================================================================
[3989]662
663END MODULE orbit
Note: See TracBrowser for help on using the repository browser.