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 branches/2017/dev_r8600_xios_read/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2017/dev_r8600_xios_read/NEMOGCM/NEMO/OPA_SRC/SBC/sbcisf.F90 @ 8612

Last change on this file since 8612 was 8612, checked in by andmirek, 6 years ago

#1953 read single file restart with XIOS

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