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

source: branches/2015/dev_r5218_CNRS17_coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 5299

Last change on this file since 5299 was 5299, checked in by smasson, 9 years ago

dev_r5218_CNRS17_coupling: bugfixes

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