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.
iceitd.F90 in NEMO/trunk/src/ICE – NEMO

source: NEMO/trunk/src/ICE/iceitd.F90

Last change on this file was 15046, checked in by clem, 3 years ago

trunk: bug found by ecmwf: do not do ice remapping if, by some mystery and in very rare occasion, hR=hL=hice.

  • Property svn:keywords set to Id
File size: 38.0 KB
Line 
1MODULE iceitd
2   !!======================================================================
3   !!                       ***  MODULE iceitd ***
4   !!   sea-ice : ice thickness distribution
5   !!======================================================================
6   !! History :  3.0  !  2005-12  (M. Vancoppenolle) original code (based on CICE)
7   !!            4.0  !  2018     (many people)      SI3 [aka Sea Ice cube]
8   !!----------------------------------------------------------------------
9#if defined key_si3
10   !!----------------------------------------------------------------------
11   !!   'key_si3'                                       SI3 sea-ice model
12   !!----------------------------------------------------------------------
13   !!   ice_itd_rem   : redistribute ice thicknesses after thermo growth and melt
14   !!   itd_glinear   : build g(h) satisfying area and volume constraints
15   !!   itd_shiftice  : shift ice across category boundaries, conserving everything
16   !!   ice_itd_reb   : rebin ice thicknesses into bounded categories
17   !!   ice_itd_init  : read ice thicknesses mean and min from namelist
18   !!----------------------------------------------------------------------
19   USE dom_oce        ! ocean domain
20   USE phycst         ! physical constants
21   USE ice1D          ! sea-ice: thermodynamic variables
22   USE ice            ! sea-ice: variables
23   USE icevar         ! sea-ice: operations
24   USE icectl         ! sea-ice: conservation tests
25   USE icetab         ! sea-ice: convert 1D<=>2D
26   !
27   USE in_out_manager ! I/O manager
28   USE lib_mpp        ! MPP library
29   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
30   USE prtctl         ! Print control
31   USE timing         ! Timing
32
33   IMPLICIT NONE
34   PRIVATE
35
36   PUBLIC   ice_itd_init  ! called in icestp
37   PUBLIC   ice_itd_rem   ! called in icethd
38   PUBLIC   ice_itd_reb   ! called in icecor
39
40   INTEGER            ::   nice_catbnd     ! choice of the type of ice category function
41   !                                       ! associated indices:
42   INTEGER, PARAMETER ::   np_cathfn = 1   ! categories defined by a function
43   INTEGER, PARAMETER ::   np_catusr = 2   ! categories defined by the user
44   !
45   !                                           !! ** namelist (namitd) **
46   LOGICAL                    ::   ln_cat_hfn   ! ice categories are defined by function like rn_himean**(-0.05)
47   REAL(wp)                   ::   rn_himean    ! mean thickness of the domain
48   LOGICAL                    ::   ln_cat_usr   ! ice categories are defined by rn_catbnd
49   REAL(wp), DIMENSION(0:100) ::   rn_catbnd    ! ice categories bounds
50   REAL(wp)                   ::   rn_himax     ! maximum ice thickness allowed
51   !
52   !! * Substitutions
53#  include "do_loop_substitute.h90"
54   !!----------------------------------------------------------------------
55   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
56   !! $Id$
57   !! Software governed by the CeCILL license (see ./LICENSE)
58   !!----------------------------------------------------------------------
59CONTAINS
60
61   SUBROUTINE ice_itd_rem( kt )
62      !!------------------------------------------------------------------
63      !!                ***  ROUTINE ice_itd_rem ***
64      !!
65      !! ** Purpose :   computes the redistribution of ice thickness
66      !!                after thermodynamic growth of ice thickness
67      !!
68      !! ** Method  :   Linear remapping
69      !!
70      !! References :   W.H. Lipscomb, JGR 2001
71      !!------------------------------------------------------------------
72      INTEGER , INTENT (in) ::   kt      ! Ocean time step
73      !
74      INTEGER  ::   ji, jj, jl, jcat     ! dummy loop index
75      INTEGER  ::   ipti                 ! local integer
76      REAL(wp) ::   zx1, zwk1, zdh0, zetamin, zdamax   ! local scalars
77      REAL(wp) ::   zx2, zwk2, zda0, zetamax           !   -      -
78      REAL(wp) ::   zx3
79      REAL(wp) ::   zslope          ! used to compute local thermodynamic "speeds"
80      !
81      INTEGER , DIMENSION(jpij)       ::   iptidx          ! compute remapping or not
82      INTEGER , DIMENSION(jpij,jpl-1) ::   jdonor          ! donor category index
83      REAL(wp), DIMENSION(jpij,jpl)   ::   zdhice          ! ice thickness increment
84      REAL(wp), DIMENSION(jpij,jpl)   ::   g0, g1          ! coefficients for fitting the line of the ITD
85      REAL(wp), DIMENSION(jpij,jpl)   ::   hL, hR          ! left and right boundary for the ITD for each thickness
86      REAL(wp), DIMENSION(jpij,jpl-1) ::   zdaice, zdvice  ! local increment of ice area and volume
87      REAL(wp), DIMENSION(jpij)       ::   zhb0, zhb1      ! category boundaries for thinnes categories
88      REAL(wp), DIMENSION(jpij,0:jpl) ::   zhbnew          ! new boundaries of ice categories
89      !!------------------------------------------------------------------
90      IF( ln_timing )   CALL timing_start('iceitd_rem')
91
92      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_itd_rem: remapping ice thickness distribution'
93
94      IF( ln_icediachk )   CALL ice_cons_hsm(0, 'iceitd_rem', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft)
95      IF( ln_icediachk )   CALL ice_cons2D  (0, 'iceitd_rem',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft)
96
97      !-----------------------------------------------------------------------------------------------
98      !  1) Identify grid cells with ice
99      !-----------------------------------------------------------------------------------------------
100      at_i(:,:) = SUM( a_i, dim=3 )
101      !
102      npti = 0   ;   nptidx(:) = 0
103      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
104         IF ( at_i(ji,jj) > epsi10 ) THEN
105            npti = npti + 1
106            nptidx( npti ) = (jj - 1) * jpi + ji
107         ENDIF
108      END_2D
109
110      !-----------------------------------------------------------------------------------------------
111      !  2) Compute new category boundaries
112      !-----------------------------------------------------------------------------------------------
113      IF( npti > 0 ) THEN
114         !
115         zdhice(:,:) = 0._wp
116         zhbnew(:,:) = 0._wp
117         !
118         CALL tab_3d_2d( npti, nptidx(1:npti), h_i_2d (1:npti,1:jpl), h_i   )
119         CALL tab_3d_2d( npti, nptidx(1:npti), h_ib_2d(1:npti,1:jpl), h_i_b )
120         CALL tab_3d_2d( npti, nptidx(1:npti), a_i_2d (1:npti,1:jpl), a_i   )
121         CALL tab_3d_2d( npti, nptidx(1:npti), a_ib_2d(1:npti,1:jpl), a_i_b )
122         !
123         DO jl = 1, jpl
124            ! Compute thickness change in each ice category
125            DO ji = 1, npti
126               IF( a_i_2d(ji,jl) > epsi10 )   zdhice(ji,jl) = h_i_2d(ji,jl) - h_ib_2d(ji,jl)
127            END DO
128         END DO
129         !
130         ! --- New boundaries for category 1:jpl-1 --- !
131         DO jl = 1, jpl - 1
132            !
133            DO ji = 1, npti
134               !
135               ! --- New boundary: Hn* = Hn + Fn*dt --- !
136               !     Fn*dt = ( fn + (fn+1 - fn)/(hn+1 - hn) * (Hn - hn) ) * dt = zdhice + zslope * (Hmax - h_i_b)
137               !
138               IF    ( a_ib_2d(ji,jl) >  epsi10 .AND. a_ib_2d(ji,jl+1) >  epsi10 ) THEN   ! a(jl+1) & a(jl) /= 0
139                  zslope        = ( zdhice(ji,jl+1) - zdhice(ji,jl) ) / ( h_ib_2d(ji,jl+1) - h_ib_2d(ji,jl) )
140                  zhbnew(ji,jl) = hi_max(jl) + zdhice(ji,jl) + zslope * ( hi_max(jl) - h_ib_2d(ji,jl) )
141               ELSEIF( a_ib_2d(ji,jl) >  epsi10 .AND. a_ib_2d(ji,jl+1) <= epsi10 ) THEN   ! a(jl+1)=0 => Hn* = Hn + fn*dt
142                  zhbnew(ji,jl) = hi_max(jl) + zdhice(ji,jl)
143               ELSEIF( a_ib_2d(ji,jl) <= epsi10 .AND. a_ib_2d(ji,jl+1) >  epsi10 ) THEN   ! a(jl)=0 => Hn* = Hn + fn+1*dt
144                  zhbnew(ji,jl) = hi_max(jl) + zdhice(ji,jl+1)
145               ELSE                                                                       ! a(jl+1) & a(jl) = 0
146                  zhbnew(ji,jl) = hi_max(jl)
147               ENDIF
148               !
149               ! --- 2 conditions for remapping --- !
150               ! 1) hn(t+1)+espi < Hn* < hn+1(t+1)-epsi
151               !    Note: hn(t+1) must not be too close to either HR or HL otherwise a division by nearly 0 is possible
152               !          in itd_glinear in the case (HR-HL) = 3(Hice - HL) or = 3(HR - Hice)
153# if defined key_single
154               IF( a_i_2d(ji,jl  ) > epsi10 .AND. h_i_2d(ji,jl  ) > ( zhbnew(ji,jl) - epsi06 ) )   nptidx(ji) = 0
155               IF( a_i_2d(ji,jl+1) > epsi10 .AND. h_i_2d(ji,jl+1) < ( zhbnew(ji,jl) + epsi06 ) )   nptidx(ji) = 0
156# else
157               IF( a_i_2d(ji,jl  ) > epsi10 .AND. h_i_2d(ji,jl  ) > ( zhbnew(ji,jl) - epsi10 ) )   nptidx(ji) = 0
158               IF( a_i_2d(ji,jl+1) > epsi10 .AND. h_i_2d(ji,jl+1) < ( zhbnew(ji,jl) + epsi10 ) )   nptidx(ji) = 0
159# endif
160               !
161               ! 2) Hn-1 < Hn* < Hn+1
162               IF( zhbnew(ji,jl) < hi_max(jl-1) )   nptidx(ji) = 0
163               IF( zhbnew(ji,jl) > hi_max(jl+1) )   nptidx(ji) = 0
164               !
165            END DO
166         END DO
167         !
168         ! --- New boundaries for category jpl --- !
169         DO ji = 1, npti
170            IF( a_i_2d(ji,jpl) > epsi10 ) THEN
171               zhbnew(ji,jpl) = MAX( hi_max(jpl-1), 3._wp * h_i_2d(ji,jpl) - 2._wp * zhbnew(ji,jpl-1) )
172            ELSE
173               zhbnew(ji,jpl) = hi_max(jpl)
174            ENDIF
175            !
176            ! --- 1 additional condition for remapping (1st category) --- !
177            ! H0+epsi < h1(t) < H1-epsi
178            !    h1(t) must not be too close to either HR or HL otherwise a division by nearly 0 is possible
179            !    in itd_glinear in the case (HR-HL) = 3(Hice - HL) or = 3(HR - Hice)
180# if defined key_single
181            IF( h_ib_2d(ji,1) < ( hi_max(0) + epsi06 ) )   nptidx(ji) = 0
182            IF( h_ib_2d(ji,1) > ( hi_max(1) - epsi06 ) )   nptidx(ji) = 0
183# else
184            IF( h_ib_2d(ji,1) < ( hi_max(0) + epsi10 ) )   nptidx(ji) = 0
185            IF( h_ib_2d(ji,1) > ( hi_max(1) - epsi10 ) )   nptidx(ji) = 0
186# endif
187         END DO
188         !
189         !-----------------------------------------------------------------------------------------------
190         !  3) Identify cells where remapping
191         !-----------------------------------------------------------------------------------------------
192         ipti = 0   ;   iptidx(:) = 0
193         DO ji = 1, npti
194            IF( nptidx(ji) /= 0 ) THEN
195               ipti = ipti + 1
196               iptidx(ipti)   = nptidx(ji)
197               zhbnew(ipti,:) = zhbnew(ji,:)  ! adjust zhbnew to new indices
198            ENDIF
199         END DO
200         nptidx(:) = iptidx(:)
201         npti      = ipti
202         !
203      ENDIF
204
205      !-----------------------------------------------------------------------------------------------
206      !  4) Compute g(h)
207      !-----------------------------------------------------------------------------------------------
208      IF( npti > 0 ) THEN
209         !
210         zhb0(:) = hi_max(0)   ;   zhb1(:) = hi_max(1)
211         g0(:,:) = 0._wp       ;   g1(:,:) = 0._wp
212         hL(:,:) = 0._wp       ;   hR(:,:) = 0._wp
213         !
214         DO jl = 1, jpl
215            !
216            CALL tab_2d_1d( npti, nptidx(1:npti), h_ib_1d(1:npti), h_i_b(:,:,jl) )
217            CALL tab_2d_1d( npti, nptidx(1:npti), h_i_1d (1:npti), h_i  (:,:,jl) )
218            CALL tab_2d_1d( npti, nptidx(1:npti), a_i_1d (1:npti), a_i  (:,:,jl) )
219            CALL tab_2d_1d( npti, nptidx(1:npti), v_i_1d (1:npti), v_i  (:,:,jl) )
220            !
221            IF( jl == 1 ) THEN
222               !
223               ! --- g(h) for category 1 --- !
224               CALL itd_glinear( zhb0(1:npti)  , zhb1(1:npti)  , h_ib_1d(1:npti)  , a_i_1d(1:npti)  ,  &   ! in
225                  &              g0  (1:npti,1), g1  (1:npti,1), hL     (1:npti,1), hR    (1:npti,1)   )   ! out
226               !
227               ! Area lost due to melting of thin ice
228               DO ji = 1, npti
229                  !
230                  IF( a_i_1d(ji) > epsi10 ) THEN
231                     !
232                     zdh0 =  h_i_1d(ji) - h_ib_1d(ji)
233                     IF( zdh0 < 0.0 ) THEN      ! remove area from category 1
234                        zdh0 = MIN( -zdh0, hi_max(1) )
235                        !Integrate g(1) from 0 to dh0 to estimate area melted
236                        zetamax = MIN( zdh0, hR(ji,1) ) - hL(ji,1)
237                        !
238                        IF( zetamax > 0.0 ) THEN
239                           zx1    = zetamax
240                           zx2    = 0.5 * zetamax * zetamax
241                           zda0   = g1(ji,1) * zx2 + g0(ji,1) * zx1                ! ice area removed
242                           zdamax = a_i_1d(ji) * (1.0 - h_i_1d(ji) / h_ib_1d(ji) ) ! Constrain new thickness <= h_i
243                           zda0   = MIN( zda0, zdamax )                            ! ice area lost due to melting of thin ice (zdamax > 0)
244                           ! Remove area, conserving volume
245                           h_i_1d(ji) = h_i_1d(ji) * a_i_1d(ji) / ( a_i_1d(ji) - zda0 )
246                           a_i_1d(ji) = a_i_1d(ji) - zda0
247                           v_i_1d(ji) = a_i_1d(ji) * h_i_1d(ji) ! useless ?
248                        ENDIF
249                        !
250                     ELSE ! if ice accretion zdh0 > 0
251                        ! zhbnew was 0, and is shifted to the right to account for thin ice growth in openwater (F0 = f1)
252                        zhbnew(ji,0) = MIN( zdh0, hi_max(1) )
253                     ENDIF
254                     !
255                  ENDIF
256                  !
257               END DO
258               !
259               CALL tab_1d_2d( npti, nptidx(1:npti), h_i_1d(1:npti), h_i(:,:,jl) )
260               CALL tab_1d_2d( npti, nptidx(1:npti), a_i_1d(1:npti), a_i(:,:,jl) )
261               CALL tab_1d_2d( npti, nptidx(1:npti), v_i_1d(1:npti), v_i(:,:,jl) )
262               !
263            ENDIF ! jl=1
264            !
265            ! --- g(h) for each thickness category --- !
266            CALL itd_glinear( zhbnew(1:npti,jl-1), zhbnew(1:npti,jl), h_i_1d(1:npti)   , a_i_1d(1:npti)   ,  &   ! in
267               &              g0    (1:npti,jl  ), g1    (1:npti,jl), hL    (1:npti,jl), hR    (1:npti,jl)   )   ! out
268            !
269         END DO
270
271         !-----------------------------------------------------------------------------------------------
272         !  5) Compute area and volume to be shifted across each boundary (Eq. 18)
273         !-----------------------------------------------------------------------------------------------
274         DO jl = 1, jpl - 1
275            !
276            DO ji = 1, npti
277               !
278               ! left and right integration limits in eta space
279               IF (zhbnew(ji,jl) > hi_max(jl)) THEN ! Hn* > Hn => transfer from jl to jl+1
280                  zetamin = MAX( hi_max(jl)   , hL(ji,jl) ) - hL(ji,jl)   ! hi_max(jl) - hL
281                  zetamax = MIN( zhbnew(ji,jl), hR(ji,jl) ) - hL(ji,jl)   ! hR - hL
282                  jdonor(ji,jl) = jl
283               ELSE                                 ! Hn* <= Hn => transfer from jl+1 to jl
284                  zetamin = 0.0
285                  zetamax = MIN( hi_max(jl), hR(ji,jl+1) ) - hL(ji,jl+1)  ! hi_max(jl) - hL
286                  jdonor(ji,jl) = jl + 1
287               ENDIF
288               zetamax = MAX( zetamax, zetamin ) ! no transfer if etamax < etamin
289               !
290               zx1  = zetamax - zetamin
291               zwk1 = zetamin * zetamin
292               zwk2 = zetamax * zetamax
293               zx2  = 0.5 * ( zwk2 - zwk1 )
294               zwk1 = zwk1 * zetamin
295               zwk2 = zwk2 * zetamax
296               zx3  = 1.0 / 3.0 * ( zwk2 - zwk1 )
297               jcat = jdonor(ji,jl)
298               zdaice(ji,jl) = g1(ji,jcat)*zx2 + g0(ji,jcat)*zx1
299               zdvice(ji,jl) = g1(ji,jcat)*zx3 + g0(ji,jcat)*zx2 + zdaice(ji,jl)*hL(ji,jcat)
300               !
301            END DO
302         END DO
303
304         !----------------------------------------------------------------------------------------------
305         ! 6) Shift ice between categories
306         !----------------------------------------------------------------------------------------------
307         CALL itd_shiftice ( jdonor(1:npti,:), zdaice(1:npti,:), zdvice(1:npti,:) )
308
309         !----------------------------------------------------------------------------------------------
310         ! 7) Make sure h_i >= minimum ice thickness hi_min
311         !----------------------------------------------------------------------------------------------
312         CALL tab_2d_1d( npti, nptidx(1:npti), h_i_1d (1:npti), h_i (:,:,1) )
313         CALL tab_2d_1d( npti, nptidx(1:npti), a_i_1d (1:npti), a_i (:,:,1) )
314         CALL tab_2d_1d( npti, nptidx(1:npti), a_ip_1d(1:npti), a_ip(:,:,1) )
315         !
316         DO ji = 1, npti
317            IF ( a_i_1d(ji) > epsi10 .AND. h_i_1d(ji) < rn_himin ) THEN
318               a_i_1d(ji) = a_i_1d(ji) * h_i_1d(ji) / rn_himin
319               IF( ln_pnd_LEV .OR. ln_pnd_TOPO )   a_ip_1d(ji) = a_ip_1d(ji) * h_i_1d(ji) / rn_himin
320               h_i_1d(ji) = rn_himin
321            ENDIF
322         END DO
323         !
324         CALL tab_1d_2d( npti, nptidx(1:npti), h_i_1d (1:npti), h_i (:,:,1) )
325         CALL tab_1d_2d( npti, nptidx(1:npti), a_i_1d (1:npti), a_i (:,:,1) )
326         CALL tab_1d_2d( npti, nptidx(1:npti), a_ip_1d(1:npti), a_ip(:,:,1) )
327         !
328      ENDIF
329      !
330      IF( ln_icediachk )   CALL ice_cons_hsm(1, 'iceitd_rem', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft)
331      IF( ln_icediachk )   CALL ice_cons2D  (1, 'iceitd_rem',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft)
332      IF( ln_timing    )   CALL timing_stop ('iceitd_rem')
333      !
334   END SUBROUTINE ice_itd_rem
335
336
337   SUBROUTINE itd_glinear( HbL, Hbr, phice, paice, pg0, pg1, phL, phR )
338      !!------------------------------------------------------------------
339      !!                ***  ROUTINE itd_glinear ***
340      !!
341      !! ** Purpose :   build g(h) satisfying area and volume constraints (Eq. 6 and 9)
342      !!
343      !! ** Method  :   g(h) is linear and written as: g(eta) = g1(eta) + g0
344      !!                with eta = h - HL
345      !!------------------------------------------------------------------
346      REAL(wp), DIMENSION(:), INTENT(in   ) ::   HbL, HbR      ! left and right category boundaries
347      REAL(wp), DIMENSION(:), INTENT(in   ) ::   phice, paice  ! ice thickness and concentration
348      REAL(wp), DIMENSION(:), INTENT(inout) ::   pg0, pg1      ! coefficients in linear equation for g(eta)
349      REAL(wp), DIMENSION(:), INTENT(inout) ::   phL, phR      ! min and max value of range over which g(h) > 0
350      !
351      INTEGER  ::   ji           ! horizontal indices
352      REAL(wp) ::   z1_3 , z2_3  ! 1/3 , 2/3
353      REAL(wp) ::   zh13         ! HbL + 1/3 * (HbR - HbL)
354      REAL(wp) ::   zh23         ! HbL + 2/3 * (HbR - HbL)
355      REAL(wp) ::   zdhr         ! 1 / (hR - hL)
356      REAL(wp) ::   zwk1, zwk2   ! temporary variables
357      !!------------------------------------------------------------------
358      !
359      z1_3 = 1._wp / 3._wp
360      z2_3 = 2._wp / 3._wp
361      !
362      DO ji = 1, npti
363         !
364         IF( paice(ji) > epsi10  .AND. phice(ji) > epsi10 )  THEN
365            !
366            ! Initialize hL and hR
367            phL(ji) = HbL(ji)
368            phR(ji) = HbR(ji)
369            !
370            ! Change hL or hR if hice falls outside central third of range,
371            ! so that hice is in the central third of the range [HL HR]
372            zh13 = z1_3 * ( 2._wp * phL(ji) +         phR(ji) )
373            zh23 = z1_3 * (         phL(ji) + 2._wp * phR(ji) )
374            !
375            IF    ( phice(ji) < zh13 ) THEN   ;   phR(ji) = 3._wp * phice(ji) - 2._wp * phL(ji) ! move HR to the left
376            ELSEIF( phice(ji) > zh23 ) THEN   ;   phL(ji) = 3._wp * phice(ji) - 2._wp * phR(ji) ! move HL to the right
377            ENDIF
378            !
379            ! Compute coefficients of g(eta) = g0 + g1*eta
380            IF( phR(ji) > phL(ji) ) THEN   ;   zdhr = 1._wp / (phR(ji) - phL(ji))
381            ELSE                           ;   zdhr = 0._wp ! if hR=hL=hice => no remapping
382            ENDIF
383            !!zdhr = 1._wp / (phR(ji) - phL(ji))
384            zwk1 = 6._wp * paice(ji) * zdhr
385            zwk2 = ( phice(ji) - phL(ji) ) * zdhr
386            pg0(ji) = zwk1 * ( z2_3 - zwk2 )                    ! Eq. 14
387            pg1(ji) = 2._wp * zdhr * zwk1 * ( zwk2 - 0.5_wp )   ! Eq. 14
388            !
389         ELSE  ! remap_flag = .false. or a_i < epsi10
390            phL(ji) = 0._wp
391            phR(ji) = 0._wp
392            pg0(ji) = 0._wp
393            pg1(ji) = 0._wp
394         ENDIF
395         !
396      END DO
397      !
398   END SUBROUTINE itd_glinear
399
400
401   SUBROUTINE itd_shiftice( kdonor, pdaice, pdvice )
402      !!------------------------------------------------------------------
403      !!                ***  ROUTINE itd_shiftice ***
404      !!
405      !! ** Purpose :   shift ice across category boundaries, conserving everything
406      !!              ( area, volume, energy, age*vol, and mass of salt )
407      !!------------------------------------------------------------------
408      INTEGER , DIMENSION(:,:), INTENT(in) ::   kdonor   ! donor category index
409      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pdaice   ! ice area transferred across boundary
410      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pdvice   ! ice volume transferred across boundary
411      !
412      INTEGER  ::   ji, jl, jk         ! dummy loop indices
413      INTEGER  ::   jl2, jl1           ! local integers
414      REAL(wp) ::   ztrans             ! ice/snow transferred
415      REAL(wp), DIMENSION(jpij)            ::   zworka, zworkv   ! workspace
416      REAL(wp), DIMENSION(jpij,jpl)        ::   zaTsfn           !  -    -
417      REAL(wp), DIMENSION(jpij,nlay_i,jpl) ::   ze_i_2d
418      REAL(wp), DIMENSION(jpij,nlay_s,jpl) ::   ze_s_2d
419      !!------------------------------------------------------------------
420
421      CALL tab_3d_2d( npti, nptidx(1:npti), h_i_2d (1:npti,1:jpl), h_i  )
422      CALL tab_3d_2d( npti, nptidx(1:npti), a_i_2d (1:npti,1:jpl), a_i  )
423      CALL tab_3d_2d( npti, nptidx(1:npti), v_i_2d (1:npti,1:jpl), v_i  )
424      CALL tab_3d_2d( npti, nptidx(1:npti), v_s_2d (1:npti,1:jpl), v_s  )
425      CALL tab_3d_2d( npti, nptidx(1:npti), oa_i_2d(1:npti,1:jpl), oa_i )
426      CALL tab_3d_2d( npti, nptidx(1:npti), sv_i_2d(1:npti,1:jpl), sv_i )
427      CALL tab_3d_2d( npti, nptidx(1:npti), a_ip_2d(1:npti,1:jpl), a_ip )
428      CALL tab_3d_2d( npti, nptidx(1:npti), v_ip_2d(1:npti,1:jpl), v_ip )
429      CALL tab_3d_2d( npti, nptidx(1:npti), v_il_2d(1:npti,1:jpl), v_il )
430      CALL tab_3d_2d( npti, nptidx(1:npti), t_su_2d(1:npti,1:jpl), t_su )
431      DO jl = 1, jpl
432         DO jk = 1, nlay_s
433            CALL tab_2d_1d( npti, nptidx(1:npti), ze_s_2d(1:npti,jk,jl), e_s(:,:,jk,jl) )
434         END DO
435         DO jk = 1, nlay_i
436            CALL tab_2d_1d( npti, nptidx(1:npti), ze_i_2d(1:npti,jk,jl), e_i(:,:,jk,jl) )
437         END DO
438      END DO
439      ! to correct roundoff errors on a_i
440      CALL tab_2d_1d( npti, nptidx(1:npti), rn_amax_1d(1:npti), rn_amax_2d )
441
442      !----------------------------------------------------------------------------------------------
443      ! 1) Define a variable equal to a_i*T_su
444      !----------------------------------------------------------------------------------------------
445      DO jl = 1, jpl
446         DO ji = 1, npti
447            zaTsfn(ji,jl) = a_i_2d(ji,jl) * t_su_2d(ji,jl)
448         END DO
449      END DO
450
451      !-------------------------------------------------------------------------------
452      ! 2) Transfer volume and energy between categories
453      !-------------------------------------------------------------------------------
454      DO jl = 1, jpl - 1
455         DO ji = 1, npti
456            !
457            jl1 = kdonor(ji,jl)
458            !
459            IF( jl1 > 0 ) THEN
460               !
461               IF ( jl1 == jl  ) THEN   ;   jl2 = jl1+1
462               ELSE                     ;   jl2 = jl
463               ENDIF
464               !
465               IF( v_i_2d(ji,jl1) >= epsi10 ) THEN   ;   zworkv(ji) = pdvice(ji,jl) / v_i_2d(ji,jl1)
466               ELSE                                  ;   zworkv(ji) = 0._wp
467               ENDIF
468               IF( a_i_2d(ji,jl1) >= epsi10 ) THEN   ;   zworka(ji) = pdaice(ji,jl) / a_i_2d(ji,jl1)
469               ELSE                                  ;   zworka(ji) = 0._wp
470               ENDIF
471               !
472               a_i_2d(ji,jl1) = a_i_2d(ji,jl1) - pdaice(ji,jl)       ! Ice areas
473               a_i_2d(ji,jl2) = a_i_2d(ji,jl2) + pdaice(ji,jl)
474               !
475               v_i_2d(ji,jl1) = v_i_2d(ji,jl1) - pdvice(ji,jl)       ! Ice volumes
476               v_i_2d(ji,jl2) = v_i_2d(ji,jl2) + pdvice(ji,jl)
477               !
478               ztrans         = v_s_2d(ji,jl1) * zworkv(ji)          ! Snow volumes
479               v_s_2d(ji,jl1) = v_s_2d(ji,jl1) - ztrans
480               v_s_2d(ji,jl2) = v_s_2d(ji,jl2) + ztrans
481               !
482               ztrans          = oa_i_2d(ji,jl1) * zworka(ji)        ! Ice age
483               oa_i_2d(ji,jl1) = oa_i_2d(ji,jl1) - ztrans
484               oa_i_2d(ji,jl2) = oa_i_2d(ji,jl2) + ztrans
485               !
486               ztrans          = sv_i_2d(ji,jl1) * zworkv(ji)        ! Ice salinity
487               sv_i_2d(ji,jl1) = sv_i_2d(ji,jl1) - ztrans
488               sv_i_2d(ji,jl2) = sv_i_2d(ji,jl2) + ztrans
489               !
490               ztrans          = zaTsfn(ji,jl1) * zworka(ji)         ! Surface temperature
491               zaTsfn(ji,jl1)  = zaTsfn(ji,jl1) - ztrans
492               zaTsfn(ji,jl2)  = zaTsfn(ji,jl2) + ztrans
493               !
494               IF ( ln_pnd_LEV .OR. ln_pnd_TOPO ) THEN
495                  ztrans          = a_ip_2d(ji,jl1) * zworka(ji)     ! Pond fraction
496                  a_ip_2d(ji,jl1) = a_ip_2d(ji,jl1) - ztrans
497                  a_ip_2d(ji,jl2) = a_ip_2d(ji,jl2) + ztrans
498                  !
499                  ztrans          = v_ip_2d(ji,jl1) * zworkv(ji)     ! Pond volume
500                  v_ip_2d(ji,jl1) = v_ip_2d(ji,jl1) - ztrans
501                  v_ip_2d(ji,jl2) = v_ip_2d(ji,jl2) + ztrans
502                  !
503                  IF ( ln_pnd_lids ) THEN                            ! Pond lid volume
504                     ztrans          = v_il_2d(ji,jl1) * zworkv(ji)
505                     v_il_2d(ji,jl1) = v_il_2d(ji,jl1) - ztrans
506                     v_il_2d(ji,jl2) = v_il_2d(ji,jl2) + ztrans
507                  ENDIF
508               ENDIF
509               !
510            ENDIF   ! jl1 >0
511         END DO
512         !
513         DO jk = 1, nlay_s         !--- Snow heat content
514            DO ji = 1, npti
515               !
516               jl1 = kdonor(ji,jl)
517               !
518               IF( jl1 > 0 ) THEN
519                  IF(jl1 == jl) THEN  ;  jl2 = jl+1
520                  ELSE                ;  jl2 = jl
521                  ENDIF
522                  ztrans             = ze_s_2d(ji,jk,jl1) * zworkv(ji)
523                  ze_s_2d(ji,jk,jl1) = ze_s_2d(ji,jk,jl1) - ztrans
524                  ze_s_2d(ji,jk,jl2) = ze_s_2d(ji,jk,jl2) + ztrans
525               ENDIF
526            END DO
527         END DO
528         !
529         DO jk = 1, nlay_i         !--- Ice heat content
530            DO ji = 1, npti
531               !
532               jl1 = kdonor(ji,jl)
533               !
534               IF( jl1 > 0 ) THEN
535                  IF(jl1 == jl) THEN  ;  jl2 = jl+1
536                  ELSE                ;  jl2 = jl
537                  ENDIF
538                  ztrans             = ze_i_2d(ji,jk,jl1) * zworkv(ji)
539                  ze_i_2d(ji,jk,jl1) = ze_i_2d(ji,jk,jl1) - ztrans
540                  ze_i_2d(ji,jk,jl2) = ze_i_2d(ji,jk,jl2) + ztrans
541               ENDIF
542            END DO
543         END DO
544         !
545      END DO                   ! boundaries, 1 to jpl-1
546
547      !-------------------
548      ! 3) roundoff errors
549      !-------------------
550      ! clem: The transfer between one category to another can lead to very small negative values (-1.e-20)
551      !       because of truncation error ( i.e. 1. - 1. /= 0 )
552      CALL ice_var_roundoff( a_i_2d, v_i_2d, v_s_2d, sv_i_2d, oa_i_2d, a_ip_2d, v_ip_2d, v_il_2d, ze_s_2d, ze_i_2d )
553
554      ! at_i must be <= rn_amax
555      zworka(1:npti) = SUM( a_i_2d(1:npti,:), dim=2 )
556      DO jl  = 1, jpl
557         WHERE( zworka(1:npti) > rn_amax_1d(1:npti) )   &
558            &   a_i_2d(1:npti,jl) = a_i_2d(1:npti,jl) * rn_amax_1d(1:npti) / zworka(1:npti)
559      END DO
560
561      !-------------------------------------------------------------------------------
562      ! 4) Update ice thickness and temperature
563      !-------------------------------------------------------------------------------
564# if defined key_single
565      WHERE( a_i_2d(1:npti,:) >= epsi06 )
566# else
567      WHERE( a_i_2d(1:npti,:) >= epsi20 )
568# endif
569         h_i_2d (1:npti,:)  =  v_i_2d(1:npti,:) / a_i_2d(1:npti,:)
570         t_su_2d(1:npti,:)  =  zaTsfn(1:npti,:) / a_i_2d(1:npti,:)
571      ELSEWHERE
572         h_i_2d (1:npti,:)  = 0._wp
573         t_su_2d(1:npti,:)  = rt0
574      END WHERE
575      !
576      CALL tab_2d_3d( npti, nptidx(1:npti), h_i_2d (1:npti,1:jpl), h_i  )
577      CALL tab_2d_3d( npti, nptidx(1:npti), a_i_2d (1:npti,1:jpl), a_i  )
578      CALL tab_2d_3d( npti, nptidx(1:npti), v_i_2d (1:npti,1:jpl), v_i  )
579      CALL tab_2d_3d( npti, nptidx(1:npti), v_s_2d (1:npti,1:jpl), v_s  )
580      CALL tab_2d_3d( npti, nptidx(1:npti), oa_i_2d(1:npti,1:jpl), oa_i )
581      CALL tab_2d_3d( npti, nptidx(1:npti), sv_i_2d(1:npti,1:jpl), sv_i )
582      CALL tab_2d_3d( npti, nptidx(1:npti), a_ip_2d(1:npti,1:jpl), a_ip )
583      CALL tab_2d_3d( npti, nptidx(1:npti), v_ip_2d(1:npti,1:jpl), v_ip )
584      CALL tab_2d_3d( npti, nptidx(1:npti), v_il_2d(1:npti,1:jpl), v_il )
585      CALL tab_2d_3d( npti, nptidx(1:npti), t_su_2d(1:npti,1:jpl), t_su )
586      DO jl = 1, jpl
587         DO jk = 1, nlay_s
588            CALL tab_1d_2d( npti, nptidx(1:npti), ze_s_2d(1:npti,jk,jl), e_s(:,:,jk,jl) )
589         END DO
590         DO jk = 1, nlay_i
591            CALL tab_1d_2d( npti, nptidx(1:npti), ze_i_2d(1:npti,jk,jl), e_i(:,:,jk,jl) )
592         END DO
593      END DO
594      !
595   END SUBROUTINE itd_shiftice
596
597
598   SUBROUTINE ice_itd_reb( kt )
599      !!------------------------------------------------------------------
600      !!                ***  ROUTINE ice_itd_reb ***
601      !!
602      !! ** Purpose : rebin - rebins thicknesses into defined categories
603      !!
604      !! ** Method  : If a category thickness is out of bounds, shift part (for down to top)
605      !!              or entire (for top to down) area, volume, and energy
606      !!              to the neighboring category
607      !!------------------------------------------------------------------
608      INTEGER , INTENT (in) ::   kt      ! Ocean time step
609      INTEGER ::   ji, jj, jl   ! dummy loop indices
610      !
611      INTEGER , DIMENSION(jpij,jpl-1) ::   jdonor           ! donor category index
612      REAL(wp), DIMENSION(jpij,jpl-1) ::   zdaice, zdvice   ! ice area and volume transferred
613      !!------------------------------------------------------------------
614      IF( ln_timing )   CALL timing_start('iceitd_reb')
615      !
616      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_itd_reb: rebining ice thickness distribution'
617      !
618      IF( ln_icediachk )   CALL ice_cons_hsm(0, 'iceitd_reb', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft)
619      IF( ln_icediachk )   CALL ice_cons2D  (0, 'iceitd_reb',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft)
620      !
621      jdonor(:,:) = 0
622      zdaice(:,:) = 0._wp
623      zdvice(:,:) = 0._wp
624      !
625      !                       !---------------------------------------
626      DO jl = 1, jpl-1        ! identify thicknesses that are too big
627         !                    !---------------------------------------
628         npti = 0   ;   nptidx(:) = 0
629         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
630            IF( a_i(ji,jj,jl) > 0._wp .AND. v_i(ji,jj,jl) > (a_i(ji,jj,jl) * hi_max(jl)) ) THEN
631               npti = npti + 1
632               nptidx( npti ) = (jj - 1) * jpi + ji
633            ENDIF
634         END_2D
635         !
636         IF( npti > 0 ) THEN
637            !!clem   CALL tab_2d_1d( npti, nptidx(1:npti), h_i_1d(1:npti), h_i(:,:,jl) )
638            CALL tab_2d_1d( npti, nptidx(1:npti), a_i_1d(1:npti), a_i(:,:,jl) )
639            CALL tab_2d_1d( npti, nptidx(1:npti), v_i_1d(1:npti), v_i(:,:,jl) )
640            !
641            DO ji = 1, npti
642               jdonor(ji,jl)  = jl
643               ! how much of a_i you send in cat sup is somewhat arbitrary
644               ! these are from CICE => transfer everything
645               !!zdaice(ji,jl)  = a_i_1d(ji)
646               !!zdvice(ji,jl)  = v_i_1d(ji)
647               ! these are from LLN => transfer only half of the category
648               zdaice(ji,jl)  =                       0.5_wp  * a_i_1d(ji)
649               zdvice(ji,jl)  = v_i_1d(ji) - (1._wp - 0.5_wp) * a_i_1d(ji) * hi_mean(jl)
650            END DO
651            !
652            CALL itd_shiftice( jdonor(1:npti,:), zdaice(1:npti,:), zdvice(1:npti,:) )  ! Shift jl=>jl+1
653            ! Reset shift parameters
654            jdonor(1:npti,jl) = 0
655            zdaice(1:npti,jl) = 0._wp
656            zdvice(1:npti,jl) = 0._wp
657         ENDIF
658         !
659      END DO
660
661      !                       !-----------------------------------------
662      DO jl = jpl-1, 1, -1    ! Identify thicknesses that are too small
663         !                    !-----------------------------------------
664         npti = 0 ; nptidx(:) = 0
665         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
666            IF( a_i(ji,jj,jl+1) > 0._wp .AND. v_i(ji,jj,jl+1) <= (a_i(ji,jj,jl+1) * hi_max(jl)) ) THEN
667               npti = npti + 1
668               nptidx( npti ) = (jj - 1) * jpi + ji
669            ENDIF
670         END_2D
671         !
672         IF( npti > 0 ) THEN
673            CALL tab_2d_1d( npti, nptidx(1:npti), a_i_1d(1:npti), a_i(:,:,jl+1) ) ! jl+1 is ok
674            CALL tab_2d_1d( npti, nptidx(1:npti), v_i_1d(1:npti), v_i(:,:,jl+1) ) ! jl+1 is ok
675            DO ji = 1, npti
676               jdonor(ji,jl) = jl + 1
677               zdaice(ji,jl) = a_i_1d(ji)
678               zdvice(ji,jl) = v_i_1d(ji)
679            END DO
680            !
681            CALL itd_shiftice( jdonor(1:npti,:), zdaice(1:npti,:), zdvice(1:npti,:) )  ! Shift jl+1=>jl
682            ! Reset shift parameters
683            jdonor(1:npti,jl) = 0
684            zdaice(1:npti,jl) = 0._wp
685            zdvice(1:npti,jl) = 0._wp
686         ENDIF
687         !
688      END DO
689      !
690      IF( ln_icediachk )   CALL ice_cons_hsm(1, 'iceitd_reb', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft)
691      IF( ln_icediachk )   CALL ice_cons2D  (1, 'iceitd_reb',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft)
692      IF( ln_timing    )   CALL timing_stop ('iceitd_reb')
693      !
694   END SUBROUTINE ice_itd_reb
695
696
697   SUBROUTINE ice_itd_init
698      !!------------------------------------------------------------------
699      !!                ***  ROUTINE ice_itd_init ***
700      !!
701      !! ** Purpose :   Initializes the ice thickness distribution
702      !! ** Method  :   ...
703      !! ** input   :   Namelist namitd
704      !!-------------------------------------------------------------------
705      INTEGER  ::   jl            ! dummy loop index
706      INTEGER  ::   ios, ioptio   ! Local integer output status for namelist read
707      REAL(wp) ::   zhmax, znum, zden, zalpha   !   -      -
708      !
709      NAMELIST/namitd/ ln_cat_hfn, rn_himean, ln_cat_usr, rn_catbnd, rn_himin, rn_himax
710      !!------------------------------------------------------------------
711      !
712      READ  ( numnam_ice_ref, namitd, IOSTAT = ios, ERR = 901)
713901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namitd in reference namelist' )
714      READ  ( numnam_ice_cfg, namitd, IOSTAT = ios, ERR = 902 )
715902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namitd in configuration namelist' )
716      IF(lwm) WRITE( numoni, namitd )
717      !
718      IF(lwp) THEN                  ! control print
719         WRITE(numout,*)
720         WRITE(numout,*) 'ice_itd_init: Initialization of ice cat distribution '
721         WRITE(numout,*) '~~~~~~~~~~~~'
722         WRITE(numout,*) '   Namelist namitd: '
723         WRITE(numout,*) '      Ice categories are defined by a function of rn_himean**(-0.05)    ln_cat_hfn = ', ln_cat_hfn
724         WRITE(numout,*) '         mean ice thickness in the domain                               rn_himean  = ', rn_himean
725         WRITE(numout,*) '      Ice categories are defined by rn_catbnd                           ln_cat_usr = ', ln_cat_usr
726         WRITE(numout,*) '      minimum ice thickness allowed                                     rn_himin   = ', rn_himin
727         WRITE(numout,*) '      maximum ice thickness allowed                                     rn_himax   = ', rn_himax
728      ENDIF
729      !
730      !-----------------------------------!
731      !  Thickness categories boundaries  !
732      !-----------------------------------!
733      !                             !== set the choice of ice categories ==!
734      ioptio = 0
735      IF( ln_cat_hfn ) THEN   ;   ioptio = ioptio + 1   ;   nice_catbnd = np_cathfn    ;   ENDIF
736      IF( ln_cat_usr ) THEN   ;   ioptio = ioptio + 1   ;   nice_catbnd = np_catusr    ;   ENDIF
737      IF( ioptio /= 1 )   CALL ctl_stop( 'ice_itd_init: choose one and only one ice categories boundaries' )
738      !
739      SELECT CASE( nice_catbnd )
740      !                                !------------------------!
741      CASE( np_cathfn )                ! h^(-alpha) function
742         !                             !------------------------!
743         zalpha = 0.05_wp
744         zhmax  = 3._wp * rn_himean
745         hi_max(0) = 0._wp
746         DO jl = 1, jpl
747            znum = jpl * ( zhmax+1 )**zalpha
748            zden = REAL( jpl-jl , wp ) * ( zhmax + 1._wp )**zalpha + REAL( jl , wp )
749            hi_max(jl) = ( znum / zden )**(1./zalpha) - 1
750         END DO
751         !                             !------------------------!
752      CASE( np_catusr )                ! user defined
753         !                             !------------------------!
754         DO jl = 0, jpl
755            hi_max(jl) = rn_catbnd(jl)
756         END DO
757         !
758      END SELECT
759      !
760      DO jl = 1, jpl                ! mean thickness by category
761         hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp
762      END DO
763      !
764      hi_max(jpl) = rn_himax        ! set to a big value to ensure that all ice is thinner than hi_max(jpl)
765      !
766      IF(lwp) WRITE(numout,*)
767      IF(lwp) WRITE(numout,*) '   ===>>>   resulting thickness category boundaries :'
768      IF(lwp) WRITE(numout,*) '            hi_max(:)= ', hi_max(0:jpl)
769      !
770      IF( hi_max(1) < rn_himin )   CALL ctl_stop('ice_itd_init: the upper bound of the 1st category must be bigger than rn_himin')
771      !
772   END SUBROUTINE ice_itd_init
773
774#else
775   !!----------------------------------------------------------------------
776   !!   Default option :         Empty module         NO SI3 sea-ice model
777   !!----------------------------------------------------------------------
778#endif
779
780   !!======================================================================
781END MODULE iceitd
Note: See TracBrowser for help on using the repository browser.