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/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/SBC/sbcisf.F90 @ 5944

Last change on this file since 5944 was 5944, checked in by mathiot, 8 years ago

ISF: change related to reviewers comments

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