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.
sbccpl.F90 in branches/UKMO/r5936_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/r5936_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 7168

Last change on this file since 7168 was 7168, checked in by jcastill, 7 years ago

First trial to changes needed for wave coupling

File size: 141.3 KB
Line 
1MODULE sbccpl
2   !!======================================================================
3   !!                       ***  MODULE  sbccpl  ***
4   !! Surface Boundary Condition :  momentum, heat and freshwater fluxes in coupled mode
5   !!======================================================================
6   !! History :  2.0  ! 2007-06  (R. Redler, N. Keenlyside, W. Park) Original code split into flxmod & taumod
7   !!            3.0  ! 2008-02  (G. Madec, C Talandier)  surface module
8   !!            3.1  ! 2009_02  (G. Madec, S. Masson, E. Maisonave, A. Caubel) generic coupled interface
9   !!            3.4  ! 2011_11  (C. Harris) more flexibility + multi-category fields
10   !!----------------------------------------------------------------------
11   !!----------------------------------------------------------------------
12   !!   namsbc_cpl      : coupled formulation namlist
13   !!   sbc_cpl_init    : initialisation of the coupled exchanges
14   !!   sbc_cpl_rcv     : receive fields from the atmosphere over the ocean (ocean only)
15   !!                     receive stress from the atmosphere over the ocean (ocean-ice case)
16   !!   sbc_cpl_ice_tau : receive stress from the atmosphere over ice
17   !!   sbc_cpl_ice_flx : receive fluxes from the atmosphere over ice
18   !!   sbc_cpl_snd     : send     fields to the atmosphere
19   !!----------------------------------------------------------------------
20   USE dom_oce         ! ocean space and time domain
21   USE sbc_oce         ! Surface boundary condition: ocean fields
22   USE sbc_ice         ! Surface boundary condition: ice fields
23   USE sbcapr
24   USE sbcdcy          ! surface boundary condition: diurnal cycle
25   USE sbcwave         ! surface boundary condition: waves
26   USE phycst          ! physical constants
27#if defined key_lim3
28   USE ice             ! ice variables
29#endif
30#if defined key_lim2
31   USE par_ice_2       ! ice parameters
32   USE ice_2           ! ice variables
33#endif
34   USE cpl_oasis3      ! OASIS3 coupling
35   USE geo2ocean       !
36   USE oce   , ONLY : tsn, un, vn, sshn, ub, vb, sshb, fraqsr_1lev
37   USE albedo          !
38   USE in_out_manager  ! I/O manager
39   USE iom             ! NetCDF library
40   USE lib_mpp         ! distribued memory computing library
41   USE wrk_nemo        ! work arrays
42   USE timing          ! Timing
43   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
44   USE eosbn2
45   USE sbcrnf   , ONLY : l_rnfcpl
46#if defined key_cpl_carbon_cycle
47   USE p4zflx, ONLY : oce_co2
48#endif
49#if defined key_cice
50   USE ice_domain_size, only: ncat
51#endif
52#if defined key_lim3
53   USE limthd_dh       ! for CALL lim_thd_snwblow
54#endif
55
56   IMPLICIT NONE
57   PRIVATE
58
59   PUBLIC   sbc_cpl_init       ! routine called by sbcmod.F90
60   PUBLIC   sbc_cpl_rcv        ! routine called by sbc_ice_lim(_2).F90
61   PUBLIC   sbc_cpl_snd        ! routine called by step.F90
62   PUBLIC   sbc_cpl_ice_tau    ! routine called by sbc_ice_lim(_2).F90
63   PUBLIC   sbc_cpl_ice_flx    ! routine called by sbc_ice_lim(_2).F90
64   PUBLIC   sbc_cpl_alloc      ! routine called in sbcice_cice.F90
65
66   INTEGER, PARAMETER ::   jpr_otx1   =  1            ! 3 atmosphere-ocean stress components on grid 1
67   INTEGER, PARAMETER ::   jpr_oty1   =  2            !
68   INTEGER, PARAMETER ::   jpr_otz1   =  3            !
69   INTEGER, PARAMETER ::   jpr_otx2   =  4            ! 3 atmosphere-ocean stress components on grid 2
70   INTEGER, PARAMETER ::   jpr_oty2   =  5            !
71   INTEGER, PARAMETER ::   jpr_otz2   =  6            !
72   INTEGER, PARAMETER ::   jpr_itx1   =  7            ! 3 atmosphere-ice   stress components on grid 1
73   INTEGER, PARAMETER ::   jpr_ity1   =  8            !
74   INTEGER, PARAMETER ::   jpr_itz1   =  9            !
75   INTEGER, PARAMETER ::   jpr_itx2   = 10            ! 3 atmosphere-ice   stress components on grid 2
76   INTEGER, PARAMETER ::   jpr_ity2   = 11            !
77   INTEGER, PARAMETER ::   jpr_itz2   = 12            !
78   INTEGER, PARAMETER ::   jpr_qsroce = 13            ! Qsr above the ocean
79   INTEGER, PARAMETER ::   jpr_qsrice = 14            ! Qsr above the ice
80   INTEGER, PARAMETER ::   jpr_qsrmix = 15 
81   INTEGER, PARAMETER ::   jpr_qnsoce = 16            ! Qns above the ocean
82   INTEGER, PARAMETER ::   jpr_qnsice = 17            ! Qns above the ice
83   INTEGER, PARAMETER ::   jpr_qnsmix = 18
84   INTEGER, PARAMETER ::   jpr_rain   = 19            ! total liquid precipitation (rain)
85   INTEGER, PARAMETER ::   jpr_snow   = 20            ! solid precipitation over the ocean (snow)
86   INTEGER, PARAMETER ::   jpr_tevp   = 21            ! total evaporation
87   INTEGER, PARAMETER ::   jpr_ievp   = 22            ! solid evaporation (sublimation)
88   INTEGER, PARAMETER ::   jpr_sbpr   = 23            ! sublimation - liquid precipitation - solid precipitation
89   INTEGER, PARAMETER ::   jpr_semp   = 24            ! solid freshwater budget (sublimation - snow)
90   INTEGER, PARAMETER ::   jpr_oemp   = 25            ! ocean freshwater budget (evap - precip)
91   INTEGER, PARAMETER ::   jpr_w10m   = 26            ! 10m wind
92   INTEGER, PARAMETER ::   jpr_dqnsdt = 27            ! d(Q non solar)/d(temperature)
93   INTEGER, PARAMETER ::   jpr_rnf    = 28            ! runoffs
94   INTEGER, PARAMETER ::   jpr_cal    = 29            ! calving
95   INTEGER, PARAMETER ::   jpr_taum   = 30            ! wind stress module
96   INTEGER, PARAMETER ::   jpr_co2    = 31
97   INTEGER, PARAMETER ::   jpr_topm   = 32            ! topmeltn
98   INTEGER, PARAMETER ::   jpr_botm   = 33            ! botmeltn
99   INTEGER, PARAMETER ::   jpr_sflx   = 34            ! salt flux
100   INTEGER, PARAMETER ::   jpr_toce   = 35            ! ocean temperature
101   INTEGER, PARAMETER ::   jpr_soce   = 36            ! ocean salinity
102   INTEGER, PARAMETER ::   jpr_ocx1   = 37            ! ocean current on grid 1
103   INTEGER, PARAMETER ::   jpr_ocy1   = 38            !
104   INTEGER, PARAMETER ::   jpr_ssh    = 39            ! sea surface height
105   INTEGER, PARAMETER ::   jpr_fice   = 40            ! ice fraction         
106   INTEGER, PARAMETER ::   jpr_e3t1st = 41            ! first T level thickness
107   INTEGER, PARAMETER ::   jpr_fraqsr = 42            ! fraction of solar net radiation absorbed in the first ocean level
108   INTEGER, PARAMETER ::   jpr_mslp   = 43            ! mean sea level pressure
109   INTEGER, PARAMETER ::   jpr_hsig   = 44            ! Hsig
110   INTEGER, PARAMETER ::   jpr_phioc  = 45            ! Wave=>ocean energy flux
111   INTEGER, PARAMETER ::   jpr_sdrftx = 46            ! Stokes drift on grid 1
112   INTEGER, PARAMETER ::   jpr_sdrfty = 47            ! Stokes drift on grid 2
113   INTEGER, PARAMETER ::   jpr_wper   = 48            ! Mean wave period
114   INTEGER, PARAMETER ::   jpr_wnum   = 49            ! Mean wavenumber
115   INTEGER, PARAMETER ::   jpr_wstrf  = 50            ! Stress fraction adsorbed by waves
116   INTEGER, PARAMETER ::   jpr_wdrag  = 51            ! Neutral surface drag coefficient
117   INTEGER, PARAMETER ::   jprcv      = 51            ! total number of fields received 
118
119   INTEGER, PARAMETER ::   jps_fice   =  1            ! ice fraction sent to the atmosphere
120   INTEGER, PARAMETER ::   jps_toce   =  2            ! ocean temperature
121   INTEGER, PARAMETER ::   jps_tice   =  3            ! ice   temperature
122   INTEGER, PARAMETER ::   jps_tmix   =  4            ! mixed temperature (ocean+ice)
123   INTEGER, PARAMETER ::   jps_albice =  5            ! ice   albedo
124   INTEGER, PARAMETER ::   jps_albmix =  6            ! mixed albedo
125   INTEGER, PARAMETER ::   jps_hice   =  7            ! ice  thickness
126   INTEGER, PARAMETER ::   jps_hsnw   =  8            ! snow thickness
127   INTEGER, PARAMETER ::   jps_ocx1   =  9            ! ocean current on grid 1
128   INTEGER, PARAMETER ::   jps_ocy1   = 10            !
129   INTEGER, PARAMETER ::   jps_ocz1   = 11            !
130   INTEGER, PARAMETER ::   jps_ivx1   = 12            ! ice   current on grid 1
131   INTEGER, PARAMETER ::   jps_ivy1   = 13            !
132   INTEGER, PARAMETER ::   jps_ivz1   = 14            !
133   INTEGER, PARAMETER ::   jps_co2    = 15
134   INTEGER, PARAMETER ::   jps_soce   = 16            ! ocean salinity
135   INTEGER, PARAMETER ::   jps_ssh    = 17            ! sea surface height
136   INTEGER, PARAMETER ::   jps_qsroce = 18            ! Qsr above the ocean
137   INTEGER, PARAMETER ::   jps_qnsoce = 19            ! Qns above the ocean
138   INTEGER, PARAMETER ::   jps_oemp   = 20            ! ocean freshwater budget (evap - precip)
139   INTEGER, PARAMETER ::   jps_sflx   = 21            ! salt flux
140   INTEGER, PARAMETER ::   jps_otx1   = 22            ! 2 atmosphere-ocean stress components on grid 1
141   INTEGER, PARAMETER ::   jps_oty1   = 23            !
142   INTEGER, PARAMETER ::   jps_rnf    = 24            ! runoffs
143   INTEGER, PARAMETER ::   jps_taum   = 25            ! wind stress module
144   INTEGER, PARAMETER ::   jps_fice2  = 26            ! ice fraction sent to OPA (by SAS when doing SAS-OPA coupling)
145   INTEGER, PARAMETER ::   jps_e3t1st = 27            ! first level depth (vvl)
146   INTEGER, PARAMETER ::   jps_fraqsr = 28            ! fraction of solar net radiation absorbed in the first ocean level
147   INTEGER, PARAMETER ::   jps_ficet  = 29            ! total ice fraction 
148   INTEGER, PARAMETER ::   jps_ocxw   = 30            ! currents on grid 1 
149   INTEGER, PARAMETER ::   jps_ocyw   = 31            ! currents on grid 2
150   INTEGER, PARAMETER ::   jps_wlev   = 32            ! water level
151   INTEGER, PARAMETER ::   jpsnd      = 32            ! total number of fields sent
152
153   !                                                         !!** namelist namsbc_cpl **
154   TYPE ::   FLD_C
155      CHARACTER(len = 32) ::   cldes                  ! desciption of the coupling strategy
156      CHARACTER(len = 32) ::   clcat                  ! multiple ice categories strategy
157      CHARACTER(len = 32) ::   clvref                 ! reference of vector ('spherical' or 'cartesian')
158      CHARACTER(len = 32) ::   clvor                  ! orientation of vector fields ('eastward-northward' or 'local grid')
159      CHARACTER(len = 32) ::   clvgrd                 ! grids on which is located the vector fields
160   END TYPE FLD_C
161   ! Send to the atmosphere                           !
162   TYPE(FLD_C) ::   sn_snd_temp, sn_snd_alb, sn_snd_thick, sn_snd_crt, sn_snd_co2                       
163   ! Received from the atmosphere                     !
164   TYPE(FLD_C) ::   sn_rcv_w10m, sn_rcv_taumod, sn_rcv_tau, sn_rcv_dqnsdt, sn_rcv_qsr, sn_rcv_qns, sn_rcv_emp, sn_rcv_rnf
165   TYPE(FLD_C) ::   sn_rcv_cal, sn_rcv_iceflx, sn_rcv_co2, sn_rcv_mslp                           
166   ! Send to waves
167   TYPE(FLD_C) ::   sn_snd_ifrac, sn_snd_crtw, sn_snd_wlev 
168   ! Received from waves
169   TYPE(FLD_C) ::   sn_rcv_hsig,sn_rcv_phioc,sn_rcv_sdrfx,sn_rcv_sdrfy,sn_rcv_wper,sn_rcv_wnum,sn_rcv_wstrf,sn_rcv_wdrag
170   ! Other namelist parameters                        !
171   INTEGER     ::   nn_cplmodel            ! Maximum number of models to/from which NEMO is potentialy sending/receiving data
172   LOGICAL     ::   ln_usecplmask          !  use a coupling mask file to merge data received from several models
173                                           !   -> file cplmask.nc with the float variable called cplmask (jpi,jpj,nn_cplmodel)
174   TYPE ::   DYNARR     
175      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   z3   
176   END TYPE DYNARR
177
178   TYPE( DYNARR ), SAVE, DIMENSION(jprcv) ::   frcv                      ! all fields recieved from the atmosphere
179
180   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   albedo_oce_mix     ! ocean albedo sent to atmosphere (mix clear/overcast sky)
181
182   REAL(wp) ::   rpref = 101000._wp   ! reference atmospheric pressure[N/m2]
183   REAL(wp) ::   r1_grau              ! = 1.e0 / (grav * rau0)
184
185   INTEGER , ALLOCATABLE, SAVE, DIMENSION(    :) ::   nrcvinfo           ! OASIS info argument
186
187   !! Substitution
188#  include "domzgr_substitute.h90"
189#  include "vectopt_loop_substitute.h90"
190   !!----------------------------------------------------------------------
191   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
192   !! $Id$
193   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
194   !!----------------------------------------------------------------------
195
196CONTAINS
197 
198   INTEGER FUNCTION sbc_cpl_alloc()
199      !!----------------------------------------------------------------------
200      !!             ***  FUNCTION sbc_cpl_alloc  ***
201      !!----------------------------------------------------------------------
202      INTEGER :: ierr(4)
203      !!----------------------------------------------------------------------
204      ierr(:) = 0
205      !
206      ALLOCATE( albedo_oce_mix(jpi,jpj), nrcvinfo(jprcv),  STAT=ierr(1) )
207     
208#if ! defined key_lim3 && ! defined key_lim2 && ! defined key_cice
209      ALLOCATE( a_i(jpi,jpj,1) , STAT=ierr(2) )  ! used in sbcice_if.F90 (done here as there is no sbc_ice_if_init)
210#endif
211      ALLOCATE( xcplmask(jpi,jpj,0:nn_cplmodel) , STAT=ierr(3) )
212      !
213      IF( .NOT. ln_apr_dyn ) ALLOCATE( ssh_ib(jpi,jpj), ssh_ibb(jpi,jpj), apr(jpi, jpj), STAT=ierr(4) ) 
214
215      sbc_cpl_alloc = MAXVAL( ierr )
216      IF( lk_mpp            )   CALL mpp_sum ( sbc_cpl_alloc )
217      IF( sbc_cpl_alloc > 0 )   CALL ctl_warn('sbc_cpl_alloc: allocation of arrays failed')
218      !
219   END FUNCTION sbc_cpl_alloc
220
221
222   SUBROUTINE sbc_cpl_init( k_ice )     
223      !!----------------------------------------------------------------------
224      !!             ***  ROUTINE sbc_cpl_init  ***
225      !!
226      !! ** Purpose :   Initialisation of send and received information from
227      !!                the atmospheric component
228      !!
229      !! ** Method  : * Read namsbc_cpl namelist
230      !!              * define the receive interface
231      !!              * define the send    interface
232      !!              * initialise the OASIS coupler
233      !!----------------------------------------------------------------------
234      INTEGER, INTENT(in) ::   k_ice       ! ice management in the sbc (=0/1/2/3)
235      !!
236      INTEGER ::   jn   ! dummy loop index
237      INTEGER ::   ios  ! Local integer output status for namelist read
238      INTEGER ::   inum 
239      REAL(wp), POINTER, DIMENSION(:,:) ::   zacs, zaos
240      !!
241      NAMELIST/namsbc_cpl/  sn_snd_temp , sn_snd_alb  , sn_snd_thick , sn_snd_crt   , sn_snd_co2,      & 
242         &                  sn_rcv_w10m, sn_rcv_taumod, sn_rcv_tau   , sn_rcv_dqnsdt, sn_rcv_qsr,      & 
243         &                  sn_snd_ifrac, sn_snd_crtw , sn_snd_wlev  , sn_rcv_hsig  , sn_rcv_phioc ,   & 
244         &                  sn_rcv_sdrfx, sn_rcv_sdrfy, sn_rcv_wper  , sn_rcv_wnum  , sn_rcv_wstrf ,   &
245         &                  sn_rcv_wdrag, sn_rcv_qns  , sn_rcv_emp   , sn_rcv_rnf   , sn_rcv_cal   ,   &
246         &                  sn_rcv_iceflx,sn_rcv_co2  , nn_cplmodel  , ln_usecplmask, sn_rcv_mslp 
247      !!---------------------------------------------------------------------
248      !
249      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_init')
250      !
251      CALL wrk_alloc( jpi,jpj, zacs, zaos )
252
253      ! ================================ !
254      !      Namelist informations       !
255      ! ================================ !
256
257      REWIND( numnam_ref )              ! Namelist namsbc_cpl in reference namelist : Variables for OASIS coupling
258      READ  ( numnam_ref, namsbc_cpl, IOSTAT = ios, ERR = 901)
259901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_cpl in reference namelist', lwp )
260
261      REWIND( numnam_cfg )              ! Namelist namsbc_cpl in configuration namelist : Variables for OASIS coupling
262      READ  ( numnam_cfg, namsbc_cpl, IOSTAT = ios, ERR = 902 )
263902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_cpl in configuration namelist', lwp )
264      IF(lwm) WRITE ( numond, namsbc_cpl )
265
266      IF(lwp) THEN                        ! control print
267         WRITE(numout,*)
268         WRITE(numout,*)'sbc_cpl_init : namsbc_cpl namelist '
269         WRITE(numout,*)'~~~~~~~~~~~~'
270      ENDIF
271      IF( lwp .AND. ln_cpl ) THEN                        ! control print
272         WRITE(numout,*)'  received fields (mutiple ice categogies)'
273         WRITE(numout,*)'      10m wind module                 = ', TRIM(sn_rcv_w10m%cldes  ), ' (', TRIM(sn_rcv_w10m%clcat  ), ')'
274         WRITE(numout,*)'      stress module                   = ', TRIM(sn_rcv_taumod%cldes), ' (', TRIM(sn_rcv_taumod%clcat), ')'
275         WRITE(numout,*)'      surface stress                  = ', TRIM(sn_rcv_tau%cldes   ), ' (', TRIM(sn_rcv_tau%clcat   ), ')'
276         WRITE(numout,*)'                     - referential    = ', sn_rcv_tau%clvref
277         WRITE(numout,*)'                     - orientation    = ', sn_rcv_tau%clvor
278         WRITE(numout,*)'                     - mesh           = ', sn_rcv_tau%clvgrd
279         WRITE(numout,*)'      non-solar heat flux sensitivity = ', TRIM(sn_rcv_dqnsdt%cldes), ' (', TRIM(sn_rcv_dqnsdt%clcat), ')'
280         WRITE(numout,*)'      solar heat flux                 = ', TRIM(sn_rcv_qsr%cldes   ), ' (', TRIM(sn_rcv_qsr%clcat   ), ')'
281         WRITE(numout,*)'      non-solar heat flux             = ', TRIM(sn_rcv_qns%cldes   ), ' (', TRIM(sn_rcv_qns%clcat   ), ')'
282         WRITE(numout,*)'      freshwater budget               = ', TRIM(sn_rcv_emp%cldes   ), ' (', TRIM(sn_rcv_emp%clcat   ), ')'
283         WRITE(numout,*)'      runoffs                         = ', TRIM(sn_rcv_rnf%cldes   ), ' (', TRIM(sn_rcv_rnf%clcat   ), ')'
284         WRITE(numout,*)'      calving                         = ', TRIM(sn_rcv_cal%cldes   ), ' (', TRIM(sn_rcv_cal%clcat   ), ')'
285         WRITE(numout,*)'      sea ice heat fluxes             = ', TRIM(sn_rcv_iceflx%cldes), ' (', TRIM(sn_rcv_iceflx%clcat), ')'
286         WRITE(numout,*)'      atm co2                         = ', TRIM(sn_rcv_co2%cldes   ), ' (', TRIM(sn_rcv_co2%clcat   ), ')'
287         WRITE(numout,*)'      significant wave heigth         = ', TRIM(sn_rcv_hsig%cldes  ), ' (', TRIM(sn_rcv_hsig%clcat  ), ')' 
288         WRITE(numout,*)'      wave to oce energy flux         = ', TRIM(sn_rcv_phioc%cldes ), ' (', TRIM(sn_rcv_phioc%clcat ), ')' 
289         WRITE(numout,*)'      Surface Stokes drift grid u     = ', TRIM(sn_rcv_sdrfx%cldes ), ' (', TRIM(sn_rcv_sdrfx%clcat ), ')' 
290         WRITE(numout,*)'      Surface Stokes drift grid v     = ', TRIM(sn_rcv_sdrfy%cldes ), ' (', TRIM(sn_rcv_sdrfy%clcat ), ')' 
291         WRITE(numout,*)'      Mean wave period                = ', TRIM(sn_rcv_wper%cldes  ), ' (', TRIM(sn_rcv_wper%clcat  ), ')' 
292         WRITE(numout,*)'      Mean wave number                = ', TRIM(sn_rcv_wnum%cldes  ), ' (', TRIM(sn_rcv_wnum%clcat  ), ')' 
293         WRITE(numout,*)'      Stress frac adsorbed by waves   = ', TRIM(sn_rcv_wstrf%cldes ), ' (', TRIM(sn_rcv_wstrf%clcat ), ')' 
294         WRITE(numout,*)'      Neutral surf drag coefficient   = ', TRIM(sn_rcv_wdrag%cldes ), ' (', TRIM(sn_rcv_wdrag%clcat ), ')' 
295         WRITE(numout,*)'  sent fields (multiple ice categories)'
296         WRITE(numout,*)'      surface temperature             = ', TRIM(sn_snd_temp%cldes  ), ' (', TRIM(sn_snd_temp%clcat  ), ')'
297         WRITE(numout,*)'      albedo                          = ', TRIM(sn_snd_alb%cldes   ), ' (', TRIM(sn_snd_alb%clcat   ), ')'
298         WRITE(numout,*)'      ice/snow thickness              = ', TRIM(sn_snd_thick%cldes ), ' (', TRIM(sn_snd_thick%clcat ), ')'
299         WRITE(numout,*)'      total ice fraction              = ', TRIM(sn_snd_ifrac%cldes ), ' (', TRIM(sn_snd_ifrac%clcat ), ')' 
300         WRITE(numout,*)'      surface current                 = ', TRIM(sn_snd_crt%cldes   ), ' (', TRIM(sn_snd_crt%clcat   ), ')'
301         WRITE(numout,*)'                      - referential   = ', sn_snd_crt%clvref 
302         WRITE(numout,*)'                      - orientation   = ', sn_snd_crt%clvor
303         WRITE(numout,*)'                      - mesh          = ', sn_snd_crt%clvgrd
304         WRITE(numout,*)'      oce co2 flux                    = ', TRIM(sn_snd_co2%cldes   ), ' (', TRIM(sn_snd_co2%clcat   ), ')'
305         WRITE(numout,*)'      water level                     = ', TRIM(sn_snd_wlev%cldes  ), ' (', TRIM(sn_snd_wlev%clcat  ), ')' 
306         WRITE(numout,*)'      mean sea level pressure         = ', TRIM(sn_rcv_mslp%cldes  ), ' (', TRIM(sn_rcv_mslp%clcat  ), ')' 
307         WRITE(numout,*)'      surface current to waves        = ', TRIM(sn_snd_crtw%cldes   ), ' (', TRIM(sn_snd_crtw%clcat   ), ')' 
308         WRITE(numout,*)'                      - referential   = ', sn_snd_crtw%clvref 
309         WRITE(numout,*)'                      - orientation   = ', sn_snd_crtw%clvor 
310         WRITE(numout,*)'                      - mesh          = ', sn_snd_crtw%clvgrd 
311         WRITE(numout,*)'  nn_cplmodel                         = ', nn_cplmodel
312         WRITE(numout,*)'  ln_usecplmask                       = ', ln_usecplmask
313      ENDIF
314
315      !                                   ! allocate sbccpl arrays
316      IF( sbc_cpl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_cpl_alloc : unable to allocate arrays' )
317     
318      ! ================================ !
319      !   Define the receive interface   !
320      ! ================================ !
321      nrcvinfo(:) = OASIS_idle   ! needed by nrcvinfo(jpr_otx1) if we do not receive ocean stress
322
323      ! for each field: define the OASIS name                              (srcv(:)%clname)
324      !                 define receive or not from the namelist parameters (srcv(:)%laction)
325      !                 define the north fold type of lbc                  (srcv(:)%nsgn)
326
327      ! default definitions of srcv
328      srcv(:)%laction = .FALSE.   ;   srcv(:)%clgrid = 'T'   ;   srcv(:)%nsgn = 1.   ;   srcv(:)%nct = 1
329
330      !                                                      ! ------------------------- !
331      !                                                      ! ice and ocean wind stress !   
332      !                                                      ! ------------------------- !
333      !                                                           ! Name
334      srcv(jpr_otx1)%clname = 'O_OTaux1'      ! 1st ocean component on grid ONE (T or U)
335      srcv(jpr_oty1)%clname = 'O_OTauy1'      ! 2nd   -      -         -     -
336      srcv(jpr_otz1)%clname = 'O_OTauz1'      ! 3rd   -      -         -     -
337      srcv(jpr_otx2)%clname = 'O_OTaux2'      ! 1st ocean component on grid TWO (V)
338      srcv(jpr_oty2)%clname = 'O_OTauy2'      ! 2nd   -      -         -     -
339      srcv(jpr_otz2)%clname = 'O_OTauz2'      ! 3rd   -      -         -     -
340      !
341      srcv(jpr_itx1)%clname = 'O_ITaux1'      ! 1st  ice  component on grid ONE (T, F, I or U)
342      srcv(jpr_ity1)%clname = 'O_ITauy1'      ! 2nd   -      -         -     -
343      srcv(jpr_itz1)%clname = 'O_ITauz1'      ! 3rd   -      -         -     -
344      srcv(jpr_itx2)%clname = 'O_ITaux2'      ! 1st  ice  component on grid TWO (V)
345      srcv(jpr_ity2)%clname = 'O_ITauy2'      ! 2nd   -      -         -     -
346      srcv(jpr_itz2)%clname = 'O_ITauz2'      ! 3rd   -      -         -     -
347      !
348      ! Vectors: change of sign at north fold ONLY if on the local grid
349      IF( TRIM( sn_rcv_tau%clvor ) == 'local grid' )   srcv(jpr_otx1:jpr_itz2)%nsgn = -1.
350     
351      !                                                           ! Set grid and action
352      SELECT CASE( TRIM( sn_rcv_tau%clvgrd ) )      !  'T', 'U,V', 'U,V,I', 'U,V,F', 'T,I', 'T,F', or 'T,U,V'
353      CASE( 'T' ) 
354         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
355         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
356         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
357      CASE( 'U,V' ) 
358         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
359         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
360         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'U'        ! ice components given at U-point
361         srcv(jpr_itx2:jpr_itz2)%clgrid  = 'V'        !           and           V-point
362         srcv(jpr_otx1:jpr_itz2)%laction = .TRUE.     ! receive oce and ice components on both grid 1 & 2
363      CASE( 'U,V,T' )
364         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
365         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
366         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'T'        ! ice components given at T-point
367         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
368         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
369      CASE( 'U,V,I' )
370         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
371         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
372         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'I'        ! ice components given at I-point
373         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
374         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
375      CASE( 'U,V,F' )
376         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
377         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
378         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'F'        ! ice components given at F-point
379         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
380         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
381      CASE( 'T,I' ) 
382         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
383         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'I'        ! ice components given at I-point
384         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
385         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
386      CASE( 'T,F' ) 
387         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
388         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'F'        ! ice components given at F-point
389         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
390         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
391      CASE( 'T,U,V' )
392         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'T'        ! oce components given at T-point
393         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'U'        ! ice components given at U-point
394         srcv(jpr_itx2:jpr_itz2)%clgrid  = 'V'        !           and           V-point
395         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1 only
396         srcv(jpr_itx1:jpr_itz2)%laction = .TRUE.     ! receive ice components on grid 1 & 2
397      CASE default   
398         CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_rcv_tau%clvgrd' )
399      END SELECT
400      !
401      IF( TRIM( sn_rcv_tau%clvref ) == 'spherical' )   &           ! spherical: 3rd component not received
402         &     srcv( (/jpr_otz1, jpr_otz2, jpr_itz1, jpr_itz2/) )%laction = .FALSE. 
403      !
404      IF( TRIM( sn_rcv_tau%clvor  ) == 'local grid' ) THEN        ! already on local grid -> no need of the second grid
405            srcv(jpr_otx2:jpr_otz2)%laction = .FALSE. 
406            srcv(jpr_itx2:jpr_itz2)%laction = .FALSE. 
407            srcv(jpr_oty1)%clgrid = srcv(jpr_oty2)%clgrid   ! not needed but cleaner...
408            srcv(jpr_ity1)%clgrid = srcv(jpr_ity2)%clgrid   ! not needed but cleaner...
409      ENDIF
410      !
411      IF( TRIM( sn_rcv_tau%cldes ) /= 'oce and ice' ) THEN        ! 'oce and ice' case ocean stress on ocean mesh used
412         srcv(jpr_itx1:jpr_itz2)%laction = .FALSE.    ! ice components not received
413         srcv(jpr_itx1)%clgrid = 'U'                  ! ocean stress used after its transformation
414         srcv(jpr_ity1)%clgrid = 'V'                  ! i.e. it is always at U- & V-points for i- & j-comp. resp.
415      ENDIF
416       
417      !                                                      ! ------------------------- !
418      !                                                      !    freshwater budget      !   E-P
419      !                                                      ! ------------------------- !
420      ! we suppose that atmosphere modele do not make the difference between precipiration (liquide or solid)
421      ! over ice of free ocean within the same atmospheric cell.cd
422      srcv(jpr_rain)%clname = 'OTotRain'      ! Rain = liquid precipitation
423      srcv(jpr_snow)%clname = 'OTotSnow'      ! Snow = solid precipitation
424      srcv(jpr_tevp)%clname = 'OTotEvap'      ! total evaporation (over oce + ice sublimation)
425      srcv(jpr_ievp)%clname = 'OIceEvap'      ! evaporation over ice = sublimation
426      srcv(jpr_sbpr)%clname = 'OSubMPre'      ! sublimation - liquid precipitation - solid precipitation
427      srcv(jpr_semp)%clname = 'OISubMSn'      ! ice solid water budget = sublimation - solid precipitation
428      srcv(jpr_oemp)%clname = 'OOEvaMPr'      ! ocean water budget = ocean Evap - ocean precip
429      SELECT CASE( TRIM( sn_rcv_emp%cldes ) )
430      CASE( 'none'          )       ! nothing to do
431      CASE( 'oce only'      )   ;   srcv(                                 jpr_oemp   )%laction = .TRUE. 
432      CASE( 'conservative'  )
433         srcv( (/jpr_rain, jpr_snow, jpr_ievp, jpr_tevp/) )%laction = .TRUE.
434         IF ( k_ice <= 1 )  srcv(jpr_ievp)%laction = .FALSE.
435      CASE( 'oce and ice'   )   ;   srcv( (/jpr_ievp, jpr_sbpr, jpr_semp, jpr_oemp/) )%laction = .TRUE.
436      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_rcv_emp%cldes' )
437      END SELECT
438
439      !                                                      ! ------------------------- !
440      !                                                      !     Runoffs & Calving     !   
441      !                                                      ! ------------------------- !
442      srcv(jpr_rnf   )%clname = 'O_Runoff'
443      IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) THEN
444         srcv(jpr_rnf)%laction = .TRUE.
445         l_rnfcpl              = .TRUE.                      ! -> no need to read runoffs in sbcrnf
446         ln_rnf                = nn_components /= jp_iam_sas ! -> force to go through sbcrnf if not sas
447         IF(lwp) WRITE(numout,*)
448         IF(lwp) WRITE(numout,*) '   runoffs received from oasis -> force ln_rnf = ', ln_rnf
449      ENDIF
450      !
451      srcv(jpr_cal   )%clname = 'OCalving'   ;   IF( TRIM( sn_rcv_cal%cldes ) == 'coupled' )   srcv(jpr_cal)%laction = .TRUE.
452
453      !                                                      ! ------------------------- !
454      !                                                      !    non solar radiation    !   Qns
455      !                                                      ! ------------------------- !
456      srcv(jpr_qnsoce)%clname = 'O_QnsOce'
457      srcv(jpr_qnsice)%clname = 'O_QnsIce'
458      srcv(jpr_qnsmix)%clname = 'O_QnsMix'
459      SELECT CASE( TRIM( sn_rcv_qns%cldes ) )
460      CASE( 'none'          )       ! nothing to do
461      CASE( 'oce only'      )   ;   srcv(               jpr_qnsoce   )%laction = .TRUE.
462      CASE( 'conservative'  )   ;   srcv( (/jpr_qnsice, jpr_qnsmix/) )%laction = .TRUE.
463      CASE( 'oce and ice'   )   ;   srcv( (/jpr_qnsice, jpr_qnsoce/) )%laction = .TRUE.
464      CASE( 'mixed oce-ice' )   ;   srcv(               jpr_qnsmix   )%laction = .TRUE. 
465      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_rcv_qns%cldes' )
466      END SELECT
467      IF( TRIM( sn_rcv_qns%cldes ) == 'mixed oce-ice' .AND. jpl > 1 ) &
468         CALL ctl_stop( 'sbc_cpl_init: sn_rcv_qns%cldes not currently allowed to be mixed oce-ice for multi-category ice' )
469      !                                                      ! ------------------------- !
470      !                                                      !    solar radiation        !   Qsr
471      !                                                      ! ------------------------- !
472      srcv(jpr_qsroce)%clname = 'O_QsrOce'
473      srcv(jpr_qsrice)%clname = 'O_QsrIce'
474      srcv(jpr_qsrmix)%clname = 'O_QsrMix'
475      SELECT CASE( TRIM( sn_rcv_qsr%cldes ) )
476      CASE( 'none'          )       ! nothing to do
477      CASE( 'oce only'      )   ;   srcv(               jpr_qsroce   )%laction = .TRUE.
478      CASE( 'conservative'  )   ;   srcv( (/jpr_qsrice, jpr_qsrmix/) )%laction = .TRUE.
479      CASE( 'oce and ice'   )   ;   srcv( (/jpr_qsrice, jpr_qsroce/) )%laction = .TRUE.
480      CASE( 'mixed oce-ice' )   ;   srcv(               jpr_qsrmix   )%laction = .TRUE. 
481      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_rcv_qsr%cldes' )
482      END SELECT
483      IF( TRIM( sn_rcv_qsr%cldes ) == 'mixed oce-ice' .AND. jpl > 1 ) &
484         CALL ctl_stop( 'sbc_cpl_init: sn_rcv_qsr%cldes not currently allowed to be mixed oce-ice for multi-category ice' )
485      !                                                      ! ------------------------- !
486      !                                                      !   non solar sensitivity   !   d(Qns)/d(T)
487      !                                                      ! ------------------------- !
488      srcv(jpr_dqnsdt)%clname = 'O_dQnsdT'   
489      IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'coupled' )   srcv(jpr_dqnsdt)%laction = .TRUE.
490      !
491      ! non solar sensitivity mandatory for LIM ice model
492      IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 .AND. nn_components /= jp_iam_sas ) &
493         CALL ctl_stop( 'sbc_cpl_init: sn_rcv_dqnsdt%cldes must be coupled in namsbc_cpl namelist' )
494      ! non solar sensitivity mandatory for mixed oce-ice solar radiation coupling technique
495      IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. TRIM( sn_rcv_qns%cldes ) == 'mixed oce-ice' ) &
496         CALL ctl_stop( 'sbc_cpl_init: namsbc_cpl namelist mismatch between sn_rcv_qns%cldes and sn_rcv_dqnsdt%cldes' )
497      !                                                      ! ------------------------- !
498      !                                                      !      10m wind module      !   
499      !                                                      ! ------------------------- !
500      srcv(jpr_w10m)%clname = 'O_Wind10'   ;   IF( TRIM(sn_rcv_w10m%cldes  ) == 'coupled' )   srcv(jpr_w10m)%laction = .TRUE. 
501      !
502      !                                                      ! ------------------------- !
503      !                                                      !   wind stress module      !   
504      !                                                      ! ------------------------- !
505      srcv(jpr_taum)%clname = 'O_TauMod'   ;   IF( TRIM(sn_rcv_taumod%cldes) == 'coupled' )   srcv(jpr_taum)%laction = .TRUE.
506      lhftau = srcv(jpr_taum)%laction
507
508      !                                                      ! ------------------------- !
509      !                                                      !      Atmospheric CO2      !
510      !                                                      ! ------------------------- !
511      srcv(jpr_co2 )%clname = 'O_AtmCO2'   ;   IF( TRIM(sn_rcv_co2%cldes   ) == 'coupled' )    srcv(jpr_co2 )%laction = .TRUE.
512
513      !                                                      ! ------------------------- !
514      !                                                      ! Mean Sea Level Pressure   !
515      !                                                      ! ------------------------- !
516      srcv(jpr_mslp)%clname = 'O_MSLP'     ;   IF( TRIM(sn_rcv_mslp%cldes  ) == 'coupled' )    srcv(jpr_mslp)%laction = .TRUE. 
517
518      !                                                      ! ------------------------- !
519      !                                                      !   topmelt and botmelt     !   
520      !                                                      ! ------------------------- !
521      srcv(jpr_topm )%clname = 'OTopMlt'
522      srcv(jpr_botm )%clname = 'OBotMlt'
523      IF( TRIM(sn_rcv_iceflx%cldes) == 'coupled' ) THEN
524         IF ( TRIM( sn_rcv_iceflx%clcat ) == 'yes' ) THEN
525            srcv(jpr_topm:jpr_botm)%nct = jpl
526         ELSE
527            CALL ctl_stop( 'sbc_cpl_init: sn_rcv_iceflx%clcat should always be set to yes currently' )
528         ENDIF
529         srcv(jpr_topm:jpr_botm)%laction = .TRUE.
530      ENDIF
531      !                                                      ! ------------------------- !
532      !                                                      !      Wave breaking        !   
533      !                                                      ! ------------------------- !
534      srcv(jpr_hsig)%clname  = 'O_Hsigwa'    ! significant wave height
535      IF( TRIM(sn_rcv_hsig%cldes  ) == 'coupled' )  THEN
536         srcv(jpr_hsig)%laction = .TRUE.
537         cpl_hsig = .TRUE.
538      ENDIF
539      srcv(jpr_phioc)%clname = 'O_PhiOce'    ! wave to ocean energy
540      IF( TRIM(sn_rcv_phioc%cldes ) == 'coupled' )  THEN
541         srcv(jpr_phioc)%laction = .TRUE.
542         cpl_phioc = .TRUE.
543      ENDIF
544      srcv(jpr_sdrftx)%clname = 'O_Sdrfx'    ! Stokes drift in the u direction
545      IF( TRIM(sn_rcv_sdrfx%cldes ) == 'coupled' )  THEN
546         srcv(jpr_sdrftx)%laction = .TRUE.
547         cpl_sdrftx = .TRUE.
548      ENDIF
549      srcv(jpr_sdrfty)%clname = 'O_Sdrfy'    ! Stokes drift in the v direction
550      IF( TRIM(sn_rcv_sdrfy%cldes ) == 'coupled' )  THEN
551         srcv(jpr_sdrfty)%laction = .TRUE.
552         cpl_sdrfty = .TRUE.
553      ENDIF
554      srcv(jpr_wper)%clname = 'O_WPer'       ! mean wave period
555      IF( TRIM(sn_rcv_wper%cldes  ) == 'coupled' )  THEN
556         srcv(jpr_wper)%laction = .TRUE.
557         cpl_wper = .TRUE.
558      ENDIF
559      srcv(jpr_wnum)%clname = 'O_WNum'       ! mean wave number
560      IF( TRIM(sn_rcv_wnum%cldes ) == 'coupled' )  THEN
561         srcv(jpr_wnum)%laction = .TRUE.
562         cpl_wnum = .TRUE.
563      ENDIF
564      srcv(jpr_wstrf)%clname = 'O_WStrf'     ! stress fraction adsorbed by the wave
565      IF( TRIM(sn_rcv_wstrf%cldes ) == 'coupled' )  THEN
566         srcv(jpr_wstrf)%laction = .TRUE.
567         cpl_wstrf = .TRUE.
568      ENDIF
569      srcv(jpr_wdrag)%clname = 'O_WDrag'     ! neutral surface drag coefficient
570      IF( TRIM(sn_rcv_wdrag%cldes ) == 'coupled' )  THEN
571         srcv(jpr_wdrag)%laction = .TRUE.
572         cpl_wdrag = .TRUE.
573      ENDIF
574      !
575      !                                                      ! ------------------------------- !
576      !                                                      !   OPA-SAS coupling - rcv by opa !   
577      !                                                      ! ------------------------------- !
578      srcv(jpr_sflx)%clname = 'O_SFLX'
579      srcv(jpr_fice)%clname = 'RIceFrc'
580      !
581      IF( nn_components == jp_iam_opa ) THEN    ! OPA coupled to SAS via OASIS: force received field by OPA (sent by SAS)
582         srcv(:)%laction = .FALSE.   ! force default definition in case of opa <-> sas coupling
583         srcv(:)%clgrid  = 'T'       ! force default definition in case of opa <-> sas coupling
584         srcv(:)%nsgn    = 1.        ! force default definition in case of opa <-> sas coupling
585         srcv( (/jpr_qsroce, jpr_qnsoce, jpr_oemp, jpr_sflx, jpr_fice, jpr_otx1, jpr_oty1, jpr_taum/) )%laction = .TRUE.
586         srcv(jpr_otx1)%clgrid = 'U'        ! oce components given at U-point
587         srcv(jpr_oty1)%clgrid = 'V'        !           and           V-point
588         ! Vectors: change of sign at north fold ONLY if on the local grid
589         srcv( (/jpr_otx1,jpr_oty1/) )%nsgn = -1.
590         sn_rcv_tau%clvgrd = 'U,V'
591         sn_rcv_tau%clvor = 'local grid'
592         sn_rcv_tau%clvref = 'spherical'
593         sn_rcv_emp%cldes = 'oce only'
594         !
595         IF(lwp) THEN                        ! control print
596            WRITE(numout,*)
597            WRITE(numout,*)'               Special conditions for SAS-OPA coupling  '
598            WRITE(numout,*)'               OPA component  '
599            WRITE(numout,*)
600            WRITE(numout,*)'  received fields from SAS component '
601            WRITE(numout,*)'                  ice cover '
602            WRITE(numout,*)'                  oce only EMP  '
603            WRITE(numout,*)'                  salt flux  '
604            WRITE(numout,*)'                  mixed oce-ice solar flux  '
605            WRITE(numout,*)'                  mixed oce-ice non solar flux  '
606            WRITE(numout,*)'                  wind stress U,V on local grid and sperical coordinates '
607            WRITE(numout,*)'                  wind stress module'
608            WRITE(numout,*)
609         ENDIF
610      ENDIF
611      !                                                      ! -------------------------------- !
612      !                                                      !   OPA-SAS coupling - rcv by sas  !   
613      !                                                      ! -------------------------------- !
614      srcv(jpr_toce  )%clname = 'I_SSTSST'
615      srcv(jpr_soce  )%clname = 'I_SSSal'
616      srcv(jpr_ocx1  )%clname = 'I_OCurx1'
617      srcv(jpr_ocy1  )%clname = 'I_OCury1'
618      srcv(jpr_ssh   )%clname = 'I_SSHght'
619      srcv(jpr_e3t1st)%clname = 'I_E3T1st'   
620      srcv(jpr_fraqsr)%clname = 'I_FraQsr'   
621      !
622      IF( nn_components == jp_iam_sas ) THEN
623         IF( .NOT. ln_cpl ) srcv(:)%laction = .FALSE.   ! force default definition in case of opa <-> sas coupling
624         IF( .NOT. ln_cpl ) srcv(:)%clgrid  = 'T'       ! force default definition in case of opa <-> sas coupling
625         IF( .NOT. ln_cpl ) srcv(:)%nsgn    = 1.        ! force default definition in case of opa <-> sas coupling
626         srcv( (/jpr_toce, jpr_soce, jpr_ssh, jpr_fraqsr, jpr_ocx1, jpr_ocy1/) )%laction = .TRUE.
627         srcv( jpr_e3t1st )%laction = lk_vvl
628         srcv(jpr_ocx1)%clgrid = 'U'        ! oce components given at U-point
629         srcv(jpr_ocy1)%clgrid = 'V'        !           and           V-point
630         ! Vectors: change of sign at north fold ONLY if on the local grid
631         srcv(jpr_ocx1:jpr_ocy1)%nsgn = -1.
632         ! Change first letter to couple with atmosphere if already coupled OPA
633         ! this is nedeed as each variable name used in the namcouple must be unique:
634         ! for example O_Runoff received by OPA from SAS and therefore O_Runoff received by SAS from the Atmosphere
635         DO jn = 1, jprcv
636            IF ( srcv(jn)%clname(1:1) == "O" ) srcv(jn)%clname = "S"//srcv(jn)%clname(2:LEN(srcv(jn)%clname))
637         END DO
638         !
639         IF(lwp) THEN                        ! control print
640            WRITE(numout,*)
641            WRITE(numout,*)'               Special conditions for SAS-OPA coupling  '
642            WRITE(numout,*)'               SAS component  '
643            WRITE(numout,*)
644            IF( .NOT. ln_cpl ) THEN
645               WRITE(numout,*)'  received fields from OPA component '
646            ELSE
647               WRITE(numout,*)'  Additional received fields from OPA component : '
648            ENDIF
649            WRITE(numout,*)'               sea surface temperature (Celcius) '
650            WRITE(numout,*)'               sea surface salinity ' 
651            WRITE(numout,*)'               surface currents ' 
652            WRITE(numout,*)'               sea surface height ' 
653            WRITE(numout,*)'               thickness of first ocean T level '       
654            WRITE(numout,*)'               fraction of solar net radiation absorbed in the first ocean level'
655            WRITE(numout,*)
656         ENDIF
657      ENDIF
658     
659      ! =================================================== !
660      ! Allocate all parts of frcv used for received fields !
661      ! =================================================== !
662      DO jn = 1, jprcv
663         IF ( srcv(jn)%laction ) ALLOCATE( frcv(jn)%z3(jpi,jpj,srcv(jn)%nct) )
664      END DO
665      ! Allocate taum part of frcv which is used even when not received as coupling field
666      IF ( .NOT. srcv(jpr_taum)%laction ) ALLOCATE( frcv(jpr_taum)%z3(jpi,jpj,srcv(jpr_taum)%nct) )
667      ! Allocate w10m part of frcv which is used even when not received as coupling field
668      IF ( .NOT. srcv(jpr_w10m)%laction ) ALLOCATE( frcv(jpr_w10m)%z3(jpi,jpj,srcv(jpr_w10m)%nct) )
669      ! Allocate jpr_otx1 part of frcv which is used even when not received as coupling field
670      IF ( .NOT. srcv(jpr_otx1)%laction ) ALLOCATE( frcv(jpr_otx1)%z3(jpi,jpj,srcv(jpr_otx1)%nct) )
671      IF ( .NOT. srcv(jpr_oty1)%laction ) ALLOCATE( frcv(jpr_oty1)%z3(jpi,jpj,srcv(jpr_oty1)%nct) )
672      ! Allocate itx1 and ity1 as they are used in sbc_cpl_ice_tau even if srcv(jpr_itx1)%laction = .FALSE.
673      IF( k_ice /= 0 ) THEN
674         IF ( .NOT. srcv(jpr_itx1)%laction ) ALLOCATE( frcv(jpr_itx1)%z3(jpi,jpj,srcv(jpr_itx1)%nct) )
675         IF ( .NOT. srcv(jpr_ity1)%laction ) ALLOCATE( frcv(jpr_ity1)%z3(jpi,jpj,srcv(jpr_ity1)%nct) )
676      END IF
677
678      ! ================================ !
679      !     Define the send interface    !
680      ! ================================ !
681      ! for each field: define the OASIS name                           (ssnd(:)%clname)
682      !                 define send or not from the namelist parameters (ssnd(:)%laction)
683      !                 define the north fold type of lbc               (ssnd(:)%nsgn)
684     
685      ! default definitions of nsnd
686      ssnd(:)%laction = .FALSE.   ;   ssnd(:)%clgrid = 'T'   ;   ssnd(:)%nsgn = 1.  ; ssnd(:)%nct = 1
687         
688      !                                                      ! ------------------------- !
689      !                                                      !    Surface temperature    !
690      !                                                      ! ------------------------- !
691      ssnd(jps_toce)%clname = 'O_SSTSST'
692      ssnd(jps_tice)%clname = 'O_TepIce'
693      ssnd(jps_tmix)%clname = 'O_TepMix'
694      SELECT CASE( TRIM( sn_snd_temp%cldes ) )
695      CASE( 'none'                                 )       ! nothing to do
696      CASE( 'oce only'                             )   ;   ssnd( jps_toce )%laction = .TRUE.
697      CASE( 'oce and ice' , 'weighted oce and ice' )
698         ssnd( (/jps_toce, jps_tice/) )%laction = .TRUE.
699         IF ( TRIM( sn_snd_temp%clcat ) == 'yes' )  ssnd(jps_tice)%nct = jpl
700      CASE( 'mixed oce-ice'                        )   ;   ssnd( jps_tmix )%laction = .TRUE.
701      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_temp%cldes' )
702      END SELECT
703           
704      !                                                      ! ------------------------- !
705      !                                                      !          Albedo           !
706      !                                                      ! ------------------------- !
707      ssnd(jps_albice)%clname = 'O_AlbIce' 
708      ssnd(jps_albmix)%clname = 'O_AlbMix'
709      SELECT CASE( TRIM( sn_snd_alb%cldes ) )
710      CASE( 'none'                 )     ! nothing to do
711      CASE( 'ice' , 'weighted ice' )   ; ssnd(jps_albice)%laction = .TRUE.
712      CASE( 'mixed oce-ice'        )   ; ssnd(jps_albmix)%laction = .TRUE.
713      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_alb%cldes' )
714      END SELECT
715      !
716      ! Need to calculate oceanic albedo if
717      !     1. sending mixed oce-ice albedo or
718      !     2. receiving mixed oce-ice solar radiation
719      IF ( TRIM ( sn_snd_alb%cldes ) == 'mixed oce-ice' .OR. TRIM ( sn_rcv_qsr%cldes ) == 'mixed oce-ice' ) THEN
720         CALL albedo_oce( zaos, zacs )
721         ! Due to lack of information on nebulosity : mean clear/overcast sky
722         albedo_oce_mix(:,:) = ( zacs(:,:) + zaos(:,:) ) * 0.5
723      ENDIF
724
725      !                                                      ! ------------------------- !
726      !                                                      !  Ice fraction & Thickness !
727      !                                                      ! ------------------------- !
728      ssnd(jps_fice)%clname = 'OIceFrc'
729      ssnd(jps_ficet)%clname = 'OIceFrcT' 
730      ssnd(jps_hice)%clname = 'OIceTck'
731      ssnd(jps_hsnw)%clname = 'OSnwTck'
732      IF( k_ice /= 0 ) THEN
733         ssnd(jps_fice)%laction = .TRUE.                  ! if ice treated in the ocean (even in climato case)
734! Currently no namelist entry to determine sending of multi-category ice fraction so use the thickness entry for now
735         IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) ssnd(jps_fice)%nct = jpl
736      ENDIF
737     
738      IF (TRIM( sn_snd_ifrac%cldes )  == 'coupled') ssnd(jps_ficet)%laction = .TRUE. 
739
740      SELECT CASE ( TRIM( sn_snd_thick%cldes ) )
741      CASE( 'none'         )       ! nothing to do
742      CASE( 'ice and snow' ) 
743         ssnd(jps_hice:jps_hsnw)%laction = .TRUE.
744         IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) THEN
745            ssnd(jps_hice:jps_hsnw)%nct = jpl
746         ENDIF
747      CASE ( 'weighted ice and snow' ) 
748         ssnd(jps_hice:jps_hsnw)%laction = .TRUE.
749         IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) ssnd(jps_hice:jps_hsnw)%nct = jpl
750      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_thick%cldes' )
751      END SELECT
752
753      !                                                      ! ------------------------- !
754      !                                                      !      Surface current      !
755      !                                                      ! ------------------------- !
756      !        ocean currents              !            ice velocities
757      ssnd(jps_ocx1)%clname = 'O_OCurx1'   ;   ssnd(jps_ivx1)%clname = 'O_IVelx1'
758      ssnd(jps_ocy1)%clname = 'O_OCury1'   ;   ssnd(jps_ivy1)%clname = 'O_IVely1'
759      ssnd(jps_ocz1)%clname = 'O_OCurz1'   ;   ssnd(jps_ivz1)%clname = 'O_IVelz1'
760      ssnd(jps_ocxw)%clname = 'O_OCurxw' 
761      ssnd(jps_ocyw)%clname = 'O_OCuryw' 
762      !
763      ssnd(jps_ocx1:jps_ivz1)%nsgn = -1.   ! vectors: change of the sign at the north fold
764
765      IF( sn_snd_crt%clvgrd == 'U,V' ) THEN
766         ssnd(jps_ocx1)%clgrid = 'U' ; ssnd(jps_ocy1)%clgrid = 'V'
767      ELSE IF( sn_snd_crt%clvgrd /= 'T' ) THEN 
768         CALL ctl_stop( 'sn_snd_crt%clvgrd must be equal to T' )
769         ssnd(jps_ocx1:jps_ivz1)%clgrid  = 'T'      ! all oce and ice components on the same unique grid
770      ENDIF
771      ssnd(jps_ocx1:jps_ivz1)%laction = .TRUE.   ! default: all are send
772      IF( TRIM( sn_snd_crt%clvref ) == 'spherical' )   ssnd( (/jps_ocz1, jps_ivz1/) )%laction = .FALSE. 
773      IF( TRIM( sn_snd_crt%clvor ) == 'eastward-northward' ) ssnd(jps_ocx1:jps_ivz1)%nsgn = 1.
774      SELECT CASE( TRIM( sn_snd_crt%cldes ) )
775      CASE( 'none'                 )   ;   ssnd(jps_ocx1:jps_ivz1)%laction = .FALSE.
776      CASE( 'oce only'             )   ;   ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE.
777      CASE( 'weighted oce and ice' )   !   nothing to do
778      CASE( 'mixed oce-ice'        )   ;   ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE.
779      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_crt%cldes' )
780      END SELECT
781
782      ssnd(jps_ocxw:jps_ocyw)%nsgn = -1.   ! vectors: change of the sign at the north fold
783       
784      IF( sn_snd_crtw%clvgrd == 'U,V' ) THEN
785         ssnd(jps_ocxw)%clgrid = 'U' ; ssnd(jps_ocyw)%clgrid = 'V' 
786      ELSE IF( sn_snd_crtw%clvgrd /= 'T' ) THEN
787         CALL ctl_stop( 'sn_snd_crtw%clvgrd must be equal to T' ) 
788      ENDIF
789      IF( TRIM( sn_snd_crtw%clvor ) == 'eastward-northward' ) ssnd(jps_ocxw:jps_ocyw)%nsgn = 1. 
790      SELECT CASE( TRIM( sn_snd_crtw%cldes ) ) 
791         CASE( 'none'                 )   ; ssnd(jps_ocxw:jps_ocyw)%laction = .FALSE. 
792         CASE( 'oce only'             )   ; ssnd(jps_ocxw:jps_ocyw)%laction = .TRUE. 
793         CASE( 'weighted oce and ice' )   !   nothing to do
794         CASE( 'mixed oce-ice'        )   ; ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE. 
795         CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_crtw%cldes' ) 
796      END SELECT 
797
798      !                                                      ! ------------------------- !
799      !                                                      !          CO2 flux         !
800      !                                                      ! ------------------------- !
801      ssnd(jps_co2)%clname = 'O_CO2FLX' ;  IF( TRIM(sn_snd_co2%cldes) == 'coupled' )    ssnd(jps_co2 )%laction = .TRUE.
802
803      !                                                      ! ------------------------- !
804      !                                                      !     Sea surface height    !
805      !                                                      ! ------------------------- !
806      ssnd(jps_wlev)%clname = 'O_Wlevel' ;  IF( TRIM(sn_snd_wlev%cldes) == 'coupled' )   ssnd(jps_wlev)%laction = .TRUE. 
807
808      !                                                      ! ------------------------------- !
809      !                                                      !   OPA-SAS coupling - snd by opa !   
810      !                                                      ! ------------------------------- !
811      ssnd(jps_ssh   )%clname = 'O_SSHght' 
812      ssnd(jps_soce  )%clname = 'O_SSSal' 
813      ssnd(jps_e3t1st)%clname = 'O_E3T1st'   
814      ssnd(jps_fraqsr)%clname = 'O_FraQsr'
815      !
816      IF( nn_components == jp_iam_opa ) THEN
817         ssnd(:)%laction = .FALSE.   ! force default definition in case of opa <-> sas coupling
818         ssnd( (/jps_toce, jps_soce, jps_ssh, jps_fraqsr, jps_ocx1, jps_ocy1/) )%laction = .TRUE.
819         ssnd( jps_e3t1st )%laction = lk_vvl
820         ! vector definition: not used but cleaner...
821         ssnd(jps_ocx1)%clgrid  = 'U'        ! oce components given at U-point
822         ssnd(jps_ocy1)%clgrid  = 'V'        !           and           V-point
823         sn_snd_crt%clvgrd = 'U,V'
824         sn_snd_crt%clvor = 'local grid'
825         sn_snd_crt%clvref = 'spherical'
826         !
827         IF(lwp) THEN                        ! control print
828            WRITE(numout,*)
829            WRITE(numout,*)'  sent fields to SAS component '
830            WRITE(numout,*)'               sea surface temperature (T before, Celcius) '
831            WRITE(numout,*)'               sea surface salinity ' 
832            WRITE(numout,*)'               surface currents U,V on local grid and spherical coordinates' 
833            WRITE(numout,*)'               sea surface height ' 
834            WRITE(numout,*)'               thickness of first ocean T level '       
835            WRITE(numout,*)'               fraction of solar net radiation absorbed in the first ocean level'
836            WRITE(numout,*)
837         ENDIF
838      ENDIF
839      !                                                      ! ------------------------------- !
840      !                                                      !   OPA-SAS coupling - snd by sas !   
841      !                                                      ! ------------------------------- !
842      ssnd(jps_sflx  )%clname = 'I_SFLX'     
843      ssnd(jps_fice2 )%clname = 'IIceFrc'
844      ssnd(jps_qsroce)%clname = 'I_QsrOce'   
845      ssnd(jps_qnsoce)%clname = 'I_QnsOce'   
846      ssnd(jps_oemp  )%clname = 'IOEvaMPr' 
847      ssnd(jps_otx1  )%clname = 'I_OTaux1'   
848      ssnd(jps_oty1  )%clname = 'I_OTauy1'   
849      ssnd(jps_rnf   )%clname = 'I_Runoff'   
850      ssnd(jps_taum  )%clname = 'I_TauMod'   
851      !
852      IF( nn_components == jp_iam_sas ) THEN
853         IF( .NOT. ln_cpl ) ssnd(:)%laction = .FALSE.   ! force default definition in case of opa <-> sas coupling
854         ssnd( (/jps_qsroce, jps_qnsoce, jps_oemp, jps_fice2, jps_sflx, jps_otx1, jps_oty1, jps_taum/) )%laction = .TRUE.
855         !
856         ! Change first letter to couple with atmosphere if already coupled with sea_ice
857         ! this is nedeed as each variable name used in the namcouple must be unique:
858         ! for example O_SSTSST sent by OPA to SAS and therefore S_SSTSST sent by SAS to the Atmosphere
859         DO jn = 1, jpsnd
860            IF ( ssnd(jn)%clname(1:1) == "O" ) ssnd(jn)%clname = "S"//ssnd(jn)%clname(2:LEN(ssnd(jn)%clname))
861         END DO
862         !
863         IF(lwp) THEN                        ! control print
864            WRITE(numout,*)
865            IF( .NOT. ln_cpl ) THEN
866               WRITE(numout,*)'  sent fields to OPA component '
867            ELSE
868               WRITE(numout,*)'  Additional sent fields to OPA component : '
869            ENDIF
870            WRITE(numout,*)'                  ice cover '
871            WRITE(numout,*)'                  oce only EMP  '
872            WRITE(numout,*)'                  salt flux  '
873            WRITE(numout,*)'                  mixed oce-ice solar flux  '
874            WRITE(numout,*)'                  mixed oce-ice non solar flux  '
875            WRITE(numout,*)'                  wind stress U,V components'
876            WRITE(numout,*)'                  wind stress module'
877         ENDIF
878      ENDIF
879
880      !
881      ! ================================ !
882      !   initialisation of the coupler  !
883      ! ================================ !
884
885      CALL cpl_define(jprcv, jpsnd, nn_cplmodel)
886     
887      IF (ln_usecplmask) THEN
888         xcplmask(:,:,:) = 0.
889         CALL iom_open( 'cplmask', inum )
890         CALL iom_get( inum, jpdom_unknown, 'cplmask', xcplmask(1:nlci,1:nlcj,1:nn_cplmodel),   &
891            &          kstart = (/ mig(1),mjg(1),1 /), kcount = (/ nlci,nlcj,nn_cplmodel /) )
892         CALL iom_close( inum )
893      ELSE
894         xcplmask(:,:,:) = 1.
895      ENDIF
896      xcplmask(:,:,0) = 1. - SUM( xcplmask(:,:,1:nn_cplmodel), dim = 3 )
897      !
898      ncpl_qsr_freq = cpl_freq( 'O_QsrOce' ) + cpl_freq( 'O_QsrMix' ) + cpl_freq( 'I_QsrOce' ) + cpl_freq( 'I_QsrMix' )
899      IF( ln_dm2dc .AND. ln_cpl .AND. ncpl_qsr_freq /= 86400 )   &
900         &   CALL ctl_stop( 'sbc_cpl_init: diurnal cycle reconstruction (ln_dm2dc) needs daily couping for solar radiation' )
901      ncpl_qsr_freq = 86400 / ncpl_qsr_freq
902
903      CALL wrk_dealloc( jpi,jpj, zacs, zaos )
904      !
905      IF( nn_timing == 1 )  CALL timing_stop('sbc_cpl_init')
906      !
907   END SUBROUTINE sbc_cpl_init
908
909
910   SUBROUTINE sbc_cpl_rcv( kt, k_fsbc, k_ice )     
911      !!----------------------------------------------------------------------
912      !!             ***  ROUTINE sbc_cpl_rcv  ***
913      !!
914      !! ** Purpose :   provide the stress over the ocean and, if no sea-ice,
915      !!                provide the ocean heat and freshwater fluxes.
916      !!
917      !! ** Method  : - Receive all the atmospheric fields (stored in frcv array). called at each time step.
918      !!                OASIS controls if there is something do receive or not. nrcvinfo contains the info
919      !!                to know if the field was really received or not
920      !!
921      !!              --> If ocean stress was really received:
922      !!
923      !!                  - transform the received ocean stress vector from the received
924      !!                 referential and grid into an atmosphere-ocean stress in
925      !!                 the (i,j) ocean referencial and at the ocean velocity point.
926      !!                    The received stress are :
927      !!                     - defined by 3 components (if cartesian coordinate)
928      !!                            or by 2 components (if spherical)
929      !!                     - oriented along geographical   coordinate (if eastward-northward)
930      !!                            or  along the local grid coordinate (if local grid)
931      !!                     - given at U- and V-point, resp.   if received on 2 grids
932      !!                            or at T-point               if received on 1 grid
933      !!                    Therefore and if necessary, they are successively
934      !!                  processed in order to obtain them
935      !!                     first  as  2 components on the sphere
936      !!                     second as  2 components oriented along the local grid
937      !!                     third  as  2 components on the U,V grid
938      !!
939      !!              -->
940      !!
941      !!              - In 'ocean only' case, non solar and solar ocean heat fluxes
942      !!             and total ocean freshwater fluxes 
943      !!
944      !! ** Method  :   receive all fields from the atmosphere and transform
945      !!              them into ocean surface boundary condition fields
946      !!
947      !! ** Action  :   update  utau, vtau   ocean stress at U,V grid
948      !!                        taum         wind stress module at T-point
949      !!                        wndm         wind speed  module at T-point over free ocean or leads in presence of sea-ice
950      !!                        qns          non solar heat fluxes including emp heat content    (ocean only case)
951      !!                                     and the latent heat flux of solid precip. melting
952      !!                        qsr          solar ocean heat fluxes   (ocean only case)
953      !!                        emp          upward mass flux [evap. - precip. (- runoffs) (- calving)] (ocean only case)
954      !!----------------------------------------------------------------------
955      USE zdf_oce,  ONLY : ln_zdfqiao
956
957      IMPLICIT NONE
958
959      INTEGER, INTENT(in)           ::   kt          ! ocean model time step index
960      INTEGER, INTENT(in)           ::   k_fsbc      ! frequency of sbc (-> ice model) computation
961      INTEGER, INTENT(in)           ::   k_ice       ! ice management in the sbc (=0/1/2/3)
962
963      !!
964      LOGICAL  ::   llnewtx, llnewtau      ! update wind stress components and module??
965      INTEGER  ::   ji, jj, jn             ! dummy loop indices
966      INTEGER  ::   isec                   ! number of seconds since nit000 (assuming rdttra did not change since nit000)
967      REAL(wp) ::   zcumulneg, zcumulpos   ! temporary scalars     
968      REAL(wp) ::   zcoef                  ! temporary scalar
969      REAL(wp) ::   zrhoa  = 1.22          ! Air density kg/m3
970      REAL(wp) ::   zcdrag = 1.5e-3        ! drag coefficient
971      REAL(wp) ::   zzx, zzy               ! temporary variables
972      REAL(wp), POINTER, DIMENSION(:,:) ::   ztx, zty, zmsk, zemp, zqns, zqsr
973      !!----------------------------------------------------------------------
974      !
975      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_rcv')
976      !
977      CALL wrk_alloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr )
978      !
979      IF( ln_mixcpl )   zmsk(:,:) = 1. - xcplmask(:,:,0)
980      !
981      !                                                      ! ======================================================= !
982      !                                                      ! Receive all the atmos. fields (including ice information)
983      !                                                      ! ======================================================= !
984      isec = ( kt - nit000 ) * NINT( rdttra(1) )                ! date of exchanges
985      DO jn = 1, jprcv                                          ! received fields sent by the atmosphere
986         IF( srcv(jn)%laction )   CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask(:,:,1:nn_cplmodel), nrcvinfo(jn) )
987      END DO
988
989      !                                                      ! ========================= !
990      IF( srcv(jpr_otx1)%laction ) THEN                      !  ocean stress components  !
991         !                                                   ! ========================= !
992         ! define frcv(jpr_otx1)%z3(:,:,1) and frcv(jpr_oty1)%z3(:,:,1): stress at U/V point along model grid
993         ! => need to be done only when we receive the field
994         IF(  nrcvinfo(jpr_otx1) == OASIS_Rcv ) THEN
995            !
996            IF( TRIM( sn_rcv_tau%clvref ) == 'cartesian' ) THEN            ! 2 components on the sphere
997               !                                                       ! (cartesian to spherical -> 3 to 2 components)
998               !
999               CALL geo2oce( frcv(jpr_otx1)%z3(:,:,1), frcv(jpr_oty1)%z3(:,:,1), frcv(jpr_otz1)%z3(:,:,1),   &
1000                  &          srcv(jpr_otx1)%clgrid, ztx, zty )
1001               frcv(jpr_otx1)%z3(:,:,1) = ztx(:,:)   ! overwrite 1st comp. on the 1st grid
1002               frcv(jpr_oty1)%z3(:,:,1) = zty(:,:)   ! overwrite 2nd comp. on the 1st grid
1003               !
1004               IF( srcv(jpr_otx2)%laction ) THEN
1005                  CALL geo2oce( frcv(jpr_otx2)%z3(:,:,1), frcv(jpr_oty2)%z3(:,:,1), frcv(jpr_otz2)%z3(:,:,1),   &
1006                     &          srcv(jpr_otx2)%clgrid, ztx, zty )
1007                  frcv(jpr_otx2)%z3(:,:,1) = ztx(:,:)   ! overwrite 1st comp. on the 2nd grid
1008                  frcv(jpr_oty2)%z3(:,:,1) = zty(:,:)   ! overwrite 2nd comp. on the 2nd grid
1009               ENDIF
1010               !
1011            ENDIF
1012            !
1013            IF( TRIM( sn_rcv_tau%clvor ) == 'eastward-northward' ) THEN   ! 2 components oriented along the local grid
1014               !                                                       ! (geographical to local grid -> rotate the components)
1015               CALL rot_rep( frcv(jpr_otx1)%z3(:,:,1), frcv(jpr_oty1)%z3(:,:,1), srcv(jpr_otx1)%clgrid, 'en->i', ztx )   
1016               IF( srcv(jpr_otx2)%laction ) THEN
1017                  CALL rot_rep( frcv(jpr_otx2)%z3(:,:,1), frcv(jpr_oty2)%z3(:,:,1), srcv(jpr_otx2)%clgrid, 'en->j', zty )   
1018               ELSE 
1019                  CALL rot_rep( frcv(jpr_otx1)%z3(:,:,1), frcv(jpr_oty1)%z3(:,:,1), srcv(jpr_otx1)%clgrid, 'en->j', zty ) 
1020               ENDIF
1021               frcv(jpr_otx1)%z3(:,:,1) = ztx(:,:)      ! overwrite 1st component on the 1st grid
1022               frcv(jpr_oty1)%z3(:,:,1) = zty(:,:)      ! overwrite 2nd component on the 2nd grid
1023            ENDIF
1024            !                             
1025            IF( srcv(jpr_otx1)%clgrid == 'T' ) THEN
1026               DO jj = 2, jpjm1                                          ! T ==> (U,V)
1027                  DO ji = fs_2, fs_jpim1   ! vector opt.
1028                     frcv(jpr_otx1)%z3(ji,jj,1) = 0.5 * ( frcv(jpr_otx1)%z3(ji+1,jj  ,1) + frcv(jpr_otx1)%z3(ji,jj,1) )
1029                     frcv(jpr_oty1)%z3(ji,jj,1) = 0.5 * ( frcv(jpr_oty1)%z3(ji  ,jj+1,1) + frcv(jpr_oty1)%z3(ji,jj,1) )
1030                  END DO
1031               END DO
1032               CALL lbc_lnk( frcv(jpr_otx1)%z3(:,:,1), 'U',  -1. )   ;   CALL lbc_lnk( frcv(jpr_oty1)%z3(:,:,1), 'V',  -1. )
1033            ENDIF
1034            llnewtx = .TRUE.
1035         ELSE
1036            llnewtx = .FALSE.
1037         ENDIF
1038         !                                                   ! ========================= !
1039      ELSE                                                   !   No dynamical coupling   !
1040         !                                                   ! ========================= !
1041         frcv(jpr_otx1)%z3(:,:,1) = 0.e0                               ! here simply set to zero
1042         frcv(jpr_oty1)%z3(:,:,1) = 0.e0                               ! an external read in a file can be added instead
1043         llnewtx = .TRUE.
1044         !
1045      ENDIF
1046      !                                                      ! ========================= !
1047      !                                                      !    wind stress module     !   (taum)
1048      !                                                      ! ========================= !
1049      !
1050      IF( .NOT. srcv(jpr_taum)%laction ) THEN                    ! compute wind stress module from its components if not received
1051         ! => need to be done only when otx1 was changed
1052         IF( llnewtx ) THEN
1053            DO jj = 2, jpjm1
1054               DO ji = fs_2, fs_jpim1   ! vect. opt.
1055                  zzx = frcv(jpr_otx1)%z3(ji-1,jj  ,1) + frcv(jpr_otx1)%z3(ji,jj,1)
1056                  zzy = frcv(jpr_oty1)%z3(ji  ,jj-1,1) + frcv(jpr_oty1)%z3(ji,jj,1)
1057                  frcv(jpr_taum)%z3(ji,jj,1) = 0.5 * SQRT( zzx * zzx + zzy * zzy )
1058               END DO
1059            END DO
1060            CALL lbc_lnk( frcv(jpr_taum)%z3(:,:,1), 'T', 1. )
1061            llnewtau = .TRUE.
1062         ELSE
1063            llnewtau = .FALSE.
1064         ENDIF
1065      ELSE
1066         llnewtau = nrcvinfo(jpr_taum) == OASIS_Rcv
1067         ! Stress module can be negative when received (interpolation problem)
1068         IF( llnewtau ) THEN
1069            frcv(jpr_taum)%z3(:,:,1) = MAX( 0._wp, frcv(jpr_taum)%z3(:,:,1) )
1070         ENDIF
1071      ENDIF
1072      !
1073      !                                                      ! ========================= !
1074      !                                                      !      10 m wind speed      !   (wndm)
1075      !                                                      ! ========================= !
1076      !
1077      IF( .NOT. srcv(jpr_w10m)%laction ) THEN                    ! compute wind spreed from wind stress module if not received 
1078         ! => need to be done only when taumod was changed
1079         IF( llnewtau ) THEN
1080            zcoef = 1. / ( zrhoa * zcdrag ) 
1081            DO jj = 1, jpj
1082               DO ji = 1, jpi 
1083                  frcv(jpr_w10m)%z3(ji,jj,1) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef )
1084               END DO
1085            END DO
1086         ENDIF
1087      ENDIF
1088
1089      ! u(v)tau and taum will be modified by ice model
1090      ! -> need to be reset before each call of the ice/fsbc     
1091      IF( MOD( kt-1, k_fsbc ) == 0 ) THEN
1092         !
1093         IF( ln_mixcpl ) THEN
1094            utau(:,:) = utau(:,:) * xcplmask(:,:,0) + frcv(jpr_otx1)%z3(:,:,1) * zmsk(:,:)
1095            vtau(:,:) = vtau(:,:) * xcplmask(:,:,0) + frcv(jpr_oty1)%z3(:,:,1) * zmsk(:,:)
1096            taum(:,:) = taum(:,:) * xcplmask(:,:,0) + frcv(jpr_taum)%z3(:,:,1) * zmsk(:,:)
1097            wndm(:,:) = wndm(:,:) * xcplmask(:,:,0) + frcv(jpr_w10m)%z3(:,:,1) * zmsk(:,:)
1098         ELSE
1099            utau(:,:) = frcv(jpr_otx1)%z3(:,:,1)
1100            vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1)
1101            taum(:,:) = frcv(jpr_taum)%z3(:,:,1)
1102            wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1)
1103         ENDIF
1104         CALL iom_put( "taum_oce", taum )   ! output wind stress module
1105         
1106      ENDIF
1107
1108#if defined key_cpl_carbon_cycle
1109      !                                                      ! ================== !
1110      !                                                      ! atmosph. CO2 (ppm) !
1111      !                                                      ! ================== !
1112      IF( srcv(jpr_co2)%laction )   atm_co2(:,:) = frcv(jpr_co2)%z3(:,:,1)
1113#endif
1114      !
1115      !                                                      ! ========================= !
1116      !                                                      ! Mean Sea Level Pressure   !   (taum)
1117      !                                                      ! ========================= !
1118      !
1119      IF( srcv(jpr_mslp)%laction ) THEN                    ! UKMO SHELF effect of atmospheric pressure on SSH
1120          IF( kt /= nit000 )   ssh_ibb(:,:) = ssh_ib(:,:)    !* Swap of ssh_ib fields
1121
1122          r1_grau = 1.e0 / (grav * rau0)               !* constant for optimization
1123          ssh_ib(:,:) = - ( frcv(jpr_mslp)%z3(:,:,1) - rpref ) * r1_grau    ! equivalent ssh (inverse barometer)
1124          apr   (:,:) =     frcv(jpr_mslp)%z3(:,:,1)                         !atmospheric pressure
1125   
1126          IF( kt == nit000 ) ssh_ibb(:,:) = ssh_ib(:,:)  ! correct this later (read from restart if possible)
1127      END IF 
1128      !
1129      IF( ln_sdw ) THEN  ! Stokes Drift correction activated
1130      !                                                      ! ========================= !
1131      !                                                      !       Stokes drift u      !
1132      !                                                      ! ========================= !
1133         IF( srcv(jpr_sdrftx)%laction ) zusd2dt(:,:) = frcv(jpr_sdrftx)%z3(:,:,1)
1134      !
1135      !                                                      ! ========================= !
1136      !                                                      !       Stokes drift v      !
1137      !                                                      ! ========================= !
1138         IF( srcv(jpr_sdrfty)%laction ) zvsd2dt(:,:) = frcv(jpr_sdrfty)%z3(:,:,1)
1139      !
1140      !                                                      ! ========================= !
1141      !                                                      !      Wave mean period     !
1142      !                                                      ! ========================= !
1143         IF( srcv(jpr_wper)%laction ) wmp(:,:) = frcv(jpr_wper)%z3(:,:,1)
1144      !
1145      !                                                      ! ========================= !
1146      !                                                      !  Significant wave height  !
1147      !                                                      ! ========================= !
1148         IF( srcv(jpr_hsig)%laction ) swh(:,:) = frcv(jpr_hsig)%z3(:,:,1)
1149      !
1150      !                                                      ! ========================= !
1151      !                                                      !    Vertical mixing Qiao   !
1152      !                                                      ! ========================= !
1153         IF( srcv(jpr_wnum)%laction .AND. ln_zdfqiao ) wnum(:,:) = frcv(jpr_wnum)%z3(:,:,1)
1154
1155         ! Calculate the 3D Stokes drift both in coupled and not fully uncoupled mode
1156         IF( srcv(jpr_sdrftx)%laction .OR. srcv(jpr_sdrfty)%laction .OR. srcv(jpr_wper)%laction &
1157                                                                    .OR. srcv(jpr_hsig)%laction ) THEN
1158            CALL sbc_stokes()
1159            IF( ln_zdfqiao .AND. .NOT. srcv(jpr_wnum)%laction ) CALL sbc_qiao()
1160         ENDIF
1161         IF( ln_zdfqiao .AND. srcv(jpr_wnum)%laction ) CALL sbc_qiao()
1162      ENDIF
1163      !                                                      ! ========================= !
1164      !                                                      ! Stress adsorbed by waves  !
1165      !                                                      ! ========================= !
1166      IF( srcv(jpr_wstrf)%laction .AND. ln_tauoc ) tauoc_wave(:,:) = frcv(jpr_wstrf)%z3(:,:,1)
1167
1168      !                                                      ! ========================= !
1169      !                                                      !   Wave drag coefficient   !
1170      !                                                      ! ========================= !
1171      IF( srcv(jpr_wdrag)%laction .AND. ln_cdgw ) cdn_wave(:,:) = frcv(jpr_wdrag)%z3(:,:,1)
1172
1173      !  Fields received by SAS when OASIS coupling
1174      !  (arrays no more filled at sbcssm stage)
1175      !                                                      ! ================== !
1176      !                                                      !        SSS         !
1177      !                                                      ! ================== !
1178      IF( srcv(jpr_soce)%laction ) THEN                      ! received by sas in case of opa <-> sas coupling
1179         sss_m(:,:) = frcv(jpr_soce)%z3(:,:,1)
1180         CALL iom_put( 'sss_m', sss_m )
1181      ENDIF
1182      !                                               
1183      !                                                      ! ================== !
1184      !                                                      !        SST         !
1185      !                                                      ! ================== !
1186      IF( srcv(jpr_toce)%laction ) THEN                      ! received by sas in case of opa <-> sas coupling
1187         sst_m(:,:) = frcv(jpr_toce)%z3(:,:,1)
1188         IF( srcv(jpr_soce)%laction .AND. ln_useCT ) THEN    ! make sure that sst_m is the potential temperature
1189            sst_m(:,:) = eos_pt_from_ct( sst_m(:,:), sss_m(:,:) )
1190         ENDIF
1191      ENDIF
1192      !                                                      ! ================== !
1193      !                                                      !        SSH         !
1194      !                                                      ! ================== !
1195      IF( srcv(jpr_ssh )%laction ) THEN                      ! received by sas in case of opa <-> sas coupling
1196         ssh_m(:,:) = frcv(jpr_ssh )%z3(:,:,1)
1197         CALL iom_put( 'ssh_m', ssh_m )
1198      ENDIF
1199      !                                                      ! ================== !
1200      !                                                      !  surface currents  !
1201      !                                                      ! ================== !
1202      IF( srcv(jpr_ocx1)%laction ) THEN                      ! received by sas in case of opa <-> sas coupling
1203         ssu_m(:,:) = frcv(jpr_ocx1)%z3(:,:,1)
1204         ub (:,:,1) = ssu_m(:,:)                             ! will be used in sbcice_lim in the call of lim_sbc_tau
1205         CALL iom_put( 'ssu_m', ssu_m )
1206      ENDIF
1207      IF( srcv(jpr_ocy1)%laction ) THEN
1208         ssv_m(:,:) = frcv(jpr_ocy1)%z3(:,:,1)
1209         vb (:,:,1) = ssv_m(:,:)                             ! will be used in sbcice_lim in the call of lim_sbc_tau
1210         CALL iom_put( 'ssv_m', ssv_m )
1211      ENDIF
1212      !                                                      ! ======================== !
1213      !                                                      !  first T level thickness !
1214      !                                                      ! ======================== !
1215      IF( srcv(jpr_e3t1st )%laction ) THEN                   ! received by sas in case of opa <-> sas coupling
1216         e3t_m(:,:) = frcv(jpr_e3t1st )%z3(:,:,1)
1217         CALL iom_put( 'e3t_m', e3t_m(:,:) )
1218      ENDIF
1219      !                                                      ! ================================ !
1220      !                                                      !  fraction of solar net radiation !
1221      !                                                      ! ================================ !
1222      IF( srcv(jpr_fraqsr)%laction ) THEN                    ! received by sas in case of opa <-> sas coupling
1223         frq_m(:,:) = frcv(jpr_fraqsr)%z3(:,:,1)
1224         CALL iom_put( 'frq_m', frq_m )
1225      ENDIF
1226     
1227      !                                                      ! ========================= !
1228      IF( k_ice <= 1 .AND. MOD( kt-1, k_fsbc ) == 0 ) THEN   !  heat & freshwater fluxes ! (Ocean only case)
1229         !                                                   ! ========================= !
1230         !
1231         !                                                       ! total freshwater fluxes over the ocean (emp)
1232         IF( srcv(jpr_oemp)%laction .OR. srcv(jpr_rain)%laction ) THEN
1233            SELECT CASE( TRIM( sn_rcv_emp%cldes ) )                                    ! evaporation - precipitation
1234            CASE( 'conservative' )
1235               zemp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) )
1236            CASE( 'oce only', 'oce and ice' )
1237               zemp(:,:) = frcv(jpr_oemp)%z3(:,:,1)
1238            CASE default
1239               CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' )
1240            END SELECT
1241         ELSE
1242            zemp(:,:) = 0._wp
1243         ENDIF
1244         !
1245         !                                                        ! runoffs and calving (added in emp)
1246         IF( srcv(jpr_rnf)%laction )     rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1)
1247         IF( srcv(jpr_cal)%laction )     zemp(:,:) = zemp(:,:) - frcv(jpr_cal)%z3(:,:,1)
1248         
1249         IF( ln_mixcpl ) THEN   ;   emp(:,:) = emp(:,:) * xcplmask(:,:,0) + zemp(:,:) * zmsk(:,:)
1250         ELSE                   ;   emp(:,:) =                              zemp(:,:)
1251         ENDIF
1252         !
1253         !                                                       ! non solar heat flux over the ocean (qns)
1254         IF(      srcv(jpr_qnsoce)%laction ) THEN   ;   zqns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1)
1255         ELSE IF( srcv(jpr_qnsmix)%laction ) THEN   ;   zqns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1)
1256         ELSE                                       ;   zqns(:,:) = 0._wp
1257         END IF
1258         ! update qns over the free ocean with:
1259         IF( nn_components /= jp_iam_opa ) THEN
1260            zqns(:,:) =  zqns(:,:) - zemp(:,:) * sst_m(:,:) * rcp         ! remove heat content due to mass flux (assumed to be at SST)
1261            IF( srcv(jpr_snow  )%laction ) THEN
1262               zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus    ! energy for melting solid precipitation over the free ocean
1263            ENDIF
1264         ENDIF
1265         IF( ln_mixcpl ) THEN   ;   qns(:,:) = qns(:,:) * xcplmask(:,:,0) + zqns(:,:) * zmsk(:,:)
1266         ELSE                   ;   qns(:,:) =                              zqns(:,:)
1267         ENDIF
1268
1269         !                                                       ! solar flux over the ocean          (qsr)
1270         IF     ( srcv(jpr_qsroce)%laction ) THEN   ;   zqsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1)
1271         ELSE IF( srcv(jpr_qsrmix)%laction ) then   ;   zqsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1)
1272         ELSE                                       ;   zqsr(:,:) = 0._wp
1273         ENDIF
1274         IF( ln_dm2dc .AND. ln_cpl )   zqsr(:,:) = sbc_dcy( zqsr )   ! modify qsr to include the diurnal cycle
1275         IF( ln_mixcpl ) THEN   ;   qsr(:,:) = qsr(:,:) * xcplmask(:,:,0) + zqsr(:,:) * zmsk(:,:)
1276         ELSE                   ;   qsr(:,:) =                              zqsr(:,:)
1277         ENDIF
1278         !
1279         ! salt flux over the ocean (received by opa in case of opa <-> sas coupling)
1280         IF( srcv(jpr_sflx )%laction )   sfx(:,:) = frcv(jpr_sflx  )%z3(:,:,1)
1281         ! Ice cover  (received by opa in case of opa <-> sas coupling)
1282         IF( srcv(jpr_fice )%laction )   fr_i(:,:) = frcv(jpr_fice )%z3(:,:,1)
1283         !
1284
1285      ENDIF
1286      !
1287      CALL wrk_dealloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr )
1288      !
1289      IF( nn_timing == 1 )  CALL timing_stop('sbc_cpl_rcv')
1290      !
1291   END SUBROUTINE sbc_cpl_rcv
1292   
1293
1294   SUBROUTINE sbc_cpl_ice_tau( p_taui, p_tauj )     
1295      !!----------------------------------------------------------------------
1296      !!             ***  ROUTINE sbc_cpl_ice_tau  ***
1297      !!
1298      !! ** Purpose :   provide the stress over sea-ice in coupled mode
1299      !!
1300      !! ** Method  :   transform the received stress from the atmosphere into
1301      !!             an atmosphere-ice stress in the (i,j) ocean referencial
1302      !!             and at the velocity point of the sea-ice model (cp_ice_msh):
1303      !!                'C'-grid : i- (j-) components given at U- (V-) point
1304      !!                'I'-grid : B-grid lower-left corner: both components given at I-point
1305      !!
1306      !!                The received stress are :
1307      !!                 - defined by 3 components (if cartesian coordinate)
1308      !!                        or by 2 components (if spherical)
1309      !!                 - oriented along geographical   coordinate (if eastward-northward)
1310      !!                        or  along the local grid coordinate (if local grid)
1311      !!                 - given at U- and V-point, resp.   if received on 2 grids
1312      !!                        or at a same point (T or I) if received on 1 grid
1313      !!                Therefore and if necessary, they are successively
1314      !!             processed in order to obtain them
1315      !!                 first  as  2 components on the sphere
1316      !!                 second as  2 components oriented along the local grid
1317      !!                 third  as  2 components on the cp_ice_msh point
1318      !!
1319      !!                Except in 'oce and ice' case, only one vector stress field
1320      !!             is received. It has already been processed in sbc_cpl_rcv
1321      !!             so that it is now defined as (i,j) components given at U-
1322      !!             and V-points, respectively. Therefore, only the third
1323      !!             transformation is done and only if the ice-grid is a 'I'-grid.
1324      !!
1325      !! ** Action  :   return ptau_i, ptau_j, the stress over the ice at cp_ice_msh point
1326      !!----------------------------------------------------------------------
1327      REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_taui   ! i- & j-components of atmos-ice stress [N/m2]
1328      REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_tauj   ! at I-point (B-grid) or U & V-point (C-grid)
1329      !!
1330      INTEGER ::   ji, jj                          ! dummy loop indices
1331      INTEGER ::   itx                             ! index of taux over ice
1332      REAL(wp), POINTER, DIMENSION(:,:) ::   ztx, zty 
1333      !!----------------------------------------------------------------------
1334      !
1335      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_ice_tau')
1336      !
1337      CALL wrk_alloc( jpi,jpj, ztx, zty )
1338
1339      IF( srcv(jpr_itx1)%laction ) THEN   ;   itx =  jpr_itx1   
1340      ELSE                                ;   itx =  jpr_otx1
1341      ENDIF
1342
1343      ! do something only if we just received the stress from atmosphere
1344      IF(  nrcvinfo(itx) == OASIS_Rcv ) THEN
1345
1346         !                                                      ! ======================= !
1347         IF( srcv(jpr_itx1)%laction ) THEN                      !   ice stress received   !
1348            !                                                   ! ======================= !
1349           
1350            IF( TRIM( sn_rcv_tau%clvref ) == 'cartesian' ) THEN            ! 2 components on the sphere
1351               !                                                       ! (cartesian to spherical -> 3 to 2 components)
1352               CALL geo2oce(  frcv(jpr_itx1)%z3(:,:,1), frcv(jpr_ity1)%z3(:,:,1), frcv(jpr_itz1)%z3(:,:,1),   &
1353                  &          srcv(jpr_itx1)%clgrid, ztx, zty )
1354               frcv(jpr_itx1)%z3(:,:,1) = ztx(:,:)   ! overwrite 1st comp. on the 1st grid
1355               frcv(jpr_ity1)%z3(:,:,1) = zty(:,:)   ! overwrite 2nd comp. on the 1st grid
1356               !
1357               IF( srcv(jpr_itx2)%laction ) THEN
1358                  CALL geo2oce( frcv(jpr_itx2)%z3(:,:,1), frcv(jpr_ity2)%z3(:,:,1), frcv(jpr_itz2)%z3(:,:,1),   &
1359                     &          srcv(jpr_itx2)%clgrid, ztx, zty )
1360                  frcv(jpr_itx2)%z3(:,:,1) = ztx(:,:)   ! overwrite 1st comp. on the 2nd grid
1361                  frcv(jpr_ity2)%z3(:,:,1) = zty(:,:)   ! overwrite 2nd comp. on the 2nd grid
1362               ENDIF
1363               !
1364            ENDIF
1365            !
1366            IF( TRIM( sn_rcv_tau%clvor ) == 'eastward-northward' ) THEN   ! 2 components oriented along the local grid
1367               !                                                       ! (geographical to local grid -> rotate the components)
1368               CALL rot_rep( frcv(jpr_itx1)%z3(:,:,1), frcv(jpr_ity1)%z3(:,:,1), srcv(jpr_itx1)%clgrid, 'en->i', ztx )   
1369               IF( srcv(jpr_itx2)%laction ) THEN
1370                  CALL rot_rep( frcv(jpr_itx2)%z3(:,:,1), frcv(jpr_ity2)%z3(:,:,1), srcv(jpr_itx2)%clgrid, 'en->j', zty )   
1371               ELSE
1372                  CALL rot_rep( frcv(jpr_itx1)%z3(:,:,1), frcv(jpr_ity1)%z3(:,:,1), srcv(jpr_itx1)%clgrid, 'en->j', zty ) 
1373               ENDIF
1374               frcv(jpr_itx1)%z3(:,:,1) = ztx(:,:)      ! overwrite 1st component on the 1st grid
1375               frcv(jpr_ity1)%z3(:,:,1) = zty(:,:)      ! overwrite 2nd component on the 1st grid
1376            ENDIF
1377            !                                                   ! ======================= !
1378         ELSE                                                   !     use ocean stress    !
1379            !                                                   ! ======================= !
1380            frcv(jpr_itx1)%z3(:,:,1) = frcv(jpr_otx1)%z3(:,:,1)
1381            frcv(jpr_ity1)%z3(:,:,1) = frcv(jpr_oty1)%z3(:,:,1)
1382            !
1383         ENDIF
1384         !                                                      ! ======================= !
1385         !                                                      !     put on ice grid     !
1386         !                                                      ! ======================= !
1387         !   
1388         !                                                  j+1   j     -----V---F
1389         ! ice stress on ice velocity point (cp_ice_msh)                 !       |
1390         ! (C-grid ==>(U,V) or B-grid ==> I or F)                 j      |   T   U
1391         !                                                               |       |
1392         !                                                   j    j-1   -I-------|
1393         !                                               (for I)         |       |
1394         !                                                              i-1  i   i
1395         !                                                               i      i+1 (for I)
1396         SELECT CASE ( cp_ice_msh )
1397            !
1398         CASE( 'I' )                                         ! B-grid ==> I
1399            SELECT CASE ( srcv(jpr_itx1)%clgrid )
1400            CASE( 'U' )
1401               DO jj = 2, jpjm1                                   ! (U,V) ==> I
1402                  DO ji = 2, jpim1   ! NO vector opt.
1403                     p_taui(ji,jj) = 0.5 * ( frcv(jpr_itx1)%z3(ji-1,jj  ,1) + frcv(jpr_itx1)%z3(ji-1,jj-1,1) )
1404                     p_tauj(ji,jj) = 0.5 * ( frcv(jpr_ity1)%z3(ji  ,jj-1,1) + frcv(jpr_ity1)%z3(ji-1,jj-1,1) )
1405                  END DO
1406               END DO
1407            CASE( 'F' )
1408               DO jj = 2, jpjm1                                   ! F ==> I
1409                  DO ji = 2, jpim1   ! NO vector opt.
1410                     p_taui(ji,jj) = frcv(jpr_itx1)%z3(ji-1,jj-1,1)
1411                     p_tauj(ji,jj) = frcv(jpr_ity1)%z3(ji-1,jj-1,1)
1412                  END DO
1413               END DO
1414            CASE( 'T' )
1415               DO jj = 2, jpjm1                                   ! T ==> I
1416                  DO ji = 2, jpim1   ! NO vector opt.
1417                     p_taui(ji,jj) = 0.25 * ( frcv(jpr_itx1)%z3(ji,jj  ,1) + frcv(jpr_itx1)%z3(ji-1,jj  ,1)   &
1418                        &                   + frcv(jpr_itx1)%z3(ji,jj-1,1) + frcv(jpr_itx1)%z3(ji-1,jj-1,1) ) 
1419                     p_tauj(ji,jj) = 0.25 * ( frcv(jpr_ity1)%z3(ji,jj  ,1) + frcv(jpr_ity1)%z3(ji-1,jj  ,1)   &
1420                        &                   + frcv(jpr_oty1)%z3(ji,jj-1,1) + frcv(jpr_ity1)%z3(ji-1,jj-1,1) )
1421                  END DO
1422               END DO
1423            CASE( 'I' )
1424               p_taui(:,:) = frcv(jpr_itx1)%z3(:,:,1)                   ! I ==> I
1425               p_tauj(:,:) = frcv(jpr_ity1)%z3(:,:,1)
1426            END SELECT
1427            IF( srcv(jpr_itx1)%clgrid /= 'I' ) THEN
1428               CALL lbc_lnk( p_taui, 'I',  -1. )   ;   CALL lbc_lnk( p_tauj, 'I',  -1. )
1429            ENDIF
1430            !
1431         CASE( 'F' )                                         ! B-grid ==> F
1432            SELECT CASE ( srcv(jpr_itx1)%clgrid )
1433            CASE( 'U' )
1434               DO jj = 2, jpjm1                                   ! (U,V) ==> F
1435                  DO ji = fs_2, fs_jpim1   ! vector opt.
1436                     p_taui(ji,jj) = 0.5 * ( frcv(jpr_itx1)%z3(ji,jj,1) + frcv(jpr_itx1)%z3(ji  ,jj+1,1) )
1437                     p_tauj(ji,jj) = 0.5 * ( frcv(jpr_ity1)%z3(ji,jj,1) + frcv(jpr_ity1)%z3(ji+1,jj  ,1) )
1438                  END DO
1439               END DO
1440            CASE( 'I' )
1441               DO jj = 2, jpjm1                                   ! I ==> F
1442                  DO ji = 2, jpim1   ! NO vector opt.
1443                     p_taui(ji,jj) = frcv(jpr_itx1)%z3(ji+1,jj+1,1)
1444                     p_tauj(ji,jj) = frcv(jpr_ity1)%z3(ji+1,jj+1,1)
1445                  END DO
1446               END DO
1447            CASE( 'T' )
1448               DO jj = 2, jpjm1                                   ! T ==> F
1449                  DO ji = 2, jpim1   ! NO vector opt.
1450                     p_taui(ji,jj) = 0.25 * ( frcv(jpr_itx1)%z3(ji,jj  ,1) + frcv(jpr_itx1)%z3(ji+1,jj  ,1)   &
1451                        &                   + frcv(jpr_itx1)%z3(ji,jj+1,1) + frcv(jpr_itx1)%z3(ji+1,jj+1,1) ) 
1452                     p_tauj(ji,jj) = 0.25 * ( frcv(jpr_ity1)%z3(ji,jj  ,1) + frcv(jpr_ity1)%z3(ji+1,jj  ,1)   &
1453                        &                   + frcv(jpr_ity1)%z3(ji,jj+1,1) + frcv(jpr_ity1)%z3(ji+1,jj+1,1) )
1454                  END DO
1455               END DO
1456            CASE( 'F' )
1457               p_taui(:,:) = frcv(jpr_itx1)%z3(:,:,1)                   ! F ==> F
1458               p_tauj(:,:) = frcv(jpr_ity1)%z3(:,:,1)
1459            END SELECT
1460            IF( srcv(jpr_itx1)%clgrid /= 'F' ) THEN
1461               CALL lbc_lnk( p_taui, 'F',  -1. )   ;   CALL lbc_lnk( p_tauj, 'F',  -1. )
1462            ENDIF
1463            !
1464         CASE( 'C' )                                         ! C-grid ==> U,V
1465            SELECT CASE ( srcv(jpr_itx1)%clgrid )
1466            CASE( 'U' )
1467               p_taui(:,:) = frcv(jpr_itx1)%z3(:,:,1)                   ! (U,V) ==> (U,V)
1468               p_tauj(:,:) = frcv(jpr_ity1)%z3(:,:,1)
1469            CASE( 'F' )
1470               DO jj = 2, jpjm1                                   ! F ==> (U,V)
1471                  DO ji = fs_2, fs_jpim1   ! vector opt.
1472                     p_taui(ji,jj) = 0.5 * ( frcv(jpr_itx1)%z3(ji,jj,1) + frcv(jpr_itx1)%z3(ji  ,jj-1,1) )
1473                     p_tauj(ji,jj) = 0.5 * ( frcv(jpr_ity1)%z3(jj,jj,1) + frcv(jpr_ity1)%z3(ji-1,jj  ,1) )
1474                  END DO
1475               END DO
1476            CASE( 'T' )
1477               DO jj = 2, jpjm1                                   ! T ==> (U,V)
1478                  DO ji = fs_2, fs_jpim1   ! vector opt.
1479                     p_taui(ji,jj) = 0.5 * ( frcv(jpr_itx1)%z3(ji+1,jj  ,1) + frcv(jpr_itx1)%z3(ji,jj,1) )
1480                     p_tauj(ji,jj) = 0.5 * ( frcv(jpr_ity1)%z3(ji  ,jj+1,1) + frcv(jpr_ity1)%z3(ji,jj,1) )
1481                  END DO
1482               END DO
1483            CASE( 'I' )
1484               DO jj = 2, jpjm1                                   ! I ==> (U,V)
1485                  DO ji = 2, jpim1   ! NO vector opt.
1486                     p_taui(ji,jj) = 0.5 * ( frcv(jpr_itx1)%z3(ji+1,jj+1,1) + frcv(jpr_itx1)%z3(ji+1,jj  ,1) )
1487                     p_tauj(ji,jj) = 0.5 * ( frcv(jpr_ity1)%z3(ji+1,jj+1,1) + frcv(jpr_ity1)%z3(ji  ,jj+1,1) )
1488                  END DO
1489               END DO
1490            END SELECT
1491            IF( srcv(jpr_itx1)%clgrid /= 'U' ) THEN
1492               CALL lbc_lnk( p_taui, 'U',  -1. )   ;   CALL lbc_lnk( p_tauj, 'V',  -1. )
1493            ENDIF
1494         END SELECT
1495
1496      ENDIF
1497      !   
1498      CALL wrk_dealloc( jpi,jpj, ztx, zty )
1499      !
1500      IF( nn_timing == 1 )  CALL timing_stop('sbc_cpl_ice_tau')
1501      !
1502   END SUBROUTINE sbc_cpl_ice_tau
1503   
1504
1505   SUBROUTINE sbc_cpl_ice_flx( p_frld, palbi, psst, pist )
1506      !!----------------------------------------------------------------------
1507      !!             ***  ROUTINE sbc_cpl_ice_flx  ***
1508      !!
1509      !! ** Purpose :   provide the heat and freshwater fluxes of the
1510      !!              ocean-ice system.
1511      !!
1512      !! ** Method  :   transform the fields received from the atmosphere into
1513      !!             surface heat and fresh water boundary condition for the
1514      !!             ice-ocean system. The following fields are provided:
1515      !!              * total non solar, solar and freshwater fluxes (qns_tot,
1516      !!             qsr_tot and emp_tot) (total means weighted ice-ocean flux)
1517      !!             NB: emp_tot include runoffs and calving.
1518      !!              * fluxes over ice (qns_ice, qsr_ice, emp_ice) where
1519      !!             emp_ice = sublimation - solid precipitation as liquid
1520      !!             precipitation are re-routed directly to the ocean and
1521      !!             runoffs and calving directly enter the ocean.
1522      !!              * solid precipitation (sprecip), used to add to qns_tot
1523      !!             the heat lost associated to melting solid precipitation
1524      !!             over the ocean fraction.
1525      !!       ===>> CAUTION here this changes the net heat flux received from
1526      !!             the atmosphere
1527      !!
1528      !!                  - the fluxes have been separated from the stress as
1529      !!                 (a) they are updated at each ice time step compare to
1530      !!                 an update at each coupled time step for the stress, and
1531      !!                 (b) the conservative computation of the fluxes over the
1532      !!                 sea-ice area requires the knowledge of the ice fraction
1533      !!                 after the ice advection and before the ice thermodynamics,
1534      !!                 so that the stress is updated before the ice dynamics
1535      !!                 while the fluxes are updated after it.
1536      !!
1537      !! ** Action  :   update at each nf_ice time step:
1538      !!                   qns_tot, qsr_tot  non-solar and solar total heat fluxes
1539      !!                   qns_ice, qsr_ice  non-solar and solar heat fluxes over the ice
1540      !!                   emp_tot            total evaporation - precipitation(liquid and solid) (-runoff)(-calving)
1541      !!                   emp_ice            ice sublimation - solid precipitation over the ice
1542      !!                   dqns_ice           d(non-solar heat flux)/d(Temperature) over the ice
1543      !!                   sprecip             solid precipitation over the ocean 
1544      !!----------------------------------------------------------------------
1545      REAL(wp), INTENT(in   ), DIMENSION(:,:)   ::   p_frld     ! lead fraction                [0 to 1]
1546      ! optional arguments, used only in 'mixed oce-ice' case
1547      REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   palbi      ! all skies ice albedo
1548      REAL(wp), INTENT(in   ), DIMENSION(:,:  ), OPTIONAL ::   psst       ! sea surface temperature     [Celsius]
1549      REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   pist       ! ice surface temperature     [Kelvin]
1550      !
1551      INTEGER ::   jl         ! dummy loop index
1552      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zcptn, ztmp, zicefr, zmsk
1553      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot
1554      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zqns_ice, zqsr_ice, zdqns_ice
1555      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zevap, zsnw, zqns_oce, zqsr_oce, zqprec_ice, zqemp_oce ! for LIM3
1556      !!----------------------------------------------------------------------
1557      !
1558      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_ice_flx')
1559      !
1560      CALL wrk_alloc( jpi,jpj,     zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot )
1561      CALL wrk_alloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice )
1562
1563      IF( ln_mixcpl )   zmsk(:,:) = 1. - xcplmask(:,:,0)
1564      zicefr(:,:) = 1.- p_frld(:,:)
1565      zcptn(:,:) = rcp * sst_m(:,:)
1566      !
1567      !                                                      ! ========================= !
1568      !                                                      !    freshwater budget      !   (emp)
1569      !                                                      ! ========================= !
1570      !
1571      !                                                           ! total Precipitation - total Evaporation (emp_tot)
1572      !                                                           ! solid precipitation - sublimation       (emp_ice)
1573      !                                                           ! solid Precipitation                     (sprecip)
1574      !                                                           ! liquid + solid Precipitation            (tprecip)
1575      SELECT CASE( TRIM( sn_rcv_emp%cldes ) )
1576      CASE( 'conservative'  )   ! received fields: jpr_rain, jpr_snow, jpr_ievp, jpr_tevp
1577         zsprecip(:,:) = frcv(jpr_snow)%z3(:,:,1)                  ! May need to ensure positive here
1578         ztprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + zsprecip(:,:)  ! May need to ensure positive here
1579         zemp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ztprecip(:,:)
1580         zemp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1)
1581            CALL iom_put( 'rain'         , frcv(jpr_rain)%z3(:,:,1)              )   ! liquid precipitation
1582         IF( iom_use('hflx_rain_cea') )   &
1583            CALL iom_put( 'hflx_rain_cea', frcv(jpr_rain)%z3(:,:,1) * zcptn(:,:) )   ! heat flux from liq. precip.
1584         IF( iom_use('evap_ao_cea') .OR. iom_use('hflx_evap_cea') )   &
1585            ztmp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:)
1586         IF( iom_use('evap_ao_cea'  ) )   &
1587            CALL iom_put( 'evap_ao_cea'  , ztmp                   )   ! ice-free oce evap (cell average)
1588         IF( iom_use('hflx_evap_cea') )   &
1589            CALL iom_put( 'hflx_evap_cea', ztmp(:,:) * zcptn(:,:) )   ! heat flux from from evap (cell average)
1590      CASE( 'oce and ice'   )   ! received fields: jpr_sbpr, jpr_semp, jpr_oemp, jpr_ievp
1591         zemp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1)
1592         zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1)
1593         zsprecip(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_semp)%z3(:,:,1)
1594         ztprecip(:,:) = frcv(jpr_semp)%z3(:,:,1) - frcv(jpr_sbpr)%z3(:,:,1) + zsprecip(:,:)
1595      END SELECT
1596
1597      IF( iom_use('subl_ai_cea') )   &
1598         CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) )   ! Sublimation over sea-ice         (cell average)
1599      !   
1600      !                                                           ! runoffs and calving (put in emp_tot)
1601      IF( srcv(jpr_rnf)%laction )   rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1)
1602      IF( srcv(jpr_cal)%laction ) THEN
1603         zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1)
1604         CALL iom_put( 'calving_cea', frcv(jpr_cal)%z3(:,:,1) )
1605      ENDIF
1606
1607      IF( ln_mixcpl ) THEN
1608         emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:)
1609         emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:)
1610         sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:)
1611         tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:)
1612      ELSE
1613         emp_tot(:,:) =                                  zemp_tot(:,:)
1614         emp_ice(:,:) =                                  zemp_ice(:,:)
1615         sprecip(:,:) =                                  zsprecip(:,:)
1616         tprecip(:,:) =                                  ztprecip(:,:)
1617      ENDIF
1618
1619         CALL iom_put( 'snowpre'    , sprecip                                )   ! Snow
1620      IF( iom_use('snow_ao_cea') )   &
1621         CALL iom_put( 'snow_ao_cea', sprecip(:,:) * p_frld(:,:)             )   ! Snow        over ice-free ocean  (cell average)
1622      IF( iom_use('snow_ai_cea') )   &
1623         CALL iom_put( 'snow_ai_cea', sprecip(:,:) * zicefr(:,:)             )   ! Snow        over sea-ice         (cell average)
1624
1625      !                                                      ! ========================= !
1626      SELECT CASE( TRIM( sn_rcv_qns%cldes ) )                !   non solar heat fluxes   !   (qns)
1627      !                                                      ! ========================= !
1628      CASE( 'oce only' )                                     ! the required field is directly provided
1629         zqns_tot(:,:  ) = frcv(jpr_qnsoce)%z3(:,:,1)
1630      CASE( 'conservative' )                                      ! the required fields are directly provided
1631         zqns_tot(:,:  ) = frcv(jpr_qnsmix)%z3(:,:,1)
1632         IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN
1633            zqns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl)
1634         ELSE
1635            ! Set all category values equal for the moment
1636            DO jl=1,jpl
1637               zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)
1638            ENDDO
1639         ENDIF
1640      CASE( 'oce and ice' )       ! the total flux is computed from ocean and ice fluxes
1641         zqns_tot(:,:  ) =  p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1)
1642         IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN
1643            DO jl=1,jpl
1644               zqns_tot(:,:   ) = zqns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl)   
1645               zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl)
1646            ENDDO
1647         ELSE
1648            qns_tot(:,:   ) = qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1)
1649            DO jl=1,jpl
1650               zqns_tot(:,:   ) = zqns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1)
1651               zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)
1652            ENDDO
1653         ENDIF
1654      CASE( 'mixed oce-ice' )     ! the ice flux is cumputed from the total flux, the SST and ice informations
1655! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED **
1656         zqns_tot(:,:  ) = frcv(jpr_qnsmix)%z3(:,:,1)
1657         zqns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1)    &
1658            &            + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,:  ) ) * p_frld(:,:)   &
1659            &                                                   +          pist(:,:,1)   * zicefr(:,:) ) )
1660      END SELECT
1661!!gm
1662!!    currently it is taken into account in leads budget but not in the zqns_tot, and thus not in
1663!!    the flux that enter the ocean....
1664!!    moreover 1 - it is not diagnose anywhere....
1665!!             2 - it is unclear for me whether this heat lost is taken into account in the atmosphere or not...
1666!!
1667!! similar job should be done for snow and precipitation temperature
1668      !                                     
1669      IF( srcv(jpr_cal)%laction ) THEN                            ! Iceberg melting
1670         ztmp(:,:) = frcv(jpr_cal)%z3(:,:,1) * lfus               ! add the latent heat of iceberg melting
1671         zqns_tot(:,:) = zqns_tot(:,:) - ztmp(:,:)
1672         IF( iom_use('hflx_cal_cea') )   &
1673            CALL iom_put( 'hflx_cal_cea', ztmp + frcv(jpr_cal)%z3(:,:,1) * zcptn(:,:) )   ! heat flux from calving
1674      ENDIF
1675
1676      ztmp(:,:) = p_frld(:,:) * zsprecip(:,:) * lfus
1677      IF( iom_use('hflx_snow_cea') )    CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) )   ! heat flux from snow (cell average)
1678
1679#if defined key_lim3
1680      CALL wrk_alloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 
1681
1682      ! --- evaporation --- !
1683      ! clem: evap_ice is set to 0 for LIM3 since we still do not know what to do with sublimation
1684      ! the problem is: the atm. imposes both mass evaporation and heat removed from the snow/ice
1685      !                 but it is incoherent WITH the ice model 
1686      DO jl=1,jpl
1687         evap_ice(:,:,jl) = 0._wp  ! should be: frcv(jpr_ievp)%z3(:,:,1)
1688      ENDDO
1689      zevap(:,:) = zemp_tot(:,:) + ztprecip(:,:) ! evaporation over ocean
1690
1691      ! --- evaporation minus precipitation --- !
1692      emp_oce(:,:) = emp_tot(:,:) - emp_ice(:,:)
1693
1694      ! --- non solar flux over ocean --- !
1695      !         note: p_frld cannot be = 0 since we limit the ice concentration to amax
1696      zqns_oce = 0._wp
1697      WHERE( p_frld /= 0._wp )  zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / p_frld(:,:)
1698
1699      ! --- heat flux associated with emp --- !
1700      zsnw(:,:) = 0._wp
1701      CALL lim_thd_snwblow( p_frld, zsnw )  ! snow distribution over ice after wind blowing
1702      zqemp_oce(:,:) = -      zevap(:,:)                   * p_frld(:,:)      *   zcptn(:,:)   &      ! evap
1703         &             + ( ztprecip(:,:) - zsprecip(:,:) )                    *   zcptn(:,:)   &      ! liquid precip
1704         &             +   zsprecip(:,:)                   * ( 1._wp - zsnw ) * ( zcptn(:,:) - lfus ) ! solid precip over ocean
1705      qemp_ice(:,:)  = -   frcv(jpr_ievp)%z3(:,:,1)        * zicefr(:,:)      *   zcptn(:,:)   &      ! ice evap
1706         &             +   zsprecip(:,:)                   * zsnw             * ( zcptn(:,:) - lfus ) ! solid precip over ice
1707
1708      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
1709      zqprec_ice(:,:) = rhosn * ( zcptn(:,:) - lfus )
1710
1711      ! --- total non solar flux --- !
1712      zqns_tot(:,:) = zqns_tot(:,:) + qemp_ice(:,:) + zqemp_oce(:,:)
1713
1714      ! --- in case both coupled/forced are active, we must mix values --- !
1715      IF( ln_mixcpl ) THEN
1716         qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:)
1717         qns_oce(:,:) = qns_oce(:,:) * xcplmask(:,:,0) + zqns_oce(:,:)* zmsk(:,:)
1718         DO jl=1,jpl
1719            qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) +  zqns_ice(:,:,jl)* zmsk(:,:)
1720         ENDDO
1721         qprec_ice(:,:) = qprec_ice(:,:) * xcplmask(:,:,0) + zqprec_ice(:,:)* zmsk(:,:)
1722         qemp_oce (:,:) =  qemp_oce(:,:) * xcplmask(:,:,0) +  zqemp_oce(:,:)* zmsk(:,:)
1723!!clem         evap_ice(:,:) = evap_ice(:,:) * xcplmask(:,:,0)
1724      ELSE
1725         qns_tot  (:,:  ) = zqns_tot  (:,:  )
1726         qns_oce  (:,:  ) = zqns_oce  (:,:  )
1727         qns_ice  (:,:,:) = zqns_ice  (:,:,:)
1728         qprec_ice(:,:)   = zqprec_ice(:,:)
1729         qemp_oce (:,:)   = zqemp_oce (:,:)
1730      ENDIF
1731
1732      CALL wrk_dealloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 
1733#else
1734
1735      ! clem: this formulation is certainly wrong... but better than it was...
1736      zqns_tot(:,:) = zqns_tot(:,:)                       &            ! zqns_tot update over free ocean with:
1737         &          - ztmp(:,:)                           &            ! remove the latent heat flux of solid precip. melting
1738         &          - (  zemp_tot(:,:)                    &            ! remove the heat content of mass flux (assumed to be at SST)
1739         &             - zemp_ice(:,:) * zicefr(:,:)  ) * zcptn(:,:) 
1740
1741     IF( ln_mixcpl ) THEN
1742         qns_tot(:,:) = qns(:,:) * p_frld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk
1743         qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) +  zqns_tot(:,:)* zmsk(:,:)
1744         DO jl=1,jpl
1745            qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) +  zqns_ice(:,:,jl)* zmsk(:,:)
1746         ENDDO
1747      ELSE
1748         qns_tot(:,:  ) = zqns_tot(:,:  )
1749         qns_ice(:,:,:) = zqns_ice(:,:,:)
1750      ENDIF
1751
1752#endif
1753
1754      !                                                      ! ========================= !
1755      SELECT CASE( TRIM( sn_rcv_qsr%cldes ) )                !      solar heat fluxes    !   (qsr)
1756      !                                                      ! ========================= !
1757      CASE( 'oce only' )
1758         zqsr_tot(:,:  ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) )
1759      CASE( 'conservative' )
1760         zqsr_tot(:,:  ) = frcv(jpr_qsrmix)%z3(:,:,1)
1761         IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN
1762            zqsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl)
1763         ELSE
1764            ! Set all category values equal for the moment
1765            DO jl=1,jpl
1766               zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)
1767            ENDDO
1768         ENDIF
1769         zqsr_tot(:,:  ) = frcv(jpr_qsrmix)%z3(:,:,1)
1770         zqsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1)
1771      CASE( 'oce and ice' )
1772         zqsr_tot(:,:  ) =  p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1)
1773         IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN
1774            DO jl=1,jpl
1775               zqsr_tot(:,:   ) = zqsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl)   
1776               zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl)
1777            ENDDO
1778         ELSE
1779            qsr_tot(:,:   ) = qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1)
1780            DO jl=1,jpl
1781               zqsr_tot(:,:   ) = zqsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1)
1782               zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)
1783            ENDDO
1784         ENDIF
1785      CASE( 'mixed oce-ice' )
1786         zqsr_tot(:,:  ) = frcv(jpr_qsrmix)%z3(:,:,1)
1787! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED **
1788!       Create solar heat flux over ice using incoming solar heat flux and albedos
1789!       ( see OASIS3 user guide, 5th edition, p39 )
1790         zqsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) )   &
1791            &            / (  1.- ( albedo_oce_mix(:,:  ) * p_frld(:,:)       &
1792            &                     + palbi         (:,:,1) * zicefr(:,:) ) )
1793      END SELECT
1794      IF( ln_dm2dc .AND. ln_cpl ) THEN   ! modify qsr to include the diurnal cycle
1795         zqsr_tot(:,:  ) = sbc_dcy( zqsr_tot(:,:  ) )
1796         DO jl=1,jpl
1797            zqsr_ice(:,:,jl) = sbc_dcy( zqsr_ice(:,:,jl) )
1798         ENDDO
1799      ENDIF
1800
1801#if defined key_lim3
1802      CALL wrk_alloc( jpi,jpj, zqsr_oce ) 
1803      ! --- solar flux over ocean --- !
1804      !         note: p_frld cannot be = 0 since we limit the ice concentration to amax
1805      zqsr_oce = 0._wp
1806      WHERE( p_frld /= 0._wp )  zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / p_frld(:,:)
1807
1808      IF( ln_mixcpl ) THEN   ;   qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) +  zqsr_oce(:,:)* zmsk(:,:)
1809      ELSE                   ;   qsr_oce(:,:) = zqsr_oce(:,:)   ;   ENDIF
1810
1811      CALL wrk_dealloc( jpi,jpj, zqsr_oce ) 
1812#endif
1813
1814      IF( ln_mixcpl ) THEN
1815         qsr_tot(:,:) = qsr(:,:) * p_frld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk
1816         qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) +  zqsr_tot(:,:)* zmsk(:,:)
1817         DO jl=1,jpl
1818            qsr_ice(:,:,jl) = qsr_ice(:,:,jl) * xcplmask(:,:,0) +  zqsr_ice(:,:,jl)* zmsk(:,:)
1819         ENDDO
1820      ELSE
1821         qsr_tot(:,:  ) = zqsr_tot(:,:  )
1822         qsr_ice(:,:,:) = zqsr_ice(:,:,:)
1823      ENDIF
1824
1825      !                                                      ! ========================= !
1826      SELECT CASE( TRIM( sn_rcv_dqnsdt%cldes ) )             !          d(qns)/dt        !
1827      !                                                      ! ========================= !
1828      CASE ('coupled')
1829         IF ( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN
1830            zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl)
1831         ELSE
1832            ! Set all category values equal for the moment
1833            DO jl=1,jpl
1834               zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1)
1835            ENDDO
1836         ENDIF
1837      END SELECT
1838     
1839      IF( ln_mixcpl ) THEN
1840         DO jl=1,jpl
1841            dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:)
1842         ENDDO
1843      ELSE
1844         dqns_ice(:,:,:) = zdqns_ice(:,:,:)
1845      ENDIF
1846     
1847      !                                                      ! ========================= !
1848      SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) )             !    topmelt and botmelt    !
1849      !                                                      ! ========================= !
1850      CASE ('coupled')
1851         topmelt(:,:,:)=frcv(jpr_topm)%z3(:,:,:)
1852         botmelt(:,:,:)=frcv(jpr_botm)%z3(:,:,:)
1853      END SELECT
1854
1855      ! Surface transimission parameter io (Maykut Untersteiner , 1971 ; Ebert and Curry, 1993 )
1856      ! Used for LIM2 and LIM3
1857      ! Coupled case: since cloud cover is not received from atmosphere
1858      !               ===> used prescribed cloud fraction representative for polar oceans in summer (0.81)
1859      fr1_i0(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice )
1860      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice )
1861
1862      CALL wrk_dealloc( jpi,jpj,     zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot )
1863      CALL wrk_dealloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice )
1864      !
1865      IF( nn_timing == 1 )  CALL timing_stop('sbc_cpl_ice_flx')
1866      !
1867   END SUBROUTINE sbc_cpl_ice_flx
1868   
1869   
1870   SUBROUTINE sbc_cpl_snd( kt )
1871      !!----------------------------------------------------------------------
1872      !!             ***  ROUTINE sbc_cpl_snd  ***
1873      !!
1874      !! ** Purpose :   provide the ocean-ice informations to the atmosphere
1875      !!
1876      !! ** Method  :   send to the atmosphere through a call to cpl_snd
1877      !!              all the needed fields (as defined in sbc_cpl_init)
1878      !!----------------------------------------------------------------------
1879      INTEGER, INTENT(in) ::   kt
1880      !
1881      INTEGER ::   ji, jj, jl   ! dummy loop indices
1882      INTEGER ::   isec, info   ! local integer
1883      REAL(wp) ::   zumax, zvmax
1884      REAL(wp), POINTER, DIMENSION(:,:)   ::   zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1
1885      REAL(wp), POINTER, DIMENSION(:,:,:) ::   ztmp3, ztmp4   
1886      !!----------------------------------------------------------------------
1887      !
1888      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_snd')
1889      !
1890      CALL wrk_alloc( jpi,jpj, zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 )
1891      CALL wrk_alloc( jpi,jpj,jpl, ztmp3, ztmp4 )
1892
1893      isec = ( kt - nit000 ) * NINT(rdttra(1))        ! date of exchanges
1894
1895      zfr_l(:,:) = 1.- fr_i(:,:)
1896      !                                                      ! ------------------------- !
1897      !                                                      !    Surface temperature    !   in Kelvin
1898      !                                                      ! ------------------------- !
1899      IF( ssnd(jps_toce)%laction .OR. ssnd(jps_tice)%laction .OR. ssnd(jps_tmix)%laction ) THEN
1900         
1901         IF ( nn_components == jp_iam_opa ) THEN
1902            ztmp1(:,:) = tsn(:,:,1,jp_tem)   ! send temperature as it is (potential or conservative) -> use of ln_useCT on the received part
1903         ELSE
1904            ! we must send the surface potential temperature
1905            IF( ln_useCT )  THEN    ;   ztmp1(:,:) = eos_pt_from_ct( tsn(:,:,1,jp_tem), tsn(:,:,1,jp_sal) )
1906            ELSE                    ;   ztmp1(:,:) = tsn(:,:,1,jp_tem)
1907            ENDIF
1908            !
1909            SELECT CASE( sn_snd_temp%cldes)
1910            CASE( 'oce only'             )   ;   ztmp1(:,:) =   ztmp1(:,:) + rt0
1911            CASE( 'oce and ice'          )   ;   ztmp1(:,:) =   ztmp1(:,:) + rt0
1912               SELECT CASE( sn_snd_temp%clcat )
1913               CASE( 'yes' )   
1914                  ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl)
1915               CASE( 'no' )
1916                  WHERE( SUM( a_i, dim=3 ) /= 0. )
1917                     ztmp3(:,:,1) = SUM( tn_ice * a_i, dim=3 ) / SUM( a_i, dim=3 )
1918                  ELSEWHERE
1919                     ztmp3(:,:,1) = rt0 ! TODO: Is freezing point a good default? (Maybe SST is better?)
1920                  END WHERE
1921               CASE default   ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' )
1922               END SELECT
1923            CASE( 'weighted oce and ice' )   ;   ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:)   
1924               SELECT CASE( sn_snd_temp%clcat )
1925               CASE( 'yes' )   
1926                  ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl)
1927               CASE( 'no' )
1928                  ztmp3(:,:,:) = 0.0
1929                  DO jl=1,jpl
1930                     ztmp3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl)
1931                  ENDDO
1932               CASE default                  ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' )
1933               END SELECT
1934            CASE( 'mixed oce-ice'        )   
1935               ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 
1936               DO jl=1,jpl
1937                  ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl)
1938               ENDDO
1939            CASE default                     ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' )
1940            END SELECT
1941         ENDIF
1942         IF( ssnd(jps_toce)%laction )   CALL cpl_snd( jps_toce, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info )
1943         IF( ssnd(jps_tice)%laction )   CALL cpl_snd( jps_tice, isec, ztmp3, info )
1944         IF( ssnd(jps_tmix)%laction )   CALL cpl_snd( jps_tmix, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info )
1945      ENDIF
1946      !                                                      ! ------------------------- !
1947      !                                                      !           Albedo          !
1948      !                                                      ! ------------------------- !
1949      IF( ssnd(jps_albice)%laction ) THEN                         ! ice
1950         SELECT CASE( sn_snd_alb%cldes )
1951         CASE( 'ice'          )   ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl)
1952         CASE( 'weighted ice' )   ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) * a_i(:,:,1:jpl)
1953         CASE default             ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_alb%cldes' )
1954         END SELECT
1955         CALL cpl_snd( jps_albice, isec, ztmp3, info )
1956      ENDIF
1957      IF( ssnd(jps_albmix)%laction ) THEN                         ! mixed ice-ocean
1958         ztmp1(:,:) = albedo_oce_mix(:,:) * zfr_l(:,:)
1959         DO jl=1,jpl
1960            ztmp1(:,:) = ztmp1(:,:) + alb_ice(:,:,jl) * a_i(:,:,jl)
1961         ENDDO
1962         CALL cpl_snd( jps_albmix, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info )
1963      ENDIF
1964      !                                                      ! ------------------------- !
1965      !                                                      !  Ice fraction & Thickness !
1966      !                                                      ! ------------------------- !
1967      ! Send ice fraction field to atmosphere
1968      IF( ssnd(jps_fice)%laction ) THEN
1969         SELECT CASE( sn_snd_thick%clcat )
1970         CASE( 'yes' )   ;   ztmp3(:,:,1:jpl) =  a_i(:,:,1:jpl)
1971         CASE( 'no'  )   ;   ztmp3(:,:,1    ) = fr_i(:,:      )
1972         CASE default    ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' )
1973         END SELECT
1974         IF( ssnd(jps_fice)%laction )   CALL cpl_snd( jps_fice, isec, ztmp3, info )
1975      ENDIF
1976     
1977      ! Send ice fraction field to OPA (sent by SAS in SAS-OPA coupling)
1978      IF( ssnd(jps_fice2)%laction ) THEN
1979         ztmp3(:,:,1) = fr_i(:,:)
1980         IF( ssnd(jps_fice2)%laction )   CALL cpl_snd( jps_fice2, isec, ztmp3, info )
1981      ENDIF
1982
1983      ! Send ice and snow thickness field
1984      IF( ssnd(jps_hice)%laction .OR. ssnd(jps_hsnw)%laction ) THEN
1985         SELECT CASE( sn_snd_thick%cldes)
1986         CASE( 'none'                  )       ! nothing to do
1987         CASE( 'weighted ice and snow' )   
1988            SELECT CASE( sn_snd_thick%clcat )
1989            CASE( 'yes' )   
1990               ztmp3(:,:,1:jpl) =  ht_i(:,:,1:jpl) * a_i(:,:,1:jpl)
1991               ztmp4(:,:,1:jpl) =  ht_s(:,:,1:jpl) * a_i(:,:,1:jpl)
1992            CASE( 'no' )
1993               ztmp3(:,:,:) = 0.0   ;  ztmp4(:,:,:) = 0.0
1994               DO jl=1,jpl
1995                  ztmp3(:,:,1) = ztmp3(:,:,1) + ht_i(:,:,jl) * a_i(:,:,jl)
1996                  ztmp4(:,:,1) = ztmp4(:,:,1) + ht_s(:,:,jl) * a_i(:,:,jl)
1997               ENDDO
1998            CASE default                  ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' )
1999            END SELECT
2000         CASE( 'ice and snow'         )   
2001            SELECT CASE( sn_snd_thick%clcat )
2002            CASE( 'yes' )
2003               ztmp3(:,:,1:jpl) = ht_i(:,:,1:jpl)
2004               ztmp4(:,:,1:jpl) = ht_s(:,:,1:jpl)
2005            CASE( 'no' )
2006               WHERE( SUM( a_i, dim=3 ) /= 0. )
2007                  ztmp3(:,:,1) = SUM( ht_i * a_i, dim=3 ) / SUM( a_i, dim=3 )
2008                  ztmp4(:,:,1) = SUM( ht_s * a_i, dim=3 ) / SUM( a_i, dim=3 )
2009               ELSEWHERE
2010                 ztmp3(:,:,1) = 0.
2011                 ztmp4(:,:,1) = 0.
2012               END WHERE
2013            CASE default                  ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' )
2014            END SELECT
2015         CASE default                     ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%cldes' )
2016         END SELECT
2017         IF( ssnd(jps_hice)%laction )   CALL cpl_snd( jps_hice, isec, ztmp3, info )
2018         IF( ssnd(jps_hsnw)%laction )   CALL cpl_snd( jps_hsnw, isec, ztmp4, info )
2019      ENDIF
2020      !
2021#if defined key_cpl_carbon_cycle
2022      !                                                      ! ------------------------- !
2023      !                                                      !  CO2 flux from PISCES     !
2024      !                                                      ! ------------------------- !
2025      IF( ssnd(jps_co2)%laction )   CALL cpl_snd( jps_co2, isec, RESHAPE ( oce_co2, (/jpi,jpj,1/) ) , info )
2026      !
2027#endif
2028      !                                                      ! ------------------------- !
2029      IF( ssnd(jps_ocx1)%laction ) THEN                      !      Surface current      !
2030         !                                                   ! ------------------------- !
2031         !   
2032         !                                                  j+1   j     -----V---F
2033         ! surface velocity always sent from T point                     !       |
2034         !                                                        j      |   T   U
2035         !                                                               |       |
2036         !                                                   j    j-1   -I-------|
2037         !                                               (for I)         |       |
2038         !                                                              i-1  i   i
2039         !                                                               i      i+1 (for I)
2040         IF( nn_components == jp_iam_opa ) THEN
2041            zotx1(:,:) = un(:,:,1) 
2042            zoty1(:,:) = vn(:,:,1) 
2043         ELSE       
2044            SELECT CASE( TRIM( sn_snd_crt%cldes ) )
2045            CASE( 'oce only'             )      ! C-grid ==> T
2046               DO jj = 2, jpjm1
2047                  DO ji = fs_2, fs_jpim1   ! vector opt.
2048                     zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj  ,1) )
2049                     zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji  ,jj-1,1) ) 
2050                  END DO
2051               END DO
2052            CASE( 'weighted oce and ice' )   
2053               SELECT CASE ( cp_ice_msh )
2054               CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
2055                  DO jj = 2, jpjm1
2056                     DO ji = fs_2, fs_jpim1   ! vector opt.
2057                        zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj) 
2058                        zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj)
2059                        zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj)
2060                        zity1(ji,jj) = 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj)
2061                     END DO
2062                  END DO
2063               CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
2064                  DO jj = 2, jpjm1
2065                     DO ji = 2, jpim1   ! NO vector opt.
2066                        zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj) 
2067                        zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj) 
2068                        zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     &
2069                           &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
2070                        zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     &
2071                           &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
2072                     END DO
2073                  END DO
2074               CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
2075                  DO jj = 2, jpjm1
2076                     DO ji = 2, jpim1   ! NO vector opt.
2077                        zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj) 
2078                        zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj) 
2079                        zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     &
2080                           &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
2081                        zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     &
2082                           &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
2083                     END DO
2084                  END DO
2085               END SELECT
2086               CALL lbc_lnk( zitx1, 'T', -1. )   ;   CALL lbc_lnk( zity1, 'T', -1. )
2087            CASE( 'mixed oce-ice'        )
2088               SELECT CASE ( cp_ice_msh )
2089               CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
2090                  DO jj = 2, jpjm1
2091                     DO ji = fs_2, fs_jpim1   ! vector opt.
2092                        zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj)   &
2093                           &         + 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj)
2094                        zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj)   &
2095                           &         + 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj)
2096                     END DO
2097                  END DO
2098               CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
2099                  DO jj = 2, jpjm1
2100                     DO ji = 2, jpim1   ! NO vector opt.
2101                        zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   &   
2102                           &         + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     &
2103                           &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
2104                        zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
2105                           &         + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     &
2106                           &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
2107                     END DO
2108                  END DO
2109               CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
2110                  DO jj = 2, jpjm1
2111                     DO ji = 2, jpim1   ! NO vector opt.
2112                        zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   &   
2113                           &         + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     &
2114                           &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
2115                        zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
2116                           &         + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     &
2117                           &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
2118                     END DO
2119                  END DO
2120               END SELECT
2121            END SELECT
2122            CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. )   ;   CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. )
2123            !
2124         ENDIF
2125         !
2126         !
2127         IF( TRIM( sn_snd_crt%clvor ) == 'eastward-northward' ) THEN             ! Rotation of the components
2128            !                                                                     ! Ocean component
2129            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocx1)%clgrid, 'ij->e', ztmp1 )       ! 1st component
2130            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocx1)%clgrid, 'ij->n', ztmp2 )       ! 2nd component
2131            zotx1(:,:) = ztmp1(:,:)                                                   ! overwrite the components
2132            zoty1(:,:) = ztmp2(:,:)
2133            IF( ssnd(jps_ivx1)%laction ) THEN                                     ! Ice component
2134               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->e', ztmp1 )    ! 1st component
2135               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->n', ztmp2 )    ! 2nd component
2136               zitx1(:,:) = ztmp1(:,:)                                                ! overwrite the components
2137               zity1(:,:) = ztmp2(:,:)
2138            ENDIF
2139         ENDIF
2140         !
2141         ! spherical coordinates to cartesian -> 2 components to 3 components
2142         IF( TRIM( sn_snd_crt%clvref ) == 'cartesian' ) THEN
2143            ztmp1(:,:) = zotx1(:,:)                     ! ocean currents
2144            ztmp2(:,:) = zoty1(:,:)
2145            CALL oce2geo ( ztmp1, ztmp2, 'T', zotx1, zoty1, zotz1 )
2146            !
2147            IF( ssnd(jps_ivx1)%laction ) THEN           ! ice velocities
2148               ztmp1(:,:) = zitx1(:,:)
2149               ztmp1(:,:) = zity1(:,:)
2150               CALL oce2geo ( ztmp1, ztmp2, 'T', zitx1, zity1, zitz1 )
2151            ENDIF
2152         ENDIF
2153         !
2154         IF( ssnd(jps_ocx1)%laction )   CALL cpl_snd( jps_ocx1, isec, RESHAPE ( zotx1, (/jpi,jpj,1/) ), info )   ! ocean x current 1st grid
2155         IF( ssnd(jps_ocy1)%laction )   CALL cpl_snd( jps_ocy1, isec, RESHAPE ( zoty1, (/jpi,jpj,1/) ), info )   ! ocean y current 1st grid
2156         IF( ssnd(jps_ocz1)%laction )   CALL cpl_snd( jps_ocz1, isec, RESHAPE ( zotz1, (/jpi,jpj,1/) ), info )   ! ocean z current 1st grid
2157         !
2158         IF( ssnd(jps_ivx1)%laction )   CALL cpl_snd( jps_ivx1, isec, RESHAPE ( zitx1, (/jpi,jpj,1/) ), info )   ! ice   x current 1st grid
2159         IF( ssnd(jps_ivy1)%laction )   CALL cpl_snd( jps_ivy1, isec, RESHAPE ( zity1, (/jpi,jpj,1/) ), info )   ! ice   y current 1st grid
2160         IF( ssnd(jps_ivz1)%laction )   CALL cpl_snd( jps_ivz1, isec, RESHAPE ( zitz1, (/jpi,jpj,1/) ), info )   ! ice   z current 1st grid
2161         !
2162      ENDIF
2163      !
2164      !                                                      ! ------------------------- !
2165      !                                                      !  Surface current to waves !
2166      !                                                      ! ------------------------- !
2167      IF( ssnd(jps_ocxw)%laction .OR. ssnd(jps_ocyw)%laction ) THEN 
2168          !     
2169          !                                                  j+1  j     -----V---F
2170          ! surface velocity always sent from T point                    !       |
2171          !                                                       j      |   T   U
2172          !                                                              |       |
2173          !                                                   j   j-1   -I-------|
2174          !                                               (for I)        |       |
2175          !                                                             i-1  i   i
2176          !                                                              i      i+1 (for I)
2177          SELECT CASE( TRIM( sn_snd_crtw%cldes ) ) 
2178          CASE( 'oce only'             )      ! C-grid ==> T
2179             DO jj = 2, jpjm1 
2180                DO ji = fs_2, fs_jpim1   ! vector opt.
2181                   zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj  ,1) ) 
2182                   zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji , jj-1,1) ) 
2183                END DO
2184             END DO
2185          CASE( 'weighted oce and ice' )   
2186             SELECT CASE ( cp_ice_msh ) 
2187             CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
2188                DO jj = 2, jpjm1 
2189                   DO ji = fs_2, fs_jpim1   ! vector opt.
2190                      zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj)   
2191                      zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj) 
2192                      zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj) 
2193                      zity1(ji,jj) = 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj) 
2194                   END DO
2195                END DO
2196             CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
2197                DO jj = 2, jpjm1 
2198                   DO ji = 2, jpim1   ! NO vector opt.
2199                      zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   
2200                      zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   
2201                      zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     & 
2202                         &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2203                      zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     & 
2204                         &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2205                   END DO
2206                END DO
2207             CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
2208                DO jj = 2, jpjm1 
2209                   DO ji = 2, jpim1   ! NO vector opt.
2210                      zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   
2211                      zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   
2212                      zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     & 
2213                         &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2214                      zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     & 
2215                         &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2216                   END DO
2217                END DO
2218             END SELECT
2219             CALL lbc_lnk( zitx1, 'T', -1. )   ;   CALL lbc_lnk( zity1, 'T', -1. ) 
2220          CASE( 'mixed oce-ice'        ) 
2221             SELECT CASE ( cp_ice_msh ) 
2222             CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
2223                DO jj = 2, jpjm1 
2224                   DO ji = fs_2, fs_jpim1   ! vector opt.
2225                      zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj)   & 
2226                         &         + 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj) 
2227                      zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
2228                         &         + 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj) 
2229                   END DO
2230                END DO
2231             CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
2232                DO jj = 2, jpjm1 
2233                   DO ji = 2, jpim1   ! NO vector opt.
2234                      zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   &   
2235                         &         + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     & 
2236                         &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2237                      zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
2238                         &         + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     & 
2239                         &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2240                   END DO
2241                END DO
2242             CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
2243                DO jj = 2, jpjm1 
2244                   DO ji = 2, jpim1   ! NO vector opt.
2245                      zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj  ,1) ) * zfr_l(ji,jj)   &   
2246                         &         + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     & 
2247                         &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2248                      zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
2249                         &         + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     & 
2250                         &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj) 
2251                   END DO
2252                END DO
2253             END SELECT
2254          END SELECT
2255         CALL lbc_lnk( zotx1, ssnd(jps_ocxw)%clgrid, -1. )   ;   CALL lbc_lnk( zoty1, ssnd(jps_ocyw)%clgrid, -1. ) 
2256         !
2257         !
2258         IF( TRIM( sn_snd_crtw%clvor ) == 'eastward-northward' ) THEN             ! Rotation of the components
2259         !                                                                        ! Ocean component
2260            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocxw)%clgrid, 'ij->e', ztmp1 )       ! 1st component 
2261            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocxw)%clgrid, 'ij->n', ztmp2 )       ! 2nd component 
2262            zotx1(:,:) = ztmp1(:,:)                                                   ! overwrite the components 
2263            zoty1(:,:) = ztmp2(:,:) 
2264            IF( ssnd(jps_ivx1)%laction ) THEN                                     ! Ice component
2265               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->e', ztmp1 )    ! 1st component 
2266               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->n', ztmp2 )    ! 2nd component 
2267               zitx1(:,:) = ztmp1(:,:)                                                ! overwrite the components 
2268               zity1(:,:) = ztmp2(:,:) 
2269            ENDIF
2270         ENDIF 
2271         !
2272!         ! spherical coordinates to cartesian -> 2 components to 3 components
2273!         IF( TRIM( sn_snd_crtw%clvref ) == 'cartesian' ) THEN
2274!            ztmp1(:,:) = zotx1(:,:)                     ! ocean currents
2275!            ztmp2(:,:) = zoty1(:,:)
2276!            CALL oce2geo ( ztmp1, ztmp2, 'T', zotx1, zoty1, zotz1 )
2277!            !
2278!            IF( ssnd(jps_ivx1)%laction ) THEN           ! ice velocities
2279!               ztmp1(:,:) = zitx1(:,:)
2280!               ztmp1(:,:) = zity1(:,:)
2281!               CALL oce2geo ( ztmp1, ztmp2, 'T', zitx1, zity1, zitz1 )
2282!            ENDIF
2283!         ENDIF
2284         !
2285         IF( ssnd(jps_ocxw)%laction )   CALL cpl_snd( jps_ocxw, isec, RESHAPE ( zotx1, (/jpi,jpj,1/) ), info )   ! ocean x current 1st grid
2286         IF( ssnd(jps_ocyw)%laction )   CALL cpl_snd( jps_ocyw, isec, RESHAPE ( zoty1, (/jpi,jpj,1/) ), info )   ! ocean y current 1st grid
2287         
2288      ENDIF 
2289      !
2290      IF( ssnd(jps_ficet)%laction ) THEN
2291         CALL cpl_snd( jps_ficet, isec, RESHAPE ( fr_i, (/jpi,jpj,1/) ), info ) 
2292      END IF 
2293      !                                                      ! ------------------------- !
2294      !                                                      !   Water levels to waves   !
2295      !                                                      ! ------------------------- !
2296      IF( ssnd(jps_wlev)%laction ) THEN
2297         IF( ln_apr_dyn ) THEN 
2298            IF( kt /= nit000 ) THEN 
2299               ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) ) 
2300            ELSE 
2301               ztmp1(:,:) = sshb(:,:) 
2302            ENDIF 
2303         ELSE 
2304            ztmp1(:,:) = sshn(:,:) 
2305         ENDIF 
2306         CALL cpl_snd( jps_wlev  , isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info ) 
2307      END IF 
2308      !
2309      !  Fields sent by OPA to SAS when doing OPA<->SAS coupling
2310      !                                                        ! SSH
2311      IF( ssnd(jps_ssh )%laction )  THEN
2312         !                          ! removed inverse barometer ssh when Patm
2313         !                          forcing is used (for sea-ice dynamics)
2314         IF( ln_apr_dyn ) THEN   ;   ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) )
2315         ELSE                    ;   ztmp1(:,:) = sshn(:,:)
2316         ENDIF
2317         CALL cpl_snd( jps_ssh   , isec, RESHAPE ( ztmp1            , (/jpi,jpj,1/) ), info )
2318
2319      ENDIF
2320      !                                                        ! SSS
2321      IF( ssnd(jps_soce  )%laction )  THEN
2322         CALL cpl_snd( jps_soce  , isec, RESHAPE ( tsn(:,:,1,jp_sal), (/jpi,jpj,1/) ), info )
2323      ENDIF
2324      !                                                        ! first T level thickness
2325      IF( ssnd(jps_e3t1st )%laction )  THEN
2326         CALL cpl_snd( jps_e3t1st, isec, RESHAPE ( fse3t_n(:,:,1)   , (/jpi,jpj,1/) ), info )
2327      ENDIF
2328      !                                                        ! Qsr fraction
2329      IF( ssnd(jps_fraqsr)%laction )  THEN
2330         CALL cpl_snd( jps_fraqsr, isec, RESHAPE ( fraqsr_1lev(:,:) , (/jpi,jpj,1/) ), info )
2331      ENDIF
2332      !
2333      !  Fields sent by SAS to OPA when OASIS coupling
2334      !                                                        ! Solar heat flux
2335      IF( ssnd(jps_qsroce)%laction )  CALL cpl_snd( jps_qsroce, isec, RESHAPE ( qsr , (/jpi,jpj,1/) ), info )
2336      IF( ssnd(jps_qnsoce)%laction )  CALL cpl_snd( jps_qnsoce, isec, RESHAPE ( qns , (/jpi,jpj,1/) ), info )
2337      IF( ssnd(jps_oemp  )%laction )  CALL cpl_snd( jps_oemp  , isec, RESHAPE ( emp , (/jpi,jpj,1/) ), info )
2338      IF( ssnd(jps_sflx  )%laction )  CALL cpl_snd( jps_sflx  , isec, RESHAPE ( sfx , (/jpi,jpj,1/) ), info )
2339      IF( ssnd(jps_otx1  )%laction )  CALL cpl_snd( jps_otx1  , isec, RESHAPE ( utau, (/jpi,jpj,1/) ), info )
2340      IF( ssnd(jps_oty1  )%laction )  CALL cpl_snd( jps_oty1  , isec, RESHAPE ( vtau, (/jpi,jpj,1/) ), info )
2341      IF( ssnd(jps_rnf   )%laction )  CALL cpl_snd( jps_rnf   , isec, RESHAPE ( rnf , (/jpi,jpj,1/) ), info )
2342      IF( ssnd(jps_taum  )%laction )  CALL cpl_snd( jps_taum  , isec, RESHAPE ( taum, (/jpi,jpj,1/) ), info )
2343
2344      CALL wrk_dealloc( jpi,jpj, zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 )
2345      CALL wrk_dealloc( jpi,jpj,jpl, ztmp3, ztmp4 )
2346      !
2347      IF( nn_timing == 1 )  CALL timing_stop('sbc_cpl_snd')
2348      !
2349   END SUBROUTINE sbc_cpl_snd
2350   
2351   !!======================================================================
2352END MODULE sbccpl
Note: See TracBrowser for help on using the repository browser.