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/2013/dev_MERGE_2013/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/LIM_SRC_2/limsbc_2.F90 @ 4317

Last change on this file since 4317 was 4306, checked in by cetlod, 10 years ago

dev_MERGE_2013 : merge in the solar mean flux branch from MERCATOR, see ticket #1187

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