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

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

source: branches/UKMO/dev_r5518_ukv_mslp/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 5801

Last change on this file since 5801 was 5801, checked in by jcastill, 9 years ago

Version as in svn://fcm3/NEMO_svn/NEMO/branches/dev/cloneill/r10320_vn3.6_ukv_mslp

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