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.
sbcisf.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcisf.F90 @ 7788

Last change on this file since 7788 was 7788, checked in by cetlod, 7 years ago

trunk : representation of ice shelf melting in coupled mode

  • Property svn:keywords set to Id
File size: 43.3 KB
Line 
1MODULE sbcisf
2   !!======================================================================
3   !!                       ***  MODULE  sbcisf  ***
4   !! Surface module :  update surface ocean boundary condition under ice
5   !!                   shelf
6   !!======================================================================
7   !! History :  3.2   !  2011-02  (C.Harris  ) Original code isf cav
8   !!            X.X   !  2006-02  (C. Wang   ) Original code bg03
9   !!            3.4   !  2013-03  (P. Mathiot) Merging + parametrization
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   sbc_isf        : update sbc under ice shelf
14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and tracers
16   USE dom_oce         ! ocean space and time domain
17   USE phycst          ! physical constants
18   USE eosbn2          ! equation of state
19   USE sbc_oce         ! surface boundary condition: ocean fields
20   USE zdfbfr          !
21   !
22   USE in_out_manager  ! I/O manager
23   USE iom             ! I/O manager library
24   USE fldread         ! read input field at current time step
25   USE lbclnk          !
26   USE wrk_nemo        ! Memory allocation
27   USE timing          ! Timing
28   USE lib_fortran     ! glob_sum
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   sbc_isf, sbc_isf_div, sbc_isf_alloc  ! routine called in sbcmod and divhor
34
35   ! public in order to be able to output then
36
37   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   risf_tsc_b, risf_tsc  !: before and now T & S isf contents [K.m/s & PSU.m/s] 
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   qisf                  !: net heat flux from ice shelf      [W/m2]
39   REAL(wp), PUBLIC ::   rn_hisf_tbl                 !: thickness of top boundary layer [m]
40   INTEGER , PUBLIC ::   nn_isf                      !: flag to choose between explicit/param/specified 
41   INTEGER , PUBLIC ::   nn_isfblk                   !: flag to choose the bulk formulation to compute the ice shelf melting
42   INTEGER , PUBLIC ::   nn_gammablk                 !: flag to choose how the exchange coefficient is computed
43   REAL(wp), PUBLIC ::   rn_gammat0                  !: temperature exchange coeficient []
44   REAL(wp), PUBLIC ::   rn_gammas0                  !: salinity    exchange coeficient []
45
46   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  rzisf_tbl              !:depth of calving front (shallowest point) nn_isf ==2/3
47   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  rhisf_tbl, rhisf_tbl_0 !:thickness of tbl  [m]
48   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  r1_hisf_tbl            !:1/thickness of tbl
49   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  ralpha                 !:proportion of bottom cell influenced by tbl
50   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  risfLeff               !:effective length (Leff) BG03 nn_isf==2
51   REAL(wp)   , PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)     ::  ttbl, stbl, utbl, vtbl !:top boundary layer variable at T point
52   INTEGER,    PUBLIC, ALLOCATABLE, SAVE, DIMENSION (:,:)      ::  misfkt, misfkb         !:Level of ice shelf base
53
54   LOGICAL, PUBLIC ::   l_isfcpl = .false.       ! isf recieved from oasis
55
56   REAL(wp), PUBLIC, SAVE ::   rcpi     = 2000.0_wp     ! specific heat of ice shelf             [J/kg/K]
57   REAL(wp), PUBLIC, SAVE ::   rkappa   = 1.54e-6_wp    ! heat diffusivity through the ice-shelf [m2/s]
58   REAL(wp), PUBLIC, SAVE ::   rhoisf   = 920.0_wp      ! volumic mass of ice shelf              [kg/m3]
59   REAL(wp), PUBLIC, SAVE ::   tsurf    = -20.0_wp      ! air temperature on top of ice shelf    [C]
60   REAL(wp), PUBLIC, SAVE ::   rlfusisf = 0.334e6_wp    ! latent heat of fusion of ice shelf     [J/kg]
61
62!: Variable used in fldread to read the forcing file (nn_isf == 4 .OR. nn_isf == 3)
63   CHARACTER(len=100), PUBLIC           :: cn_dirisf  = './' !: Root directory for location of ssr files
64   TYPE(FLD_N)       , PUBLIC           :: sn_fwfisf         !: information about the isf melting file to be read
65   TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_fwfisf
66   TYPE(FLD_N)       , PUBLIC           :: sn_rnfisf         !: information about the isf melting param.   file to be read
67   TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_rnfisf           
68   TYPE(FLD_N)       , PUBLIC           :: sn_depmax_isf     !: information about the grounding line depth file to be read
69   TYPE(FLD_N)       , PUBLIC           :: sn_depmin_isf     !: information about the calving   line depth file to be read
70   TYPE(FLD_N)       , PUBLIC           :: sn_Leff_isf       !: information about the effective length     file to be read
71   
72   !!----------------------------------------------------------------------
73   !! NEMO/OPA 3.7 , LOCEAN-IPSL (2015)
74   !! $Id$
75   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
76   !!----------------------------------------------------------------------
77CONTAINS
78 
79  SUBROUTINE sbc_isf(kt)
80      !!---------------------------------------------------------------------
81      !!                  ***  ROUTINE sbc_isf  ***
82      !!
83      !! ** Purpose : Compute Salt and Heat fluxes related to ice_shelf
84      !!              melting and freezing
85      !!
86      !! ** Method  :  4 parameterizations are available according to nn_isf
87      !!               nn_isf = 1 : Realistic ice_shelf formulation
88      !!                        2 : Beckmann & Goose parameterization
89      !!                        3 : Specified runoff in deptht (Mathiot & al. )
90      !!                        4 : specified fwf and heat flux forcing beneath the ice shelf
91      !!----------------------------------------------------------------------
92      INTEGER, INTENT( in ) :: kt                   ! ocean time step
93      !
94      INTEGER               :: ji, jj, jk           ! loop index
95      INTEGER               :: ikt, ikb             ! loop index
96      REAL(wp), DIMENSION (:,:), POINTER :: zt_frz, zdep ! freezing temperature (zt_frz) at depth (zdep)
97      REAL(wp), DIMENSION(:,:,:), POINTER :: zfwfisf3d, zqhcisf3d, zqlatisf3d
98      REAL(wp), DIMENSION(:,:  ), POINTER :: zqhcisf2d
99      !!---------------------------------------------------------------------
100      !
101      IF( MOD( kt-1, nn_fsbc) == 0 ) THEN
102         ! allocation
103         CALL wrk_alloc( jpi,jpj, zt_frz, zdep  )
104
105         ! compute salt and heat flux
106         SELECT CASE ( nn_isf )
107         CASE ( 1 )    ! realistic ice shelf formulation
108            ! compute T/S/U/V for the top boundary layer
109            CALL sbc_isf_tbl(tsn(:,:,:,jp_tem),ttbl(:,:),'T')
110            CALL sbc_isf_tbl(tsn(:,:,:,jp_sal),stbl(:,:),'T')
111            CALL sbc_isf_tbl(un(:,:,:)        ,utbl(:,:),'U')
112            CALL sbc_isf_tbl(vn(:,:,:)        ,vtbl(:,:),'V')
113            ! iom print
114            CALL iom_put('ttbl',ttbl(:,:))
115            CALL iom_put('stbl',stbl(:,:))
116            CALL iom_put('utbl',utbl(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:))
117            CALL iom_put('vtbl',vtbl(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:))
118            ! compute fwf and heat flux
119            ! compute fwf and heat flux
120            IF( .NOT.l_isfcpl ) THEN    ;   CALL sbc_isf_cav (kt)
121            ELSE                        ;   qisf(:,:)  = fwfisf(:,:) * rlfusisf  ! heat        flux
122            ENDIF
123
124         CASE ( 2 )    ! Beckmann and Goosse parametrisation
125            stbl(:,:)   = soce
126            CALL sbc_isf_bg03(kt)
127
128         CASE ( 3 )    ! specified runoff in depth (Mathiot et al., XXXX in preparation)
129            ! specified runoff in depth (Mathiot et al., XXXX in preparation)
130            IF( .NOT.l_isfcpl ) THEN
131               CALL fld_read ( kt, nn_fsbc, sf_rnfisf   )
132               fwfisf(:,:) = - sf_rnfisf(1)%fnow(:,:,1)         ! fresh water flux from the isf (fwfisf <0 mean melting)
133            ENDIF
134            qisf(:,:)   = fwfisf(:,:) * rlfusisf             ! heat flux
135            stbl(:,:)   = soce
136
137         CASE ( 4 )    ! specified fwf and heat flux forcing beneath the ice shelf
138           ! specified fwf and heat flux forcing beneath the ice shelf
139            IF( .NOT.l_isfcpl ) THEN
140               CALL fld_read ( kt, nn_fsbc, sf_fwfisf   )
141               !CALL fld_read ( kt, nn_fsbc, sf_qisf   )
142               fwfisf(:,:) = -sf_fwfisf(1)%fnow(:,:,1)            ! fwf
143            ENDIF
144            qisf(:,:)   = fwfisf(:,:) * rlfusisf               ! heat flux
145            stbl(:,:)   = soce
146
147         END SELECT
148
149         ! compute tsc due to isf
150         ! isf melting implemented as a volume flux and we assume that melt water is at 0 PSU.
151         ! WARNING water add at temp = 0C, need to add a correction term (fwfisf * tfreez / rau0).
152         ! compute freezing point beneath ice shelf (or top cell if nn_isf = 3)
153         DO jj = 1,jpj
154            DO ji = 1,jpi
155               zdep(ji,jj)=gdepw_n(ji,jj,misfkt(ji,jj))
156            END DO
157         END DO
158         CALL eos_fzp( stbl(:,:), zt_frz(:,:), zdep(:,:) )
159         
160         risf_tsc(:,:,jp_tem) = qisf(:,:) * r1_rau0_rcp - fwfisf(:,:) * zt_frz(:,:) * r1_rau0 !
161         risf_tsc(:,:,jp_sal) = 0.0_wp
162
163         ! lbclnk
164         CALL lbc_lnk(risf_tsc(:,:,jp_tem),'T',1.)
165         CALL lbc_lnk(risf_tsc(:,:,jp_sal),'T',1.)
166         CALL lbc_lnk(fwfisf(:,:)   ,'T',1.)
167         CALL lbc_lnk(qisf(:,:)     ,'T',1.)
168
169         ! output
170         CALL iom_put('qlatisf' , qisf)
171         CALL iom_put('fwfisf'  , fwfisf)
172
173        ! Diagnostics
174        IF ( iom_use('fwfisf3d') .OR. iom_use('qlatisf3d') .OR. iom_use('qhcisf3d') .OR. iom_use('qhcisf')) THEN
175            CALL wrk_alloc( jpi,jpj,jpk, zfwfisf3d, zqhcisf3d, zqlatisf3d )
176            CALL wrk_alloc( jpi,jpj,     zqhcisf2d                        )
177
178            zfwfisf3d(:,:,:) = 0.0_wp                         ! 3d ice shelf melting (kg/m2/s)
179            zqhcisf3d(:,:,:) = 0.0_wp                         ! 3d heat content flux (W/m2)
180            zqlatisf3d(:,:,:)= 0.0_wp                         ! 3d ice shelf melting latent heat flux (W/m2)
181            zqhcisf2d(:,:)   = fwfisf(:,:) * zt_frz * rcp     ! 2d heat content flux (W/m2)
182
183            DO jj = 1,jpj
184               DO ji = 1,jpi
185                  ikt = misfkt(ji,jj)
186                  ikb = misfkb(ji,jj)
187                  DO jk = ikt, ikb - 1
188                     zfwfisf3d (ji,jj,jk) = zfwfisf3d (ji,jj,jk) + fwfisf   (ji,jj) * r1_hisf_tbl(ji,jj) * e3t_n(ji,jj,jk)
189                     zqhcisf3d (ji,jj,jk) = zqhcisf3d (ji,jj,jk) + zqhcisf2d(ji,jj) * r1_hisf_tbl(ji,jj) * e3t_n(ji,jj,jk)
190                     zqlatisf3d(ji,jj,jk) = zqlatisf3d(ji,jj,jk) + qisf     (ji,jj) * r1_hisf_tbl(ji,jj) * e3t_n(ji,jj,jk)
191                  END DO
192                  zfwfisf3d (ji,jj,jk) = zfwfisf3d (ji,jj,jk) + fwfisf   (ji,jj) * r1_hisf_tbl(ji,jj) * ralpha(ji,jj) * e3t_n(ji,jj,jk)
193                  zqhcisf3d (ji,jj,jk) = zqhcisf3d (ji,jj,jk) + zqhcisf2d(ji,jj) * r1_hisf_tbl(ji,jj) * ralpha(ji,jj) * e3t_n(ji,jj,jk)
194                  zqlatisf3d(ji,jj,jk) = zqlatisf3d(ji,jj,jk) + qisf     (ji,jj) * r1_hisf_tbl(ji,jj) * ralpha(ji,jj) * e3t_n(ji,jj,jk)
195               END DO
196            END DO
197
198            CALL iom_put('fwfisf3d' , zfwfisf3d (:,:,:))
199            CALL iom_put('qlatisf3d', zqlatisf3d(:,:,:))
200            CALL iom_put('qhcisf3d' , zqhcisf3d (:,:,:))
201            CALL iom_put('qhcisf'   , zqhcisf2d (:,:  ))
202
203            CALL wrk_dealloc( jpi,jpj,jpk, zfwfisf3d, zqhcisf3d, zqlatisf3d )
204            CALL wrk_dealloc( jpi,jpj,     zqhcisf2d                        )
205          END IF
206          ! deallocation
207          CALL wrk_dealloc( jpi,jpj, zt_frz, zdep  )
208          !
209        END IF
210
211        IF( kt == nit000 ) THEN                         !   set the forcing field at nit000 - 1    !
212           IF( ln_rstart .AND.    &                     ! Restart: read in restart file
213                 & iom_varid( numror, 'fwf_isf_b', ldstop = .FALSE. ) > 0 ) THEN
214               IF(lwp) WRITE(numout,*) '          nit000-1 isf tracer content forcing fields read in the restart file'
215               CALL iom_get( numror, jpdom_autoglo, 'fwf_isf_b', fwfisf_b(:,:) )   ! before salt content isf_tsc trend
216               CALL iom_get( numror, jpdom_autoglo, 'isf_sc_b', risf_tsc_b(:,:,jp_sal) )   ! before salt content isf_tsc trend
217               CALL iom_get( numror, jpdom_autoglo, 'isf_hc_b', risf_tsc_b(:,:,jp_tem) )   ! before salt content isf_tsc trend
218           ELSE
219               fwfisf_b(:,:)    = fwfisf(:,:)
220               risf_tsc_b(:,:,:)= risf_tsc(:,:,:)
221           END IF
222         END IF
223         !
224         IF( lrst_oce ) THEN
225            IF(lwp) WRITE(numout,*)
226            IF(lwp) WRITE(numout,*) 'sbc : isf surface tracer content forcing fields written in ocean restart file ',   &
227               &                    'at it= ', kt,' date= ', ndastp
228            IF(lwp) WRITE(numout,*) '~~~~'
229            CALL iom_rstput( kt, nitrst, numrow, 'fwf_isf_b', fwfisf(:,:) )
230            CALL iom_rstput( kt, nitrst, numrow, 'isf_hc_b' , risf_tsc(:,:,jp_tem) )
231            CALL iom_rstput( kt, nitrst, numrow, 'isf_sc_b' , risf_tsc(:,:,jp_sal) )
232         ENDIF
233         !
234  END SUBROUTINE sbc_isf
235
236
237  INTEGER FUNCTION sbc_isf_alloc()
238      !!----------------------------------------------------------------------
239      !!               ***  FUNCTION sbc_isf_rnf_alloc  ***
240      !!----------------------------------------------------------------------
241      sbc_isf_alloc = 0       ! set to zero if no array to be allocated
242      IF( .NOT. ALLOCATED( qisf ) ) THEN
243         ALLOCATE(  risf_tsc(jpi,jpj,jpts), risf_tsc_b(jpi,jpj,jpts), qisf(jpi,jpj)   , &
244               &    rhisf_tbl(jpi,jpj)    , r1_hisf_tbl(jpi,jpj), rzisf_tbl(jpi,jpj)  , &
245               &    ttbl(jpi,jpj)         , stbl(jpi,jpj)       , utbl(jpi,jpj)       , &
246               &    vtbl(jpi, jpj)        , risfLeff(jpi,jpj)   , rhisf_tbl_0(jpi,jpj), &
247               &    ralpha(jpi,jpj)       , misfkt(jpi,jpj)     , misfkb(jpi,jpj)     , &
248               &    STAT= sbc_isf_alloc )
249         !
250         IF( lk_mpp             )   CALL mpp_sum ( sbc_isf_alloc )
251         IF( sbc_isf_alloc /= 0 )   CALL ctl_warn('sbc_isf_alloc: failed to allocate arrays.')
252         !
253      END IF
254  END FUNCTION
255
256  SUBROUTINE sbc_isf_init
257      !!---------------------------------------------------------------------
258      !!                  ***  ROUTINE sbc_isf_init  ***
259      !!
260      !! ** Purpose : Initialisation of variables for iceshelf fluxes formulation
261      !!
262      !! ** Method  :  4 parameterizations are available according to nn_isf
263      !!               nn_isf = 1 : Realistic ice_shelf formulation
264      !!                        2 : Beckmann & Goose parameterization
265      !!                        3 : Specified runoff in deptht (Mathiot & al. )
266      !!                        4 : specified fwf and heat flux forcing beneath the ice shelf
267      !!----------------------------------------------------------------------
268      INTEGER               :: ji, jj, jk           ! loop index
269      INTEGER               :: ik                   ! current level index
270      INTEGER               :: ikt, ikb             ! top and bottom level of the isf boundary layer
271      INTEGER               :: inum, ierror
272      INTEGER               :: ios                  ! Local integer output status for namelist read
273      REAL(wp)              :: zhk
274      CHARACTER(len=256)    :: cvarzisf, cvarhisf   ! name for isf file
275      CHARACTER(LEN=32 )    :: cvarLeff             ! variable name for efficient Length scale
276      !!----------------------------------------------------------------------
277      NAMELIST/namsbc_isf/ nn_isfblk, rn_hisf_tbl, rn_gammat0, rn_gammas0, nn_gammablk, nn_isf, &
278                         & sn_fwfisf, sn_rnfisf, sn_depmax_isf, sn_depmin_isf, sn_Leff_isf
279      !!----------------------------------------------------------------------
280
281      REWIND( numnam_ref )              ! Namelist namsbc_rnf in reference namelist : Runoffs
282      READ  ( numnam_ref, namsbc_isf, IOSTAT = ios, ERR = 901)
283901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_isf in reference namelist', lwp )
284
285      REWIND( numnam_cfg )              ! Namelist namsbc_rnf in configuration namelist : Runoffs
286      READ  ( numnam_cfg, namsbc_isf, IOSTAT = ios, ERR = 902 )
287902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_isf in configuration namelist', lwp )
288      IF(lwm) WRITE ( numond, namsbc_isf )
289
290      IF ( lwp ) WRITE(numout,*)
291      IF ( lwp ) WRITE(numout,*) 'sbc_isf: heat flux of the ice shelf'
292      IF ( lwp ) WRITE(numout,*) '~~~~~~~~~'
293      IF ( lwp ) WRITE(numout,*) 'sbcisf :' 
294      IF ( lwp ) WRITE(numout,*) '~~~~~~~~'
295      IF ( lwp ) WRITE(numout,*) '        nn_isf      = ', nn_isf
296      IF ( lwp ) WRITE(numout,*) '        nn_isfblk   = ', nn_isfblk
297      IF ( lwp ) WRITE(numout,*) '        rn_hisf_tbl = ', rn_hisf_tbl
298      IF ( lwp ) WRITE(numout,*) '        nn_gammablk = ', nn_gammablk 
299      IF ( lwp ) WRITE(numout,*) '        rn_gammat0  = ', rn_gammat0 
300      IF ( lwp ) WRITE(numout,*) '        rn_gammas0  = ', rn_gammas0 
301      IF ( lwp ) WRITE(numout,*) '        rn_tfri2    = ', rn_tfri2 
302      !
303      ! Allocate public variable
304      IF ( sbc_isf_alloc()  /= 0 )         CALL ctl_stop( 'STOP', 'sbc_isf : unable to allocate arrays' )
305      !
306      ! initialisation
307      qisf(:,:)        = 0._wp  ; fwfisf  (:,:) = 0._wp
308      risf_tsc(:,:,:)  = 0._wp  ; fwfisf_b(:,:) = 0._wp
309      !
310      ! define isf tbl tickness, top and bottom indice
311      SELECT CASE ( nn_isf )
312      CASE ( 1 ) 
313         rhisf_tbl(:,:) = rn_hisf_tbl
314         misfkt(:,:)    = mikt(:,:)         ! same indice for bg03 et cav => used in isfdiv
315
316      CASE ( 2 , 3 )
317         IF( .NOT.l_isfcpl ) THEN
318             ALLOCATE( sf_rnfisf(1), STAT=ierror )
319             ALLOCATE( sf_rnfisf(1)%fnow(jpi,jpj,1), sf_rnfisf(1)%fdta(jpi,jpj,1,2) )
320             CALL fld_fill( sf_rnfisf, (/ sn_rnfisf /), cn_dirisf, 'sbc_isf_init', 'read fresh water flux isf data', 'namsbc_isf' )
321          ENDIF
322          !  read effective lenght (BG03)
323          IF (nn_isf == 2) THEN
324            CALL iom_open( sn_Leff_isf%clname, inum )
325            cvarLeff = TRIM(sn_Leff_isf%clvar)
326            CALL iom_get( inum, jpdom_data, cvarLeff, risfLeff , 1)
327            CALL iom_close(inum)
328            !
329            risfLeff = risfLeff*1000.0_wp           !: convertion in m
330          END IF
331         ! read depth of the top and bottom of the isf top boundary layer (in this case, isf front depth and grounding line depth)
332         CALL iom_open( sn_depmax_isf%clname, inum )
333         cvarhisf = TRIM(sn_depmax_isf%clvar)
334         CALL iom_get( inum, jpdom_data, cvarhisf, rhisf_tbl, 1) !: depth of deepest point of the ice shelf base
335         CALL iom_close(inum)
336         !
337         CALL iom_open( sn_depmin_isf%clname, inum )
338         cvarzisf = TRIM(sn_depmin_isf%clvar)
339         CALL iom_get( inum, jpdom_data, cvarzisf, rzisf_tbl, 1) !: depth of shallowest point of the ice shelves base
340         CALL iom_close(inum)
341         !
342         rhisf_tbl(:,:) = rhisf_tbl(:,:) - rzisf_tbl(:,:)        !: tickness isf boundary layer
343
344         !! compute first level of the top boundary layer
345         DO ji = 1, jpi
346            DO jj = 1, jpj
347                ik = 2
348                DO WHILE ( ik <= mbkt(ji,jj) .AND. gdepw_n(ji,jj,ik) < rzisf_tbl(ji,jj) ) ;  ik = ik + 1 ;  END DO
349                misfkt(ji,jj) = ik-1
350            END DO
351         END DO
352
353      CASE ( 4 ) 
354         ! as in nn_isf == 1
355         rhisf_tbl(:,:) = rn_hisf_tbl
356         misfkt(:,:)    = mikt(:,:)         ! same indice for bg03 et cav => used in isfdiv
357         
358         ! load variable used in fldread (use for temporal interpolation of isf fwf forcing)
359         IF( .NOT.l_isfcpl ) THEN
360           ALLOCATE( sf_fwfisf(1), STAT=ierror )
361           ALLOCATE( sf_fwfisf(1)%fnow(jpi,jpj,1), sf_fwfisf(1)%fdta(jpi,jpj,1,2) )
362           CALL fld_fill( sf_fwfisf, (/ sn_fwfisf /), cn_dirisf, 'sbc_isf_init', 'read fresh water flux isf data', 'namsbc_isf' )
363         ENDIF
364
365      END SELECT
366         
367      rhisf_tbl_0(:,:) = rhisf_tbl(:,:)
368
369      ! compute bottom level of isf tbl and thickness of tbl below the ice shelf
370      DO jj = 1,jpj
371         DO ji = 1,jpi
372            ikt = misfkt(ji,jj)
373            ikb = misfkt(ji,jj)
374            ! thickness of boundary layer at least the top level thickness
375            rhisf_tbl(ji,jj) = MAX(rhisf_tbl_0(ji,jj), e3t_n(ji,jj,ikt))
376
377            ! determine the deepest level influenced by the boundary layer
378            DO jk = ikt+1, mbkt(ji,jj)
379               IF ( (SUM(e3t_n(ji,jj,ikt:jk-1)) < rhisf_tbl(ji,jj)) .AND. (tmask(ji,jj,jk) == 1) ) ikb = jk
380            END DO
381            rhisf_tbl(ji,jj) = MIN(rhisf_tbl(ji,jj), SUM(e3t_n(ji,jj,ikt:ikb))) ! limit the tbl to water thickness.
382            misfkb(ji,jj) = ikb                                                   ! last wet level of the tbl
383            r1_hisf_tbl(ji,jj) = 1._wp / rhisf_tbl(ji,jj)
384
385            zhk           = SUM( e3t_n(ji, jj, ikt:ikb - 1)) * r1_hisf_tbl(ji,jj) ! proportion of tbl cover by cell from ikt to ikb - 1
386            ralpha(ji,jj) = rhisf_tbl(ji,jj) * (1._wp - zhk ) / e3t_n(ji,jj,ikb)  ! proportion of bottom cell influenced by boundary layer
387         END DO
388      END DO
389
390  END SUBROUTINE sbc_isf_init
391
392  SUBROUTINE sbc_isf_bg03(kt)
393      !!---------------------------------------------------------------------
394      !!                  ***  ROUTINE sbc_isf_bg03  ***
395      !!
396      !! ** Purpose : add net heat and fresh water flux from ice shelf melting
397      !!          into the adjacent ocean
398      !!
399      !! ** Method  :   See reference
400      !!
401      !! ** Reference : Beckmann and Goosse (2003), "A parameterization of ice shelf-ocean
402      !!         interaction for climate models", Ocean Modelling 5(2003) 157-170.
403      !!         (hereafter BG)
404      !! History :
405      !!         06-02  (C. Wang) Original code
406      !!----------------------------------------------------------------------
407      INTEGER, INTENT ( in ) :: kt
408      !
409      INTEGER  :: ji, jj, jk ! dummy loop index
410      INTEGER  :: ik         ! current level
411      REAL(wp) :: zt_sum     ! sum of the temperature between 200m and 600m
412      REAL(wp) :: zt_ave     ! averaged temperature between 200m and 600m
413      REAL(wp) :: zt_frz     ! freezing point temperature at depth z
414      REAL(wp) :: zpress     ! pressure to compute the freezing point in depth
415      !!----------------------------------------------------------------------
416
417      IF ( nn_timing == 1 ) CALL timing_start('sbc_isf_bg03')
418      !
419      DO ji = 1, jpi
420         DO jj = 1, jpj
421            ik = misfkt(ji,jj)
422            !! Initialize arrays to 0 (each step)
423            zt_sum = 0.e0_wp
424            IF ( ik > 1 ) THEN
425               ! 1. -----------the average temperature between 200m and 600m ---------------------
426               DO jk = misfkt(ji,jj),misfkb(ji,jj)
427                  ! freezing point temperature  at ice shelf base BG eq. 2 (JMM sign pb ??? +7.64e-4 !!!)
428                  ! after verif with UNESCO, wrong sign in BG eq. 2
429                  ! Calculate freezing temperature
430                  CALL eos_fzp(stbl(ji,jj), zt_frz, zpress) 
431                  zt_sum = zt_sum + (tsn(ji,jj,jk,jp_tem)-zt_frz) * e3t_n(ji,jj,jk) * tmask(ji,jj,jk)  ! sum temp
432               END DO
433               zt_ave = zt_sum/rhisf_tbl(ji,jj) ! calcul mean value
434               ! 2. ------------Net heat flux and fresh water flux due to the ice shelf
435               ! For those corresponding to zonal boundary   
436               qisf(ji,jj) = - rau0 * rcp * rn_gammat0 * risfLeff(ji,jj) * e1t(ji,jj) * zt_ave  &
437                           & * r1_e1e2t(ji,jj) * tmask(ji,jj,jk)
438             
439               fwfisf(ji,jj) = qisf(ji,jj) / rlfusisf          !fresh water flux kg/(m2s)                 
440               fwfisf(ji,jj) = fwfisf(ji,jj) * ( soce / stbl(ji,jj) )
441               !add to salinity trend
442            ELSE
443               qisf(ji,jj) = 0._wp ; fwfisf(ji,jj) = 0._wp
444            END IF
445         END DO
446      END DO
447      !
448      IF( nn_timing == 1 )  CALL timing_stop('sbc_isf_bg03')
449      !
450  END SUBROUTINE sbc_isf_bg03
451
452  SUBROUTINE sbc_isf_cav( kt )
453      !!---------------------------------------------------------------------
454      !!                     ***  ROUTINE sbc_isf_cav  ***
455      !!
456      !! ** Purpose :   handle surface boundary condition under ice shelf
457      !!
458      !! ** Method  : -
459      !!
460      !! ** Action  :   utau, vtau : remain unchanged
461      !!                taum, wndm : remain unchanged
462      !!                qns        : update heat flux below ice shelf
463      !!                emp, emps  : update freshwater flux below ice shelf
464      !!---------------------------------------------------------------------
465      INTEGER, INTENT(in)          ::   kt         ! ocean time step
466      !
467      INTEGER  ::   ji, jj     ! dummy loop indices
468      INTEGER  ::   nit
469      REAL(wp) ::   zlamb1, zlamb2, zlamb3
470      REAL(wp) ::   zeps1,zeps2,zeps3,zeps4,zeps6,zeps7
471      REAL(wp) ::   zaqe,zbqe,zcqe,zaqer,zdis,zsfrz,zcfac
472      REAL(wp) ::   zeps = 1.e-20_wp       
473      REAL(wp) ::   zerr
474      REAL(wp), DIMENSION(:,:), POINTER ::   zfrz
475      REAL(wp), DIMENSION(:,:), POINTER ::   zgammat, zgammas 
476      REAL(wp), DIMENSION(:,:), POINTER ::   zfwflx, zhtflx, zhtflx_b
477      LOGICAL  ::   lit
478      !!---------------------------------------------------------------------
479      ! coeficient for linearisation of potential tfreez
480      ! Crude approximation for pressure (but commonly used)
481      zlamb1 =-0.0573_wp
482      zlamb2 = 0.0832_wp
483      zlamb3 =-7.53e-08_wp * grav * rau0
484      IF( nn_timing == 1 )  CALL timing_start('sbc_isf_cav')
485      !
486      CALL wrk_alloc( jpi,jpj, zfrz  , zgammat, zgammas  )
487      CALL wrk_alloc( jpi,jpj, zfwflx, zhtflx , zhtflx_b )
488
489      ! initialisation
490      zgammat(:,:) = rn_gammat0 ; zgammas (:,:) = rn_gammas0
491      zhtflx (:,:) = 0.0_wp     ; zhtflx_b(:,:) = 0.0_wp   
492      zfwflx (:,:) = 0.0_wp
493
494      ! compute ice shelf melting
495      nit = 1 ; lit = .TRUE.
496      DO WHILE ( lit )    ! maybe just a constant number of iteration as in blk_core is fine
497         SELECT CASE ( nn_isfblk )
498         CASE ( 1 )   !  ISOMIP formulation (2 equations) for volume flux (Hunter et al., 2006)
499            ! Calculate freezing temperature
500            CALL eos_fzp( stbl(:,:), zfrz(:,:), risfdep(:,:) )
501
502            ! compute gammat every where (2d)
503            CALL sbc_isf_gammats(zgammat, zgammas, zhtflx, zfwflx)
504           
505            ! compute upward heat flux zhtflx and upward water flux zwflx
506            DO jj = 1, jpj
507               DO ji = 1, jpi
508                  zhtflx(ji,jj) =   zgammat(ji,jj)*rcp*rau0*(ttbl(ji,jj)-zfrz(ji,jj))
509                  zfwflx(ji,jj) = - zhtflx(ji,jj)/rlfusisf
510               END DO
511            END DO
512
513            ! Compute heat flux and upward fresh water flux
514            qisf  (:,:) = - zhtflx(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:)
515            fwfisf(:,:) =   zfwflx(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:)
516
517         CASE ( 2 )  ! ISOMIP+ formulation (3 equations) for volume flux (Asay-Davis et al., 2015)
518            ! compute gammat every where (2d)
519            CALL sbc_isf_gammats(zgammat, zgammas, zhtflx, zfwflx)
520
521            ! compute upward heat flux zhtflx and upward water flux zwflx
522            ! Resolution of a 2d equation from equation 21, 22 and 23 to find Sb (Asay-Davis et al., 2015)
523            DO jj = 1, jpj
524               DO ji = 1, jpi
525                  ! compute coeficient to solve the 2nd order equation
526                  zeps1 = rcp*rau0*zgammat(ji,jj)
527                  zeps2 = rlfusisf*rau0*zgammas(ji,jj)
528                  zeps3 = rhoisf*rcpi*rkappa/MAX(risfdep(ji,jj),zeps)
529                  zeps4 = zlamb2+zlamb3*risfdep(ji,jj)
530                  zeps6 = zeps4-ttbl(ji,jj)
531                  zeps7 = zeps4-tsurf
532                  zaqe  = zlamb1 * (zeps1 + zeps3)
533                  zaqer = 0.5_wp/MIN(zaqe,-zeps)
534                  zbqe  = zeps1*zeps6+zeps3*zeps7-zeps2
535                  zcqe  = zeps2*stbl(ji,jj)
536                  zdis  = zbqe*zbqe-4.0_wp*zaqe*zcqe               
537
538                  ! Presumably zdis can never be negative because gammas is very small compared to gammat
539                  ! compute s freeze
540                  zsfrz=(-zbqe-SQRT(zdis))*zaqer
541                  IF ( zsfrz < 0.0_wp ) zsfrz=(-zbqe+SQRT(zdis))*zaqer
542
543                  ! compute t freeze (eq. 22)
544                  zfrz(ji,jj)=zeps4+zlamb1*zsfrz
545 
546                  ! zfwflx is upward water flux
547                  ! zhtflx is upward heat flux (out of ocean)
548                  ! compute the upward water and heat flux (eq. 28 and eq. 29)
549                  zfwflx(ji,jj) = rau0 * zgammas(ji,jj) * (zsfrz-stbl(ji,jj)) / MAX(zsfrz,zeps)
550                  zhtflx(ji,jj) = zgammat(ji,jj) * rau0 * rcp * (ttbl(ji,jj) - zfrz(ji,jj) ) 
551               END DO
552            END DO
553
554            ! compute heat and water flux
555            qisf  (:,:) = - zhtflx(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:)
556            fwfisf(:,:) =   zfwflx(:,:) * (1._wp - tmask(:,:,1)) * ssmask(:,:)
557
558         END SELECT
559
560         ! define if we need to iterate (nn_gammablk 0/1 do not need iteration)
561         IF ( nn_gammablk <  2 ) THEN ; lit = .FALSE.
562         ELSE                           
563            ! check total number of iteration
564            IF (nit >= 100) THEN ; CALL ctl_stop( 'STOP', 'sbc_isf_hol99 : too many iteration ...' )
565            ELSE                 ; nit = nit + 1
566            END IF
567
568            ! compute error between 2 iterations
569            ! if needed save gammat and compute zhtflx_b for next iteration
570            zerr = MAXVAL(ABS(zhtflx-zhtflx_b))
571            IF ( zerr <= 0.01_wp ) THEN ; lit = .FALSE.
572            ELSE                        ; zhtflx_b(:,:) = zhtflx(:,:)
573            END IF
574         END IF
575      END DO
576      !
577      CALL iom_put('isfgammat', zgammat)
578      CALL iom_put('isfgammas', zgammas)
579      !
580      CALL wrk_dealloc( jpi,jpj, zfrz  , zgammat, zgammas  )
581      CALL wrk_dealloc( jpi,jpj, zfwflx, zhtflx , zhtflx_b )
582      !
583      IF( nn_timing == 1 )  CALL timing_stop('sbc_isf_cav')
584      !
585   END SUBROUTINE sbc_isf_cav
586
587
588   SUBROUTINE sbc_isf_gammats(pgt, pgs, pqhisf, pqwisf )
589      !!----------------------------------------------------------------------
590      !! ** Purpose    : compute the coefficient echange for heat flux
591      !!
592      !! ** Method     : gamma assume constant or depends of u* and stability
593      !!
594      !! ** References : Holland and Jenkins, 1999, JPO, p1787-1800, eq 14
595      !!                Jenkins et al., 2010, JPO, p2298-2312
596      !!---------------------------------------------------------------------
597      REAL(wp), DIMENSION(:,:), INTENT(out) :: pgt, pgs
598      REAL(wp), DIMENSION(:,:), INTENT(in ) :: pqhisf, pqwisf
599      !
600      INTEGER  :: ikt                       
601      INTEGER  :: ji, jj                     ! loop index
602      REAL(wp), DIMENSION(:,:), POINTER :: zustar           ! U, V at T point and friction velocity
603      REAL(wp) :: zdku, zdkv                 ! U, V shear
604      REAL(wp) :: zPr, zSc, zRc              ! Prandtl, Scmidth and Richardson number
605      REAL(wp) :: zmob, zmols                ! Monin Obukov length, coriolis factor at T point
606      REAL(wp) :: zbuofdep, zhnu             ! Bouyancy length scale, sublayer tickness
607      REAL(wp) :: zhmax                      ! limitation of mol
608      REAL(wp) :: zetastar                   ! stability parameter
609      REAL(wp) :: zgmolet, zgmoles, zgturb   ! contribution of modelecular sublayer and turbulence
610      REAL(wp) :: zcoef                      ! temporary coef
611      REAL(wp) :: zdep
612      REAL(wp) :: zeps = 1.0e-20_wp   
613      REAL(wp), PARAMETER :: zxsiN = 0.052_wp   ! dimensionless constant
614      REAL(wp), PARAMETER :: znu   = 1.95e-6_wp ! kinamatic viscosity of sea water (m2.s-1)
615      REAL(wp), DIMENSION(2) :: zts, zab
616      !!---------------------------------------------------------------------
617      CALL wrk_alloc( jpi,jpj, zustar )
618      !
619      SELECT CASE ( nn_gammablk )
620      CASE ( 0 ) ! gamma is constant (specified in namelist)
621         !! ISOMIP formulation (Hunter et al, 2006)
622         pgt(:,:) = rn_gammat0
623         pgs(:,:) = rn_gammas0
624
625      CASE ( 1 ) ! gamma is assume to be proportional to u*
626         !! Jenkins et al., 2010, JPO, p2298-2312
627         !! Adopted by Asay-Davis et al. (2015)
628
629         !! compute ustar (eq. 24)
630         zustar(:,:) = SQRT( rn_tfri2 * (utbl(:,:) * utbl(:,:) + vtbl(:,:) * vtbl(:,:) + rn_tfeb2) )
631
632         !! Compute gammats
633         pgt(:,:) = zustar(:,:) * rn_gammat0
634         pgs(:,:) = zustar(:,:) * rn_gammas0
635     
636      CASE ( 2 ) ! gamma depends of stability of boundary layer
637         !! Holland and Jenkins, 1999, JPO, p1787-1800, eq 14
638         !! as MOL depends of flux and flux depends of MOL, best will be iteration (TO DO)
639         !! compute ustar
640         zustar(:,:) = SQRT( rn_tfri2 * (utbl(:,:) * utbl(:,:) + vtbl(:,:) * vtbl(:,:) + rn_tfeb2) )
641
642         !! compute Pr and Sc number (can be improved)
643         zPr =   13.8_wp
644         zSc = 2432.0_wp
645
646         !! compute gamma mole
647         zgmolet = 12.5_wp * zPr ** (2.0/3.0) - 6.0_wp
648         zgmoles = 12.5_wp * zSc ** (2.0/3.0) - 6.0_wp
649
650         !! compute gamma
651         DO ji=2,jpi
652            DO jj=2,jpj
653               ikt = mikt(ji,jj)
654
655               IF (zustar(ji,jj) == 0._wp) THEN           ! only for kt = 1 I think
656                  pgt = rn_gammat0
657                  pgs = rn_gammas0
658               ELSE
659                  !! compute Rc number (as done in zdfric.F90)
660                  zcoef = 0.5_wp / e3w_n(ji,jj,ikt)
661                  !                                            ! shear of horizontal velocity
662                  zdku = zcoef * (  un(ji-1,jj  ,ikt  ) + un(ji,jj,ikt  )  &
663                     &             -un(ji-1,jj  ,ikt+1) - un(ji,jj,ikt+1)  )
664                  zdkv = zcoef * (  vn(ji  ,jj-1,ikt  ) + vn(ji,jj,ikt  )  &
665                     &             -vn(ji  ,jj-1,ikt+1) - vn(ji,jj,ikt+1)  )
666                  !                                            ! richardson number (minimum value set to zero)
667                  zRc = rn2(ji,jj,ikt+1) / MAX( zdku*zdku + zdkv*zdkv, zeps )
668
669                  !! compute bouyancy
670                  zts(jp_tem) = ttbl(ji,jj)
671                  zts(jp_sal) = stbl(ji,jj)
672                  zdep        = gdepw_n(ji,jj,ikt)
673                  !
674                  CALL eos_rab( zts, zdep, zab )
675                  !
676                  !! compute length scale
677                  zbuofdep = grav * ( zab(jp_tem) * pqhisf(ji,jj) - zab(jp_sal) * pqwisf(ji,jj) )  !!!!!!!!!!!!!!!!!!!!!!!!!!!!
678
679                  !! compute Monin Obukov Length
680                  ! Maximum boundary layer depth
681                  zhmax = gdept_n(ji,jj,mbkt(ji,jj)) - gdepw_n(ji,jj,mikt(ji,jj)) -0.001_wp
682                  ! Compute Monin obukhov length scale at the surface and Ekman depth:
683                  zmob   = zustar(ji,jj) ** 3 / (vkarmn * (zbuofdep + zeps))
684                  zmols  = SIGN(1._wp, zmob) * MIN(ABS(zmob), zhmax) * tmask(ji,jj,ikt)
685
686                  !! compute eta* (stability parameter)
687                  zetastar = 1._wp / ( SQRT(1._wp + MAX(zxsiN * zustar(ji,jj) / ( ABS(ff_f(ji,jj)) * zmols * zRc ), 0._wp)))
688
689                  !! compute the sublayer thickness
690                  zhnu = 5 * znu / zustar(ji,jj)
691
692                  !! compute gamma turb
693                  zgturb = 1._wp / vkarmn * LOG(zustar(ji,jj) * zxsiN * zetastar * zetastar / ( ABS(ff_f(ji,jj)) * zhnu )) &
694                  &      + 1._wp / ( 2 * zxsiN * zetastar ) - 1._wp / vkarmn
695
696                  !! compute gammats
697                  pgt(ji,jj) = zustar(ji,jj) / (zgturb + zgmolet)
698                  pgs(ji,jj) = zustar(ji,jj) / (zgturb + zgmoles)
699               END IF
700            END DO
701         END DO
702         CALL lbc_lnk(pgt(:,:),'T',1.)
703         CALL lbc_lnk(pgs(:,:),'T',1.)
704      END SELECT
705      CALL wrk_dealloc( jpi,jpj, zustar )
706      !
707   END SUBROUTINE sbc_isf_gammats
708
709   SUBROUTINE sbc_isf_tbl( pvarin, pvarout, cd_ptin )
710      !!----------------------------------------------------------------------
711      !!                  ***  SUBROUTINE sbc_isf_tbl  ***
712      !!
713      !! ** Purpose : compute mean T/S/U/V in the boundary layer at T- point
714      !!
715      !!----------------------------------------------------------------------
716      REAL(wp), DIMENSION(:,:,:), INTENT( in  ) :: pvarin
717      REAL(wp), DIMENSION(:,:)  , INTENT( out ) :: pvarout
718      CHARACTER(len=1),           INTENT( in  ) :: cd_ptin ! point of variable in/out
719      !
720      REAL(wp) :: ze3, zhk
721      REAL(wp), DIMENSION(:,:), POINTER :: zhisf_tbl ! thickness of the tbl
722
723      INTEGER :: ji, jj, jk                  ! loop index
724      INTEGER :: ikt, ikb                    ! top and bottom index of the tbl
725      !!----------------------------------------------------------------------
726      ! allocation
727      CALL wrk_alloc( jpi,jpj, zhisf_tbl)
728     
729      ! initialisation
730      pvarout(:,:)=0._wp
731   
732      SELECT CASE ( cd_ptin )
733      CASE ( 'U' ) ! compute U in the top boundary layer at T- point
734         DO jj = 1,jpj
735            DO ji = 1,jpi
736               ikt = miku(ji,jj) ; ikb = miku(ji,jj)
737               ! thickness of boundary layer at least the top level thickness
738               zhisf_tbl(ji,jj) = MAX(rhisf_tbl_0(ji,jj), e3u_n(ji,jj,ikt))
739
740               ! determine the deepest level influenced by the boundary layer
741               DO jk = ikt+1, mbku(ji,jj)
742                  IF ( (SUM(e3u_n(ji,jj,ikt:jk-1)) < zhisf_tbl(ji,jj)) .AND. (umask(ji,jj,jk) == 1) ) ikb = jk
743               END DO
744               zhisf_tbl(ji,jj) = MIN(zhisf_tbl(ji,jj), SUM(e3u_n(ji,jj,ikt:ikb)))  ! limit the tbl to water thickness.
745
746               ! level fully include in the ice shelf boundary layer
747               DO jk = ikt, ikb - 1
748                  ze3 = e3u_n(ji,jj,jk)
749                  pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,jk) / zhisf_tbl(ji,jj) * ze3
750               END DO
751
752               ! level partially include in ice shelf boundary layer
753               zhk = SUM( e3u_n(ji, jj, ikt:ikb - 1)) / zhisf_tbl(ji,jj)
754               pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,ikb) * (1._wp - zhk)
755            END DO
756         END DO
757         DO jj = 2,jpj
758            DO ji = 2,jpi
759               pvarout(ji,jj) = 0.5_wp * (pvarout(ji,jj) + pvarout(ji-1,jj))
760            END DO
761         END DO
762         CALL lbc_lnk(pvarout,'T',-1.)
763     
764      CASE ( 'V' ) ! compute V in the top boundary layer at T- point
765         DO jj = 1,jpj
766            DO ji = 1,jpi
767               ikt = mikv(ji,jj) ; ikb = mikv(ji,jj)
768               ! thickness of boundary layer at least the top level thickness
769               zhisf_tbl(ji,jj) = MAX(rhisf_tbl_0(ji,jj), e3v_n(ji,jj,ikt))
770
771               ! determine the deepest level influenced by the boundary layer
772               DO jk = ikt+1, mbkv(ji,jj)
773                  IF ( (SUM(e3v_n(ji,jj,ikt:jk-1)) < zhisf_tbl(ji,jj)) .AND. (vmask(ji,jj,jk) == 1) ) ikb = jk
774               END DO
775               zhisf_tbl(ji,jj) = MIN(zhisf_tbl(ji,jj), SUM(e3v_n(ji,jj,ikt:ikb)))  ! limit the tbl to water thickness.
776
777               ! level fully include in the ice shelf boundary layer
778               DO jk = ikt, ikb - 1
779                  ze3 = e3v_n(ji,jj,jk)
780                  pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,jk) / zhisf_tbl(ji,jj) * ze3
781               END DO
782
783               ! level partially include in ice shelf boundary layer
784               zhk = SUM( e3v_n(ji, jj, ikt:ikb - 1)) / zhisf_tbl(ji,jj)
785               pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,ikb) * (1._wp - zhk)
786            END DO
787         END DO
788         DO jj = 2,jpj
789            DO ji = 2,jpi
790               pvarout(ji,jj) = 0.5_wp * (pvarout(ji,jj) + pvarout(ji,jj-1))
791            END DO
792         END DO
793         CALL lbc_lnk(pvarout,'T',-1.)
794
795      CASE ( 'T' ) ! compute T in the top boundary layer at T- point
796         DO jj = 1,jpj
797            DO ji = 1,jpi
798               ikt = misfkt(ji,jj)
799               ikb = misfkb(ji,jj)
800
801               ! level fully include in the ice shelf boundary layer
802               DO jk = ikt, ikb - 1
803                  ze3 = e3t_n(ji,jj,jk)
804                  pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,jk) * r1_hisf_tbl(ji,jj) * ze3
805               END DO
806
807               ! level partially include in ice shelf boundary layer
808               zhk = SUM( e3t_n(ji, jj, ikt:ikb - 1)) * r1_hisf_tbl(ji,jj)
809               pvarout(ji,jj) = pvarout(ji,jj) + pvarin(ji,jj,ikb) * (1._wp - zhk)
810            END DO
811         END DO
812      END SELECT
813
814      ! mask mean tbl value
815      pvarout(:,:) = pvarout(:,:) * ssmask(:,:)
816
817      ! deallocation
818      CALL wrk_dealloc( jpi,jpj, zhisf_tbl )     
819      !
820   END SUBROUTINE sbc_isf_tbl
821     
822
823   SUBROUTINE sbc_isf_div( phdivn )
824      !!----------------------------------------------------------------------
825      !!                  ***  SUBROUTINE sbc_isf_div  ***
826      !!       
827      !! ** Purpose :   update the horizontal divergence with the runoff inflow
828      !!
829      !! ** Method  :   
830      !!                CAUTION : risf_tsc(:,:,jp_sal) is negative (outflow) increase the
831      !!                          divergence and expressed in m/s
832      !!
833      !! ** Action  :   phdivn   decreased by the runoff inflow
834      !!----------------------------------------------------------------------
835      REAL(wp), DIMENSION(:,:,:), INTENT( inout ) ::   phdivn   ! horizontal divergence
836      !
837      INTEGER  ::   ji, jj, jk   ! dummy loop indices
838      INTEGER  ::   ikt, ikb 
839      REAL(wp) ::   zhk
840      REAL(wp) ::   zfact     ! local scalar
841      !!----------------------------------------------------------------------
842      !
843      zfact   = 0.5_wp
844      !
845      IF(.NOT.ln_linssh ) THEN     ! need to re compute level distribution of isf fresh water
846         DO jj = 1,jpj
847            DO ji = 1,jpi
848               ikt = misfkt(ji,jj)
849               ikb = misfkt(ji,jj)
850               ! thickness of boundary layer at least the top level thickness
851               rhisf_tbl(ji,jj) = MAX(rhisf_tbl_0(ji,jj), e3t_n(ji,jj,ikt))
852
853               ! determine the deepest level influenced by the boundary layer
854               DO jk = ikt, mbkt(ji,jj)
855                  IF ( (SUM(e3t_n(ji,jj,ikt:jk-1)) .LT. rhisf_tbl(ji,jj)) .AND. (tmask(ji,jj,jk) == 1) ) ikb = jk
856               END DO
857               rhisf_tbl(ji,jj) = MIN(rhisf_tbl(ji,jj), SUM(e3t_n(ji,jj,ikt:ikb)))  ! limit the tbl to water thickness.
858               misfkb(ji,jj) = ikb                                                  ! last wet level of the tbl
859               r1_hisf_tbl(ji,jj) = 1._wp / rhisf_tbl(ji,jj)
860
861               zhk           = SUM( e3t_n(ji, jj, ikt:ikb - 1)) * r1_hisf_tbl(ji,jj)  ! proportion of tbl cover by cell from ikt to ikb - 1
862               ralpha(ji,jj) = rhisf_tbl(ji,jj) * (1._wp - zhk ) / e3t_n(ji,jj,ikb)  ! proportion of bottom cell influenced by boundary layer
863            END DO
864         END DO
865      END IF 
866      !
867      !==   ice shelf melting distributed over several levels   ==!
868      DO jj = 1,jpj
869         DO ji = 1,jpi
870               ikt = misfkt(ji,jj)
871               ikb = misfkb(ji,jj)
872               ! level fully include in the ice shelf boundary layer
873               DO jk = ikt, ikb - 1
874                  phdivn(ji,jj,jk) = phdivn(ji,jj,jk) + ( fwfisf(ji,jj) + fwfisf_b(ji,jj) ) &
875                    &              * r1_hisf_tbl(ji,jj) * r1_rau0 * zfact
876               END DO
877               ! level partially include in ice shelf boundary layer
878               phdivn(ji,jj,ikb) = phdivn(ji,jj,ikb) + ( fwfisf(ji,jj) &
879                    &            + fwfisf_b(ji,jj) ) * r1_hisf_tbl(ji,jj) * r1_rau0 * zfact * ralpha(ji,jj) 
880         END DO
881      END DO
882      !
883   END SUBROUTINE sbc_isf_div
884   !!======================================================================
885END MODULE sbcisf
Note: See TracBrowser for help on using the repository browser.