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/UKMO/dev_isf_flx_UKESM_r9321/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/dev_isf_flx_UKESM_r9321/NEMOGCM/NEMO/OPA_SRC/SBC/sbcisf.F90 @ 9855

Last change on this file since 9855 was 9855, checked in by mathiot, 6 years ago

update sbcisf from trunk

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