source: trunk/LMDZ.COMMON/libf/evolution/maths.F90 @ 4076

Last change on this file since 4076 was 4065, checked in by jbclement, 2 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: 10.3 KB
Line 
1module maths
2!-----------------------------------------------------------------------
3! NAME
4!     maths
5!
6! DESCRIPTION
7!     The module contains all the mathematical subroutines used in the PEM.
8!
9! AUTHORS & DATE
10!     L. Lange
11!     JB Clement, 2023-2025
12!
13! NOTES
14!     Adapted from Schorghofer MSIM (N.S, Icarus 2010)
15!-----------------------------------------------------------------------
16
17! DEPENDENCIES
18! ------------
19use numerics, only: dp, di, k4, minieps
20
21! DECLARATION
22! -----------
23implicit none
24
25! PARAMETERS
26! ----------
27real(dp),    parameter :: pi = 4._dp*atan(1._dp) ! PI = 3.14159...
28integer(di), parameter :: limiter_minmod  = 1
29integer(di), parameter :: limiter_MC      = 2
30integer(di), parameter :: limiter_vanLeer = 3
31
32contains
33!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34
35!=======================================================================
36SUBROUTINE deriv1(z,nz,y,y0,ybot,dzY)
37!-----------------------------------------------------------------------
38! NAME
39!     deriv1
40!
41! DESCRIPTION
42!     Compute the first derivative of a function y(z) on an irregular grid.
43!
44! AUTHORS & DATE
45!     N.S (Icarus 2010)
46!     L. Lange
47!     JB Clement, 2023-2025
48!
49! NOTES
50!     Upper boundary conditions: y(0) = y0
51!     Lower boundary condition: yp = ybottom
52!-----------------------------------------------------------------------
53
54! DECLARATION
55! -----------
56implicit none
57
58! ARGUMENTS
59! ---------
60integer(di),             intent(in)  :: nz       ! number of layer
61real(dp), dimension(nz), intent(in)  :: z        ! depth layer
62real(dp), dimension(nz), intent(in)  :: y        ! function which needs to be derived
63real(dp),                intent(in)  :: y0, ybot ! boundary conditions
64real(dp), dimension(nz), intent(out) :: dzY      ! derivative of y w.r.t depth
65
66! LOCAL VARIABLES
67! ---------------
68integer(di) :: j
69real(dp)    :: hm, hp, c1, c2, c3
70
71! CODE
72! ----
73hp = z(2) - z(1)
74c1 = z(1)/(hp*z(2))
75c2 = 1/z(1) - 1/(z(2) - z(1))
76c3 = -hp/(z(1)*z(2))
77dzY(1) = c1*y(2) + c2*y(1) + c3*y0
78do j = 2,nz - 1
79    hp = z(j + 1) - z(j)
80    hm = z(j) - z(j - 1)
81    c1 = +hm/(hp*(z(j + 1) - z(j - 1)))
82    c2 = 1._dp/hm - 1._dp/hp
83    c3 = -hp/(hm*(z(j + 1) - z(j - 1)))
84    dzY(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1)
85end do
86dzY(nz) = (ybot - y(nz - 1))/(2._dp*(z(nz) - z(nz - 1)))
87
88END SUBROUTINE deriv1
89!=======================================================================
90
91!=======================================================================
92SUBROUTINE deriv2_simple(z,nz,y,y0,yNp1,yp2)
93!-----------------------------------------------------------------------
94! NAME
95!     deriv2_simple
96!
97! DESCRIPTION
98!     Compute the second derivative of a function y(z) on an irregular grid.
99!
100! AUTHORS & DATE
101!     N.S (Icarus 2010)
102!     JB Clement, 2023-2025
103!
104! NOTES
105!     Boundary conditions: y(0) = y0, y(nz + 1) = yNp1
106!-----------------------------------------------------------------------
107
108! DECLARATION
109! -----------
110implicit none
111
112! ARGUMENTS
113! ---------
114integer(di),             intent(in)  :: nz
115real(dp), dimension(nz), intent(in)  :: z, y
116real(dp),                intent(in)  :: y0, yNp1
117real(dp), dimension(nz), intent(out) :: yp2
118
119! LOCAL VARIABLES
120! ---------------
121integer(di) :: j
122real(dp)    :: hm, hp, c1, c2, c3
123
124! CODE
125! ----
126c1 = +2._dp/((z(2) - z(1))*z(2))
127c2 = -2._dp/((z(2) - z(1))*z(1))
128c3 = +2._dp/(z(1)*z(2))
129yp2(1) = c1*y(2) + c2*y(1) + c3*y0
130do j = 2,nz - 1
131    hp = z(j + 1) - z(j)
132    hm = z(j) - z(j - 1)
133    c1 = +2._dp/(hp*(z(j + 1) - z(j - 1)))
134    c2 = -2._dp/(hp*hm)
135    c3 = +2._dp/(hm*(z(j + 1) - z(j - 1)))
136    yp2(j) = c1*y(j + 1) + c2*y(j) + c3*y(j - 1)
137end do
138yp2(nz) = (yNp1 - 2._dp*y(nz) + y(nz - 1))/(z(nz) - z(nz - 1))**2
139
140END SUBROUTINE deriv2_simple
141!=======================================================================
142
143!=======================================================================
144SUBROUTINE  deriv1_onesided(j,z,nz,y,dy_zj)
145!-----------------------------------------------------------------------
146! NAME
147!     deriv1_onesided
148!
149! DESCRIPTION
150!     First derivative of function y(z) at z(j) one-sided derivative
151!     on irregular grid.
152!
153! AUTHORS & DATE
154!     N.S (Icarus 2010)
155!     JB Clement, 2023-2025
156!
157! NOTES
158!
159!-----------------------------------------------------------------------
160
161! DECLARATION
162! -----------
163implicit none
164
165! ARGUMENTS
166! ---------
167integer(di),             intent(in)  :: nz, j
168real(dp), dimension(nz), intent(in)  :: z, y
169real(dp),                intent(out) :: dy_zj
170
171! LOCAL VARIABLES
172! ---------------
173real(dp) :: h1, h2, c1, c2, c3
174
175! CODE
176! ----
177if (j < 1 .or. j > nz - 2) then
178    dy_zj = -1._dp
179else
180    h1 = z(j + 1) - z(j)
181    h2 = z(j + 2)- z(j + 1)
182    c1 = -(2._dp*h1 + h2)/(h1*(h1 + h2))
183    c2 = (h1 + h2)/(h1*h2)
184    c3 = -h1/(h2*(h1 + h2))
185    dy_zj = c1*y(j) + c2*y(j + 1) + c3*y(j + 2)
186end if
187
188END SUBROUTINE deriv1_onesided
189!=======================================================================
190
191!=======================================================================
192PURE SUBROUTINE colint(y,z,nz,i1,i2,integral)
193!-----------------------------------------------------------------------
194! NAME
195!     colint
196!
197! DESCRIPTION
198!     Column integrates y on irregular grid.
199!
200! AUTHORS & DATE
201!     N.S (Icarus 2010)
202!     JB Clement, 2023-2025
203!
204! NOTES
205!
206!-----------------------------------------------------------------------
207
208! DECLARATION
209! -----------
210implicit none
211
212! ARGUMENTS
213! ---------
214integer(di),             intent(in)  :: nz, i1, i2
215real(dp), dimension(nz), intent(in)  :: y, z
216real(dp),                intent(out) :: integral
217
218! LOCAL VARIABLES
219! ---------------
220integer(di)             :: i
221real(dp), dimension(nz) :: dz
222
223! CODE
224! ----
225dz(1) = (z(2) - 0.)/2
226do i = 2,nz - 1
227    dz(i) = (z(i + 1) - z(i - 1))/2.
228end do
229dz(nz) = z(nz) - z(nz - 1)
230integral = sum(y(i1:i2)*dz(i1:i2))
231
232END SUBROUTINE colint
233!=======================================================================
234
235!=======================================================================
236SUBROUTINE findroot(y1,y2,z1,z2,zr)
237!-----------------------------------------------------------------------
238! NAME
239!     findroot
240!
241! DESCRIPTION
242!     Compute the root zr, between two values y1 and y2 at depth z1,z2.
243!
244! AUTHORS & DATE
245!     L. Lange
246!     JB Clement, 2023-2025
247!
248! NOTES
249!
250!-----------------------------------------------------------------------
251
252! DECLARATION
253! -----------
254implicit none
255
256! ARGUMENTS
257! ---------
258real(dp), intent(in)  :: y1, y2
259real(dp), intent(in)  :: z1, z2
260real(dp), intent(out) :: zr
261
262! CODE
263! ----
264zr = (y1*z2 - y2*z1)/(y1 - y2)
265
266END SUBROUTINE findroot
267!=======================================================================
268
269!=======================================================================
270SUBROUTINE solve_tridiag(a,b,c,d,n,x,error)
271!-----------------------------------------------------------------------
272! NAME
273!     solve_tridiag
274!
275! DESCRIPTION
276!     Solve a tridiagonal system Ax = d using the Thomas' algorithm
277!     a: sub-diagonal
278!     b: main diagonal
279!     c: super-diagonal
280!     d: right-hand side
281!     x: solution
282!
283! AUTHORS & DATE
284!     JB Clement, 2025
285!
286! NOTES
287!
288!-----------------------------------------------------------------------
289
290! DECLARATION
291! -----------
292implicit none
293
294! ARGUMENTS
295! ---------
296integer(di),                intent(in)  :: n
297real(dp), dimension(n),     intent(in)  :: b, d
298real(dp), dimension(n - 1), intent(in)  :: a, c
299real(dp), dimension(n),     intent(out) :: x
300integer(di),                intent(out) :: error
301
302! LOCAL VARIABLES
303! ---------------
304integer(di)            :: i
305real(dp)               :: m
306real(dp), dimension(n) :: cp, dp
307
308! CODE
309! ----
310! Check stability: diagonally dominant condition
311error = 0
312if (abs(b(1)) < abs(c(1))) then
313    error = 1
314    return
315end if
316do i = 2,n - 1
317    if (abs(b(i)) < abs(a(i - 1)) + abs(c(i))) then
318        error = 1
319        return
320    end if
321end do
322if (abs(b(n)) < abs(a(n - 1))) then
323    error = 1
324    return
325end if
326
327! Initialization
328cp(1) = c(1)/b(1)
329dp(1) = d(1)/b(1)
330
331! Forward sweep
332do i = 2,n - 1
333    m = b(i) - a(i - 1)*cp(i - 1)
334    cp(i) = c(i)/m
335    dp(i) = (d(i) - a(i - 1)*dp(i - 1))/m
336end do
337m = b(n) - a(n - 1)*cp(n - 1)
338dp(n) = (d(n) - a(n - 1)*dp(n - 1))/m
339
340! Backward substitution
341x(n) = dp(n)
342do i = n - 1,1,-1
343    x(i) = dp(i) - cp(i)*x(i + 1)
344end do
345
346END SUBROUTINE solve_tridiag
347!=======================================================================
348
349!=======================================================================
350SUBROUTINE solve_steady_heat(n,z,mz,kappa,mkappa,T_left,q_right,T)
351!-----------------------------------------------------------------------
352! NAME
353!     solve_steady_heat
354!
355! DESCRIPTION
356!     Solve 1D steady-state heat equation with space-dependent thermal diffusivity.
357!     Uses Thomas algorithm to solve tridiagonal system.
358!     Left boundary: prescribed temperature (Dirichlet).
359!     Right boundary: prescribed thermal flux (Neumann).
360!
361! AUTHORS & DATE
362!     JB Clement, 2025
363!
364! NOTES
365!     Grid points at z, mid-grid points at mz.
366!     Thermal diffusivity specified at both grids (kappa, mkappa).
367!-----------------------------------------------------------------------
368
369! DEPENDENCIES
370! ------------
371use stoppage, only: stop_clean
372
373! DECLARATION
374! -----------
375implicit none
376
377! ARGUMENTS
378! ---------
379integer(di),                intent(in)  :: n
380real(dp), dimension(n),     intent(in)  :: z, kappa
381real(dp), dimension(n - 1), intent(in)  :: mz, mkappa
382real(dp),                   intent(in)  :: T_left, q_right
383real(dp), dimension(n),     intent(out) :: T
384
385! LOCAL VARIABLES
386! ---------------
387integer(di)                :: i, error
388real(dp), dimension(n)     :: b, d
389real(dp), dimension(n - 1) :: a, c
390
391! CODE
392! ----
393! Initialization
394a = 0._dp; b = 0._dp; c = 0._dp; d = 0._dp
395
396! Left boundary condition (Dirichlet: prescribed temperature)
397b(1) = 1._dp
398d(1) = T_left
399
400! Internal points
401do i = 2,n - 1
402    a(i - 1) = -mkappa(i - 1)/((mz(i) - mz(i - 1))*(z(i) - z(i - 1)))
403    c(i) = -mkappa(i)/((mz(i) - mz(i - 1))*(z(i + 1) - z(i)))
404    b(i) = -(a(i - 1) + c(i))
405end do
406
407! Right boundary condition (Neumann: prescribed temperature)
408a(n - 1) = kappa(n - 1)/(z(n) - z(n - 1))
409b(n) = -kappa(n)/(z(n) - z(n - 1))
410d(n) = q_right
411
412! Solve the tridiagonal linear system with the Thomas' algorithm
413call solve_tridiag(a,b,c,d,n,T,error)
414if (error /= 0) call stop_clean(__FILE__,__LINE__,"unstable solving!",1)
415
416END SUBROUTINE solve_steady_heat
417!=======================================================================
418
419end module maths
Note: See TracBrowser for help on using the repository browser.