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

Last change on this file since 4074 was 4074, checked in by jbclement, 11 days ago

PEM:

  • Correct management of H2O ice tendency in 1D when there is not enough ice anymore.
  • Clean initialization of allocatable module arrays (especially needed when no slope)
  • One more renaming for consistency + few small updates thoughout the code.

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