New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
limthd_dif.F90 in branches/2014/dev_4728_CNRS04_coupled_interface/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2014/dev_4728_CNRS04_coupled_interface/NEMOGCM/NEMO/LIM_SRC_3/limthd_dif.F90 @ 4733

Last change on this file since 4733 was 4733, checked in by vancop, 10 years ago

Fix energy budget in coupled case

  • Property svn:keywords set to Id
File size: 40.7 KB
Line 
1MODULE limthd_dif
2   !!======================================================================
3   !!                       ***  MODULE limthd_dif ***
4   !!                       heat diffusion in sea ice
5   !!                   computation of surface and inner T 
6   !!======================================================================
7   !! History :  LIM  ! 02-2003 (M. Vancoppenolle) original 1D code
8   !!                 ! 06-2005 (M. Vancoppenolle) 3d version
9   !!                 ! 11-2006 (X Fettweis) Vectorization by Xavier
10   !!                 ! 04-2007 (M. Vancoppenolle) Energy conservation
11   !!            4.0  ! 2011-02 (G. Madec) dynamical allocation
12   !!             -   ! 2012-05 (C. Rousset) add penetration solar flux
13   !!----------------------------------------------------------------------
14#if defined key_lim3
15   !!----------------------------------------------------------------------
16   !!   'key_lim3'                                      LIM3 sea-ice model
17   !!----------------------------------------------------------------------
18   USE par_oce        ! ocean parameters
19   USE phycst         ! physical constants (ocean directory)
20   USE ice            ! LIM-3 variables
21   USE par_ice        ! LIM-3 parameters
22   USE thd_ice        ! LIM-3: thermodynamics
23   USE in_out_manager ! I/O manager
24   USE lib_mpp        ! MPP library
25   USE wrk_nemo       ! work arrays
26   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
27   USE cpl_oasis3, ONLY : lk_cpl
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   lim_thd_dif   ! called by lim_thd
33
34   REAL(wp) ::   epsi10 = 1.e-10_wp    !
35   !!----------------------------------------------------------------------
36   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
37   !! $Id$
38   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
39   !!----------------------------------------------------------------------
40CONTAINS
41
42   SUBROUTINE lim_thd_dif( kideb , kiut )
43      !!------------------------------------------------------------------
44      !!                ***  ROUTINE lim_thd_dif  ***
45      !! ** Purpose :
46      !!           This routine determines the time evolution of snow and sea-ice
47      !!           temperature profiles.
48      !! ** Method  :
49      !!           This is done by solving the heat equation diffusion with
50      !!           a Neumann boundary condition at the surface and a Dirichlet one
51      !!           at the bottom. Solar radiation is partially absorbed into the ice.
52      !!           The specific heat and thermal conductivities depend on ice salinity
53      !!           and temperature to take into account brine pocket melting. The
54      !!           numerical
55      !!           scheme is an iterative Crank-Nicolson on a non-uniform multilayer grid
56      !!           in the ice and snow system.
57      !!
58      !!           The successive steps of this routine are
59      !!           1.  Thermal conductivity at the interfaces of the ice layers
60      !!           2.  Internal absorbed radiation
61      !!           3.  Scale factors due to non-uniform grid
62      !!           4.  Kappa factors
63      !!           Then iterative procedure begins
64      !!           5.  specific heat in the ice
65      !!           6.  eta factors
66      !!           7.  surface flux computation
67      !!           8.  tridiagonal system terms
68      !!           9.  solving the tridiagonal system with Gauss elimination
69      !!           Iterative procedure ends according to a criterion on evolution
70      !!           of temperature
71      !!
72      !! ** Arguments :
73      !!           kideb , kiut : Starting and ending points on which the
74      !!                         the computation is applied
75      !!
76      !! ** Inputs / Ouputs : (global commons)
77      !!           surface temperature : t_su_b
78      !!           ice/snow temperatures   : t_i_b, t_s_b
79      !!           ice salinities          : s_i_b
80      !!           number of layers in the ice/snow: nlay_i, nlay_s
81      !!           profile of the ice/snow layers : z_i, z_s
82      !!           total ice/snow thickness : ht_i_b, ht_s_b
83      !!
84      !! ** External :
85      !!
86      !! ** References :
87      !!
88      !! ** History :
89      !!           (02-2003) Martin Vancoppenolle, Louvain-la-Neuve, Belgium
90      !!           (06-2005) Martin Vancoppenolle, 3d version
91      !!           (11-2006) Vectorized by Xavier Fettweis (UCL-ASTR)
92      !!           (04-2007) Energy conservation tested by M. Vancoppenolle
93      !!------------------------------------------------------------------
94      INTEGER , INTENT(in) ::   kideb, kiut   ! Start/End point on which the  the computation is applied
95
96      !! * Local variables
97      INTEGER ::   ji          ! spatial loop index
98      INTEGER ::   ii, ij      ! temporary dummy loop index
99      INTEGER ::   numeq       ! current reference number of equation
100      INTEGER ::   layer       ! vertical dummy loop index
101      INTEGER ::   nconv       ! number of iterations in iterative procedure
102      INTEGER ::   minnumeqmin, maxnumeqmax
103      INTEGER, POINTER, DIMENSION(:) ::   numeqmin   ! reference number of top equation
104      INTEGER, POINTER, DIMENSION(:) ::   numeqmax   ! reference number of bottom equation
105      INTEGER, POINTER, DIMENSION(:) ::   isnow      ! switch for presence (1) or absence (0) of snow
106      REAL(wp) ::   zg1s      =  2._wp        ! for the tridiagonal system
107      REAL(wp) ::   zg1       =  2._wp        !
108      REAL(wp) ::   zgamma    =  18009._wp    ! for specific heat
109      REAL(wp) ::   zbeta     =  0.117_wp     ! for thermal conductivity (could be 0.13)
110      REAL(wp) ::   zraext_s  =  1.e+8_wp     ! extinction coefficient of radiation in the snow
111      REAL(wp) ::   zkimin    =  0.10_wp      ! minimum ice thermal conductivity
112      REAL(wp) ::   ztsu_err  =  1.e-5_wp     ! range around which t_su is considered as 0°C
113      REAL(wp) ::   ztmelt_i    ! ice melting temperature
114      REAL(wp) ::   zerritmax   ! current maximal error on temperature
115      REAL(wp), POINTER, DIMENSION(:) ::   ztfs        ! ice melting point
116      REAL(wp), POINTER, DIMENSION(:) ::   ztsuold     ! old surface temperature (before the iterative procedure )
117      REAL(wp), POINTER, DIMENSION(:) ::   ztsuoldit   ! surface temperature at previous iteration
118      REAL(wp), POINTER, DIMENSION(:) ::   zh_i        ! ice layer thickness
119      REAL(wp), POINTER, DIMENSION(:) ::   zh_s        ! snow layer thickness
120      REAL(wp), POINTER, DIMENSION(:) ::   zfsw        ! solar radiation absorbed at the surface
121      REAL(wp), POINTER, DIMENSION(:) ::   zf          ! surface flux function
122      REAL(wp), POINTER, DIMENSION(:) ::   dzf         ! derivative of the surface flux function
123      REAL(wp), POINTER, DIMENSION(:) ::   zerrit      ! current error on temperature
124      REAL(wp), POINTER, DIMENSION(:) ::   zdifcase    ! case of the equation resolution (1->4)
125      REAL(wp), POINTER, DIMENSION(:) ::   zftrice     ! solar radiation transmitted through the ice
126      REAL(wp), POINTER, DIMENSION(:) ::   zihic, zhsu
127      REAL(wp), POINTER, DIMENSION(:,:) ::   ztcond_i    ! Ice thermal conductivity
128      REAL(wp), POINTER, DIMENSION(:,:) ::   zradtr_i    ! Radiation transmitted through the ice
129      REAL(wp), POINTER, DIMENSION(:,:) ::   zradab_i    ! Radiation absorbed in the ice
130      REAL(wp), POINTER, DIMENSION(:,:) ::   zkappa_i    ! Kappa factor in the ice
131      REAL(wp), POINTER, DIMENSION(:,:) ::   ztiold      ! Old temperature in the ice
132      REAL(wp), POINTER, DIMENSION(:,:) ::   zeta_i      ! Eta factor in the ice
133      REAL(wp), POINTER, DIMENSION(:,:) ::   ztitemp     ! Temporary temperature in the ice to check the convergence
134      REAL(wp), POINTER, DIMENSION(:,:) ::   zspeche_i   ! Ice specific heat
135      REAL(wp), POINTER, DIMENSION(:,:) ::   z_i         ! Vertical cotes of the layers in the ice
136      REAL(wp), POINTER, DIMENSION(:,:) ::   zradtr_s    ! Radiation transmited through the snow
137      REAL(wp), POINTER, DIMENSION(:,:) ::   zradab_s    ! Radiation absorbed in the snow
138      REAL(wp), POINTER, DIMENSION(:,:) ::   zkappa_s    ! Kappa factor in the snow
139      REAL(wp), POINTER, DIMENSION(:,:) ::   zeta_s       ! Eta factor in the snow
140      REAL(wp), POINTER, DIMENSION(:,:) ::   ztstemp      ! Temporary temperature in the snow to check the convergence
141      REAL(wp), POINTER, DIMENSION(:,:) ::   ztsold       ! Temporary temperature in the snow
142      REAL(wp), POINTER, DIMENSION(:,:) ::   z_s          ! Vertical cotes of the layers in the snow
143      REAL(wp), POINTER, DIMENSION(:,:) ::   zindterm   ! Independent term
144      REAL(wp), POINTER, DIMENSION(:,:) ::   zindtbis   ! temporary independent term
145      REAL(wp), POINTER, DIMENSION(:,:) ::   zdiagbis
146      REAL(wp), POINTER, DIMENSION(:,:,:) ::   ztrid   ! tridiagonal system terms
147      ! diag errors on heat
148      REAL(wp), POINTER, DIMENSION(:) :: zdq, zq_ini, zhfx_err
149      !!------------------------------------------------------------------     
150      !
151      CALL wrk_alloc( jpij, numeqmin, numeqmax, isnow )
152      CALL wrk_alloc( jpij, ztfs, ztsuold, ztsuoldit, zh_i, zh_s, zfsw )
153      CALL wrk_alloc( jpij, zf, dzf, zerrit, zdifcase, zftrice, zihic, zhsu )
154      CALL wrk_alloc( jpij, nlay_i+1, ztcond_i, zradtr_i, zradab_i, zkappa_i, ztiold, zeta_i, ztitemp, z_i, zspeche_i, kjstart=0)
155      CALL wrk_alloc( jpij, nlay_s+1,           zradtr_s, zradab_s, zkappa_s, ztsold, zeta_s, ztstemp, z_s, kjstart=0)
156      CALL wrk_alloc( jpij, jkmax+2, zindterm, zindtbis, zdiagbis  )
157      CALL wrk_alloc( jpij, jkmax+2, 3, ztrid )
158
159      CALL wrk_alloc( jpij, zdq, zq_ini, zhfx_err )
160
161      ! --- diag error on heat diffusion - PART 1 --- !
162      zdq(:) = 0._wp ; zq_ini(:) = 0._wp     
163      DO ji = kideb, kiut
164         zq_ini(ji) = ( SUM( q_i_b(ji,1:nlay_i) ) * ht_i_b(ji) / REAL( nlay_i ) +  &
165            &           SUM( q_s_b(ji,1:nlay_s) ) * ht_s_b(ji) / REAL( nlay_s ) ) 
166      END DO
167
168      !------------------------------------------------------------------------------!
169      ! 1) Initialization                                                            !
170      !------------------------------------------------------------------------------!
171      ! clem clean: replace just ztfs by rtt
172      DO ji = kideb , kiut
173         ! is there snow or not
174         isnow(ji)= NINT(  1._wp - MAX( 0._wp , SIGN(1._wp, - ht_s_b(ji) ) )  )
175         ! surface temperature of fusion
176         ztfs(ji) = REAL( isnow(ji) ) * rtt + REAL( 1 - isnow(ji) ) * rtt
177         ! layer thickness
178         zh_i(ji) = ht_i_b(ji) / REAL( nlay_i )
179         zh_s(ji) = ht_s_b(ji) / REAL( nlay_s )
180      END DO
181
182      !--------------------
183      ! Ice / snow layers
184      !--------------------
185
186      z_s(:,0) = 0._wp   ! vert. coord. of the up. lim. of the 1st snow layer
187      z_i(:,0) = 0._wp   ! vert. coord. of the up. lim. of the 1st ice layer
188
189      DO layer = 1, nlay_s            ! vert. coord of the up. lim. of the layer-th snow layer
190         DO ji = kideb , kiut
191            z_s(ji,layer) = z_s(ji,layer-1) + ht_s_b(ji) / REAL( nlay_s )
192         END DO
193      END DO
194
195      DO layer = 1, nlay_i            ! vert. coord of the up. lim. of the layer-th ice layer
196         DO ji = kideb , kiut
197            z_i(ji,layer) = z_i(ji,layer-1) + ht_i_b(ji) / REAL( nlay_i )
198         END DO
199      END DO
200      !
201      !------------------------------------------------------------------------------|
202      ! 2) Radiations                                                                |
203      !------------------------------------------------------------------------------|
204      !
205      !-------------------
206      ! Computation of i0
207      !-------------------
208      ! i0 describes the fraction of solar radiation which does not contribute
209      ! to the surface energy budget but rather penetrates inside the ice.
210      ! We assume that no radiation is transmitted through the snow
211      ! If there is no no snow
212      ! zfsw    = (1-i0).qsr_ice   is absorbed at the surface
213      ! zftrice = io.qsr_ice       is below the surface
214      ! ftr_ice = io.qsr_ice.exp(-k(h_i)) transmitted below the ice
215
216      DO ji = kideb , kiut
217         ! switches
218         isnow(ji) = NINT(  1._wp - MAX( 0._wp , SIGN( 1._wp , - ht_s_b(ji) ) )  ) 
219         ! hs > 0, isnow = 1
220         zhsu (ji) = hnzst  ! threshold for the computation of i0
221         zihic(ji) = MAX( 0._wp , 1._wp - ( ht_i_b(ji) / zhsu(ji) ) )     
222
223         i0(ji)    = REAL( 1 - isnow(ji) ) * ( fr1_i0_1d(ji) + zihic(ji) * fr2_i0_1d(ji) )
224         !fr1_i0_1d = i0 for a thin ice surface
225         !fr1_i0_2d = i0 for a thick ice surface
226         !            a function of the cloud cover
227         !
228         !i0(ji)     =  (1.0-FLOAT(isnow(ji)))*3.0/(100*ht_s_b(ji)+10.0)
229         !formula used in Cice
230      END DO
231
232      !-------------------------------------------------------
233      ! Solar radiation absorbed / transmitted at the surface
234      ! Derivative of the non solar flux
235      !-------------------------------------------------------
236      DO ji = kideb , kiut
237         zfsw   (ji) =  qsr_ice_1d(ji) * ( 1 - i0(ji) )   ! Shortwave radiation absorbed at surface
238         zftrice(ji) =  qsr_ice_1d(ji) *       i0(ji)     ! Solar radiation transmitted below the surface layer
239         dzf    (ji) = dqns_ice_1d(ji)                    ! derivative of incoming nonsolar flux
240      END DO
241
242      !---------------------------------------------------------
243      ! Transmission - absorption of solar radiation in the ice
244      !---------------------------------------------------------
245
246      DO ji = kideb, kiut           ! snow initialization
247         zradtr_s(ji,0) = zftrice(ji)     ! radiation penetrating through snow
248      END DO
249
250      DO layer = 1, nlay_s          ! Radiation through snow
251         DO ji = kideb, kiut
252            !                             ! radiation transmitted below the layer-th snow layer
253            zradtr_s(ji,layer) = zradtr_s(ji,0) * EXP( - zraext_s * ( MAX ( 0._wp , z_s(ji,layer) ) ) )
254            !                             ! radiation absorbed by the layer-th snow layer
255            zradab_s(ji,layer) = zradtr_s(ji,layer-1) - zradtr_s(ji,layer)
256         END DO
257      END DO
258
259      DO ji = kideb, kiut           ! ice initialization
260         zradtr_i(ji,0) = zradtr_s(ji,nlay_s) * REAL( isnow(ji) ) + zftrice(ji) * REAL( 1 - isnow(ji) )
261      END DO
262
263      DO layer = 1, nlay_i          ! Radiation through ice
264         DO ji = kideb, kiut
265            !                             ! radiation transmitted below the layer-th ice layer
266            zradtr_i(ji,layer) = zradtr_i(ji,0) * EXP( - kappa_i * ( MAX ( 0._wp , z_i(ji,layer) ) ) )
267            !                             ! radiation absorbed by the layer-th ice layer
268            zradab_i(ji,layer) = zradtr_i(ji,layer-1) - zradtr_i(ji,layer)
269         END DO
270      END DO
271
272      DO ji = kideb, kiut           ! Radiation transmitted below the ice
273         !!!ftr_ice_1d(ji) = ftr_ice_1d(ji) + iatte_1d(ji) * zradtr_i(ji,nlay_i) * a_i_b(ji) / at_i_b(ji) ! clem modif
274         ftr_ice_1d(ji) = zradtr_i(ji,nlay_i) 
275      END DO
276
277      !
278      !------------------------------------------------------------------------------|
279      !  3) Iterative procedure begins                                               |
280      !------------------------------------------------------------------------------|
281      !
282      DO ji = kideb, kiut        ! Old surface temperature
283         ztsuold  (ji) =  t_su_b(ji)                              ! temperature at the beg of iter pr.
284         ztsuoldit(ji) =  t_su_b(ji)                              ! temperature at the previous iter
285         t_su_b   (ji) =  MIN( t_su_b(ji), ztfs(ji) - ztsu_err )  ! necessary
286         zerrit   (ji) =  1000._wp                                ! initial value of error
287      END DO
288
289      DO layer = 1, nlay_s       ! Old snow temperature
290         DO ji = kideb , kiut
291            ztsold(ji,layer) =  t_s_b(ji,layer)
292         END DO
293      END DO
294
295      DO layer = 1, nlay_i       ! Old ice temperature
296         DO ji = kideb , kiut
297            ztiold(ji,layer) =  t_i_b(ji,layer)
298         END DO
299      END DO
300
301      nconv     =  0           ! number of iterations
302      zerritmax =  1000._wp    ! maximal value of error on all points
303
304      DO WHILE ( zerritmax > maxer_i_thd .AND. nconv < nconv_i_thd )
305         !
306         nconv = nconv + 1
307         !
308         !------------------------------------------------------------------------------|
309         ! 4) Sea ice thermal conductivity                                              |
310         !------------------------------------------------------------------------------|
311         !
312         IF( thcon_i_swi == 0 ) THEN      ! Untersteiner (1964) formula
313            DO ji = kideb , kiut
314               ztcond_i(ji,0)        = rcdic + zbeta*s_i_b(ji,1) / MIN(-epsi10,t_i_b(ji,1)-rtt)
315               ztcond_i(ji,0)        = MAX(ztcond_i(ji,0),zkimin)
316            END DO
317            DO layer = 1, nlay_i-1
318               DO ji = kideb , kiut
319                  ztcond_i(ji,layer) = rcdic + zbeta*( s_i_b(ji,layer) + s_i_b(ji,layer+1) ) /  &
320                     MIN(-2.0_wp * epsi10, t_i_b(ji,layer)+t_i_b(ji,layer+1) - 2.0_wp * rtt)
321                  ztcond_i(ji,layer) = MAX(ztcond_i(ji,layer),zkimin)
322               END DO
323            END DO
324         ENDIF
325
326         IF( thcon_i_swi == 1 ) THEN      ! Pringle et al formula included: 2.11 + 0.09 S/T - 0.011.T
327            DO ji = kideb , kiut
328               ztcond_i(ji,0) = rcdic + 0.090_wp * s_i_b(ji,1) / MIN( -epsi10, t_i_b(ji,1)-rtt )   &
329                  &                   - 0.011_wp * ( t_i_b(ji,1) - rtt ) 
330               ztcond_i(ji,0) = MAX( ztcond_i(ji,0), zkimin )
331            END DO
332            DO layer = 1, nlay_i-1
333               DO ji = kideb , kiut
334                  ztcond_i(ji,layer) = rcdic + 0.090_wp * ( s_i_b(ji,layer) + s_i_b(ji,layer+1) )   &
335                     &                                  / MIN(-2.0_wp * epsi10, t_i_b(ji,layer)+t_i_b(ji,layer+1) - 2.0_wp * rtt)   &
336                     &                       - 0.0055_wp* ( t_i_b(ji,layer) + t_i_b(ji,layer+1) - 2.0*rtt ) 
337                  ztcond_i(ji,layer) = MAX( ztcond_i(ji,layer), zkimin )
338               END DO
339            END DO
340            DO ji = kideb , kiut
341               ztcond_i(ji,nlay_i) = rcdic + 0.090_wp * s_i_b(ji,nlay_i) / MIN(-epsi10,t_bo_b(ji)-rtt)   &
342                  &                        - 0.011_wp * ( t_bo_b(ji) - rtt ) 
343               ztcond_i(ji,nlay_i) = MAX( ztcond_i(ji,nlay_i), zkimin )
344            END DO
345         ENDIF
346         !
347         !------------------------------------------------------------------------------|
348         !  5) kappa factors                                                            |
349         !------------------------------------------------------------------------------|
350         !
351         DO ji = kideb, kiut
352
353            !-- Snow kappa factors
354            zkappa_s(ji,0)         = rcdsn / MAX(epsi10,zh_s(ji))
355            zkappa_s(ji,nlay_s)    = rcdsn / MAX(epsi10,zh_s(ji))
356         END DO
357
358         DO layer = 1, nlay_s-1
359            DO ji = kideb , kiut
360               zkappa_s(ji,layer)  = 2.0 * rcdsn / &
361                  MAX(epsi10,2.0*zh_s(ji))
362            END DO
363         END DO
364
365         DO layer = 1, nlay_i-1
366            DO ji = kideb , kiut
367               !-- Ice kappa factors
368               zkappa_i(ji,layer)  = 2.0*ztcond_i(ji,layer)/ &
369                  MAX(epsi10,2.0*zh_i(ji)) 
370            END DO
371         END DO
372
373         DO ji = kideb , kiut
374            zkappa_i(ji,0)        = ztcond_i(ji,0)/MAX(epsi10,zh_i(ji))
375            zkappa_i(ji,nlay_i)   = ztcond_i(ji,nlay_i) / MAX(epsi10,zh_i(ji))
376            !-- Interface
377            zkappa_s(ji,nlay_s)   = 2.0*rcdsn*ztcond_i(ji,0)/MAX(epsi10, &
378               (ztcond_i(ji,0)*zh_s(ji) + rcdsn*zh_i(ji)))
379            zkappa_i(ji,0)        = zkappa_s(ji,nlay_s)*REAL( isnow(ji) ) &
380               + zkappa_i(ji,0)*REAL( 1 - isnow(ji) )
381         END DO
382         !
383         !------------------------------------------------------------------------------|
384         ! 6) Sea ice specific heat, eta factors                                        |
385         !------------------------------------------------------------------------------|
386         !
387         DO layer = 1, nlay_i
388            DO ji = kideb , kiut
389               ztitemp(ji,layer)   = t_i_b(ji,layer)
390               zspeche_i(ji,layer) = cpic + zgamma*s_i_b(ji,layer)/ &
391                  MAX((t_i_b(ji,layer)-rtt)*(ztiold(ji,layer)-rtt),epsi10)
392               zeta_i(ji,layer)    = rdt_ice / MAX(rhoic*zspeche_i(ji,layer)*zh_i(ji), &
393                  epsi10)
394            END DO
395         END DO
396
397         DO layer = 1, nlay_s
398            DO ji = kideb , kiut
399               ztstemp(ji,layer) = t_s_b(ji,layer)
400               zeta_s(ji,layer)  = rdt_ice / MAX(rhosn*cpic*zh_s(ji),epsi10)
401            END DO
402         END DO
403         !
404         !------------------------------------------------------------------------------|
405         ! 7) surface flux computation                                                  |
406         !------------------------------------------------------------------------------|
407         !
408         IF( .NOT. lk_cpl ) THEN   !--- forced atmosphere case
409            DO ji = kideb , kiut
410               ! update of the non solar flux according to the update in T_su
411               qns_ice_1d(ji) = qns_ice_1d(ji) + dqns_ice_1d(ji) * ( t_su_b(ji) - ztsuoldit(ji) )
412            END DO
413         ENDIF
414
415         ! Update incoming flux
416         DO ji = kideb , kiut
417            ! update incoming flux
418            zf(ji)    =   zfsw(ji)              & ! net absorbed solar radiation
419               + qns_ice_1d(ji)                   ! non solar total flux
420            ! (LWup, LWdw, SH, LH)
421         END DO
422
423         !
424         !------------------------------------------------------------------------------|
425         ! 8) tridiagonal system terms                                                  |
426         !------------------------------------------------------------------------------|
427         !
428         !!layer denotes the number of the layer in the snow or in the ice
429         !!numeq denotes the reference number of the equation in the tridiagonal
430         !!system, terms of tridiagonal system are indexed as following :
431         !!1 is subdiagonal term, 2 is diagonal and 3 is superdiagonal one
432
433         !!ice interior terms (top equation has the same form as the others)
434
435         DO numeq=1,jkmax+2
436            DO ji = kideb , kiut
437               ztrid(ji,numeq,1) = 0.
438               ztrid(ji,numeq,2) = 0.
439               ztrid(ji,numeq,3) = 0.
440               zindterm(ji,numeq)= 0.
441               zindtbis(ji,numeq)= 0.
442               zdiagbis(ji,numeq)= 0.
443            ENDDO
444         ENDDO
445
446         DO numeq = nlay_s + 2, nlay_s + nlay_i 
447            DO ji = kideb , kiut
448               layer              = numeq - nlay_s - 1
449               ztrid(ji,numeq,1)  =  - zeta_i(ji,layer)*zkappa_i(ji,layer-1)
450               ztrid(ji,numeq,2)  =  1.0 + zeta_i(ji,layer)*(zkappa_i(ji,layer-1) + &
451                  zkappa_i(ji,layer))
452               ztrid(ji,numeq,3)  =  - zeta_i(ji,layer)*zkappa_i(ji,layer)
453               zindterm(ji,numeq) =  ztiold(ji,layer) + zeta_i(ji,layer)* &
454                  zradab_i(ji,layer)
455            END DO
456         ENDDO
457
458         numeq =  nlay_s + nlay_i + 1
459         DO ji = kideb , kiut
460            !!ice bottom term
461            ztrid(ji,numeq,1)  =  - zeta_i(ji,nlay_i)*zkappa_i(ji,nlay_i-1)   
462            ztrid(ji,numeq,2)  =  1.0 + zeta_i(ji,nlay_i)*( zkappa_i(ji,nlay_i)*zg1 &
463               +  zkappa_i(ji,nlay_i-1) )
464            ztrid(ji,numeq,3)  =  0.0
465            zindterm(ji,numeq) =  ztiold(ji,nlay_i) + zeta_i(ji,nlay_i)* &
466               ( zradab_i(ji,nlay_i) + zkappa_i(ji,nlay_i)*zg1 &
467               *  t_bo_b(ji) ) 
468         ENDDO
469
470
471         DO ji = kideb , kiut
472            IF ( ht_s_b(ji).gt.0.0 ) THEN
473               !
474               !------------------------------------------------------------------------------|
475               !  snow-covered cells                                                          |
476               !------------------------------------------------------------------------------|
477               !
478               !!snow interior terms (bottom equation has the same form as the others)
479               DO numeq = 3, nlay_s + 1
480                  layer =  numeq - 1
481                  ztrid(ji,numeq,1)   =  - zeta_s(ji,layer)*zkappa_s(ji,layer-1)
482                  ztrid(ji,numeq,2)   =  1.0 + zeta_s(ji,layer)*( zkappa_s(ji,layer-1) + &
483                     zkappa_s(ji,layer) )
484                  ztrid(ji,numeq,3)   =  - zeta_s(ji,layer)*zkappa_s(ji,layer)
485                  zindterm(ji,numeq)  =  ztsold(ji,layer) + zeta_s(ji,layer)* &
486                     zradab_s(ji,layer)
487               END DO
488
489               !!case of only one layer in the ice (ice equation is altered)
490               IF ( nlay_i.eq.1 ) THEN
491                  ztrid(ji,nlay_s+2,3)    =  0.0
492                  zindterm(ji,nlay_s+2)   =  zindterm(ji,nlay_s+2) + zkappa_i(ji,1)* &
493                     t_bo_b(ji) 
494               ENDIF
495
496               IF ( t_su_b(ji) .LT. rtt ) THEN
497
498                  !------------------------------------------------------------------------------|
499                  !  case 1 : no surface melting - snow present                                  |
500                  !------------------------------------------------------------------------------|
501                  zdifcase(ji)    =  1.0
502                  numeqmin(ji)    =  1
503                  numeqmax(ji)    =  nlay_i + nlay_s + 1
504
505                  !!surface equation
506                  ztrid(ji,1,1) = 0.0
507                  ztrid(ji,1,2) = dzf(ji) - zg1s*zkappa_s(ji,0)
508                  ztrid(ji,1,3) = zg1s*zkappa_s(ji,0)
509                  zindterm(ji,1) = dzf(ji)*t_su_b(ji)   - zf(ji)
510
511                  !!first layer of snow equation
512                  ztrid(ji,2,1)  =  - zkappa_s(ji,0)*zg1s*zeta_s(ji,1)
513                  ztrid(ji,2,2)  =  1.0 + zeta_s(ji,1)*(zkappa_s(ji,1) + zkappa_s(ji,0)*zg1s)
514                  ztrid(ji,2,3)  =  - zeta_s(ji,1)* zkappa_s(ji,1)
515                  zindterm(ji,2) =  ztsold(ji,1) + zeta_s(ji,1)*zradab_s(ji,1)
516
517               ELSE 
518                  !
519                  !------------------------------------------------------------------------------|
520                  !  case 2 : surface is melting - snow present                                  |
521                  !------------------------------------------------------------------------------|
522                  !
523                  zdifcase(ji)    =  2.0
524                  numeqmin(ji)    =  2
525                  numeqmax(ji)    =  nlay_i + nlay_s + 1
526
527                  !!first layer of snow equation
528                  ztrid(ji,2,1)  =  0.0
529                  ztrid(ji,2,2)  =  1.0 + zeta_s(ji,1) * ( zkappa_s(ji,1) + &
530                     zkappa_s(ji,0) * zg1s )
531                  ztrid(ji,2,3)  =  - zeta_s(ji,1)*zkappa_s(ji,1) 
532                  zindterm(ji,2) = ztsold(ji,1) + zeta_s(ji,1) *            &
533                     ( zradab_s(ji,1) +                         &
534                     zkappa_s(ji,0) * zg1s * t_su_b(ji) ) 
535               ENDIF
536            ELSE
537               !
538               !------------------------------------------------------------------------------|
539               !  cells without snow                                                          |
540               !------------------------------------------------------------------------------|
541               !
542               IF (t_su_b(ji) .LT. rtt) THEN
543                  !
544                  !------------------------------------------------------------------------------|
545                  !  case 3 : no surface melting - no snow                                       |
546                  !------------------------------------------------------------------------------|
547                  !
548                  zdifcase(ji)      =  3.0
549                  numeqmin(ji)      =  nlay_s + 1
550                  numeqmax(ji)      =  nlay_i + nlay_s + 1
551
552                  !!surface equation   
553                  ztrid(ji,numeqmin(ji),1)   =  0.0
554                  ztrid(ji,numeqmin(ji),2)   =  dzf(ji) - zkappa_i(ji,0)*zg1   
555                  ztrid(ji,numeqmin(ji),3)   =  zkappa_i(ji,0)*zg1
556                  zindterm(ji,numeqmin(ji))  =  dzf(ji)*t_su_b(ji) - zf(ji)
557
558                  !!first layer of ice equation
559                  ztrid(ji,numeqmin(ji)+1,1) =  - zkappa_i(ji,0) * zg1 * zeta_i(ji,1)
560                  ztrid(ji,numeqmin(ji)+1,2) =  1.0 + zeta_i(ji,1) * ( zkappa_i(ji,1) & 
561                     + zkappa_i(ji,0) * zg1 )
562                  ztrid(ji,numeqmin(ji)+1,3) =  - zeta_i(ji,1)*zkappa_i(ji,1) 
563                  zindterm(ji,numeqmin(ji)+1)=  ztiold(ji,1) + zeta_i(ji,1)*zradab_i(ji,1) 
564
565                  !!case of only one layer in the ice (surface & ice equations are altered)
566
567                  IF (nlay_i.eq.1) THEN
568                     ztrid(ji,numeqmin(ji),1)    =  0.0
569                     ztrid(ji,numeqmin(ji),2)    =  dzf(ji) - zkappa_i(ji,0)*2.0
570                     ztrid(ji,numeqmin(ji),3)    =  zkappa_i(ji,0)*2.0
571                     ztrid(ji,numeqmin(ji)+1,1)  =  -zkappa_i(ji,0)*2.0*zeta_i(ji,1)
572                     ztrid(ji,numeqmin(ji)+1,2)  =  1.0 + zeta_i(ji,1)*(zkappa_i(ji,0)*2.0 + &
573                        zkappa_i(ji,1))
574                     ztrid(ji,numeqmin(ji)+1,3)  =  0.0
575
576                     zindterm(ji,numeqmin(ji)+1) =  ztiold(ji,1) + zeta_i(ji,1)* &
577                        ( zradab_i(ji,1) + zkappa_i(ji,1)*t_bo_b(ji) )
578                  ENDIF
579
580               ELSE
581
582                  !
583                  !------------------------------------------------------------------------------|
584                  ! case 4 : surface is melting - no snow                                        |
585                  !------------------------------------------------------------------------------|
586                  !
587                  zdifcase(ji)    =  4.0
588                  numeqmin(ji)    =  nlay_s + 2
589                  numeqmax(ji)    =  nlay_i + nlay_s + 1
590
591                  !!first layer of ice equation
592                  ztrid(ji,numeqmin(ji),1)      =  0.0
593                  ztrid(ji,numeqmin(ji),2)      =  1.0 + zeta_i(ji,1)*(zkappa_i(ji,1) + zkappa_i(ji,0)* &
594                     zg1) 
595                  ztrid(ji,numeqmin(ji),3)      =  - zeta_i(ji,1) * zkappa_i(ji,1)
596                  zindterm(ji,numeqmin(ji))     =  ztiold(ji,1) + zeta_i(ji,1)*( zradab_i(ji,1) + &
597                     zkappa_i(ji,0) * zg1 * t_su_b(ji) ) 
598
599                  !!case of only one layer in the ice (surface & ice equations are altered)
600                  IF (nlay_i.eq.1) THEN
601                     ztrid(ji,numeqmin(ji),1)  =  0.0
602                     ztrid(ji,numeqmin(ji),2)  =  1.0 + zeta_i(ji,1)*(zkappa_i(ji,0)*2.0 + &
603                        zkappa_i(ji,1))
604                     ztrid(ji,numeqmin(ji),3)  =  0.0
605                     zindterm(ji,numeqmin(ji)) =  ztiold(ji,1) + zeta_i(ji,1)* &
606                        (zradab_i(ji,1) + zkappa_i(ji,1)*t_bo_b(ji)) &
607                        + t_su_b(ji)*zeta_i(ji,1)*zkappa_i(ji,0)*2.0
608                  ENDIF
609
610               ENDIF
611            ENDIF
612
613         END DO
614
615         !
616         !------------------------------------------------------------------------------|
617         ! 9) tridiagonal system solving                                                |
618         !------------------------------------------------------------------------------|
619         !
620
621         ! Solve the tridiagonal system with Gauss elimination method.
622         ! Thomas algorithm, from Computational fluid Dynamics, J.D. ANDERSON,
623         ! McGraw-Hill 1984. 
624
625         maxnumeqmax = 0
626         minnumeqmin = jkmax+4
627
628         DO ji = kideb , kiut
629            zindtbis(ji,numeqmin(ji)) =  zindterm(ji,numeqmin(ji))
630            zdiagbis(ji,numeqmin(ji)) =  ztrid(ji,numeqmin(ji),2)
631            minnumeqmin               =  MIN(numeqmin(ji),minnumeqmin)
632            maxnumeqmax               =  MAX(numeqmax(ji),maxnumeqmax)
633         END DO
634
635         DO layer = minnumeqmin+1, maxnumeqmax
636            DO ji = kideb , kiut
637               numeq               =  min(max(numeqmin(ji)+1,layer),numeqmax(ji))
638               zdiagbis(ji,numeq)  =  ztrid(ji,numeq,2) - ztrid(ji,numeq,1)* &
639                  ztrid(ji,numeq-1,3)/zdiagbis(ji,numeq-1)
640               zindtbis(ji,numeq)  =  zindterm(ji,numeq) - ztrid(ji,numeq,1)* &
641                  zindtbis(ji,numeq-1)/zdiagbis(ji,numeq-1)
642            END DO
643         END DO
644
645         DO ji = kideb , kiut
646            ! ice temperatures
647            t_i_b(ji,nlay_i)    =  zindtbis(ji,numeqmax(ji))/zdiagbis(ji,numeqmax(ji))
648         END DO
649
650         DO numeq = nlay_i + nlay_s + 1, nlay_s + 2, -1
651            DO ji = kideb , kiut
652               layer    =  numeq - nlay_s - 1
653               t_i_b(ji,layer)  =  (zindtbis(ji,numeq) - ztrid(ji,numeq,3)* &
654                  t_i_b(ji,layer+1))/zdiagbis(ji,numeq)
655            END DO
656         END DO
657
658         DO ji = kideb , kiut
659            ! snow temperatures     
660            IF (ht_s_b(ji).GT.0._wp) &
661               t_s_b(ji,nlay_s)     =  (zindtbis(ji,nlay_s+1) - ztrid(ji,nlay_s+1,3) &
662               *  t_i_b(ji,1))/zdiagbis(ji,nlay_s+1) &
663               *        MAX(0.0,SIGN(1.0,ht_s_b(ji))) 
664
665            ! surface temperature
666            isnow(ji)     = NINT(  1.0 - MAX( 0.0 , SIGN( 1.0 , -ht_s_b(ji) )  )  )
667            ztsuoldit(ji) = t_su_b(ji)
668            IF( t_su_b(ji) < ztfs(ji) ) &
669               t_su_b(ji) = ( zindtbis(ji,numeqmin(ji)) - ztrid(ji,numeqmin(ji),3)* ( REAL( isnow(ji) )*t_s_b(ji,1)   &
670               &          + REAL( 1 - isnow(ji) )*t_i_b(ji,1) ) ) / zdiagbis(ji,numeqmin(ji)) 
671         END DO
672         !
673         !--------------------------------------------------------------------------
674         !  10) Has the scheme converged ?, end of the iterative procedure         |
675         !--------------------------------------------------------------------------
676         !
677         ! check that nowhere it has started to melt
678         ! zerrit(ji) is a measure of error, it has to be under maxer_i_thd
679         DO ji = kideb , kiut
680            t_su_b(ji) =  MAX(  MIN( t_su_b(ji) , ztfs(ji) ) , 190._wp  )
681            zerrit(ji) =  ABS( t_su_b(ji) - ztsuoldit(ji) )     
682         END DO
683
684         DO layer  =  1, nlay_s
685            DO ji = kideb , kiut
686               t_s_b(ji,layer) = MAX(  MIN( t_s_b(ji,layer), rtt ), 190._wp  )
687               zerrit(ji)      = MAX(zerrit(ji),ABS(t_s_b(ji,layer) - ztstemp(ji,layer)))
688            END DO
689         END DO
690
691         DO layer  =  1, nlay_i
692            DO ji = kideb , kiut
693               ztmelt_i        = -tmut * s_i_b(ji,layer) + rtt 
694               t_i_b(ji,layer) =  MAX(MIN(t_i_b(ji,layer),ztmelt_i), 190._wp)
695               zerrit(ji)      =  MAX(zerrit(ji),ABS(t_i_b(ji,layer) - ztitemp(ji,layer)))
696            END DO
697         END DO
698
699         ! Compute spatial maximum over all errors
700         ! note that this could be optimized substantially by iterating only the non-converging points
701         zerritmax = 0._wp
702         DO ji = kideb, kiut
703            zerritmax = MAX( zerritmax, zerrit(ji) )   
704         END DO
705         IF( lk_mpp ) CALL mpp_max( zerritmax, kcom=ncomm_ice )
706
707      END DO  ! End of the do while iterative procedure
708
709      IF( ln_nicep .AND. lwp ) THEN
710         WRITE(numout,*) ' zerritmax : ', zerritmax
711         WRITE(numout,*) ' nconv     : ', nconv
712      ENDIF
713
714      !
715      !-------------------------------------------------------------------------!
716      !   11) Fluxes at the interfaces                                          !
717      !-------------------------------------------------------------------------!
718      DO ji = kideb, kiut
719         ! forced mode only : update of latent heat fluxes (sublimation) (always >=0, upward flux)
720         IF( .NOT. lk_cpl) qla_ice_1d (ji) = MAX( 0._wp, qla_ice_1d (ji) + dqla_ice_1d(ji) * ( t_su_b(ji) - ztsuold(ji) ) )
721         !                                ! surface ice conduction flux
722         isnow(ji)       = NINT(  1._wp - MAX( 0._wp, SIGN( 1._wp, -ht_s_b(ji) ) )  )
723         fc_su(ji)       =  -     REAL( isnow(ji) ) * zkappa_s(ji,0) * zg1s * (t_s_b(ji,1) - t_su_b(ji))   &
724            &               - REAL( 1 - isnow(ji) ) * zkappa_i(ji,0) * zg1  * (t_i_b(ji,1) - t_su_b(ji))
725         !                                ! bottom ice conduction flux
726         fc_bo_i(ji)     =  - zkappa_i(ji,nlay_i) * ( zg1*(t_bo_b(ji) - t_i_b(ji,nlay_i)) )
727      END DO
728
729      !-----------------------------------------
730      ! Heat flux used to warm/cool ice in W.m-2
731      !-----------------------------------------
732      DO ji = kideb, kiut
733         IF( t_su_b(ji) < rtt ) THEN  ! case T_su < 0degC
734            hfx_dif_1d(ji) = hfx_dif_1d(ji) + ( qns_ice_1d(ji) + qsr_ice_1d(ji) - zradtr_i(ji,nlay_i) - fc_bo_i(ji) ) * a_i_b(ji)
735         ELSE                         ! case T_su = 0degC
736            hfx_dif_1d(ji) = hfx_dif_1d(ji) + ( fc_su(ji) + i0(ji) * qsr_ice_1d(ji) - zradtr_i(ji,nlay_i) - fc_bo_i(ji) ) * a_i_b(ji)
737         ENDIF
738      END DO
739
740      ! --- computes sea ice energy of melting compulsory for limthd_dh --- !
741      CALL lim_thd_enmelt( kideb, kiut )
742
743      ! --- diag conservation imbalance on heat diffusion - PART 2 --- !
744      DO ji = kideb, kiut
745         zdq(ji)        = - zq_ini(ji) + ( SUM( q_i_b(ji,1:nlay_i) ) * ht_i_b(ji) / REAL( nlay_i ) +  &
746            &                              SUM( q_s_b(ji,1:nlay_s) ) * ht_s_b(ji) / REAL( nlay_s ) )
747         zhfx_err(ji)   = ( fc_su(ji) + i0(ji) * qsr_ice_1d(ji) - zradtr_i(ji,nlay_i) - fc_bo_i(ji) + zdq(ji) * r1_rdtice ) 
748         hfx_err_1d(ji) = hfx_err_1d(ji) + zhfx_err(ji) * a_i_b(ji)
749      END DO 
750
751      ! diagnose external surface (forced case) or bottom (forced case) from heat conservation
752      IF( .NOT. lk_cpl ) THEN   ! --- forced case: qns_ice and fc_su are diagnosed
753         !
754         DO ji = kideb, kiut
755            qns_ice_1d(ji) = qns_ice_1d(ji) - zhfx_err(ji)
756            fc_su     (ji) = fc_su(ji)      - zhfx_err(ji)
757         END DO
758         !
759      ELSE                      ! --- coupled case: ocean turbulent heat flux is diagnosed
760         !
761         DO ji = kideb, kiut
762            fhtur_1d  (ji) = fhtur_1d(ji)   - zhfx_err(ji)
763         END DO
764         !
765      ENDIF
766
767      ! --- compute diagnostic net heat flux at the surface of the snow-ice system (W.m2)
768      DO ji = kideb, kiut
769         ii = MOD( npb(ji) - 1, jpi ) + 1 ; ij = ( npb(ji) - 1 ) / jpi + 1
770         hfx_in (ii,ij) = hfx_in (ii,ij) + a_i_b(ji) * ( qsr_ice_1d(ji) + qns_ice_1d(ji) )
771      END DO
772   
773      !
774      CALL wrk_dealloc( jpij, numeqmin, numeqmax, isnow )
775      CALL wrk_dealloc( jpij, ztfs, ztsuold, ztsuoldit, zh_i, zh_s, zfsw )
776      CALL wrk_dealloc( jpij, zf, dzf, zerrit, zdifcase, zftrice, zihic, zhsu )
777      CALL wrk_dealloc( jpij, nlay_i+1, ztcond_i, zradtr_i, zradab_i, zkappa_i, ztiold, zeta_i, ztitemp, z_i, zspeche_i, kjstart = 0 )
778      CALL wrk_dealloc( jpij, nlay_s+1,           zradtr_s, zradab_s, zkappa_s, ztsold, zeta_s, ztstemp, z_s, kjstart = 0 )
779      CALL wrk_dealloc( jpij, jkmax+2, zindterm, zindtbis, zdiagbis )
780      CALL wrk_dealloc( jpij, jkmax+2, 3, ztrid )
781      CALL wrk_dealloc( jpij, zdq, zq_ini, zhfx_err )
782
783   END SUBROUTINE lim_thd_dif
784
785   SUBROUTINE lim_thd_enmelt( kideb, kiut )
786      !!-----------------------------------------------------------------------
787      !!                   ***  ROUTINE lim_thd_enmelt ***
788      !!                 
789      !! ** Purpose :   Computes sea ice energy of melting q_i (J.m-3) from temperature
790      !!
791      !! ** Method  :   Formula (Bitz and Lipscomb, 1999)
792      !!-------------------------------------------------------------------
793      INTEGER, INTENT(in) ::   kideb, kiut   ! bounds for the spatial loop
794      !
795      INTEGER  ::   ji, jk   ! dummy loop indices
796      REAL(wp) ::   ztmelts, zindb  ! local scalar
797      !!-------------------------------------------------------------------
798      !
799      DO jk = 1, nlay_i             ! Sea ice energy of melting
800         DO ji = kideb, kiut
801            ztmelts      = - tmut  * s_i_b(ji,jk) + rtt 
802            zindb        = MAX( 0._wp , SIGN( 1._wp , -(t_i_b(ji,jk) - rtt) - epsi10 ) )
803            q_i_b(ji,jk) = rhoic * ( cpic * ( ztmelts - t_i_b(ji,jk) )                                             &
804               &                   + lfus * ( 1.0 - zindb * ( ztmelts-rtt ) / MIN( t_i_b(ji,jk)-rtt, -epsi10 ) )   &
805               &                   - rcp  *                 ( ztmelts-rtt )  ) 
806         END DO
807      END DO
808      DO jk = 1, nlay_s             ! Snow energy of melting
809         DO ji = kideb, kiut
810            q_s_b(ji,jk) = rhosn * ( cpic * ( rtt - t_s_b(ji,jk) ) + lfus )
811         END DO
812      END DO
813      !
814   END SUBROUTINE lim_thd_enmelt
815
816#else
817   !!----------------------------------------------------------------------
818   !!                   Dummy Module                 No LIM-3 sea-ice model
819   !!----------------------------------------------------------------------
820CONTAINS
821   SUBROUTINE lim_thd_dif          ! Empty routine
822   END SUBROUTINE lim_thd_dif
823#endif
824   !!======================================================================
825END MODULE limthd_dif
Note: See TracBrowser for help on using the repository browser.