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

Last change on this file since 4066 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: 20.0 KB
Line 
1MODULE soil_temp
2!-----------------------------------------------------------------------
3! NAME
4!     soil_temp
5!
6! DESCRIPTION
7!     Routines to compute soil temperature evolution and initialization.
8!
9! AUTHORS & DATE
10!     L. Lange, 2023
11!     JB Clement, 2023-2025
12!
13! NOTES
14!
15!-----------------------------------------------------------------------
16
17! DEPENDENCIES
18! ------------
19use numerics, only: dp, di, k4
20
21! DECLARATION
22! -----------
23implicit none
24
25! PARAMETERS
26! ----------
27real(dp), dimension(:,:,:), allocatable, protected :: tsoil_PCM    ! Soil temperature in the PCM at the beginning [K]
28real(dp), dimension(:,:),   allocatable, protected :: flux_geo_PCM ! Geothermal flux in the PCM at the beginning [K]
29
30contains
31!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32
33!=======================================================================
34SUBROUTINE ini_soil_temp()
35!-----------------------------------------------------------------------
36! NAME
37!     ini_soil_temp
38!
39! DESCRIPTION
40!     Initialize the parameters of module 'soil_temp'.
41!
42! AUTHORS & DATE
43!     JB Clement, 12/2025
44!
45! NOTES
46!
47!-----------------------------------------------------------------------
48
49! DEPENDENCIES
50! ------------
51use geometry, only: ngrid, nsoil_PCM, nslope
52
53! DECLARATION
54! -----------
55implicit none
56
57! CODE
58! ----
59if (.not. allocated(tsoil_PCM)) allocate(tsoil_PCM(ngrid,nsoil_PCM,nslope))
60if (.not. allocated(flux_geo_PCM)) allocate(flux_geo_PCM(ngrid,nslope))
61
62END SUBROUTINE ini_soil_temp
63!=======================================================================
64
65!=======================================================================
66SUBROUTINE end_soil_temp()
67!-----------------------------------------------------------------------
68! NAME
69!     end_soil_temp
70!
71! DESCRIPTION
72!     Deallocate soil_temp arrays.
73!
74! AUTHORS & DATE
75!     JB Clement, 12/2025
76!
77! NOTES
78!
79!-----------------------------------------------------------------------
80
81! DECLARATION
82! -----------
83implicit none
84
85! CODE
86! ----
87if (allocated(tsoil_PCM)) deallocate(tsoil_PCM)
88if (allocated(flux_geo_PCM)) deallocate(flux_geo_PCM)
89
90END SUBROUTINE end_soil_temp
91!=======================================================================
92
93!=======================================================================
94SUBROUTINE set_tsoil_PCM(tsoil_PCM_in)
95!-----------------------------------------------------------------------
96! NAME
97!     set_tsoil_PCM
98!
99! DESCRIPTION
100!     Setter for 'tsoil_PCM'.
101!
102! AUTHORS & DATE
103!     JB Clement, 12/2025
104!
105! NOTES
106!
107!-----------------------------------------------------------------------
108
109! DECLARATION
110! -----------
111implicit none
112
113! ARGUMENTS
114! ---------
115real(dp), dimension(:,:,:), intent(in) :: tsoil_PCM_in
116
117! CODE
118! ----
119tsoil_PCM(:,:,:) = tsoil_PCM_in(:,:,:)
120
121END SUBROUTINE set_tsoil_PCM
122!=======================================================================
123
124!=======================================================================
125SUBROUTINE set_flux_geo_PCM(flux_geo_PCM_in)
126!-----------------------------------------------------------------------
127! NAME
128!     set_flux_geo_PCM
129!
130! DESCRIPTION
131!     Setter for 'flux_geo_PCM'.
132!
133! AUTHORS & DATE
134!     JB Clement, 12/2025
135!
136! NOTES
137!
138!-----------------------------------------------------------------------
139
140! DECLARATION
141! -----------
142implicit none
143
144! ARGUMENTS
145! ---------
146real(dp), dimension(:,:), intent(in) :: flux_geo_PCM_in
147
148! CODE
149! ----
150flux_geo_PCM(:,:) = flux_geo_PCM_in(:,:)
151
152END SUBROUTINE set_flux_geo_PCM
153!=======================================================================
154
155!=======================================================================
156SUBROUTINE compute_tsoil(ngrid,nsoil,firstcall,therm_i,timestep,tsurf,tsoil)
157!-----------------------------------------------------------------------
158! NAME
159!     compute_tsoil
160!
161! DESCRIPTION
162!     Compute soil temperature using an implicit 1st order scheme.
163!
164! AUTHORS & DATE
165!     L. Lange, 2023
166!     JB Clement, 2023-2025
167!
168! NOTES
169!
170!-----------------------------------------------------------------------
171
172! DEPENDENCIES
173! ------------
174use soil, only: layer, mlayer, mthermdiff, thermdiff, coefq, coefd, mu, alph, beta, flux_geo, volcapa
175
176! DECLARATION
177! -----------
178implicit none
179
180#include "dimensions.h"
181
182! ARGUMENTS
183! ---------
184integer(di),              intent(in)    :: ngrid     ! number of (horizontal) grid-points
185integer(di),              intent(in)    :: nsoil     ! number of soil layers
186logical(k4),              intent(in)    :: firstcall ! identifier for initialization call
187real(dp), dimension(:,:), intent(in)    :: therm_i   ! thermal inertia [SI]
188real(dp),                 intent(in)    :: timestep  ! time step [s]
189real(dp), dimension(:),   intent(in)    :: tsurf     ! surface temperature [K]
190real(dp), dimension(:,:), intent(inout) :: tsoil     ! soil (mid-layer) temperature [K]
191
192! LOCAL VARIABLES
193! ---------------
194integer(di) :: ig, ik
195
196! CODE
197! ----
198! 0. Initialisations and preprocessing step
199 if (firstcall) then
200! 0.1 Build mthermdiff(:), the mid-layer thermal diffusivities
201    do ig = 1,ngrid
202        do ik = 0,nsoil - 1
203            mthermdiff(ig,ik) = therm_i(ig,ik + 1)*therm_i(ig,ik + 1)/volcapa
204        end do
205    end do
206
207! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities
208    do ig = 1,ngrid
209        do ik = 1,nsoil - 1
210            thermdiff(ig,ik) = ((layer(ik) - mlayer(ik - 1))*mthermdiff(ig,ik) &
211                                   + (mlayer(ik) - layer(ik))*mthermdiff(ig,ik - 1))/(mlayer(ik) - mlayer(ik - 1))
212        end do
213    end do
214
215! 0.3 Build coefficients mu, q_{k+1/2}, d_k, alph
216    ! mu
217    mu = mlayer(0)/(mlayer(1) - mlayer(0))
218
219    ! q_{1/2}
220    coefq(0) = volcapa*layer(1)/timestep
221    ! q_{k+1/2}
222    do ik = 1,nsoil - 1
223        coefq(ik) = volcapa*(layer(ik + 1) - layer(ik))/timestep
224    end do
225
226    do ig = 1,ngrid
227        ! d_k
228        do ik = 1,nsoil - 1
229            coefd(ig,ik) = thermdiff(ig,ik)/(mlayer(ik)-mlayer(ik - 1))
230        end do
231
232        ! alph_PEM_{N-1}
233        alph(ig,nsoil - 1) = coefd(ig,nsoil - 1)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
234        ! alph_PEM_k
235        do ik = nsoil - 2,1,-1
236            alph(ig,ik) = coefd(ig,ik)/(coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
237        end do
238      end do ! of do ig=1,ngrid
239end if ! of if (firstcall)
240
241if (.not. firstcall) THEN
242! 2. Compute soil temperatures
243! First layer:
244    do ig = 1,ngrid
245        tsoil(ig,1) = (tsurf(ig) + mu*beta(ig,1)*thermdiff(ig,1)/mthermdiff(ig,0))/ &
246                      (1._dp + mu*(1._dp - alph(ig,1))*thermdiff(ig,1)/mthermdiff(ig,0))
247
248! Other layers:
249        do ik = 1,nsoil - 1
250                tsoil(ig,ik + 1) = alph(ig,ik)*tsoil(ig,ik) + beta(ig,ik)
251        end do
252    end do
253end if
254
255!  2. Compute beta coefficients (preprocessing for next time step)
256! Bottom layer, beta_PEM_{N-1}
257do ig = 1,ngrid
258    beta(ig,nsoil - 1) = coefq(nsoil - 1)*tsoil(ig,nsoil)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1)) &
259                             + flux_geo/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
260end do
261! Other layers
262do ik = nsoil - 2,1,-1
263    do ig = 1,ngrid
264        beta(ig,ik) = (coefq(ik)*tsoil(ig,ik + 1) + coefd(ig,ik + 1)*beta(ig,ik + 1))/ &
265                          (coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
266    end do
267end do
268
269END SUBROUTINE compute_tsoil
270!=======================================================================
271
272!=======================================================================
273SUBROUTINE ini_tsoil_pem(ngrid,nsoil,therm_i,tsurf,tsoil)
274!-----------------------------------------------------------------------
275! NAME
276!     ini_tsoil_pem
277!
278! DESCRIPTION
279!     Initialize soil with the solution of the stationary heat conduction problem.
280!     Boundary conditions: Tsurf averaged from PCM; geothermal flux at bottom.
281!
282! AUTHORS & DATE
283!     L. Lang, 2023
284!     JB Clement, 2023-2025
285!
286! NOTES
287!
288!-----------------------------------------------------------------------
289
290! DEPENDENCIES
291! ------------
292use soil, only: layer, mlayer, mthermdiff, thermdiff, coefq, coefd, mu, alph, beta, flux_geo, volcapa
293
294! DECLARATION
295! -----------
296implicit none
297
298#include "dimensions.h"
299
300! ARGUMENTS
301! ---------
302integer(di),              intent(in)    :: ngrid   ! number of (horizontal) grid-points
303integer(di),              intent(in)    :: nsoil   ! number of soil layers
304real(dp), dimension(:,:), intent(in)    :: therm_i ! thermal inertia [SI]
305real(dp), dimension(:),   intent(in)    :: tsurf   ! surface temperature [K]
306real(dp), dimension(:,:), intent(inout) :: tsoil   ! soil (mid-layer) temperature [K]
307
308! LOCAL VARIABLES
309! ---------------
310integer(di) :: ig, ik, iloop
311
312! CODE
313! ----
314! 0. Initialisations and preprocessing step
315! 0.1 Build mthermdiff(:), the mid-layer thermal diffusivities
316do ig = 1,ngrid
317    do ik = 0,nsoil - 1
318        mthermdiff(ig,ik) = therm_i(ig,ik + 1)*therm_i(ig,ik + 1)/volcapa
319    end do
320end do
321
322! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities
323do ig = 1,ngrid
324    do ik = 1,nsoil - 1
325        thermdiff(ig,ik) = ((layer(ik) - mlayer(ik - 1))*mthermdiff(ig,ik) &
326                             + (mlayer(ik) - layer(ik))*mthermdiff(ig,ik - 1))/(mlayer(ik) - mlayer(ik - 1))
327    end do
328end do
329
330! 0.3 Build coefficients mu, q_{k+1/2}, d_k, alph
331! mu
332mu = mlayer(0)/(mlayer(1) - mlayer(0))
333
334! q_{1/2}
335coefq(:) = 0.
336! q_{k+1/2}
337
338do ig = 1,ngrid
339    ! d_k
340    do ik = 1,nsoil - 1
341        coefd(ig,ik) = thermdiff(ig,ik)/(mlayer(ik) - mlayer(ik - 1))
342    end do
343
344    ! alph_PEM_{N-1}
345    alph(ig,nsoil - 1) = coefd(ig,nsoil - 1)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
346    ! alph_PEM_k
347    do ik = nsoil - 2,1,-1
348        alph(ig,ik) = coefd(ig,ik)/(coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
349    end do
350end do ! of do ig=1,ngrid
351
352!  1. Compute beta coefficients
353! Bottom layer, beta_PEM_{N-1}
354do ig = 1,ngrid
355        beta(ig,nsoil - 1) = coefq(nsoil - 1)*tsoil(ig,nsoil)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1)) &
356                                 + flux_geo/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
357end do
358! Other layers
359do ik = nsoil - 2,1,-1
360    do ig = 1,ngrid
361        beta(ig,ik) = (coefq(ik)*tsoil(ig,ik + 1) + coefd(ig,ik + 1)*beta(ig,ik + 1))/ &
362                          (coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
363    end do
364end do
365
366!  2. Compute soil temperatures
367do iloop = 1,10 !just convergence
368    do ig = 1,ngrid
369        ! First layer:
370        tsoil(ig,1) = (tsurf(ig) + mu*beta(ig,1)*thermdiff(ig,1)/mthermdiff(ig,0))/ &
371                      (1._dp + mu*(1._dp - alph(ig,1))*thermdiff(ig,1)/mthermdiff(ig,0))
372        ! Other layers:
373        do ik = 1,nsoil - 1
374            tsoil(ig,ik + 1) = alph(ig,ik)*tsoil(ig,ik) + beta(ig,ik)
375        end do
376    end do
377end do ! iloop
378
379END SUBROUTINE ini_tsoil_pem
380!=======================================================================
381
382!=======================================================================
383SUBROUTINE shift_tsoil2surf(ngrid,nsoil,nslope,zshift_surf,zlag,tsurf,tsoil)
384!-----------------------------------------------------------------------
385! NAME
386!     shift_tsoil2surf
387!
388! DESCRIPTION
389!     Shift soil temperature profile to follow surface evolution due to ice
390!     condensation/sublimation.
391!
392! AUTHORS & DATE
393!     JB Clement, 2025
394!
395! NOTES
396!
397!-----------------------------------------------------------------------
398
399! DEPENDENCIES
400! ------------
401use soil,    only: layer, mlayer, flux_geo, thermdiff, mthermdiff
402use maths,   only: solve_steady_heat
403use display, only: print_msg
404
405! DECLARATION
406! -----------
407implicit none
408
409! ARGUMENTS
410! ---------
411integer(di),                intent(in)    :: ngrid       ! number of (horizontal) grid-points
412integer(di),                intent(in)    :: nsoil       ! number of soil layers
413integer(di),                intent(in)    :: nslope      ! number of sub-slopes
414real(dp), dimension(:,:),   intent(in)    :: zshift_surf ! elevation shift for the surface [m]
415real(dp), dimension(:,:),   intent(in)    :: zlag        ! newly built lag thickness [m]
416real(dp), dimension(:,:),   intent(in)    :: tsurf       ! surface temperature [K]
417real(dp), dimension(:,:,:), intent(inout) :: tsoil       ! soil (mid-layer) temperature [K]
418
419! LOCAL VARIABLES
420! ---------------
421integer(di)                             :: ig, isoil, islope, ishift, ilag
422real(dp)                                :: z, zshift_surfloc
423real(dp), dimension(ngrid,nsoil,nslope) :: tsoil_old
424
425! CODE
426! ----
427call print_msg("> Shifting soil temperature profile to match surface evolution")
428tsoil_old = tsoil
429
430do ig = 1,ngrid
431    do islope = 1,nslope
432        zshift_surfloc = zshift_surf(ig,islope)
433
434        if (zshift_surfloc >= 0.) then ! In case of the surface is higher than initially
435            if (zshift_surfloc < mlayer(0)) then ! Surface change is too small to be taken into account
436                ! Nothing to do; we keep the soil temperature profile
437            else if (zshift_surfloc >= mlayer(nsoil - 1)) then ! Surface change is much larger than the discretization
438                tsoil(ig,:,islope) = tsurf(ig,islope)
439            else
440                ishift = 0
441                do while (mlayer(ishift) <= zshift_surfloc)
442                    ishift = ishift + 1 ! mlayer indices begin at 0 so this the good index for tsoil!
443                end do
444                ! The "new soil" temperature is set to tsurf
445                tsoil(ig,:ishift,islope) = tsurf(ig,islope)
446                do isoil = ishift + 1,nsoil
447                    ! Position in the old discretization of the depth
448                    z = mlayer(isoil - 1) - zshift_surfloc
449                    ! Interpolation of the temperature profile from the old discretization
450                    tsoil(ig,isoil,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
451                end do
452            end if
453        else ! In case of the surface is lower than initially
454            if (abs(zshift_surfloc) < mlayer(0)) then ! Surface change is too small to be taken into account
455                ! Nothing to do; we keep the soil temperature profile
456            else if (abs(zshift_surfloc) >= mlayer(nsoil - 1)) then ! Surface change is much larger than the discretization
457                call solve_steady_heat(nsoil,mlayer,layer,mthermdiff(ig,:),thermdiff(ig,:),tsurf(ig,islope),flux_geo,tsoil(ig,:,islope))
458            else
459                if (zlag(ig,islope) < mlayer(0)) then ! The lag is too thin to be taken into account
460                    ilag = 0
461                else
462                    ilag = 0
463                    do while (mlayer(ilag) <= zlag(ig,islope))
464                        ilag = ilag + 1 ! mlayer indices begin at 0 so this the good index for tsoil!
465                    end do
466                    ! Position of the lag bottom in the old discretization of the depth
467                    z = zlag(ig,islope) - zshift_surfloc
468                    ! The "new lag" temperature is set to the ice temperature just below
469                    tsoil(ig,:ilag,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
470                end if
471
472                ishift = nsoil - 1
473                z = mlayer(nsoil - 1) + zshift_surfloc
474                do while (mlayer(ishift) >= z)
475                    ishift = ishift - 1
476                end do
477                ishift = ishift + 1 ! Adding 1 is needed to match the good index for tsoil!
478                do isoil = ilag + 1,ishift
479                    ! Position in the old discretization of the depth
480                    z = mlayer(isoil - 1) - zshift_surfloc
481                    ! Interpolation of the temperature profile from the old discretization
482                    tsoil(ig,isoil,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
483                end do
484
485                ! The "new deepest layers" temperature is set by solving the steady heat equation
486                call solve_steady_heat(nsoil - ishift + 1,mlayer(ishift - 1:),layer(ishift:),mthermdiff(ig,ishift - 1:),thermdiff(ig,ishift:),tsoil(ig,ishift,islope),flux_geo,tsoil(ig,ishift:,islope))
487            end if
488        end if
489    end do
490end do
491
492END SUBROUTINE shift_tsoil2surf
493!=======================================================================
494
495!=======================================================================
496FUNCTION itp_tsoil(tsoil,tsurf,z) RESULT(tsoil_z)
497!-----------------------------------------------------------------------
498! NAME
499!     itp_tsoil
500!
501! DESCRIPTION
502!     Interpolate soil temperature profile.
503!
504! AUTHORS & DATE
505!     JB Clement, 2025
506!
507! NOTES
508!
509!-----------------------------------------------------------------------
510
511! DEPENDENCIES
512! ------------
513use soil, only: mlayer
514
515! DECLARATION
516! -----------
517implicit none
518
519! ARGUMENTS
520! ---------
521real(dp), dimension(:), intent(in) :: tsoil
522real(dp),               intent(in) :: z, tsurf
523
524! LOCAL VARIABLES
525! ---------------
526real(dp)    :: tsoil_z, tsoil_minus, mlayer_minus, a
527integer(di) :: iz
528
529! CODE
530! ----
531! Find the interval [mlayer(iz - 1),mlayer(iz)[ where the position z belongs
532iz = 0
533do while (mlayer(iz) <= z)
534    iz = iz + 1
535end do
536if (iz == 0) then
537    tsoil_minus = tsurf
538    mlayer_minus = 0._dp
539else
540    tsoil_minus = tsoil(iz)
541    mlayer_minus = mlayer(iz - 1)
542end if
543! Interpolation of the temperature profile from the old discretization
544a = (tsoil(iz + 1) - tsoil_minus)/(mlayer(iz) - mlayer_minus)
545tsoil_z = a*(z - mlayer_minus) + tsoil_minus
546
547END FUNCTION itp_tsoil
548!=======================================================================
549
550!=======================================================================
551SUBROUTINE evolve_soil_temp(tsoil_avg, tsurf_avg, tsoil_ts, tsoil_ts_old, h2o_soildensity_avg)
552!-----------------------------------------------------------------------
553! NAME
554!     evolve_soil_temp
555!
556! DESCRIPTION
557!     Update soil temperature profile and calculate water soil density
558!     time series and averages.
559!
560! AUTHORS & DATE
561!     JB Clement, 02/2026
562!
563! NOTES
564!
565!-----------------------------------------------------------------------
566
567! DEPENDENCIES
568! ------------
569use geometry,  only: ngrid, nsoil, nslope, nday
570use planet,    only: alpha_clap_h2o, beta_clap_h2o
571use tracers,   only: mmol, iPCM_qh2o
572use physics,   only: mugaz, r
573use evolution, only: dt
574use soil,      only: TI
575use stoppage,  only: stop_clean
576use display,   only: print_msg
577
578! DECLARATION
579! -----------
580implicit none
581
582! ARGUMENTS
583! ---------
584real(dp), dimension(:,:),     intent(in)    :: tsurf_avg
585real(dp), dimension(:,:,:),   intent(inout) :: tsoil_avg
586real(dp), dimension(:,:,:,:), intent(inout) :: tsoil_ts
587real(dp), dimension(:,:,:,:), intent(out)   :: tsoil_ts_old
588real(dp), dimension(:,:,:),   intent(out)   :: h2o_soildensity_avg
589
590! LOCAL VARIABLES
591! ---------------
592real(dp), dimension(ngrid,nsoil)             :: tsoil_avg_old
593real(dp), dimension(ngrid,nsoil,nslope,nday) :: h2o_soildensity_ts
594integer(di)                                  :: i, isoil, islope, iday
595
596! CODE
597! ----
598call print_msg("> Update the soil temperature")
599
600! Store current state
601tsoil_ts_old(:,:,:,:) = tsoil_ts(:,:,:,:)
602
603do islope = 1,nslope
604    ! Store the average profile before the update
605    tsoil_avg_old(:,:) = tsoil_avg(:,:,islope)
606   
607    ! Compute the new soil temperature
608    call compute_tsoil(ngrid,nsoil,.true., TI(:,:,islope),dt,tsurf_avg(:,islope),tsoil_avg(:,:,islope))
609    call compute_tsoil(ngrid,nsoil,.false.,TI(:,:,islope),dt,tsurf_avg(:,islope),tsoil_avg(:,:,islope))
610
611    do iday = 1,nday
612        do i = 1,ngrid
613            do isoil = 1,nsoil
614                ! Update soil temperature timeseries which is needed to compute the water soil density timeseries
615                tsoil_ts(i,isoil,islope,iday) = tsoil_ts(i,isoil,islope,iday)*tsoil_avg(i,isoil,islope)/tsoil_avg_old(i,isoil)
616                ! Update water soil density
617                h2o_soildensity_ts(i,isoil,islope,iday) = exp(beta_clap_h2o/tsoil_ts(i,isoil,islope,iday) + alpha_clap_h2o)/tsoil_ts(i,isoil,islope,iday)**mmol(iPCM_qh2o)/(mugaz*r)
618                ! Safety check
619                if (isnan(tsoil_avg(i,isoil,islope))) call stop_clean(__FILE__,__LINE__,"NaN detected in tsoil_avg",1)
620            end do
621        end do
622    end do
623end do
624
625! Calculate the time-averaged water soil density
626h2o_soildensity_avg = sum(h2o_soildensity_ts,4)/real(nday,dp)
627
628END SUBROUTINE evolve_soil_temp
629!=======================================================================
630
631END MODULE soil_temp
Note: See TracBrowser for help on using the repository browser.