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 @ 7263

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

Add missing variable definitions for isf diagnostics

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