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

Last change on this file since 4147 was 4138, checked in by jbclement, 2 weeks ago

PEM:

  • Move ice table variables from "ice_table" to the main program.
  • Merge "job_id_mod" and "job_timelimit_mod" into "job" which is relocated to the PEM folder.
  • Rename local variables in procedures to avoid masking variables in parent scope.
  • Few cleanings to delete remaining PEM-external "include" and "use".

JBC

File size: 19.9 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! ----
59allocate(tsoil_PCM(ngrid,nsoil_PCM,nslope))
60allocate(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! 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! ARGUMENTS
299! ---------
300integer(di),              intent(in)    :: ngrid   ! number of (horizontal) grid-points
301integer(di),              intent(in)    :: nsoil   ! number of soil layers
302real(dp), dimension(:,:), intent(in)    :: therm_i ! thermal inertia [SI]
303real(dp), dimension(:),   intent(in)    :: tsurf   ! surface temperature [K]
304real(dp), dimension(:,:), intent(inout) :: tsoil   ! soil (mid-layer) temperature [K]
305
306! LOCAL VARIABLES
307! ---------------
308integer(di) :: ig, ik, iloop
309
310! CODE
311! ----
312! 0. Initialisations and preprocessing step
313! 0.1 Build mthermdiff(:), the mid-layer thermal diffusivities
314do ig = 1,ngrid
315    do ik = 0,nsoil - 1
316        mthermdiff(ig,ik) = therm_i(ig,ik + 1)*therm_i(ig,ik + 1)/volcapa
317    end do
318end do
319
320! 0.2 Build thermdiff(:), the "interlayer" thermal diffusivities
321do ig = 1,ngrid
322    do ik = 1,nsoil - 1
323        thermdiff(ig,ik) = ((layer(ik) - mlayer(ik - 1))*mthermdiff(ig,ik) &
324                             + (mlayer(ik) - layer(ik))*mthermdiff(ig,ik - 1))/(mlayer(ik) - mlayer(ik - 1))
325    end do
326end do
327
328! 0.3 Build coefficients mu, q_{k+1/2}, d_k, alph
329! mu
330mu = mlayer(0)/(mlayer(1) - mlayer(0))
331
332! q_{1/2}
333coefq(:) = 0.
334! q_{k+1/2}
335
336do ig = 1,ngrid
337    ! d_k
338    do ik = 1,nsoil - 1
339        coefd(ig,ik) = thermdiff(ig,ik)/(mlayer(ik) - mlayer(ik - 1))
340    end do
341
342    ! alph_PEM_{N-1}
343    alph(ig,nsoil - 1) = coefd(ig,nsoil - 1)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
344    ! alph_PEM_k
345    do ik = nsoil - 2,1,-1
346        alph(ig,ik) = coefd(ig,ik)/(coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
347    end do
348end do ! of do ig=1,ngrid
349
350!  1. Compute beta coefficients
351! Bottom layer, beta_PEM_{N-1}
352do ig = 1,ngrid
353        beta(ig,nsoil - 1) = coefq(nsoil - 1)*tsoil(ig,nsoil)/(coefq(nsoil - 1) + coefd(ig,nsoil - 1)) &
354                                 + flux_geo/(coefq(nsoil - 1) + coefd(ig,nsoil - 1))
355end do
356! Other layers
357do ik = nsoil - 2,1,-1
358    do ig = 1,ngrid
359        beta(ig,ik) = (coefq(ik)*tsoil(ig,ik + 1) + coefd(ig,ik + 1)*beta(ig,ik + 1))/ &
360                          (coefq(ik) + coefd(ig,ik + 1)*(1._dp - alph(ig,ik + 1)) + coefd(ig,ik))
361    end do
362end do
363
364!  2. Compute soil temperatures
365do iloop = 1,10 !just convergence
366    do ig = 1,ngrid
367        ! First layer:
368        tsoil(ig,1) = (tsurf(ig) + mu*beta(ig,1)*thermdiff(ig,1)/mthermdiff(ig,0))/ &
369                      (1._dp + mu*(1._dp - alph(ig,1))*thermdiff(ig,1)/mthermdiff(ig,0))
370        ! Other layers:
371        do ik = 1,nsoil - 1
372            tsoil(ig,ik + 1) = alph(ig,ik)*tsoil(ig,ik) + beta(ig,ik)
373        end do
374    end do
375end do ! iloop
376
377END SUBROUTINE ini_tsoil_pem
378!=======================================================================
379
380!=======================================================================
381SUBROUTINE shift_tsoil2surf(ngrid,nsoil,nslope,zshift_surf,zlag,tsurf,tsoil)
382!-----------------------------------------------------------------------
383! NAME
384!     shift_tsoil2surf
385!
386! DESCRIPTION
387!     Shift soil temperature profile to follow surface evolution due to ice
388!     condensation/sublimation.
389!
390! AUTHORS & DATE
391!     JB Clement, 2025
392!
393! NOTES
394!
395!-----------------------------------------------------------------------
396
397! DEPENDENCIES
398! ------------
399use soil,    only: layer, mlayer, flux_geo, thermdiff, mthermdiff
400use maths,   only: solve_steady_heat
401use display, only: print_msg, LVL_NFO
402
403! DECLARATION
404! -----------
405implicit none
406
407! ARGUMENTS
408! ---------
409integer(di),                intent(in)    :: ngrid       ! number of (horizontal) grid-points
410integer(di),                intent(in)    :: nsoil       ! number of soil layers
411integer(di),                intent(in)    :: nslope      ! number of sub-slopes
412real(dp), dimension(:,:),   intent(in)    :: zshift_surf ! elevation shift for the surface [m]
413real(dp), dimension(:,:),   intent(in)    :: zlag        ! newly built lag thickness [m]
414real(dp), dimension(:,:),   intent(in)    :: tsurf       ! surface temperature [K]
415real(dp), dimension(:,:,:), intent(inout) :: tsoil       ! soil (mid-layer) temperature [K]
416
417! LOCAL VARIABLES
418! ---------------
419integer(di)                             :: ig, isoil, islope, ishift, ilag
420real(dp)                                :: z, zshift_surfloc
421real(dp), dimension(ngrid,nsoil,nslope) :: tsoil_old
422
423! CODE
424! ----
425call print_msg("> Shifting soil temperature profile to match surface evolution",LVL_NFO)
426tsoil_old = tsoil
427
428do ig = 1,ngrid
429    do islope = 1,nslope
430        zshift_surfloc = zshift_surf(ig,islope)
431
432        if (zshift_surfloc >= 0.) then ! In case of the surface is higher than initially
433            if (zshift_surfloc < mlayer(0)) then ! Surface change is too small to be taken into account
434                ! Nothing to do; we keep the soil temperature profile
435            else if (zshift_surfloc >= mlayer(nsoil - 1)) then ! Surface change is much larger than the discretization
436                tsoil(ig,:,islope) = tsurf(ig,islope)
437            else
438                ishift = 0
439                do while (mlayer(ishift) <= zshift_surfloc)
440                    ishift = ishift + 1 ! mlayer indices begin at 0 so this the good index for tsoil!
441                end do
442                ! The "new soil" temperature is set to tsurf
443                tsoil(ig,:ishift,islope) = tsurf(ig,islope)
444                do isoil = ishift + 1,nsoil
445                    ! Position in the old discretization of the depth
446                    z = mlayer(isoil - 1) - zshift_surfloc
447                    ! Interpolation of the temperature profile from the old discretization
448                    tsoil(ig,isoil,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
449                end do
450            end if
451        else ! In case of the surface is lower than initially
452            if (abs(zshift_surfloc) < mlayer(0)) then ! Surface change is too small to be taken into account
453                ! Nothing to do; we keep the soil temperature profile
454            else if (abs(zshift_surfloc) >= mlayer(nsoil - 1)) then ! Surface change is much larger than the discretization
455                call solve_steady_heat(nsoil,mlayer,layer,mthermdiff(ig,:),thermdiff(ig,:),tsurf(ig,islope),flux_geo,tsoil(ig,:,islope))
456            else
457                if (zlag(ig,islope) < mlayer(0)) then ! The lag is too thin to be taken into account
458                    ilag = 0
459                else
460                    ilag = 0
461                    do while (mlayer(ilag) <= zlag(ig,islope))
462                        ilag = ilag + 1 ! mlayer indices begin at 0 so this the good index for tsoil!
463                    end do
464                    ! Position of the lag bottom in the old discretization of the depth
465                    z = zlag(ig,islope) - zshift_surfloc
466                    ! The "new lag" temperature is set to the ice temperature just below
467                    tsoil(ig,:ilag,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
468                end if
469
470                ishift = nsoil - 1
471                z = mlayer(nsoil - 1) + zshift_surfloc
472                do while (mlayer(ishift) >= z)
473                    ishift = ishift - 1
474                end do
475                ishift = ishift + 1 ! Adding 1 is needed to match the good index for tsoil!
476                do isoil = ilag + 1,ishift
477                    ! Position in the old discretization of the depth
478                    z = mlayer(isoil - 1) - zshift_surfloc
479                    ! Interpolation of the temperature profile from the old discretization
480                    tsoil(ig,isoil,islope) = itp_tsoil(tsoil_old(ig,:,islope),tsurf(ig,islope),z)
481                end do
482
483                ! The "new deepest layers" temperature is set by solving the steady heat equation
484                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))
485            end if
486        end if
487    end do
488end do
489
490END SUBROUTINE shift_tsoil2surf
491!=======================================================================
492
493!=======================================================================
494FUNCTION itp_tsoil(tsoil,tsurf,z) RESULT(tsoil_z)
495!-----------------------------------------------------------------------
496! NAME
497!     itp_tsoil
498!
499! DESCRIPTION
500!     Interpolate soil temperature profile.
501!
502! AUTHORS & DATE
503!     JB Clement, 2025
504!
505! NOTES
506!
507!-----------------------------------------------------------------------
508
509! DEPENDENCIES
510! ------------
511use soil, only: mlayer
512
513! DECLARATION
514! -----------
515implicit none
516
517! ARGUMENTS
518! ---------
519real(dp), dimension(:), intent(in) :: tsoil
520real(dp),               intent(in) :: z, tsurf
521
522! LOCAL VARIABLES
523! ---------------
524real(dp)    :: tsoil_z, tsoil_minus, mlayer_minus, a
525integer(di) :: iz
526
527! CODE
528! ----
529! Find the interval [mlayer(iz - 1),mlayer(iz)[ where the position z belongs
530iz = 0
531do while (mlayer(iz) <= z)
532    iz = iz + 1
533end do
534if (iz == 0) then
535    tsoil_minus = tsurf
536    mlayer_minus = 0._dp
537else
538    tsoil_minus = tsoil(iz)
539    mlayer_minus = mlayer(iz - 1)
540end if
541! Interpolation of the temperature profile from the old discretization
542a = (tsoil(iz + 1) - tsoil_minus)/(mlayer(iz) - mlayer_minus)
543tsoil_z = a*(z - mlayer_minus) + tsoil_minus
544
545END FUNCTION itp_tsoil
546!=======================================================================
547
548!=======================================================================
549SUBROUTINE evolve_soil_temp(tsoil_avg,tsurf_avg,tsoil_ts,tsoil_ts_old,h2o_soildensity_avg)
550!-----------------------------------------------------------------------
551! NAME
552!     evolve_soil_temp
553!
554! DESCRIPTION
555!     Update soil temperature profile and calculate water soil density
556!     time series and averages.
557!
558! AUTHORS & DATE
559!     JB Clement, 02/2026
560!
561! NOTES
562!
563!-----------------------------------------------------------------------
564
565! DEPENDENCIES
566! ------------
567use geometry,   only: ngrid, nsoil, nslope, nday
568use tracers,    only: mmol, iPCM_qh2o
569use physics,    only: mugaz, r, alpha_clap_h2o, beta_clap_h2o
570use evolution,  only: dt
571use soil,       only: TI
572use stoppage,   only: stop_clean
573use display,    only: print_msg, LVL_NFO
574
575! DECLARATION
576! -----------
577implicit none
578
579! ARGUMENTS
580! ---------
581real(dp), dimension(:,:),     intent(in)    :: tsurf_avg
582real(dp), dimension(:,:,:),   intent(inout) :: tsoil_avg
583real(dp), dimension(:,:,:,:), intent(inout) :: tsoil_ts
584real(dp), dimension(:,:,:,:), intent(out)   :: tsoil_ts_old
585real(dp), dimension(:,:,:),   intent(out)   :: h2o_soildensity_avg
586
587! LOCAL VARIABLES
588! ---------------
589real(dp), dimension(ngrid,nsoil)             :: tsoil_avg_old
590real(dp), dimension(ngrid,nsoil,nslope,nday) :: h2o_soildensity_ts
591integer(di)                                  :: i, isoil, islope, iday
592
593! CODE
594! ----
595call print_msg("> Update the soil temperature",LVL_NFO)
596
597! Store current state
598tsoil_ts_old(:,:,:,:) = tsoil_ts(:,:,:,:)
599
600do islope = 1,nslope
601    ! Store the average profile before the update
602    tsoil_avg_old(:,:) = tsoil_avg(:,:,islope)
603   
604    ! Compute the new soil temperature
605    call compute_tsoil(ngrid,nsoil,.true., TI(:,:,islope),dt,tsurf_avg(:,islope),tsoil_avg(:,:,islope))
606    call compute_tsoil(ngrid,nsoil,.false.,TI(:,:,islope),dt,tsurf_avg(:,islope),tsoil_avg(:,:,islope))
607
608    do iday = 1,nday
609        do i = 1,ngrid
610            do isoil = 1,nsoil
611                ! Update soil temperature timeseries which is needed to compute the water soil density timeseries
612                tsoil_ts(i,isoil,islope,iday) = tsoil_ts(i,isoil,islope,iday)*tsoil_avg(i,isoil,islope)/tsoil_avg_old(i,isoil)
613                ! Update water soil density
614                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)
615                ! Safety check
616                if (isnan(tsoil_avg(i,isoil,islope))) call stop_clean(__FILE__,__LINE__,"NaN detected in tsoil_avg",1)
617            end do
618        end do
619    end do
620end do
621
622! Calculate the time-averaged water soil density
623h2o_soildensity_avg = sum(h2o_soildensity_ts,4)/real(nday,dp)
624
625END SUBROUTINE evolve_soil_temp
626!=======================================================================
627
628END MODULE soil_temp
Note: See TracBrowser for help on using the repository browser.