New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbccpl.F90 in branches/UKMO/r5936_hadgem3_cplfld/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/r5936_hadgem3_cplfld/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 7150

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

Move one line of code to allow merging with another branch

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