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/r6232_HZG_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/r6232_HZG_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 7878

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

Add Phillips vertical Stokes drift parameterization as in the HZG wave branch

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