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 @ 5817

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

Bug fix in variable name, variable declarations

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