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.F90 in trunk/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: trunk/NEMOGCM/NEMO/LIM_SRC_3/limsbc.F90 @ 4765

Last change on this file since 4765 was 4765, checked in by rblod, 10 years ago

Compilation issue, see ticket #1379

  • Property svn:keywords set to Id
File size: 22.3 KB
Line 
1MODULE limsbc
2   !!======================================================================
3   !!                       ***  MODULE limsbc   ***
4   !!           computation of the flux at the sea ice/ocean interface
5   !!======================================================================
6   !! History :   -   ! 2006-07 (M. Vancoppelle)  LIM3 original code
7   !!            3.0  ! 2008-03 (C. Tallandier)  surface module
8   !!             -   ! 2008-04 (C. Tallandier)  split in 2 + new ice-ocean coupling
9   !!            3.3  ! 2010-05 (G. Madec) decrease ocean & ice reference salinities in the Baltic sea
10   !!                 !                  + simplification of the ice-ocean stress calculation
11   !!            3.4  ! 2011-02 (G. Madec) dynamical allocation
12   !!             -   ! 2012    (D. Iovino) salt flux change
13   !!             -   ! 2012-05 (C. Rousset) add penetration solar flux
14   !!            3.5  ! 2012-10 (A. Coward, G. Madec) salt fluxes ; ice+snow mass
15   !!----------------------------------------------------------------------
16#if defined key_lim3
17   !!----------------------------------------------------------------------
18   !!   'key_lim3'                                    LIM 3.0 sea-ice model
19   !!----------------------------------------------------------------------
20   !!   lim_sbc_alloc : allocate the limsbc arrays
21   !!   lim_sbc_init  : initialisation
22   !!   lim_sbc_flx   : updates mass, heat and salt fluxes at the ocean surface
23   !!   lim_sbc_tau   : update i- and j-stresses, and its modulus at the ocean surface
24   !!----------------------------------------------------------------------
25   USE par_oce          ! ocean parameters
26   USE phycst           ! physical constants
27   USE par_ice          ! ice parameters
28   USE dom_oce          ! ocean domain
29   USE dom_ice,    ONLY : tms, area
30   USE ice              ! LIM sea-ice variables
31   USE sbc_ice          ! Surface boundary condition: sea-ice fields
32   USE sbc_oce          ! Surface boundary condition: ocean fields
33   USE sbccpl
34   USE cpl_oasis3, ONLY : lk_cpl
35   USE oce       , ONLY : iatte, oatte, sshn, sshb, snwice_mass, snwice_mass_b, snwice_fmass
36   USE albedo           ! albedo parameters
37   USE lbclnk           ! ocean lateral boundary condition - MPP exchanges
38   USE lib_mpp          ! MPP library
39   USE wrk_nemo         ! work arrays
40   USE in_out_manager   ! I/O manager
41   USE prtctl           ! Print control
42   USE lib_fortran      ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
43   USE traqsr           ! clem: add penetration of solar flux into the calculation of heat budget
44   USE iom
45   USE domvvl           ! Variable volume
46
47   IMPLICIT NONE
48   PRIVATE
49
50   PUBLIC   lim_sbc_init   ! called by ice_init
51   PUBLIC   lim_sbc_flx    ! called by sbc_ice_lim
52   PUBLIC   lim_sbc_tau    ! called by sbc_ice_lim
53
54   REAL(wp)  ::   epsi10 = 1.e-10   !
55   REAL(wp)  ::   epsi20 = 1.e-20   !
56
57   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   utau_oce, vtau_oce   ! air-ocean surface i- & j-stress     [N/m2]
58   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   tmod_io              ! modulus of the ice-ocean velocity   [m/s]
59   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   soce_0  , sice_0     ! cst SSS and ice salinity (levitating sea-ice)
60
61   !! * Substitutions
62#  include "vectopt_loop_substitute.h90"
63#  include "domzgr_substitute.h90"
64   !!----------------------------------------------------------------------
65   !! NEMO/LIM3 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()
72      !!-------------------------------------------------------------------
73      !!             ***  ROUTINE lim_sbc_alloc ***
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)
77         !
78      IF( lk_mpp             )   CALL mpp_sum( lim_sbc_alloc )
79      IF( lim_sbc_alloc /= 0 )   CALL ctl_warn('lim_sbc_alloc: failed to allocate arrays')
80   END FUNCTION lim_sbc_alloc
81
82
83   SUBROUTINE lim_sbc_flx( kt )
84      !!-------------------------------------------------------------------
85      !!                ***  ROUTINE lim_sbc_flx ***
86      !! 
87      !! ** Purpose :   Update the surface ocean boundary condition for heat
88      !!              salt and mass over areas where sea-ice is non-zero
89      !!         
90      !! ** Action  : - computes the heat and freshwater/salt fluxes
91      !!              at the ice-ocean interface.
92      !!              - Update the ocean sbc
93      !!     
94      !! ** Outputs : - qsr     : sea heat flux:     solar
95      !!              - qns     : sea heat flux: non solar
96      !!              - emp     : freshwater budget: volume flux
97      !!              - sfx     : salt flux
98      !!              - fr_i    : ice fraction
99      !!              - tn_ice  : sea-ice surface temperature
100      !!              - alb_ice : sea-ice alberdo (lk_cpl=T)
101      !!
102      !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90.
103      !!              Tartinville et al. 2001 Ocean Modelling, 3, 95-108.
104      !!---------------------------------------------------------------------
105      INTEGER, INTENT(in) ::   kt    ! number of iteration
106      !
107      INTEGER  ::   ji, jj, jl, jk           ! dummy loop indices
108      REAL(wp) ::   zinda, zemp      ! local scalars
109      REAL(wp) ::   zf_mass         ! Heat flux associated with mass exchange ice->ocean (W.m-2)
110      REAL(wp) ::   zfcm1           ! New solar flux received by the ocean
111      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalb, zalbp     ! 2D/3D workspace
112      !!---------------------------------------------------------------------
113     
114      IF( lk_cpl )   CALL wrk_alloc( jpi, jpj, jpl, zalb, zalbp )
115
116      ! make calls for heat fluxes before it is modified
117      CALL iom_put( "qsr_oce" , qsr(:,:) * pfrld(:,:) )   !     solar flux at ocean surface
118      CALL iom_put( "qns_oce" , qns(:,:) * pfrld(:,:) )   ! non-solar flux at ocean surface
119      CALL iom_put( "qsr_ice" , SUM( qsr_ice(:,:,:) * old_a_i(:,:,:), dim=3 ) )  !     solar flux at ice surface
120      CALL iom_put( "qns_ice" , SUM( qns_ice(:,:,:) * old_a_i(:,:,:), dim=3 ) )  ! non-solar flux at ice surface
121      CALL iom_put( "qtr_ice" , SUM( ftr_ice(:,:,:) * old_a_i(:,:,:), dim=3 ) )  !     solar flux transmitted thru ice
122      CALL iom_put( "qt_oce"  , ( qsr(:,:) + qns(:,:) ) * pfrld(:,:) ) 
123      CALL iom_put( "qt_ice"  , SUM( ( qns_ice(:,:,:) + qsr_ice(:,:,:) ) * old_a_i(:,:,:), dim=3 ) )
124
125      ! pfrld is the lead fraction at the previous time step (actually between TRP and THD)
126      DO jj = 1, jpj
127         DO ji = 1, jpi
128
129            !------------------------------------------!
130            !      heat flux at the ocean surface      !
131            !------------------------------------------!
132            zinda   = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ( 1._wp - pfrld(ji,jj) ) ) ) ! 1 if ice
133
134            ! Solar heat flux reaching the ocean = zfcm1 (W.m-2)
135            !---------------------------------------------------
136            IF( lk_cpl ) THEN ! be carfeful: not been tested yet
137               ! original line
138               zfcm1 = qsr_tot(ji,jj)
139               !!!zfcm1 = qsr_tot(ji,jj) + ftr_ice(ji,jj) * ( 1._wp - pfrld(ji,jj) ) / ( 1._wp - zinda + zinda * iatte(ji,jj) )
140               DO jl = 1, jpl
141                  zfcm1 = zfcm1 - ( qsr_ice(ji,jj,jl) - ftr_ice(ji,jj,jl) ) * old_a_i(ji,jj,jl)
142               END DO
143            ELSE
144               !!!zfcm1   = pfrld(ji,jj) * qsr(ji,jj)  + &
145               !!!     &    ( 1._wp - pfrld(ji,jj) ) * ftr_ice(ji,jj) / ( 1._wp - zinda + zinda * iatte(ji,jj) )
146               zfcm1   = pfrld(ji,jj) * qsr(ji,jj)
147               DO jl = 1, jpl
148                  zfcm1   = zfcm1 + old_a_i(ji,jj,jl) * ftr_ice(ji,jj,jl)
149               END DO
150            ENDIF
151
152            ! Total heat flux reaching the ocean = hfx_out (W.m-2)
153            !---------------------------------------------------
154            zf_mass        = hfx_thd(ji,jj) + hfx_dyn(ji,jj) + hfx_res(ji,jj) ! heat flux from snow is 0 (T=0 degC)
155            hfx_out(ji,jj) = hfx_out(ji,jj) + zf_mass + zfcm1
156
157            ! New qsr and qns used to compute the oceanic heat flux at the next time step
158            !---------------------------------------------------
159            qsr(ji,jj) = zfcm1                                     
160            qns(ji,jj) = hfx_out(ji,jj) - zfcm1             
161
162            !------------------------------------------!
163            !      mass flux at the ocean surface      !
164            !------------------------------------------!
165            !  case of realistic freshwater flux (Tartinville et al., 2001) (presently ACTIVATED)
166            !  -------------------------------------------------------------------------------------
167            !  The idea of this approach is that the system that we consider is the ICE-OCEAN system
168            !  Thus  FW  flux  =  External ( E-P+snow melt)
169            !       Salt flux  =  Exchanges in the ice-ocean system then converted into FW
170            !                     Associated to Ice formation AND Ice melting
171            !                     Even if i see Ice melting as a FW and SALT flux
172            !       
173            !  computing freshwater exchanges at the ice/ocean interface
174            IF( lk_cpl ) THEN
175               zemp = - emp_tot(ji,jj) + emp_ice(ji,jj) * ( 1. - pfrld(ji,jj) )    &   !
176                  &   + wfx_snw(ji,jj)
177            ELSE
178               zemp =   emp(ji,jj)     *           pfrld(ji,jj)            &   ! evaporation over oceanic fraction
179                  &   - tprecip(ji,jj) * ( 1._wp - pfrld(ji,jj) )          &   ! all precipitation reach the ocean
180                  &   + sprecip(ji,jj) * ( 1._wp - pfrld(ji,jj)**betas )       ! except solid precip intercepted by sea-ice
181            ENDIF
182
183            ! mass flux from ice/ocean
184            wfx_ice(ji,jj) = wfx_bog(ji,jj) + wfx_bom(ji,jj) + wfx_sum(ji,jj) + wfx_sni(ji,jj)   &
185                           + wfx_opw(ji,jj) + wfx_dyn(ji,jj) + wfx_res(ji,jj)
186
187            ! mass flux at the ocean/ice interface
188            fmmflx(ji,jj) = - wfx_ice(ji,jj) * rdt_ice                   ! F/M mass flux save at least for biogeochemical model
189            emp(ji,jj)    = zemp - wfx_ice(ji,jj) - wfx_snw(ji,jj) - wfx_sub(ji,jj)   ! mass flux + F/M mass flux (always ice/ocean mass exchange)
190           
191         END DO
192      END DO
193
194      !------------------------------------------!
195      !      salt flux at the ocean surface      !
196      !------------------------------------------!
197      sfx(:,:) = sfx_bog(:,:) + sfx_bom(:,:) + sfx_sum(:,:) + sfx_sni(:,:) + sfx_opw(:,:)   &
198         &     + sfx_res(:,:) + sfx_dyn(:,:) + sfx_bri(:,:)
199
200      !-------------------------------------------------------------!
201      !   mass of snow and ice per unit area for embedded sea-ice   !
202      !-------------------------------------------------------------!
203      IF( nn_ice_embd /= 0 ) THEN
204         ! save mass from the previous ice time step
205         snwice_mass_b(:,:) = snwice_mass(:,:)                 
206         ! new mass per unit area
207         snwice_mass  (:,:) = tms(:,:) * ( rhosn * vt_s(:,:) + rhoic * vt_i(:,:)  ) 
208         ! time evolution of snow+ice mass
209         snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) * r1_rdtice
210      ENDIF
211
212      !-----------------------------------------------!
213      !   Storing the transmitted variables           !
214      !-----------------------------------------------!
215      fr_i  (:,:)   = at_i(:,:)             ! Sea-ice fraction           
216      tn_ice(:,:,:) = t_su(:,:,:)           ! Ice surface temperature                     
217
218      !------------------------------------------------!
219      !    Computation of snow/ice and ocean albedo    !
220      !------------------------------------------------!
221      IF( lk_cpl ) THEN          ! coupled case
222         CALL albedo_ice( t_su, ht_i, ht_s, zalbp, zalb )                  ! snow/ice albedo
223         alb_ice(:,:,:) =  0.5_wp * zalbp(:,:,:) + 0.5_wp * zalb (:,:,:)   ! Ice albedo (mean clear and overcast skys)
224      ENDIF
225
226
227      IF(ln_ctl) THEN
228         CALL prt_ctl( tab2d_1=qsr   , clinfo1=' lim_sbc: qsr    : ', tab2d_2=qns , clinfo2=' qns     : ' )
229         CALL prt_ctl( tab2d_1=emp   , clinfo1=' lim_sbc: emp    : ', tab2d_2=sfx , clinfo2=' sfx     : ' )
230         CALL prt_ctl( tab2d_1=fr_i  , clinfo1=' lim_sbc: fr_i   : ' )
231         CALL prt_ctl( tab3d_1=tn_ice, clinfo1=' lim_sbc: tn_ice : ', kdim=jpl )
232      ENDIF
233      !
234      IF( lk_cpl )   CALL wrk_dealloc( jpi, jpj, jpl, zalb, zalbp )
235      !
236   END SUBROUTINE lim_sbc_flx
237
238
239   SUBROUTINE lim_sbc_tau( kt , pu_oce, pv_oce )
240      !!-------------------------------------------------------------------
241      !!                ***  ROUTINE lim_sbc_tau ***
242      !! 
243      !! ** Purpose : Update the ocean surface stresses due to the ice
244      !!         
245      !! ** Action  : * at each ice time step (every nn_fsbc time step):
246      !!                - compute the modulus of ice-ocean relative velocity
247      !!                  (*rho*Cd) at T-point (C-grid) or I-point (B-grid)
248      !!                      tmod_io = rhoco * | U_ice-U_oce |
249      !!                - update the modulus of stress at ocean surface
250      !!                      taum = frld * taum + (1-frld) * tmod_io * | U_ice-U_oce |
251      !!              * at each ocean time step (every kt):
252      !!                  compute linearized ice-ocean stresses as
253      !!                      Utau = tmod_io * | U_ice - pU_oce |
254      !!                using instantaneous current ocean velocity (usually before)
255      !!
256      !!    NB: - ice-ocean rotation angle no more allowed
257      !!        - here we make an approximation: taum is only computed every ice time step
258      !!          This avoids mutiple average to pass from T -> U,V grids and next from U,V grids
259      !!          to T grid. taum is used in TKE and GLS, which should not be too sensitive to this approximaton...
260      !!
261      !! ** Outputs : - utau, vtau   : surface ocean i- and j-stress (u- & v-pts) updated with ice-ocean fluxes
262      !!              - taum         : modulus of the surface ocean stress (T-point) updated with ice-ocean fluxes
263      !!---------------------------------------------------------------------
264      INTEGER ,                     INTENT(in) ::   kt               ! ocean time-step index
265      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pu_oce, pv_oce   ! surface ocean currents
266      !!
267      INTEGER  ::   ji, jj   ! dummy loop indices
268      REAL(wp) ::   zat_u, zutau_ice, zu_t, zmodt   ! local scalar
269      REAL(wp) ::   zat_v, zvtau_ice, zv_t          !   -      -
270      !!---------------------------------------------------------------------
271      !
272      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !==  Ice time-step only  ==!   (i.e. surface module time-step)
273!CDIR NOVERRCHK
274         DO jj = 2, jpjm1                             !* update the modulus of stress at ocean surface (T-point)
275!CDIR NOVERRCHK
276            DO ji = fs_2, fs_jpim1
277               !                                               ! 2*(U_ice-U_oce) at T-point
278               zu_t = u_ice(ji,jj) + u_ice(ji-1,jj) - u_oce(ji,jj) - u_oce(ji-1,jj)   
279               zv_t = v_ice(ji,jj) + v_ice(ji,jj-1) - v_oce(ji,jj) - v_oce(ji,jj-1) 
280               !                                              ! |U_ice-U_oce|^2
281               zmodt =  0.25_wp * (  zu_t * zu_t + zv_t * zv_t  )
282               !                                               ! update the ocean stress modulus
283               taum(ji,jj) = ( 1._wp - at_i(ji,jj) ) * taum(ji,jj) + at_i(ji,jj) * rhoco * zmodt
284               tmod_io(ji,jj) = rhoco * SQRT( zmodt )          ! rhoco * |U_ice-U_oce| at T-point
285            END DO
286         END DO
287         CALL lbc_lnk( taum, 'T', 1. )   ;   CALL lbc_lnk( tmod_io, 'T', 1. )
288         !
289         utau_oce(:,:) = utau(:,:)                    !* save the air-ocean stresses at ice time-step
290         vtau_oce(:,:) = vtau(:,:)
291         !
292      ENDIF
293      !
294      !                                      !==  every ocean time-step  ==!
295      !
296      DO jj = 2, jpjm1                                !* update the stress WITHOUT a ice-ocean rotation angle
297         DO ji = fs_2, fs_jpim1   ! Vect. Opt.
298            zat_u  = ( at_i(ji,jj) + at_i(ji+1,jj) ) * 0.5_wp   ! ice area at u and V-points
299            zat_v  = ( at_i(ji,jj) + at_i(ji,jj+1) ) * 0.5_wp
300            !                                                   ! linearized quadratic drag formulation
301            zutau_ice   = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji+1,jj) ) * ( u_ice(ji,jj) - pu_oce(ji,jj) )
302            zvtau_ice   = 0.5_wp * ( tmod_io(ji,jj) + tmod_io(ji,jj+1) ) * ( v_ice(ji,jj) - pv_oce(ji,jj) )
303            !                                                   ! stresses at the ocean surface
304            utau(ji,jj) = ( 1._wp - zat_u ) * utau_oce(ji,jj) + zat_u * zutau_ice
305            vtau(ji,jj) = ( 1._wp - zat_v ) * vtau_oce(ji,jj) + zat_v * zvtau_ice
306         END DO
307      END DO
308      CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )   ! lateral boundary condition
309      !
310      IF(ln_ctl)   CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau   : ', mask1=umask,   &
311         &                       tab2d_2=vtau, clinfo2=' vtau    : '        , mask2=vmask )
312     
313   END SUBROUTINE lim_sbc_tau
314
315
316   SUBROUTINE lim_sbc_init
317      !!-------------------------------------------------------------------
318      !!                  ***  ROUTINE lim_sbc_init  ***
319      !!             
320      !! ** Purpose : Preparation of the file ice_evolu for the output of
321      !!      the temporal evolution of key variables
322      !!
323      !! ** input   : Namelist namicedia
324      !!-------------------------------------------------------------------
325      REAL(wp) :: zsum, zarea
326      !
327      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
328      REAL(wp) ::   zcoefu, zcoefv, zcoeff          ! local scalar
329      IF(lwp) WRITE(numout,*)
330      IF(lwp) WRITE(numout,*) 'lim_sbc_init : LIM-3 sea-ice - surface boundary condition'
331      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   '
332
333      !                                      ! allocate lim_sbc array
334      IF( lim_sbc_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'lim_sbc_init : unable to allocate standard arrays' )
335      !
336      soce_0(:,:) = soce                     ! constant SSS and ice salinity used in levitating sea-ice case
337      sice_0(:,:) = sice
338      !
339      IF( cp_cfg == "orca" ) THEN            ! decrease ocean & ice reference salinities in the Baltic sea
340         WHERE( 14._wp <= glamt(:,:) .AND. glamt(:,:) <= 32._wp .AND.   &
341            &   54._wp <= gphit(:,:) .AND. gphit(:,:) <= 66._wp         ) 
342            soce_0(:,:) = 4._wp
343            sice_0(:,:) = 2._wp
344         END WHERE
345      ENDIF
346      ! clem modif
347      IF( .NOT. ln_rstart ) THEN
348         iatte(:,:) = 1._wp
349         oatte(:,:) = 1._wp
350      ENDIF
351      !
352      ! clem: snwice_mass in the restart file now
353      IF( .NOT. ln_rstart ) THEN
354         !                                      ! embedded sea ice
355         IF( nn_ice_embd /= 0 ) THEN            ! mass exchanges between ice and ocean (case 1 or 2) set the snow+ice mass
356            snwice_mass  (:,:) = tms(:,:) * ( rhosn * vt_s(:,:) + rhoic * vt_i(:,:)  )
357            snwice_mass_b(:,:) = snwice_mass(:,:)
358         ELSE
359            snwice_mass  (:,:) = 0.0_wp         ! no mass exchanges
360            snwice_mass_b(:,:) = 0.0_wp         ! no mass exchanges
361         ENDIF
362         IF( nn_ice_embd == 2 ) THEN            ! full embedment (case 2) deplete the initial ssh below sea-ice area
363            sshn(:,:) = sshn(:,:) - snwice_mass(:,:) * r1_rau0
364            sshb(:,:) = sshb(:,:) - snwice_mass(:,:) * r1_rau0
365#if defined key_vvl           
366           ! key_vvl necessary? clem: yes for compilation purpose
367            DO jk = 1,jpkm1                     ! adjust initial vertical scale factors
368               fse3t_n(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshn(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
369               fse3t_b(:,:,jk) = e3t_0(:,:,jk)*( 1._wp + sshb(:,:)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
370            ENDDO
371            fse3t_a(:,:,:) = fse3t_b(:,:,:)
372            ! Reconstruction of all vertical scale factors at now and before time
373            ! steps
374            ! =============================================================================
375            ! Horizontal scale factor interpolations
376            ! --------------------------------------
377            CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3u_b(:,:,:), 'U' )
378            CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3v_b(:,:,:), 'V' )
379            CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' )
380            CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' )
381            CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n(:,:,:), 'F' )
382            ! Vertical scale factor interpolations
383            ! ------------------------------------
384            CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W'  )
385            CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' )
386            CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' )
387            CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' )
388            CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' )
389            ! t- and w- points depth
390            ! ----------------------
391            fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1)
392            fsdepw_n(:,:,1) = 0.0_wp
393            fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:)
394            DO jk = 2, jpk
395               fsdept_n(:,:,jk) = fsdept_n(:,:,jk-1) + fse3w_n(:,:,jk)
396               fsdepw_n(:,:,jk) = fsdepw_n(:,:,jk-1) + fse3t_n(:,:,jk-1)
397               fsde3w_n(:,:,jk) = fsdept_n(:,:,jk  ) - sshn   (:,:)
398            END DO
399#endif
400         ENDIF
401      ENDIF ! .NOT. ln_rstart
402      !
403
404   END SUBROUTINE lim_sbc_init
405
406#else
407   !!----------------------------------------------------------------------
408   !!   Default option :        Dummy module       NO LIM 3.0 sea-ice model
409   !!----------------------------------------------------------------------
410CONTAINS
411   SUBROUTINE lim_sbc           ! Dummy routine
412   END SUBROUTINE lim_sbc
413#endif 
414
415   !!======================================================================
416END MODULE limsbc
Note: See TracBrowser for help on using the repository browser.