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/branches/2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles/src/ICE – NEMO

source: NEMO/branches/2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles/src/ICE/iceitd.F90 @ 11648

Last change on this file since 11648 was 11536, checked in by smasson, 5 years ago

trunk: merge dev_r10984_HPC-13 into the trunk

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