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

source: branches/UKMO/dev_r5518_medusa_cpl_rh/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 6171

Last change on this file since 6171 was 6171, checked in by frrh, 8 years ago

Add MEDUSA coupling interface code. Much of this is temporary and needs further
attention, e.g. for development purposes we allocate temporarry receiving
arrays and we also currently feature the original "magic number"
approach to identifying output fields which must not be allowed to
be a permanent solution and in any case should not get past a code review.

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