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/2016/dev_r7233_CMIP6_diags_trunk_version/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2016/dev_r7233_CMIP6_diags_trunk_version/NEMOGCM/NEMO/OPA_SRC/SBC/sbcisf.F90 @ 7240

Last change on this file since 7240 was 7240, checked in by timgraham, 7 years ago

Added vertically integrated zonal mass transport, temperature and salinity
Added 3D ice shelf budget terms

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