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 branches/dev_003_CPL/NEMO/LIM_SRC_3 – NEMO

source: branches/dev_003_CPL/NEMO/LIM_SRC_3/limsbc.F90 @ 991

Last change on this file since 991 was 991, checked in by smasson, 16 years ago

dev_003_CPL: preliminary draft (not working), see ticket #155

File size: 24.2 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   !!----------------------------------------------------------------------
10#if defined key_lim3
11   !!----------------------------------------------------------------------
12   !!   'key_lim3'                                    LIM 3.0 sea-ice model
13   !!----------------------------------------------------------------------
14   !!   lim_sbc  : flux at the ice / ocean interface
15   !!----------------------------------------------------------------------
16   USE oce
17   USE par_oce          ! ocean parameters
18   USE par_ice          ! ice parameters
19   USE dom_oce          ! ocean domain
20   USE sbc_ice          ! Surface boundary condition: sea-ice fields
21   USE sbc_oce          ! Surface boundary condition: ocean fields
22   USE phycst           ! physical constants
23   USE ocfzpt           ! surface ocean freezing point
24   USE ice_oce          ! sea-ice variable
25   USE ice              ! LIM sea-ice variables
26   USE iceini           ! ???
27   USE dynspg_oce       ! choice of the surface pressure gradient scheme
28
29   USE lbclnk           ! ocean lateral boundary condition
30   USE in_out_manager   ! I/O manager
31   USE albedo           ! albedo parameters
32   USE prtctl           ! Print control
33
34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC   lim_sbc_flx   ! called by sbc_ice_lim
38   PUBLIC   lim_sbc_tau   ! called by sbc_ice_lim
39
40   REAL(wp)  ::   epsi16 = 1.e-16  ! constant values
41   REAL(wp)  ::   rzero  = 0.e0   
42   REAL(wp)  ::   rone   = 1.e0
43
44   REAL(wp), DIMENSION(jpi,jpj) ::   utau_oce   !: atm.-ocean surface i-stress (ocean referential)     [N/m2]
45   REAL(wp), DIMENSION(jpi,jpj) ::   vtau_oce   !: atm.-ocean surface j-stress (ocean referential)     [N/m2]
46   !! * Substitutions
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/LIM  3.0 ,  UCL-LOCEAN-IPSL (2008)
50   !! $Id: $
51   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
52   !!----------------------------------------------------------------------
53
54CONTAINS
55
56   SUBROUTINE lim_sbc_tau( kt, kcpl )
57      !!-------------------------------------------------------------------
58      !!                ***  ROUTINE lim_sbc_tau ***
59      !! 
60      !! ** Purpose : Update the ocean surface stresses due to the ice
61      !!         
62      !! ** Action  : - compute the ice-ocean stress depending on kcpl:
63      !!              fluxes at the ice-ocean interface.
64      !!     Case 0  :  old LIM-3 way, computed at ice time-step only
65      !!     Case 1  :  at each ocean time step the stresses are computed
66      !!                using the current ocean velocity (now)
67      !!     Case 2  :  combination of half case 0 + half case 1
68      !!     
69      !! ** Outputs : - utau   : sea surface i-stress (ocean referential)
70      !!              - vtau   : sea surface j-stress (ocean referential)
71      !!
72      !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90.
73      !!              Tartinville et al. 2001 Ocean Modelling, 3, 95-108.
74      !!---------------------------------------------------------------------
75      INTEGER ::   kt     ! number of ocean iteration
76      INTEGER ::   kcpl   ! ice-ocean coupling indicator: =0  LIM-3 old case
77      !                   !                               =1  stresses computed using now ocean velocity
78      !                   !                               =2  combination of 0 and 1 cases
79      !!
80      INTEGER  ::   ji, jj           ! dummy loop indices
81      REAL(wp) ::   zfrldu, zfrldv   ! lead fraction at U- & V-points
82      REAL(wp) ::   ztglx , ztgly
83      REAL(wp) ::   zat_u, zu_ico, zutaui, zu_u, zv_u, zmodu, zmod
84      REAL(wp) ::   zat_v, zv_ico, zvtaui, zu_v, zv_v, zmodv, zsang
85
86#if defined key_coupled   
87      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   zalb     ! albedo of ice under overcast sky
88      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   zalbp    ! albedo of ice under clear sky
89#endif
90      REAL(wp), DIMENSION(jpi,jpj) ::   ztio_u, ztio_v   ! ocean stress below sea-ice
91      !!---------------------------------------------------------------------
92
93      IF( kt == nit000 ) THEN
94         IF(lwp) WRITE(numout,*)
95         IF(lwp) WRITE(numout,*) 'lim_sbc_tau : LIM 3.0 sea-ice - surface ocean momentum fluxes'
96         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ '
97      ENDIF
98
99      SELECT CASE( kcpl )
100         !                                           !--------------------------------!
101      CASE( 0 )                                   !  LIM 3 old stress computation  !  (at ice timestep only)
102         !                                        !--------------------------------!
103         ! ... ice stress over ocean with a ice-ocean rotation angle
104         DO jj = 1, jpj
105            ! ... change the sinus angle sign in the south hemisphere
106            zsang  = SIGN(1.e0, gphif(1,jj) ) * sangvg
107            DO ji = 1, jpi
108               ! ... ice velocity relative to the ocean
109               zu_ico = u_ice(ji,jj) - u_oce(ji,jj)
110               zv_ico = v_ice(ji,jj) - v_oce(ji,jj)
111               zmod  = SQRT( zu_ico * zu_ico + zv_ico * zv_ico ) 
112               ! quadratic drag formulation
113               ztglx = rhoco * zmod * ( cangvg * zu_ico - zsang * zv_ico )
114               ztgly = rhoco * zmod * ( cangvg * zv_ico + zsang * zu_ico )
115               ! IMPORTANT
116               ! these lines are bound to prevent numerical oscillations
117               ! in the ice-ocean stress
118               ! They are physically ill-based. There is a cleaner solution
119               ! to try (remember discussion in Paris Gurvan)
120               ztio_u(ji,jj) = ztglx * exp( - zmod / 0.5 ) 
121               ztio_v(ji,jj) = ztgly * exp( - zmod / 0.5 ) 
122               !
123            END DO
124         END DO
125
126         DO jj = 2, jpjm1
127            DO ji = fs_2, fs_jpim1   ! vertor opt.
128               ! ... open-ocean (lead) fraction at U- & V-points (from T-point values)
129               zfrldu = 0.5 * ( ato_i(ji,jj) + ato_i(ji+1,jj  ) )
130               zfrldv = 0.5 * ( ato_i(ji,jj) + ato_i(ji  ,jj+1) )
131               ! update surface ocean stress
132               utau(ji,jj) = zfrldu * utau(ji,jj) + ( 1. - zfrldu ) * ztio_u(ji,jj)
133               vtau(ji,jj) = zfrldv * vtau(ji,jj) + ( 1. - zfrldv ) * ztio_v(ji,jj)
134               !
135            END DO
136         END DO
137
138         ! boundary condition on the stress (utau,vtau)
139         CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )
140         !
141         !                                        !--------------------------------!
142      CASE( 1 )                                   !  Use the now velocity          !  (at each ocean timestep)
143         !                                        !--------------------------------!
144         ! ... save the air-ocean stresses at ice time-step
145         IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN
146            utau_oce(:,:) = utau(:,:)
147            vtau_oce(:,:) = vtau(:,:)
148         ENDIF
149         ! ... ice stress over ocean with a ice-ocean rotation angle
150         DO jj = 2, jpjm1
151            zsang  = SIGN(1.e0, gphif(1,jj-1) ) * sangvg
152            DO ji = fs_2, fs_jpim1
153               ! computation of wind stress over ocean in X and Y direction
154               zat_u  = ( at_i(ji,jj) + at_i(ji+1,jj) ) * 0.5      ! ice area at u and V-points
155               zat_v  = ( at_i(ji,jj) + at_i(ji,jj+1) ) * 0.5
156
157               zu_u   = u_ice(ji,jj) - un(ji,jj,1)                  ! u ice-ocean velocity at U-point
158               zv_v   = v_ice(ji,jj) - vn(ji,jj,1)                  ! v ice-ocean velocity at V-point
159               !                                                    ! u ice-ocean velocity at V-point
160               zu_v   = 0.25 * (  u_ice(ji,jj  ) - un(ji,jj  ,1)  +  u_ice(ji+1,jj  ) - un(ji+1,jj  ,1)    &
161                  &             + u_ice(ji,jj-1) - un(ji,jj-1,1)  +  u_ice(ji+1,jj-1) - un(ji+1,jj-1,1)  ) 
162               !                                                    ! v ice-ocean velocity at U-point
163               zv_u   = 0.25 * (  v_ice(ji-1,jj+1) - vn(ji-1,jj+1,1)  +  v_ice(ji,jj+1) - vn(ji,jj+1,1)    &
164                  &             + v_ice(ji-1,jj  ) - vn(ji-1,jj  ,1)  +  v_ice(ji,jj  ) - vn(ji,jj  ,1)  )
165               zmodu    = SQRT( zu_u * zu_u + zv_u * zv_u )         ! module of the ice-ocean velocity at U-point
166               zmodv    = SQRT( zu_v * zu_v + zv_v * zv_v )         ! module of the ice-ocean velocity at V-point
167               zutaui   = rhoco * zmodu * ( cangvg * zu_u - zsang * zv_u ) 
168               zvtaui   = rhoco * zmodv * ( cangvg * zv_v + zsang * zu_v ) 
169
170               utau(ji,jj) = ( 1.-zat_u ) * utau_oce(ji,jj) + zat_u * zutaui    ! stress at the ocean surface
171               vtau(ji,jj) = ( 1.-zat_v ) * vtau_oce(ji,jj) + zat_v * zvtaui
172            END DO
173         END DO
174         ! boundary condition on the stress (utau,vtau)
175         CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )
176         !
177         !                                        !--------------------------------!
178      CASE( 2 )                                   !  mixed 0 and 2 cases           !  (at each ocean timestep)
179         !                                        !--------------------------------!
180         ! ... save the air-ocean stresses at ice time-step
181         IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN
182            utau_oce(:,:) = utau(:,:)
183            vtau_oce(:,:) = vtau(:,:)
184         ENDIF
185         ! ... ice stress over ocean with a ice-ocean rotation angle
186         DO jj = 2, jpjm1
187            zsang  = SIGN(1.e0, gphif(1,jj-1) ) * sangvg
188            DO ji = fs_2, fs_jpim1
189               ! computation of wind stress over ocean in X and Y direction
190               zat_u = at_i(ji,jj) + at_i(ji+1,jj) * 0.5     ! ice area at u and V-points
191               zat_v = at_i(ji,jj) + at_i(ji,jj+1) * 0.5 
192
193               !!gm bug mixing U and V points value below     ====>>> to be corrected
194               zu_ico   = u_ice(ji,jj) - 0.5 * ( un(ji,jj,1) - ssu_m(ji,jj) )   ! ice-oce velocity using un and ssu_m
195               zv_ico   = v_ice(ji,jj) - 0.5 * ( vn(ji,jj,1) - ssu_m(ji,jj) )
196               !                                        ! module of the ice-ocean velocity
197               zmod     = SQRT( zu_ico * zu_ico + zv_ico * zv_ico )
198               !                                        ! quadratic drag formulation with rotation
199               zutaui   = rhoco * zmod * ( cangvg * zu_ico - zsang * zv_ico )
200               zvtaui   = rhoco * zmod * ( cangvg * zv_ico + zsang * zu_ico )
201               !
202               utau(ji,jj) = ( 1.-zat_u ) * utau_oce(ji,jj) + zat_u * zutaui    ! stress at the ocean surface
203               vtau(ji,jj) = ( 1.-zat_v ) * vtau_oce(ji,jj) + zat_v * zvtaui
204            END DO
205         END DO
206         ! boundary condition on the stress (utau,vtau)
207         CALL lbc_lnk( utau, 'U', -1. )   ;   CALL lbc_lnk( vtau, 'V', -1. )
208         !
209      END SELECT
210      IF(ln_ctl)   CALL prt_ctl( tab2d_1=utau, clinfo1=' lim_sbc: utau   : ', mask1=umask,   &
211         &                       tab2d_2=vtau, clinfo2=' vtau    : '        , mask2=vmask )
212     
213   END SUBROUTINE lim_sbc_tau
214
215
216   SUBROUTINE lim_sbc_flx( kt )
217      !!-------------------------------------------------------------------
218      !!                ***  ROUTINE lim_sbc_flx ***
219      !! 
220      !! ** Purpose :   Update the surface ocean boundary condition for heat
221      !!              salt and mass over areas where sea-ice is non-zero
222      !!         
223      !! ** Action  : - computes the heat and freshwater/salt fluxes
224      !!              at the ice-ocean interface.
225      !!              - Update the ocean sbc
226      !!     
227      !! ** Outputs : - qsr    : sea heat flux:     solar
228      !!              - qns    : sea heat flux: non solar
229      !!              - emp    : freshwater budget: volume flux
230      !!              - emps   : freshwater budget: concentration/dillution
231      !!
232      !! References : Goosse, H. et al. 1996, Bul. Soc. Roy. Sc. Liege, 65, 87-90.
233      !!              Tartinville et al. 2001 Ocean Modelling, 3, 95-108.
234      !!---------------------------------------------------------------------
235      INTEGER ::   kt    ! number of iteration
236      !!
237      INTEGER  ::   ji, jj           ! dummy loop indices
238      INTEGER  ::   ifvt, i1mfr, idfr               ! some switches
239      INTEGER  ::   iflt, ial, iadv, ifral, ifrdv
240      REAL(wp) ::   zinda            ! switch for testing the values of ice concentration
241      REAL(wp) ::   zfons            ! salt exchanges at the ice/ocean interface
242      REAL(wp) ::   zpme             ! freshwater exchanges at the ice/ocean interface
243      REAL(wp), DIMENSION(jpi,jpj) ::   zfcm1 , zfcm2    ! solar/non solar heat fluxes
244      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   zalb     ! albedo of ice under overcast sky
245      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   zalbp    ! albedo of ice under clear sky
246      !!---------------------------------------------------------------------
247
248      IF( kt == nit000 ) THEN
249         IF(lwp) WRITE(numout,*)
250         IF(lwp) WRITE(numout,*) 'lim_sbc_flx : LIM 3.0 sea-ice - heat salt and mass ocean surface fluxes'
251         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ '
252      ENDIF
253
254      !------------------------------------------!
255      !      heat flux at the ocean surface      !
256      !------------------------------------------!
257      ! pfrld is the lead fraction at the previous time step (actually between TRP and THD)
258      ! changed to old_frld and old ht_i
259
260      DO jj = 1, jpj
261         DO ji = 1, jpi
262            zinda   = 1.0 - MAX( rzero , SIGN( rone , - ( 1.0 - pfrld(ji,jj) ) ) )
263            ifvt    = zinda  *  MAX( rzero , SIGN( rone, -phicif  (ji,jj) ) )  !subscripts are bad here
264            i1mfr   = 1.0 - MAX( rzero , SIGN( rone ,  - ( at_i(ji,jj)       ) ) )
265            idfr    = 1.0 - MAX( rzero , SIGN( rone , ( 1.0 - at_i(ji,jj) ) - pfrld(ji,jj) ) )
266            iflt    = zinda  * (1 - i1mfr) * (1 - ifvt )
267            ial     = ifvt   * i1mfr + ( 1 - ifvt ) * idfr
268            iadv    = ( 1  - i1mfr ) * zinda
269            ifral   = ( 1  - i1mfr * ( 1 - ial ) )   
270            ifrdv   = ( 1  - ifral * ( 1 - ial ) ) * iadv 
271
272            ! switch --- 1.0 ---------------- 0.0 --------------------
273            ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
274            ! zinda   | if pfrld = 1       | if pfrld < 1            |
275            !  -> ifvt| if pfrld old_ht_i
276            ! i1mfr   | if frld = 1        | if frld  < 1            |
277            ! idfr    | if frld <= pfrld    | if frld > pfrld        |
278            ! iflt    |
279            ! ial     |
280            ! iadv    |
281            ! ifral
282            ! ifrdv
283
284            !   computation the solar flux at ocean surface
285            zfcm1(ji,jj)   = pfrld(ji,jj) * qsr(ji,jj)  + ( 1. - pfrld(ji,jj) ) * fstric(ji,jj)
286            ! fstric     Solar flux transmitted trough the ice
287            ! qsr        Net short wave heat flux on free ocean
288            ! new line
289            fscmbq(ji,jj) = ( 1.0 - pfrld(ji,jj) ) * fstric(ji,jj)
290
291            !  computation the non solar heat flux at ocean surface
292            zfcm2(ji,jj) = - zfcm1(ji,jj)                  &
293               &           + iflt    * ( fscmbq(ji,jj) )   & ! total abl -> fscmbq is given to the ocean
294               ! fscmbq and ffltbif are obsolete
295               !              &           + iflt * ffltbif(ji,jj) !!! only if one category is used
296               &           + ifral   * ( ial * qcmif(ji,jj) + (1 - ial) * qldif(ji,jj) ) / rdt_ice   &
297               &           + ifrdv   * ( qfvbq(ji,jj) + qdtcn(ji,jj) ) / rdt_ice                     &
298               &           + fhmec(ji,jj)     & ! new contribution due to snow melt in ridging!!
299               &           + fheat_rpo(ji,jj) & ! contribution from ridge formation
300               &           + fheat_res(ji,jj)
301            ! fscmbq  Part of the solar radiation transmitted through the ice and going to the ocean
302            !         computed in limthd_zdf.F90
303            ! ffltbif Total heat content of the ice (brine pockets+ice) / delta_t
304            ! qcmif   Energy needed to bring the ocean surface layer until its freezing (ok)
305            ! qldif   heat balance of the lead (or of the open ocean)
306            ! qfvbq   i think this is wrong!
307            ! ---> Array used to store energy in case of total lateral ablation
308            ! qfvbq latent heat uptake/release after accretion/ablation
309            ! qdtcn Energy from the turbulent oceanic heat flux heat flux coming in the lead
310
311            IF ( num_sal .EQ. 2 ) zfcm2(ji,jj) = zfcm2(ji,jj) + &
312               fhbri(ji,jj) ! new contribution due to brine drainage
313
314            ! bottom radiative component is sent to the computation of the
315            ! oceanic heat flux
316            fsbbq(ji,jj) = ( 1.0 - ( ifvt + iflt ) ) * fscmbq(ji,jj)     
317
318            ! used to compute the oceanic heat flux at the next time step
319            qsr(ji,jj) = zfcm1(ji,jj)                                       ! solar heat flux
320            qns(ji,jj) = zfcm2(ji,jj) - fdtcn(ji,jj)                        ! non solar heat flux
321            !                           ! fdtcn : turbulent oceanic heat flux
322
323            !!gm   this IF prevents the vertorisation of the whole loop
324            IF ( ( ji .EQ. jiindx ) .AND. ( jj .EQ. jjindx) ) THEN
325               WRITE(numout,*) ' lim_sbc : heat fluxes '
326               WRITE(numout,*) ' qsr       : ', qsr(jiindx,jjindx)
327               WRITE(numout,*) ' zfcm1     : ', zfcm1(jiindx,jjindx)
328               WRITE(numout,*) ' pfrld     : ', pfrld(jiindx,jjindx)
329               WRITE(numout,*) ' fstric    : ', fstric (jiindx,jjindx)
330               WRITE(numout,*)
331               WRITE(numout,*) ' qns       : ', qns(jiindx,jjindx)
332               WRITE(numout,*) ' zfcm2     : ', zfcm2(jiindx,jjindx)
333               WRITE(numout,*) ' zfcm1     : ', zfcm1(jiindx,jjindx)
334               WRITE(numout,*) ' ifral     : ', ifral
335               WRITE(numout,*) ' ial       : ', ial 
336               WRITE(numout,*) ' qcmif     : ', qcmif(jiindx,jjindx)
337               WRITE(numout,*) ' qldif     : ', qldif(jiindx,jjindx)
338               WRITE(numout,*) ' qcmif / dt: ', qcmif(jiindx,jjindx) / rdt_ice
339               WRITE(numout,*) ' qldif / dt: ', qldif(jiindx,jjindx) / rdt_ice
340               WRITE(numout,*) ' ifrdv     : ', ifrdv
341               WRITE(numout,*) ' qfvbq     : ', qfvbq(jiindx,jjindx)
342               WRITE(numout,*) ' qdtcn     : ', qdtcn(jiindx,jjindx)
343               WRITE(numout,*) ' qfvbq / dt: ', qfvbq(jiindx,jjindx) / rdt_ice
344               WRITE(numout,*) ' qdtcn / dt: ', qdtcn(jiindx,jjindx) / rdt_ice
345               WRITE(numout,*) ' '
346               WRITE(numout,*) ' fdtcn     : ', fdtcn(jiindx,jjindx)
347               WRITE(numout,*) ' fhmec     : ', fhmec(jiindx,jjindx)
348               WRITE(numout,*) ' fheat_rpo : ', fheat_rpo(jiindx,jjindx)
349               WRITE(numout,*) ' fhbri     : ', fhbri(jiindx,jjindx)
350               WRITE(numout,*) ' fheat_res : ', fheat_res(jiindx,jjindx)
351            ENDIF
352            !!gm   end
353         END DO
354      END DO
355
356      !------------------------------------------!
357      !      mass flux at the ocean surface      !
358      !------------------------------------------!
359
360      !!gm   optimisation: this loop have to be merged with the previous one
361      DO jj = 1, jpj
362         DO ji = 1, jpi
363            !  case of realistic freshwater flux (Tartinville et al., 2001) (presently ACTIVATED)
364            !  -------------------------------------------------------------------------------------
365            !  The idea of this approach is that the system that we consider is the ICE-OCEAN system
366            !  Thus  FW  flux  =  External ( E-P+snow melt)
367            !       Salt flux  =  Exchanges in the ice-ocean system then converted into FW
368            !                     Associated to Ice formation AND Ice melting
369            !                     Even if i see Ice melting as a FW and SALT flux
370            !       
371
372            !  computing freshwater exchanges at the ice/ocean interface
373            zpme = - emp(ji,jj)     * ( 1.0 - at_i(ji,jj) )  &   !  evaporation over oceanic fraction
374               &   + tprecip(ji,jj) *         at_i(ji,jj)    &   !  total precipitation
375               ! old fashioned way               
376               !              &   - sprecip(ji,jj) * ( 1. - pfrld(ji,jj) )  &   !  remov. snow precip over ice
377               &   - sprecip(ji,jj) * ( 1. - (pfrld(ji,jj)**betas) )  &   !  remov. snow precip over ice
378               &   - rdmsnif(ji,jj) / rdt_ice                &   !  freshwaterflux due to snow melting
379               ! new contribution from snow falling when ridging
380               &   + fmmec(ji,jj)
381
382            !  computing salt exchanges at the ice/ocean interface
383            !  sice should be the same as computed with the ice model
384            zfons =  ( soce - sice ) * ( rdmicif(ji,jj) / rdt_ice ) 
385            ! SOCE
386            zfons =  ( sss_m(ji,jj) - sice ) * ( rdmicif(ji,jj) / rdt_ice ) 
387
388            !CT useless            !  salt flux for constant salinity
389            !CT useless            fsalt(ji,jj)      =  zfons / ( sss_m(ji,jj) + epsi16 ) + fsalt_res(ji,jj)
390            !  salt flux for variable salinity
391            zinda             = 1.0 - MAX( rzero , SIGN( rone , - ( 1.0 - pfrld(ji,jj) ) ) )
392            !  correcting brine and salt fluxes
393            fsbri(ji,jj)      =  zinda*fsbri(ji,jj)
394            !  converting the salt fluxes from ice to a freshwater flux from ocean
395            fsalt_res(ji,jj)  =  fsalt_res(ji,jj) / ( sss_m(ji,jj) + epsi16 )
396            fseqv(ji,jj)      =  fseqv(ji,jj)     / ( sss_m(ji,jj) + epsi16 )
397            fsbri(ji,jj)      =  fsbri(ji,jj)     / ( sss_m(ji,jj) + epsi16 )
398            fsalt_rpo(ji,jj)  =  fsalt_rpo(ji,jj) / ( sss_m(ji,jj) + epsi16 )
399
400            !  freshwater mass exchange (positive to the ice, negative for the ocean ?)
401            !  actually it's a salt flux (so it's minus freshwater flux)
402            !  if sea ice grows, zfons is positive, fsalt also
403            !  POSITIVE SALT FLUX FROM THE ICE TO THE OCEAN
404            !  POSITIVE FRESHWATER FLUX FROM THE OCEAN TO THE ICE [kg.m-2.s-1]
405
406            emp(ji,jj) = - zpme 
407         END DO
408      END DO
409
410      IF( num_sal == 2 ) THEN      ! variable ice salinity: brine drainage included in the salt flux
411         emps(:,:) = fsbri(:,:) + fseqv(:,:) + fsalt_res(:,:) + fsalt_rpo(:,:) + emp(:,:)
412      ELSE                         ! constant ice salinity:
413         emps(:,:) =              fseqv(:,:) + fsalt_res(:,:) + fsalt_rpo(:,:) + emp(:,:)
414      ENDIF
415
416      IF( lk_dynspg_rl )    emp (:,:) = emps(:,:)      ! rigid-lid formulation : emp = emps
417
418      !-----------------------------------------------!
419      !   Storing the transmitted variables           !
420      !-----------------------------------------------!
421
422      freeze(:,:)   = at_i(:,:)             ! Sea ice cover 
423
424      IF ( lk_cpl ) THEN               
425         ! Ice surface temperature
426         tn_ice(:,:,:) = t_su(:,:,:)           
427         ! Computation of snow/ice and ocean albedo
428         CALL albedo_ice( t_su, ht_i, ht_s, zalbp, zalb )
429         alb_ice(:,:,:) =  0.5 * ( zalbp(:,:,:) + zalb (:,:,:) )   ! Ice albedo (mean clear and overcast skys)
430      ENDIF
431         
432      IF(ln_ctl) THEN
433         CALL prt_ctl( tab2d_1=qsr   , clinfo1=' lim_sbc: qsr    : ', tab2d_2=qns , clinfo2=' qns     : ' )
434         CALL prt_ctl( tab2d_1=emp   , clinfo1=' lim_sbc: emp    : ', tab2d_2=emps, clinfo2=' emps    : ' )
435         CALL prt_ctl( tab2d_1=freeze, clinfo1=' lim_sbc: freeze : ' )
436         CALL prt_ctl( tab3d_1=tn_ice, clinfo1=' lim_sbc: tn_ice : ', kdim=jpl )
437      ENDIF
438      !
439   END SUBROUTINE lim_sbc_flx
440
441
442#else
443   !!----------------------------------------------------------------------
444   !!   Default option :        Dummy module       NO LIM 3.0 sea-ice model
445   !!----------------------------------------------------------------------
446CONTAINS
447   SUBROUTINE lim_sbc           ! Dummy routine
448   END SUBROUTINE lim_sbc
449#endif 
450
451   !!======================================================================
452END MODULE limsbc
Note: See TracBrowser for help on using the repository browser.