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.
limsbc_2.F90 in branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/LIM_SRC_2/limsbc_2.F90 @ 5845

Last change on this file since 5845 was 5845, checked in by gm, 8 years ago

#1613: vvl by default: suppression of domzgr_substitute.h90

  • Property svn:keywords set to Id
File size: 29.6 KB
Line 
1MODULE limsbc_2
2   !!======================================================================
3   !!                       ***  MODULE limsbc_2   ***
4   !! LIM-2 :   updates the fluxes at the ocean surface with ice-ocean fluxes
5   !!======================================================================
6   !! History :  LIM  ! 2000-01 (H. Goosse) Original code
7   !!            1.0  ! 2002-07 (C. Ethe, G. Madec) re-writing F90
8   !!            3.0  ! 2006-07 (G. Madec) surface module
9   !!            3.3  ! 2009-05 (G. Garric, C. Bricaud) addition of the lim2_evp case
10   !!             -   ! 2010-11 (G. Madec) ice-ocean stress computed at each ocean time-step
11   !!           3.3.1 ! 2011-01 (A. R. Porter, STFC Daresbury) dynamical allocation
12   !!            3.5  ! 2012-11 ((G. Madec, Y. Aksenov, A. Coward) salt and heat fluxes associated with e-p
13   !!----------------------------------------------------------------------
14#if defined key_lim2
15   !!----------------------------------------------------------------------
16   !!   'key_lim2'                                    LIM 2.0 sea-ice model
17   !!----------------------------------------------------------------------
18   !!   lim_sbc_alloc_2 : allocate the limsbc arrays
19   !!   lim_sbc_init    : initialisation
20   !!   lim_sbc_flx_2   : update mass, heat and salt fluxes at the ocean surface
21   !!   lim_sbc_tau_2   : update i- and j-stresses, and its modulus at the ocean surface
22   !!----------------------------------------------------------------------
23   USE par_oce          ! ocean parameters
24   USE phycst           ! physical constants
25   USE dom_oce          ! ocean domain
26   USE domvvl           ! ocean vertical scale factors
27   USE dom_ice_2        ! LIM-2: ice domain
28   USE ice_2            ! LIM-2: ice variables
29   USE sbc_ice          ! surface boundary condition: ice
30   USE sbc_oce          ! surface boundary condition: ocean
31   USE sbccpl
32   USE oce       , ONLY : sshn, sshb, snwice_mass, snwice_mass_b, snwice_fmass 
33   USE albedo           ! albedo parameters
34   USE lbclnk           ! ocean lateral boundary condition - MPP exchanges
35   USE lib_mpp          ! MPP library
36   USE wrk_nemo         ! work arrays
37   USE in_out_manager   ! I/O manager
38   USE iom              ! I/O library
39   USE prtctl           ! Print control
40   USE lib_fortran      ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC   lim_sbc_init_2     ! called by ice_init_2
46   PUBLIC   lim_sbc_flx_2      ! called by sbc_ice_lim_2
47   PUBLIC   lim_sbc_tau_2      ! called by sbc_ice_lim_2
48
49   REAL(wp)  ::   r1_rdtice            ! = 1. / rdt_ice
50   REAL(wp)  ::   epsi16 = 1.e-16_wp   ! constant values
51   REAL(wp)  ::   rzero  = 0._wp       !     -      -
52   REAL(wp)  ::   rone   = 1._wp       !     -      -
53   !
54   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   soce_0, sice_0   ! constant SSS and ice salinity used in levitating sea-ice case
55   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   utau_oce, vtau_oce   ! air-ocean surface i- & j-stress              [N/m2]
56   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   tmod_io              ! modulus of the ice-ocean relative velocity   [m/s]
57
58   !! * Substitutions
59#  include "vectopt_loop_substitute.h90"
60   !!----------------------------------------------------------------------
61   !! NEMO/LIM2 4.0 , UCL - NEMO Consortium (2011)
62   !! $Id$
63   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
64   !!----------------------------------------------------------------------
65CONTAINS
66
67   INTEGER FUNCTION lim_sbc_alloc_2()
68      !!-------------------------------------------------------------------
69      !!             ***  ROUTINE lim_sbc_alloc_2 ***
70      !!-------------------------------------------------------------------
71      ALLOCATE( soce_0(jpi,jpj) , utau_oce(jpi,jpj) ,                       &
72         &      sice_0(jpi,jpj) , vtau_oce(jpi,jpj) , tmod_io(jpi,jpj), STAT=lim_sbc_alloc_2)
73         !
74      IF( lk_mpp               )   CALL mpp_sum( lim_sbc_alloc_2 )
75      IF( lim_sbc_alloc_2 /= 0 )   CALL ctl_warn('lim_sbc_alloc_2: failed to allocate arrays.')
76      !
77   END FUNCTION lim_sbc_alloc_2
78
79
80   SUBROUTINE lim_sbc_flx_2( kt )
81      !!-------------------------------------------------------------------
82      !!                ***  ROUTINE lim_sbc_2 ***
83      !! 
84      !! ** Purpose :   Update surface ocean boundary condition over areas
85      !!              that are at least partially covered by sea-ice
86      !!         
87      !! ** Action  : - comput. of the momentum, heat and freshwater/salt
88      !!                fluxes at the ice-ocean interface.
89      !!              - Update the fluxes provided to the ocean
90      !!     
91      !! ** Outputs : - qsr     : sea heat flux    : solar
92      !!              - qns     : sea heat flux    : non solar (including heat content of the mass flux)
93      !!              - emp     : freshwater budget: mass flux
94      !!              - sfx     : freshwater budget: salt flux due to Freezing/Melting
95      !!              - fr_i    : ice fraction
96      !!              - tn_ice  : sea-ice surface temperature
97      !!              - alb_ice : sea-ice albedo (ln_cpl=T)
98      !!
99      !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90.
100      !!              Tartinville et al. 2001 Ocean Modelling, 3, 95-108.
101      !!---------------------------------------------------------------------
102      INTEGER, INTENT(in) ::   kt    ! number of iteration
103      !!
104      INTEGER  ::   ji, jj   ! dummy loop indices
105      INTEGER  ::   ii0, ii1, ij0, ij1         ! local integers
106      INTEGER  ::   ifvt, i1mfr, idfr, iflt    !   -       -
107      INTEGER  ::   ial, iadv, ifral, ifrdv    !   -       -
108      REAL(wp) ::   zqsr,     zqns,   zfmm     ! local scalars
109      REAL(wp) ::   zinda,    zfsalt, zemp     !   -      -
110      REAL(wp) ::   zemp_snw, zqhc,   zcd      !   -      -
111      REAL(wp) ::   zswitch                    !   -      -
112      REAL(wp), POINTER, DIMENSION(:,:)   ::   zqnsoce       ! 2D workspace
113      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalb, zalbp   ! 2D/3D workspace
114      !!---------------------------------------------------------------------
115     
116      CALL wrk_alloc( jpi, jpj, zqnsoce )
117      CALL wrk_alloc( jpi, jpj, 1, zalb, zalbp )
118
119      SELECT CASE( nn_ice_embd )                 ! levitating or embedded sea-ice option
120        CASE( 0    )   ;   zswitch = 1           ! (0) standard levitating sea-ice : salt exchange only
121        CASE( 1, 2 )   ;   zswitch = 0           ! (1) levitating sea-ice: salt and volume exchange but no pressure effect
122                                                 ! (2) embedded sea-ice : salt and volume fluxes and pressure
123      END SELECT                                 !   
124
125      !------------------------------------------!
126      !      heat flux at the ocean surface      !
127      !------------------------------------------!
128
129      zqnsoce(:,:) = qns(:,:)
130      DO jj = 1, jpj
131         DO ji = 1, jpi
132            zinda   = 1.0   - MAX( rzero , SIGN( rone, - ( 1.0 - pfrld(ji,jj) )   ) )
133            ifvt    = zinda * MAX( rzero , SIGN( rone,  - phicif(ji,jj)           ) )
134            i1mfr   = 1.0   - MAX( rzero , SIGN( rone, - ( 1.0 - frld(ji,jj) )    ) )
135            idfr    = 1.0   - MAX( rzero , SIGN( rone, frld(ji,jj) - pfrld(ji,jj) ) )
136            iflt    = zinda  * (1 - i1mfr) * (1 - ifvt )
137            ial     = ifvt   * i1mfr + ( 1 - ifvt ) * idfr
138            iadv    = ( 1  - i1mfr ) * zinda
139            ifral   = ( 1  - i1mfr * ( 1 - ial ) )   
140            ifrdv   = ( 1  - ifral * ( 1 - ial ) ) * iadv 
141
142!!$            attempt to explain the tricky flags set above....
143!!$            zinda   = 1.0 - AINT( pfrld(ji,jj) )                   ! = 0. if ice-free ocean else 1. (after ice adv, but before ice thermo)
144!!$            i1mfr   = 1.0 - AINT(  frld(ji,jj) )                   ! = 0. if ice-free ocean else 1. (after ice thermo)
145!!$
146!!$            IF( phicif(ji,jj) <= 0. ) THEN   ;   ifvt = zinda      ! = zinda if previous thermodynamic step overmelted the ice???
147!!$            ELSE                             ;   ifvt = 0.         !
148!!$            ENDIF
149!!$
150!!$            IF( frld(ji,jj) >= pfrld(ji,jj) ) THEN   ;   idfr = 0.  !   = 0. if lead fraction increases due to ice thermodynamics
151!!$            ELSE                                     ;   idfr = 1.   
152!!$            ENDIF
153!!$
154!!$            iflt    = zinda  * (1 - i1mfr) * (1 - ifvt )    !   = 1. if ice (not only snow) at previous time and ice-free ocean currently
155!!$
156!!$            ial     = ifvt   * i1mfr    +    ( 1 - ifvt ) * idfr
157!!$                    = i1mfr if ifvt = 1 i.e.
158!!$                    = idfr  if ifvt = 0
159!!$!                 snow no ice   ice         ice or nothing  lead fraction increases
160!!$!                 at previous   now           at previous
161!!$!                -> ice area increases  ???         -> ice area decreases ???
162!!$
163!!$            iadv    = ( 1  - i1mfr ) * zinda
164!!$!                     pure ocean      ice at
165!!$!                     at current      previous
166!!$!                        -> = 1. if ice disapear between previous and current
167!!$
168!!$            ifral   = ( 1  - i1mfr * ( 1 - ial ) ) 
169!!$!                            ice at     ???
170!!$!                            current         
171!!$!                         -> ???
172!!$
173!!$            ifrdv   = ( 1  - ifral * ( 1 - ial ) ) * iadv
174!!$!                                                    ice disapear
175!!$
176!!$
177
178            !   computation the solar flux at ocean surface
179            IF( ln_cpl ) THEN
180               zqsr = qsr_tot(ji,jj) + ( fstric(ji,jj) - qsr_ice(ji,jj,1) ) * ( 1.0 - pfrld(ji,jj) )
181            ELSE
182               zqsr = pfrld(ji,jj) * qsr(ji,jj)  + ( 1.  - pfrld(ji,jj) ) * fstric(ji,jj)
183            ENDIF
184            !  computation the non solar heat flux at ocean surface
185            zqns    =  - ( 1. - thcm(ji,jj) ) * zqsr                                              &   ! part of the solar energy used in leads
186               &       + iflt    * ( fscmbq(ji,jj) + ffltbif(ji,jj) )                             &
187               &       + ifral   * ( ial * qcmif(ji,jj) + (1 - ial) * qldif(ji,jj) ) * r1_rdtice  &
188               &       + ifrdv   * (       qfvbq(ji,jj) +             qdtcn(ji,jj) ) * r1_rdtice 
189
190            fsbbq(ji,jj) = ( 1.0 - ( ifvt + iflt ) ) * fscmbq(ji,jj)     ! store residual heat flux (to put into the ocean at the next time-step)
191            zqhc = ( rdq_snw(ji,jj)                                     &
192                 & + rdq_ice(ji,jj) * ( 1.- zswitch) ) * r1_rdtice       ! heat flux due to snow ( & ice heat content,
193            !                                                            !           if ice/ocean mass exchange active)
194            qsr  (ji,jj) = zqsr                                          ! solar heat flux
195            qns  (ji,jj) = zqns - fdtcn(ji,jj) + zqhc                    ! non solar heat flux
196            !
197            !                          !------------------------------------------!
198            !                          !  mass and salt flux at the ocean surface !
199            !                          !------------------------------------------!
200            !
201            ! mass flux at the ocean-atmosphere interface (open ocean fraction = leads area)
202            !                                                  ! coupled mode:
203            IF( ln_cpl ) THEN
204               zemp = + emp_tot(ji,jj)                            &     ! net mass flux over the grid cell (ice+ocean area)
205                  &   - emp_ice(ji,jj) * ( 1. - pfrld(ji,jj) )          ! minus the mass flux intercepted by sea-ice
206            ELSE
207               !                                                  ! forced  mode:
208               zemp = + emp(ji,jj)     *         frld(ji,jj)      &     ! mass flux over open ocean fraction
209                  &   - tprecip(ji,jj) * ( 1. -  frld(ji,jj) )    &     ! liquid precip. over ice reaches directly the ocean
210                  &   + sprecip(ji,jj) * ( 1. - pfrld(ji,jj) )          ! snow is intercepted by sea-ice (previous frld)
211            ENDIF
212            !
213            ! mass flux at the ocean/ice interface (sea ice fraction)
214            zemp_snw = rdm_snw(ji,jj) * r1_rdtice                    ! snow melting = pure water that enters the ocean
215            zfmm     = rdm_ice(ji,jj) * r1_rdtice                    ! Freezing minus Melting (F-M)
216
217            fmmflx(ji,jj) = zfmm                                     ! F/M mass flux save at least for biogeochemical model
218
219            ! salt flux at the ice/ocean interface (sea ice fraction) [PSU*kg/m2/s]
220            zfsalt = - sice_0(ji,jj) * zfmm                          ! F-M salt exchange
221            zcd    =   soce_0(ji,jj) * zfmm                          ! concentration/dilution term due to F-M
222            !
223            ! salt flux only       : add concentration dilution term in salt flux  and no  F-M term in volume flux
224            ! salt and mass fluxes : non concentration dilution term in salt flux  and add F-M term in volume flux
225            sfx (ji,jj) = zfsalt +                  zswitch  * zcd   ! salt flux (+ C/D if no ice/ocean mass exchange)
226            emp (ji,jj) = zemp   + zemp_snw + ( 1.- zswitch) * zfmm  ! mass flux (+ F/M mass flux if ice/ocean mass exchange)
227            !
228         END DO
229      END DO
230      !                                !------------------------------------------!
231      !                                !    mass of snow and ice per unit area    !
232      !                                !------------------------------------------!
233      IF( nn_ice_embd /= 0 ) THEN      ! embedded sea-ice (mass required)
234         snwice_mass_b(:,:) = snwice_mass(:,:)                  ! save mass from the previous ice time step
235         !                                                      ! new mass per unit area
236         snwice_mass  (:,:) = tms(:,:) * ( rhosn * hsnif(:,:) + rhoic * hicif(:,:)  ) * ( 1.0 - frld(:,:) )
237         !                                                      ! time evolution of snow+ice mass
238         snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) / rdt_ice
239      ENDIF
240
241      IF( iom_use('hflx_ice_cea' ) )   CALL iom_put( 'hflx_ice_cea', - fdtcn(:,:) )     
242      IF( iom_use('qns_io_cea'   ) )   CALL iom_put( 'qns_io_cea', qns(:,:) - zqnsoce(:,:) * pfrld(:,:) )     
243      IF( iom_use('qsr_io_cea'   ) )   CALL iom_put( 'qsr_io_cea', fstric(:,:) * (1.e0 - pfrld(:,:)) )
244
245      IF( iom_use('isnwmlt_cea'  ) )   CALL iom_put( 'isnwmlt_cea'  ,                 rdm_snw(:,:) * r1_rdtice )
246      IF( iom_use('fsal_virt_cea') )   CALL iom_put( 'fsal_virt_cea',   soce_0(:,:) * rdm_ice(:,:) * r1_rdtice )
247      IF( iom_use('fsal_real_cea') )   CALL iom_put( 'fsal_real_cea', - sice_0(:,:) * rdm_ice(:,:) * r1_rdtice )
248
249      !-----------------------------------------------!
250      !   Coupling variables                          !
251      !-----------------------------------------------!
252
253      IF( ln_cpl) THEN
254         tn_ice(:,:,1) = sist(:,:)          ! sea-ice surface temperature       
255         ht_i(:,:,1) = hicif(:,:)
256         ht_s(:,:,1) = hsnif(:,:)
257         a_i(:,:,1) = fr_i(:,:)
258         !                                  ! Computation of snow/ice and ocean albedo
259         CALL albedo_ice( tn_ice, ht_i, ht_s, zalbp, zalb )
260         alb_ice(:,:,1) =  0.5 * ( zalbp(:,:,1) + zalb (:,:,1) )   ! Ice albedo (mean clear and overcast skys)
261         IF( iom_use('icealb_cea' ) )   CALL iom_put( 'icealb_cea', alb_ice(:,:,1) * fr_i(:,:) )  ! ice albedo
262      ENDIF
263
264      IF(ln_ctl) THEN            ! control print
265         CALL prt_ctl(tab2d_1=qsr   , clinfo1=' lim_sbc: qsr    : ', tab2d_2=qns   , clinfo2=' qns     : ')
266         CALL prt_ctl(tab2d_1=emp   , clinfo1=' lim_sbc: emp    : ', tab2d_2=sfx   , clinfo2=' sfx     : ')
267         CALL prt_ctl(tab2d_1=fr_i  , clinfo1=' lim_sbc: fr_i   : ', tab2d_2=tn_ice(:,:,1), clinfo2=' tn_ice  : ')
268      ENDIF 
269      !
270      CALL wrk_dealloc( jpi, jpj, zqnsoce )
271      CALL wrk_dealloc( jpi, jpj, 1, zalb, zalbp )
272      !
273   END SUBROUTINE lim_sbc_flx_2
274
275
276   SUBROUTINE lim_sbc_tau_2( kt , pu_oce, pv_oce )
277      !!-------------------------------------------------------------------
278      !!                ***  ROUTINE lim_sbc_tau_2 ***
279      !! 
280      !! ** Purpose : Update the ocean surface stresses due to the ice
281      !!         
282      !! ** Action  : * at each ice time step (every nn_fsbc time step):
283      !!                - compute the modulus of ice-ocean relative velocity
284      !!                  at T-point (C-grid) or I-point (B-grid)
285      !!                      tmod_io = rhoco * | U_ice-U_oce |
286      !!                - update the modulus of stress at ocean surface
287      !!                      taum = frld * taum + (1-frld) * tmod_io * | U_ice-U_oce |
288      !!              * at each ocean time step (each kt):
289      !!                  compute linearized ice-ocean stresses as
290      !!                      Utau = tmod_io * | U_ice - pU_oce |
291      !!                using instantaneous current ocean velocity (usually before)
292      !!
293      !!    NB: - the averaging operator used depends on the ice dynamics grid (cp_ice_msh='I' or 'C')
294      !!        - ice-ocean rotation angle only allowed in cp_ice_msh='I' case
295      !!        - here we make an approximation: taum is only computed every ice time step
296      !!          This avoids mutiple average to pass from T -> U,V grids and next from U,V grids
297      !!          to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximation...
298      !!
299      !! ** Outputs : - utau, vtau : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes
300      !!              - taum       : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes
301      !!---------------------------------------------------------------------
302      INTEGER ,                     INTENT(in) ::   kt               ! ocean time-step index
303      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pu_oce, pv_oce   ! surface ocean currents
304      !!
305      INTEGER  ::   ji, jj   ! dummy loop indices
306      REAL(wp) ::   zfrldu, zat_u, zu_i, zutau_ice, zu_t, zmodt   ! local scalar
307      REAL(wp) ::   zfrldv, zat_v, zv_i, zvtau_ice, zv_t, zmodi   !   -      -
308      REAL(wp) ::   zsang, zumt                                   !    -         -
309      REAL(wp), POINTER, DIMENSION(:,:) ::   ztio_u, ztio_v   ! ocean stress below sea-ice
310      !!---------------------------------------------------------------------
311      !
312      CALL wrk_alloc( jpi, jpj, ztio_u, ztio_v )
313      !
314      SELECT CASE( cp_ice_msh )     
315      !                             !-----------------------!
316      CASE( 'I' )                   !  B-grid ice dynamics  !   I-point (i.e. F-point with sea-ice indexation)
317         !                          !--=--------------------!
318         !
319         IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !==  Ice time-step only  ==! (i.e. surface module time-step)
320            !
321            DO jj = 1, jpj                               !* modulus of ice-ocean relative velocity at I-point
322               DO ji = 1, jpi
323                  zu_i  = u_ice(ji,jj) - u_oce(ji,jj)                   ! ice-ocean relative velocity at I-point
324                  zv_i  = v_ice(ji,jj) - v_oce(ji,jj)
325                  tmod_io(ji,jj) = SQRT( zu_i * zu_i + zv_i * zv_i )    ! modulus of this velocity (at I-point)
326               END DO
327            END DO
328            DO jj = 1, jpjm1                             !* update the modulus of stress at ocean surface (T-point)
329               DO ji = 1, jpim1   ! NO vector opt.
330                  !                                               ! modulus of U_ice-U_oce at T-point
331                  zumt  = 0.25_wp * (  tmod_io(ji+1,jj) + tmod_io(ji  ,jj  )    &   
332                     &               + tmod_io(ji,jj+1) + tmod_io(ji+1,jj+1)  )
333                  !                                               ! update the modulus of stress at ocean surface
334                  taum(ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zumt * zumt
335               END DO
336            END DO
337            CALL lbc_lnk( taum, 'T', 1. )
338            !
339            utau_oce(:,:) = utau(:,:)                    !* save the air-ocean stresses at ice time-step
340            vtau_oce(:,:) = vtau(:,:)
341            !
342         ENDIF
343         !
344         !                                        !==  at each ocean time-step  ==!
345         !
346         !                                               !* ice/ocean stress WITH a ice-ocean rotation angle at I-point
347         DO jj = 2, jpj
348            zsang  = SIGN( 1._wp, gphif(1,jj) ) * sangvg          ! change the cosine angle sign in the SH
349            DO ji = 2, jpi    ! NO vect. opt. possible
350               ! ... ice-ocean relative velocity at I-point using instantaneous surface ocean current at u- & v-pts
351               zu_i = u_ice(ji,jj) - 0.5_wp * ( pu_oce(ji-1,jj  ) + pu_oce(ji-1,jj-1) ) * tmu(ji,jj)
352               zv_i = v_ice(ji,jj) - 0.5_wp * ( pv_oce(ji  ,jj-1) + pv_oce(ji-1,jj-1) ) * tmu(ji,jj)
353               ! ... components of stress with a ice-ocean rotation angle
354               zmodi = rhoco * tmod_io(ji,jj)                     
355               ztio_u(ji,jj) = zmodi * ( cangvg * zu_i - zsang * zv_i )
356               ztio_v(ji,jj) = zmodi * ( cangvg * zv_i + zsang * zu_i )
357            END DO
358         END DO
359         !                                               !* surface ocean stresses at u- and v-points
360         DO jj = 2, jpjm1
361            DO ji = 2, jpim1   ! NO vector opt.
362               !                                   ! ice-ocean stress at U and V-points  (from I-point values)
363               zutau_ice  = 0.5_wp * ( ztio_u(ji+1,jj) + ztio_u(ji+1,jj+1) )
364               zvtau_ice  = 0.5_wp * ( ztio_v(ji,jj+1) + ztio_v(ji+1,jj+1) )
365               !                                   ! open-ocean (lead) fraction at U- & V-points (from T-point values)
366               zfrldu = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) )
367               zfrldv = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) )
368               !                                   ! update the surface ocean stress (ice-cover wheighted)
369               utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice
370               vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice
371            END DO
372         END DO
373         CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )     ! lateral boundary condition
374         !
375         !
376         !                          !-----------------------!
377      CASE( 'C' )                   !  C-grid ice dynamics  !   U & V-points (same as in the ocean)
378         !                          !--=--------------------!
379         !
380         IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !==  Ice time-step only  ==! (i.e. surface module time-step)
381            !
382            DO jj = 2, jpjm1                          !* modulus of the ice-ocean velocity at T-point
383               DO ji = fs_2, fs_jpim1
384                  zu_t  = u_ice(ji,jj) + u_ice(ji-1,jj) - u_oce(ji,jj) - u_oce(ji-1,jj)   ! 2*(U_ice-U_oce) at T-point
385                  zv_t  = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1)     
386                  zmodt =  0.25_wp * (  zu_t * zu_t + zv_t * zv_t  )                      ! |U_ice-U_oce|^2
387                  !                                               ! update the modulus of stress at ocean surface
388                  taum   (ji,jj) = frld(ji,jj) * taum(ji,jj) + ( 1._wp - frld(ji,jj) ) * rhoco * zmodt
389                  tmod_io(ji,jj) = SQRT( zmodt ) * rhoco          ! rhoco*|Uice-Uoce|
390               END DO
391            END DO
392            CALL lbc_lnk( taum, 'T', 1. )   ;   CALL lbc_lnk( tmod_io, 'T', 1. )
393            !
394            utau_oce(:,:) = utau(:,:)                 !* save the air-ocean stresses at ice time-step
395            vtau_oce(:,:) = vtau(:,:)
396            !
397         ENDIF
398         !
399         !                                        !==  at each ocean time-step  ==!
400         !
401         DO jj = 2, jpjm1                             !* ice stress over ocean WITHOUT a ice-ocean rotation angle
402            DO ji = fs_2, fs_jpim1
403               !                                            ! ocean area at u- & v-points
404               zfrldu  = 0.5_wp * ( frld(ji,jj) + frld(ji+1,jj) )
405               zfrldv  = 0.5_wp * ( frld(ji,jj) + frld(ji,jj+1) )
406               !                                            ! quadratic drag formulation without rotation
407               !                                            ! using instantaneous surface ocean current
408               zutau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - pu_oce(ji,jj) )
409               zvtau_ice = 0.5 * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - pv_oce(ji,jj) )
410               !                                            ! update the surface ocean stress (ice-cover wheighted)
411               utau(ji,jj) = zfrldu * utau_oce(ji,jj) + ( 1._wp - zfrldu ) * zutau_ice
412               vtau(ji,jj) = zfrldv * vtau_oce(ji,jj) + ( 1._wp - zfrldv ) * zvtau_ice
413            END DO
414         END DO
415         CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )   ! lateral boundary condition
416         !
417      END SELECT
418
419      IF(ln_ctl)   CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau   : ', mask1=umask,   &
420         &                       tab2d_2=vtau, clinfo2=' vtau    : '        , mask2=vmask )
421     
422      CALL wrk_dealloc( jpi, jpj, ztio_u, ztio_v )
423      !
424   END SUBROUTINE lim_sbc_tau_2
425
426
427   SUBROUTINE lim_sbc_init_2
428      !!-------------------------------------------------------------------
429      !!                  ***  ROUTINE lim_sbc_init  ***
430      !!             
431      !! ** Purpose : Preparation of the file ice_evolu for the output of
432      !!      the temporal evolution of key variables
433      !!
434      !! ** input   : Namelist namicedia
435      !!-------------------------------------------------------------------
436      !
437      INTEGER :: jk           ! local integer
438      !
439      IF(lwp) WRITE(numout,*)
440      IF(lwp) WRITE(numout,*) 'lim_sbc_init_2 : LIM-2 sea-ice - surface boundary condition'
441      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~   '
442
443      !                                      ! allocate lim_sbc arrays
444      IF( lim_sbc_alloc_2() /= 0 )   CALL ctl_stop( 'STOP', 'lim_sbc_flx_2 : unable to allocate arrays' )
445      !
446      r1_rdtice = 1._wp / rdt_ice
447      !
448      soce_0(:,:) = soce                     ! constant SSS and ice salinity used in levitating sea-ice case
449      sice_0(:,:) = sice
450      !
451      IF( cp_cfg == "orca" ) THEN            ! decrease ocean & ice reference salinities in the Baltic sea
452         WHERE( 14._wp <= glamt(:,:) .AND. glamt(:,:) <= 32._wp .AND.   &
453            &   54._wp <= gphit(:,:) .AND. gphit(:,:) <= 66._wp         ) 
454            soce_0(:,:) = 4._wp
455            sice_0(:,:) = 2._wp
456         END WHERE
457      ENDIF
458      !                                      ! embedded sea ice
459      IF( nn_ice_embd /= 0 ) THEN            ! mass exchanges between ice and ocean (case 1 or 2) set the snow+ice mass
460         snwice_mass  (:,:) = tms(:,:) * ( rhosn * hsnif(:,:) + rhoic * hicif(:,:)  ) * ( 1.0 - frld(:,:) )
461         snwice_mass_b(:,:) = snwice_mass(:,:)
462      ELSE
463         snwice_mass  (:,:) = 0.e0           ! no mass exchanges
464         snwice_mass_b(:,:) = 0.e0           ! no mass exchanges
465         snwice_fmass (:,:) = 0.e0           ! no mass exchanges
466      ENDIF
467      IF( nn_ice_embd == 2 .AND.          &  ! full embedment (case 2) & no restart :
468         &   .NOT.ln_rstart ) THEN           ! deplete the initial ssh below sea-ice area
469         sshn(:,:) = sshn(:,:) - snwice_mass(:,:) * r1_rau0
470         sshb(:,:) = sshb(:,:) - snwice_mass(:,:) * r1_rau0
471         do jk = 1,jpkm1                     ! adjust initial vertical scale factors
472          e3t_n(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshn(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
473          e3t_b(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshb(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
474         end do
475         e3t_a(:,:,:) = e3t_b(:,:,:)
476         ! Reconstruction of all vertical scale factors at now and before time steps
477         ! =============================================================================
478         ! Horizontal scale factor interpolations
479         ! --------------------------------------
480         CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' )
481         CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' )
482         CALL dom_vvl_interpol( e3t_n(:,:,:), e3u_n(:,:,:), 'U' )
483         CALL dom_vvl_interpol( e3t_n(:,:,:), e3v_n(:,:,:), 'V' )
484         CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F' )
485         ! Vertical scale factor interpolations
486         ! ------------------------------------
487         CALL dom_vvl_interpol( e3t_n(:,:,:), e3w_n (:,:,:), 'W'  )
488         CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' )
489         CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' )
490         CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' )
491         CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' )
492         ! t- and w- points depth
493         ! ----------------------
494         gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1)
495         gdepw_n(:,:,1) = 0.0_wp
496         gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:)
497         DO jk = 2, jpk
498            gdept_n(:,:,jk) = gdept_n(:,:,jk-1) + e3w_n(:,:,jk)
499            gdepw_n(:,:,jk) = gdepw_n(:,:,jk-1) + e3t_n(:,:,jk-1)
500            gde3w_n(:,:,jk) = gdept_n(:,:,jk  ) - sshn   (:,:)
501         END DO
502      ENDIF
503      !
504   END SUBROUTINE lim_sbc_init_2
505
506#else
507   !!----------------------------------------------------------------------
508   !!   Default option         Empty module        NO LIM 2.0 sea-ice model
509   !!----------------------------------------------------------------------
510#endif 
511
512   !!======================================================================
513END MODULE limsbc_2
Note: See TracBrowser for help on using the repository browser.