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/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_2 – NEMO

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_2/limsbc_2.F90 @ 7910

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

All wrk_alloc removed

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