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

source: branches/2015/nemo_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 7951

Last change on this file since 7951 was 7951, checked in by clem, 7 years ago

correct output, see ticket #1888

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