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

source: trunk/NEMO/OPA_SRC/SBC/sbccpl.F90 @ 1481

Last change on this file since 1481 was 1472, checked in by smasson, 15 years ago

small bugfix in the coupling interface when sending surface current, see ticket:455

  • Property svn:keywords set to Id
File size: 78.9 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  !  06-2007  (R. Redler, N. Keenlyside, W. Park) Original code split into flxmod & taumod
7   !!            3.0  !  02-2008  (G. Madec, C Talandier)  surface module
8   !!            3.1  !  02-2009  (S. Masson, E. Maisonave, A. Caubel) generic coupled interface
9   !!----------------------------------------------------------------------
10#if defined key_oasis3 || defined key_oasis4
11   !!----------------------------------------------------------------------
12   !!   'key_oasis3' or 'key_oasis4'   Coupled Ocean/Atmosphere formulation
13   !!----------------------------------------------------------------------
14   !!   namsbc_cpl      : coupled formulation namlist
15   !!   sbc_cpl_init    : initialisation of the coupled exchanges
16   !!   sbc_cpl_rcv     : receive fields from the atmosphere over the ocean (ocean only)
17   !!                     receive stress from the atmosphere over the ocean (ocean-ice case)
18   !!   sbc_cpl_ice_tau : receive stress from the atmosphere over ice
19   !!   sbc_cpl_ice_flx : receive fluxes from the atmosphere over ice
20   !!   sbc_cpl_snd     : send     fields to the atmosphere
21   !!----------------------------------------------------------------------
22   USE dom_oce         ! ocean space and time domain
23   USE sbc_oce         ! Surface boundary condition: ocean fields
24   USE sbc_ice         ! Surface boundary condition: ice fields
25#if defined key_lim3
26   USE par_ice         ! ice parameters
27#endif
28#if defined key_lim2
29   USE ice_2, ONLY : hicif, hsnif          ! Ice and Snow thickness
30#endif
31   USE cpl_oasis3      ! OASIS3 coupling
32   USE geo2ocean       !
33   USE restart         !
34   USE oce   , ONLY : tn, un, vn
35   USE phycst, ONLY : rt0
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 lbclnk          ! ocean lateral boundary conditions (or mpp link)
41   USE mod_prism_proto ! OASIS3 prism module: PRISM_* variables...
42   USE phycst, ONLY : xlsn, rhosn
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC   sbc_cpl_rcv       ! routine called by sbc_ice_lim(_2).F90
47   PUBLIC   sbc_cpl_snd       ! routine called by step.F90
48   PUBLIC   sbc_cpl_ice_tau   ! routine called by sbc_ice_lim(_2).F90
49   PUBLIC   sbc_cpl_ice_flx   ! routine called by sbc_ice_lim(_2).F90
50   
51   INTEGER, PARAMETER ::   jpr_otx1   =  1            ! 3 atmosphere-ocean stress components on grid 1
52   INTEGER, PARAMETER ::   jpr_oty1   =  2            !
53   INTEGER, PARAMETER ::   jpr_otz1   =  3            !
54   INTEGER, PARAMETER ::   jpr_otx2   =  4            ! 3 atmosphere-ocean stress components on grid 2
55   INTEGER, PARAMETER ::   jpr_oty2   =  5            !
56   INTEGER, PARAMETER ::   jpr_otz2   =  6            !
57   INTEGER, PARAMETER ::   jpr_itx1   =  7            ! 3 atmosphere-ice   stress components on grid 1
58   INTEGER, PARAMETER ::   jpr_ity1   =  8            !
59   INTEGER, PARAMETER ::   jpr_itz1   =  9            !
60   INTEGER, PARAMETER ::   jpr_itx2   = 10            ! 3 atmosphere-ice   stress components on grid 2
61   INTEGER, PARAMETER ::   jpr_ity2   = 11            !
62   INTEGER, PARAMETER ::   jpr_itz2   = 12            !
63   INTEGER, PARAMETER ::   jpr_qsroce = 13            ! Qsr above the ocean
64   INTEGER, PARAMETER ::   jpr_qsrice = 14            ! Qsr above the ice
65   INTEGER, PARAMETER ::   jpr_qsrmix = 15 
66   INTEGER, PARAMETER ::   jpr_qnsoce = 16            ! Qns above the ocean
67   INTEGER, PARAMETER ::   jpr_qnsice = 17            ! Qns above the ice
68   INTEGER, PARAMETER ::   jpr_qnsmix = 18
69   INTEGER, PARAMETER ::   jpr_rain   = 19            ! total liquid precipitation (rain)
70   INTEGER, PARAMETER ::   jpr_snow   = 20            ! solid precipitation over the ocean (snow)
71   INTEGER, PARAMETER ::   jpr_tevp   = 21            ! total evaporation
72   INTEGER, PARAMETER ::   jpr_ievp   = 22            ! solid evaporation (sublimation)
73   INTEGER, PARAMETER ::   jpr_sbpr   = 23            ! sublimation - liquid precipitation - solid precipitation
74   INTEGER, PARAMETER ::   jpr_semp   = 24            ! solid freshwater budget (sublimation - snow)
75   INTEGER, PARAMETER ::   jpr_oemp   = 25            ! ocean freshwater budget (evap - precip)
76   INTEGER, PARAMETER ::   jpr_w10m   = 26            !
77   INTEGER, PARAMETER ::   jpr_dqnsdt = 27            !
78   INTEGER, PARAMETER ::   jpr_rnf    = 28            !
79   INTEGER, PARAMETER ::   jpr_cal    = 29            !
80   INTEGER, PARAMETER ::   jprcv      = 29            ! total number of fields recieved
81   
82   INTEGER, PARAMETER ::   jps_fice   =  1            ! ice fraction
83   INTEGER, PARAMETER ::   jps_toce   =  2            ! ocean temperature
84   INTEGER, PARAMETER ::   jps_tice   =  3            ! ice   temperature
85   INTEGER, PARAMETER ::   jps_tmix   =  4            ! mixed temperature (ocean+ice)
86   INTEGER, PARAMETER ::   jps_albice =  5            ! ice   albedo
87   INTEGER, PARAMETER ::   jps_albmix =  6            ! mixed albedo
88   INTEGER, PARAMETER ::   jps_hice   =  7            ! ice  thickness
89   INTEGER, PARAMETER ::   jps_hsnw   =  8            ! snow thickness
90   INTEGER, PARAMETER ::   jps_ocx1   =  9            ! ocean current on grid 1
91   INTEGER, PARAMETER ::   jps_ocy1   = 10            !
92   INTEGER, PARAMETER ::   jps_ocz1   = 11            !
93   INTEGER, PARAMETER ::   jps_ivx1   = 12            ! ice   current on grid 1
94   INTEGER, PARAMETER ::   jps_ivy1   = 13            !
95   INTEGER, PARAMETER ::   jps_ivz1   = 14            !
96   INTEGER, PARAMETER ::   jpsnd      = 14            ! total number of fields sended
97   
98   !                                                         !!** namelist namsbc_cpl **
99   ! Send to the atmosphere                                   !
100   CHARACTER(len=100) ::   cn_snd_temperature = 'oce only'    ! 'oce only' 'weighted oce and ice' or 'mixed oce-ice'
101   CHARACTER(len=100) ::   cn_snd_albedo      = 'none'        ! 'none' 'weighted ice' or 'mixed oce-ice'
102   CHARACTER(len=100) ::   cn_snd_thickness   = 'none'        ! 'none' or 'weighted ice and snow'
103   CHARACTER(len=100) ::   cn_snd_crt_nature  = 'none'        ! 'none' 'oce only' 'weighted oce and ice' or 'mixed oce-ice'   
104   CHARACTER(len=100) ::   cn_snd_crt_refere  = 'spherical'   ! 'spherical' or 'cartesian'
105   CHARACTER(len=100) ::   cn_snd_crt_orient  = 'local grid'  ! 'eastward-northward' or 'local grid'
106   CHARACTER(len=100) ::   cn_snd_crt_grid    = 'T'           ! always at 'T' point
107   
108   ! Recieved from the atmosphere                             !
109   CHARACTER(len=100) ::   cn_rcv_tau_nature  = 'oce only'    ! 'oce only' 'oce and ice' or 'mixed oce-ice'
110   CHARACTER(len=100) ::   cn_rcv_tau_refere  = 'spherical'   ! 'spherical' or 'cartesian'
111   CHARACTER(len=100) ::   cn_rcv_tau_orient  = 'local grid'  ! 'eastward-northward' or 'local grid'
112   CHARACTER(len=100) ::   cn_rcv_tau_grid    = 'T'           ! 'T', 'U,V', 'U,V,I', 'T,I', or 'T,U,V'
113   CHARACTER(len=100) ::   cn_rcv_w10m        = 'none'        ! 'none' or 'coupled'
114   CHARACTER(len=100) ::   cn_rcv_dqnsdt      = 'none'        ! 'none' or 'coupled'
115   CHARACTER(len=100) ::   cn_rcv_qsr         = 'oce only'    ! 'oce only' 'conservative' 'oce and ice' or 'mixed oce-ice'
116   CHARACTER(len=100) ::   cn_rcv_qns         = 'oce only'    ! 'oce only' 'conservative' 'oce and ice' or 'mixed oce-ice'
117   CHARACTER(len=100) ::   cn_rcv_emp         = 'oce only'    ! 'oce only' 'conservative' or 'oce and ice'
118   CHARACTER(len=100) ::   cn_rcv_rnf         = 'coupled'     ! 'coupled' 'climato' or 'mixed'
119   CHARACTER(len=100) ::   cn_rcv_cal         = 'none'        ! 'none' or 'coupled'
120
121!!   CHARACTER(len=100), PUBLIC ::   cn_rcv_rnf   !: ???             ==>>  !!gm   treat this case in a different maner
122   
123   CHARACTER(len=100), DIMENSION(4) ::   cn_snd_crt           ! array combining cn_snd_crt_*
124   CHARACTER(len=100), DIMENSION(4) ::   cn_rcv_tau           ! array combining cn_rcv_tau_*
125
126   REAL(wp), DIMENSION(jpi,jpj)       ::   albedo_oce_mix     ! ocean albedo sent to atmosphere (mix clear/overcast sky)
127
128   REAL(wp), DIMENSION(jpi,jpj,jprcv) ::   frcv               ! all fields recieved from the atmosphere
129   INTEGER , DIMENSION(        jprcv) ::   nrcvinfo           ! OASIS info argument
130
131   !! Substitution
132#  include "vectopt_loop_substitute.h90"
133   !!----------------------------------------------------------------------
134   !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008)
135   !! $Id$
136   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
137   !!----------------------------------------------------------------------
138
139CONTAINS
140 
141   SUBROUTINE sbc_cpl_init( k_ice )     
142      !!----------------------------------------------------------------------
143      !!             ***  ROUTINE sbc_cpl_init  ***
144      !!
145      !! ** Purpose :   Initialisation of send and recieved information from
146      !!                the atmospheric component
147      !!
148      !! ** Method  : * Read namsbc_cpl namelist
149      !!              * define the receive interface
150      !!              * define the send    interface
151      !!              * initialise the OASIS coupler
152      !!----------------------------------------------------------------------
153      INTEGER, INTENT(in) ::   k_ice    ! ice management in the sbc (=0/1/2/3)
154      !!
155      INTEGER                      ::   jn           ! dummy loop index
156      REAL(wp), DIMENSION(jpi,jpj) ::   zacs, zaos   ! 2D workspace (clear & overcast sky albedos)
157      !!
158      NAMELIST/namsbc_cpl/  cn_snd_temperature, cn_snd_albedo    , cn_snd_thickness,                 &         
159         cn_snd_crt_nature, cn_snd_crt_refere , cn_snd_crt_orient, cn_snd_crt_grid , cn_rcv_w10m,    &
160         cn_rcv_tau_nature, cn_rcv_tau_refere , cn_rcv_tau_orient, cn_rcv_tau_grid ,                 &
161         cn_rcv_dqnsdt    , cn_rcv_qsr        , cn_rcv_qns       , cn_rcv_emp      , cn_rcv_rnf , cn_rcv_cal
162      !!---------------------------------------------------------------------
163
164      ! ================================ !
165      !      Namelist informations       !
166      ! ================================ !
167
168      REWIND( numnam )                    ! ... read namlist namsbc_cpl
169      READ  ( numnam, namsbc_cpl )
170
171      IF(lwp) THEN                        ! control print
172         WRITE(numout,*)
173         WRITE(numout,*)'sbc_cpl_init : namsbc_cpl namelist '
174         WRITE(numout,*)'~~~~~~~~~~~~'
175         WRITE(numout,*)'   received fields'
176         WRITE(numout,*)'       10m wind module                    cn_rcv_w10m        = ', cn_rcv_w10m 
177         WRITE(numout,*)'       surface stress - nature            cn_rcv_tau_nature  = ', cn_rcv_tau_nature
178         WRITE(numout,*)'                      - referential       cn_rcv_tau_refere  = ', cn_rcv_tau_refere
179         WRITE(numout,*)'                      - orientation       cn_rcv_tau_orient  = ', cn_rcv_tau_orient
180         WRITE(numout,*)'                      - mesh              cn_rcv_tau_grid    = ', cn_rcv_tau_grid
181         WRITE(numout,*)'       non-solar heat flux sensitivity    cn_rcv_dqnsdt      = ', cn_rcv_dqnsdt
182         WRITE(numout,*)'       solar heat flux                    cn_rcv_qsr         = ', cn_rcv_qsr 
183         WRITE(numout,*)'       non-solar heat flux                cn_rcv_qns         = ', cn_rcv_qns
184         WRITE(numout,*)'       freshwater budget                  cn_rcv_emp         = ', cn_rcv_emp
185         WRITE(numout,*)'       runoffs                            cn_rcv_rnf         = ', cn_rcv_rnf
186         WRITE(numout,*)'       calving                            cn_rcv_cal         = ', cn_rcv_cal 
187         WRITE(numout,*)'   sent fields'
188         WRITE(numout,*)'       surface temperature                cn_snd_temperature = ', cn_snd_temperature
189         WRITE(numout,*)'       albedo                             cn_snd_albedo      = ', cn_snd_albedo
190         WRITE(numout,*)'       ice/snow thickness                 cn_snd_thickness   = ', cn_snd_thickness 
191         WRITE(numout,*)'       surface current - nature           cn_snd_crt_nature  = ', cn_snd_crt_nature 
192         WRITE(numout,*)'                       - referential      cn_snd_crt_refere  = ', cn_snd_crt_refere 
193         WRITE(numout,*)'                       - orientation      cn_snd_crt_orient  = ', cn_snd_crt_orient
194         WRITE(numout,*)'                       - mesh             cn_snd_crt_grid    = ', cn_snd_crt_grid 
195      ENDIF
196
197      ! save current & stress in an array and suppress possible blank in the name
198      cn_snd_crt(1) = TRIM( cn_snd_crt_nature )   ;   cn_snd_crt(2) = TRIM( cn_snd_crt_refere )
199      cn_snd_crt(3) = TRIM( cn_snd_crt_orient )   ;   cn_snd_crt(4) = TRIM( cn_snd_crt_grid   )
200      cn_rcv_tau(1) = TRIM( cn_rcv_tau_nature )   ;   cn_rcv_tau(2) = TRIM( cn_rcv_tau_refere )
201      cn_rcv_tau(3) = TRIM( cn_rcv_tau_orient )   ;   cn_rcv_tau(4) = TRIM( cn_rcv_tau_grid   )
202     
203      ! ================================ !
204      !   Define the receive interface   !
205      ! ================================ !
206      nrcvinfo(:) = PRISM_NotDef   ! needed by nrcvinfo(jpr_otx1) if we do not receive ocea stress
207
208      ! for each field: define the OASIS name                              (srcv(:)%clname)
209      !                 define receive or not from the namelist parameters (srcv(:)%laction)
210      !                 define the north fold type of lbc                  (srcv(:)%nsgn)
211
212      ! default definitions of srcv
213      srcv(:)%laction = .FALSE.   ;   srcv(:)%clgrid = 'T'   ;   srcv(:)%nsgn = 1
214
215      !                                                      ! ------------------------- !
216      !                                                      ! ice and ocean wind stress !   
217      !                                                      ! ------------------------- !
218      !                                                           ! Name
219      srcv(jpr_otx1)%clname = 'O_OTaux1'      ! 1st ocean component on grid ONE (T or U)
220      srcv(jpr_oty1)%clname = 'O_OTauy1'      ! 2nd   -      -         -     -
221      srcv(jpr_otz1)%clname = 'O_OTauz1'      ! 3rd   -      -         -     -
222      srcv(jpr_otx2)%clname = 'O_OTaux2'      ! 1st ocean component on grid TWO (V)
223      srcv(jpr_oty2)%clname = 'O_OTauy2'      ! 2nd   -      -         -     -
224      srcv(jpr_otz2)%clname = 'O_OTauz2'      ! 3rd   -      -         -     -
225      !
226      srcv(jpr_itx1)%clname = 'O_ITaux1'      ! 1st  ice  component on grid ONE (T, F, I or U)
227      srcv(jpr_ity1)%clname = 'O_ITauy1'      ! 2nd   -      -         -     -
228      srcv(jpr_itz1)%clname = 'O_ITauz1'      ! 3rd   -      -         -     -
229      srcv(jpr_itx2)%clname = 'O_ITaux2'      ! 1st  ice  component on grid TWO (V)
230      srcv(jpr_ity2)%clname = 'O_ITauy2'      ! 2nd   -      -         -     -
231      srcv(jpr_itz2)%clname = 'O_ITauz2'      ! 3rd   -      -         -     -
232      !
233      srcv(jpr_otx1:jpr_itz2)%nsgn = -1                           ! Vectors: change of sign at north fold
234     
235      !                                                           ! Set grid and action
236      SELECT CASE( TRIM( cn_rcv_tau(4) ) )      !  'T', 'U,V', 'U,V,I', 'U,V,F', 'T,I', 'T,F', or 'T,U,V'
237      CASE( 'T' ) 
238         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
239         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
240         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
241      CASE( 'U,V' ) 
242         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
243         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
244         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'U'        ! ice components given at U-point
245         srcv(jpr_itx2:jpr_itz2)%clgrid  = 'V'        !           and           V-point
246         srcv(jpr_otx1:jpr_itz2)%laction = .TRUE.     ! receive oce and ice components on both grid 1 & 2
247      CASE( 'U,V,T' )
248         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
249         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
250         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'T'        ! ice components given at T-point
251         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
252         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
253      CASE( 'U,V,I' )
254         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
255         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
256         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'I'        ! ice components given at I-point
257         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
258         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
259      CASE( 'U,V,F' )
260         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'U'        ! oce components given at U-point
261         srcv(jpr_otx2:jpr_otz2)%clgrid  = 'V'        !           and           V-point
262         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'F'        ! ice components given at F-point
263         srcv(jpr_otx1:jpr_otz2)%laction = .TRUE.     ! receive oce components on grid 1 & 2
264         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1 only
265      CASE( 'T,I' ) 
266         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
267         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'I'        ! ice components given at I-point
268         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
269         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
270      CASE( 'T,F' ) 
271         srcv(jpr_otx1:jpr_itz2)%clgrid  = 'T'        ! oce and ice components given at T-point
272         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'F'        ! ice components given at F-point
273         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1
274         srcv(jpr_itx1:jpr_itz1)%laction = .TRUE.     ! receive ice components on grid 1
275      CASE( 'T,U,V' )
276         srcv(jpr_otx1:jpr_otz1)%clgrid  = 'T'        ! oce components given at T-point
277         srcv(jpr_itx1:jpr_itz1)%clgrid  = 'U'        ! ice components given at U-point
278         srcv(jpr_itx2:jpr_itz2)%clgrid  = 'V'        !           and           V-point
279         srcv(jpr_otx1:jpr_otz1)%laction = .TRUE.     ! receive oce components on grid 1 only
280         srcv(jpr_itx1:jpr_itz2)%laction = .TRUE.     ! receive ice components on grid 1 & 2
281      CASE default   
282         CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_rcv_tau(4)' )
283      END SELECT
284      !
285      IF( TRIM( cn_rcv_tau(2) ) == 'spherical' )   &           ! spherical: 3rd component not received
286         &     srcv( (/jpr_otz1, jpr_otz2, jpr_itz1, jpr_itz2/) )%laction = .FALSE. 
287      !
288      IF( TRIM( cn_rcv_tau(1) ) /= 'oce and ice' ) THEN        ! 'oce and ice' case ocean stress on ocean mesh used
289         srcv(jpr_itx1:jpr_itz2)%laction = .FALSE.    ! ice components not received
290         srcv(jpr_itx1)%clgrid = 'U'                  ! ocean stress used after its transformation
291         srcv(jpr_ity1)%clgrid = 'V'                  ! i.e. it is always at U- & V-points for i- & j-comp. resp.
292      ENDIF
293       
294      !                                                      ! ------------------------- !
295      !                                                      !    freshwater budget      !   E-P
296      !                                                      ! ------------------------- !
297      ! we suppose that atmosphere modele do not make the difference between precipiration (liquide or solid)
298      ! over ice of free ocean within the same atmospheric cell.cd
299      srcv(jpr_rain)%clname = 'OTotRain'      ! Rain = liquid precipitation
300      srcv(jpr_snow)%clname = 'OTotSnow'      ! Snow = solid precipitation
301      srcv(jpr_tevp)%clname = 'OTotEvap'      ! total evaporation (over oce + ice sublimation)
302      srcv(jpr_ievp)%clname = 'OIceEvap'      ! evaporation over ice = sublimation
303      srcv(jpr_sbpr)%clname = 'OSubMPre'      ! sublimation - liquid precipitation - solid precipitation
304      srcv(jpr_semp)%clname = 'OISubMSn'      ! ice solid water budget = sublimation - solid precipitation
305      srcv(jpr_oemp)%clname = 'OOEvaMPr'      ! ocean water budget = ocean Evap - ocean precip
306      SELECT CASE( TRIM( cn_rcv_emp ) )
307      CASE( 'oce only'      )   ;   srcv(                                 jpr_oemp   )%laction = .TRUE. 
308      CASE( 'conservative'  )   ;   srcv( (/jpr_rain, jpr_snow, jpr_ievp, jpr_tevp/) )%laction = .TRUE.
309      CASE( 'oce and ice'   )   ;   srcv( (/jpr_ievp, jpr_sbpr, jpr_semp, jpr_oemp/) )%laction = .TRUE.
310      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_rcv_emp' )
311      END SELECT
312
313      !                                                      ! ------------------------- !
314      !                                                      !     Runoffs & Calving     !   
315      !                                                      ! ------------------------- !
316      srcv(jpr_rnf   )%clname = 'O_Runoff'   ;   IF( TRIM( cn_rcv_rnf ) == 'coupled' )   srcv(jpr_rnf)%laction = .TRUE.
317                                                 IF( TRIM( cn_rcv_rnf ) == 'climato' )   THEN   ;   ln_rnf = .TRUE.
318                                                 ELSE                                           ;   ln_rnf = .FALSE.
319                                                 ENDIF
320      srcv(jpr_cal   )%clname = 'OCalving'   ;   IF( TRIM( cn_rcv_cal ) == 'coupled' )   srcv(jpr_cal)%laction = .TRUE.
321
322      !                                                      ! ------------------------- !
323      !                                                      !    non solar radiation    !   Qns
324      !                                                      ! ------------------------- !
325      srcv(jpr_qnsoce)%clname = 'O_QnsOce'
326      srcv(jpr_qnsice)%clname = 'O_QnsIce'
327      srcv(jpr_qnsmix)%clname = 'O_QnsMix'
328      SELECT CASE( TRIM( cn_rcv_qns ) )
329      CASE( 'oce only'      )   ;   srcv(               jpr_qnsoce   )%laction = .TRUE.
330      CASE( 'conservative'  )   ;   srcv( (/jpr_qnsice, jpr_qnsmix/) )%laction = .TRUE.
331      CASE( 'oce and ice'   )   ;   srcv( (/jpr_qnsice, jpr_qnsoce/) )%laction = .TRUE.
332      CASE( 'mixed oce-ice' )   ;   srcv(               jpr_qnsmix   )%laction = .TRUE. 
333      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_rcv_qns' )
334      END SELECT
335
336      !                                                      ! ------------------------- !
337      !                                                      !    solar radiation        !   Qsr
338      !                                                      ! ------------------------- !
339      srcv(jpr_qsroce)%clname = 'O_QsrOce'
340      srcv(jpr_qsrice)%clname = 'O_QsrIce'
341      srcv(jpr_qsrmix)%clname = 'O_QsrMix'
342      SELECT CASE( TRIM( cn_rcv_qsr ) )
343      CASE( 'oce only'      )   ;   srcv(               jpr_qsroce   )%laction = .TRUE.
344      CASE( 'conservative'  )   ;   srcv( (/jpr_qsrice, jpr_qsrmix/) )%laction = .TRUE.
345      CASE( 'oce and ice'   )   ;   srcv( (/jpr_qsrice, jpr_qsroce/) )%laction = .TRUE.
346      CASE( 'mixed oce-ice' )   ;   srcv(               jpr_qsrmix   )%laction = .TRUE. 
347      CASE default              ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_rcv_qsr' )
348      END SELECT
349
350      !                                                      ! ------------------------- !
351      !                                                      !   non solar sensitivity   !   d(Qns)/d(T)
352      !                                                      ! ------------------------- !
353      srcv(jpr_dqnsdt)%clname = 'O_dQnsdT'   
354      IF( TRIM( cn_rcv_dqnsdt ) == 'coupled' )   srcv(jpr_dqnsdt)%laction = .TRUE.
355      !
356      ! non solar sensitivity mandatory for ice model
357      IF( TRIM( cn_rcv_dqnsdt ) == 'none' .AND. k_ice /= 0 ) &
358         CALL ctl_stop( 'sbc_cpl_init: cn_rcv_dqnsdt must be coupled in namsbc_cpl namelist' )
359      ! non solar sensitivity mandatory for mixed oce-ice solar radiation coupling technique
360      IF( TRIM( cn_rcv_dqnsdt ) == 'none' .AND. TRIM( cn_rcv_qns ) == 'mixed oce-ice' ) &
361         CALL ctl_stop( 'sbc_cpl_init: namsbc_cpl namelist mismatch between cn_rcv_qns and cn_rcv_dqnsdt' )
362      !                                                      ! ------------------------- !
363      !                                                      !    Ice Qsr penetration    !   
364      !                                                      ! ------------------------- !
365      ! fraction of net shortwave radiation which is not absorbed in the thin surface layer
366      ! and penetrates inside the ice cover ( Maykut and Untersteiner, 1971 ; Elbert anbd Curry, 1993 )
367      ! Coupled case: since cloud cover is not received from atmosphere
368      !               ===> defined as constant value -> definition done in sbc_cpl_init
369      fr1_i0(:,:) = 0.18
370      fr2_i0(:,:) = 0.82
371      !                                                      ! ------------------------- !
372      !                                                      !      10m wind module      !   
373      !                                                      ! ------------------------- !
374      srcv(jpr_w10m  )%clname = 'O_Wind10'   ;   IF( TRIM(cn_rcv_w10m) == 'coupled' )   srcv(jpr_w10m)%laction = .TRUE.     
375     
376      ! ================================ !
377      !     Define the send interface    !
378      ! ================================ !
379      ! for each field: define the OASIS name                           (srcv(:)%clname)
380      !                 define send or not from the namelist parameters (srcv(:)%laction)
381      !                 define the north fold type of lbc               (srcv(:)%nsgn)
382     
383      ! default definitions of nsnd
384      ssnd(:)%laction = .FALSE.   ;   ssnd(:)%clgrid = 'T'   ;   ssnd(:)%nsgn = 1
385         
386      !                                                      ! ------------------------- !
387      !                                                      !    Surface temperature    !
388      !                                                      ! ------------------------- !
389      ssnd(jps_toce)%clname = 'O_SSTSST'
390      ssnd(jps_tice)%clname = 'O_TepIce'
391      ssnd(jps_tmix)%clname = 'O_TepMix'
392      SELECT CASE( TRIM( cn_snd_temperature ) )
393      CASE( 'oce only'             )   ;   ssnd(   jps_toce             )%laction = .TRUE.
394      CASE( 'weighted oce and ice' )   ;   ssnd( (/jps_toce, jps_tice/) )%laction = .TRUE.
395      CASE( 'mixed oce-ice'        )   ;   ssnd(   jps_tmix             )%laction = .TRUE.
396      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_snd_temperature' )
397      END SELECT
398     
399      !                                                      ! ------------------------- !
400      !                                                      !          Albedo           !
401      !                                                      ! ------------------------- !
402      ssnd(jps_albice)%clname = 'O_AlbIce' 
403      ssnd(jps_albmix)%clname = 'O_AlbMix'
404      SELECT CASE( TRIM( cn_snd_albedo ) )
405      CASE( 'none'          )       ! nothing to do
406      CASE( 'weighted ice'  )   ;   ssnd(jps_albice)%laction = .TRUE.
407      CASE( 'mixed oce-ice' )   ;   ssnd(jps_albmix)%laction = .TRUE.
408      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_snd_albedo' )
409      END SELECT
410      !
411      ! Need to calculate oceanic albedo if
412      !     1. sending mixed oce-ice albedo or
413      !     2. receiving mixed oce-ice solar radiation
414      IF ( TRIM ( cn_snd_albedo ) == 'mixed oce-ice' .OR. TRIM ( cn_rcv_qsr ) == 'mixed oce-ice' ) THEN
415         CALL albedo_oce( zaos, zacs )
416         ! Due to lack of information on nebulosity : mean clear/overcast sky
417         albedo_oce_mix(:,:) = ( zacs(:,:) + zaos(:,:) ) * 0.5
418      ENDIF
419
420      !                                                      ! ------------------------- !
421      !                                                      !  Ice fraction & Thickness !
422      !                                                      ! ------------------------- !
423      ssnd(jps_fice)%clname = 'OIceFrac'   
424      ssnd(jps_hice)%clname = 'O_IceTck'
425      ssnd(jps_hsnw)%clname = 'O_SnwTck'
426      IF( k_ice /= 0 )   ssnd(jps_fice)%laction = .TRUE.       ! if ice treated in the ocean (even in climato case)
427      IF( TRIM( cn_snd_thickness ) == 'weighted ice and snow' )   ssnd( (/jps_hice, jps_hsnw/) )%laction = .TRUE.
428         
429      !                                                      ! ------------------------- !
430      !                                                      !      Surface current      !
431      !                                                      ! ------------------------- !
432      !        ocean currents              !            ice velocities
433      ssnd(jps_ocx1)%clname = 'O_OCurx1'   ;   ssnd(jps_ivx1)%clname = 'O_IVelx1'
434      ssnd(jps_ocy1)%clname = 'O_OCury1'   ;   ssnd(jps_ivy1)%clname = 'O_IVely1'
435      ssnd(jps_ocz1)%clname = 'O_OCurz1'   ;   ssnd(jps_ivz1)%clname = 'O_IVelz1'
436      !
437      ssnd(jps_ocx1:jps_ivz1)%nsgn = -1    ! vectors: change of the sign at the north fold
438
439      IF( cn_snd_crt(4) /= 'T' )   CALL ctl_stop( 'cn_snd_crt(4) must be equal to T' )
440      ssnd(jps_ocx1:jps_ivz1)%clgrid  = 'T'      ! all oce and ice components on the same unique grid
441
442      ssnd(jps_ocx1:jps_ivz1)%laction = .TRUE.   ! default: all are send
443      IF( TRIM( cn_snd_crt(2) ) == 'spherical' )   ssnd( (/jps_ocz1, jps_ivz1/) )%laction = .FALSE. 
444      SELECT CASE( TRIM( cn_snd_crt(1) ) )
445      CASE( 'none'                 )   ;   ssnd(jps_ocx1:jps_ivz1)%laction = .FALSE.
446      CASE( 'oce only'             )   ;   ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE.
447      CASE( 'weighted oce and ice' )   !   nothing to do
448      CASE( 'mixed oce-ice'        )   ;   ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE.
449      CASE default   ;   CALL ctl_stop( 'sbc_cpl_init: wrong definition of cn_snd_crt(1)' )
450      END SELECT
451
452      ! ================================ !
453      !   initialisation of the coupler  !
454      ! ================================ !
455
456      CALL cpl_prism_define(jprcv, jpsnd)           
457      !
458   END SUBROUTINE sbc_cpl_init
459
460
461   SUBROUTINE sbc_cpl_rcv( kt, k_fsbc, k_ice )     
462      !!----------------------------------------------------------------------
463      !!             ***  ROUTINE sbc_cpl_rcv  ***
464      !!
465      !! ** Purpose :   provide the stress over the ocean and, if no sea-ice,
466      !!                provide the ocean heat and freshwater fluxes.
467      !!
468      !! ** Method  : - Receive all the atmospheric fields (stored in frcv array). called at each time step.
469      !!                OASIS controls if there is something do receive or not. nrcvinfo contains the info
470      !!                to know if the field was really received or not
471      !!
472      !!              --> If ocean stress was really received:
473      !!
474      !!                  - transform the received ocean stress vector from the received
475      !!                 referential and grid into an atmosphere-ocean stress in
476      !!                 the (i,j) ocean referencial and at the ocean velocity point.
477      !!                    The received stress are :
478      !!                     - defined by 3 components (if cartesian coordinate)
479      !!                            or by 2 components (if spherical)
480      !!                     - oriented along geographical   coordinate (if eastward-northward)
481      !!                            or  along the local grid coordinate (if local grid)
482      !!                     - given at U- and V-point, resp.   if received on 2 grids
483      !!                            or at T-point               if received on 1 grid
484      !!                    Therefore and if necessary, they are successively
485      !!                  processed in order to obtain them
486      !!                     first  as  2 components on the sphere
487      !!                     second as  2 components oriented along the local grid
488      !!                     third  as  2 components on the U,V grid
489      !!
490      !!              -->
491      !!
492      !!              - In 'ocean only' case, non solar and solar ocean heat fluxes
493      !!             and total ocean freshwater fluxes 
494      !!
495      !! ** Method  :   receive all fields from the atmosphere and transform
496      !!              them into ocean surface boundary condition fields
497      !!
498      !! ** Action  :   update  utau, vtau   ocean stress at U,V grid
499      !!                        qns , qsr    non solar and solar ocean heat fluxes   ('ocean only case)
500      !!                        emp = emps   evap. - precip. (- runoffs) (- calving) ('ocean only case)
501      !!                        wndm         10m wind speed  !!!!gm  to be checked
502      !!----------------------------------------------------------------------
503      INTEGER, INTENT(in) ::   kt       ! ocean model time step index
504      INTEGER, INTENT(in) ::   k_fsbc   ! frequency of sbc (-> ice model) computation
505      INTEGER, INTENT(in) ::   k_ice    ! ice management in the sbc (=0/1/2/3)
506      !!
507      INTEGER  ::   ji, jj, jn             ! dummy loop indices
508      INTEGER  ::   isec                   ! number of seconds since nit000 (assuming rdttra did not change since nit000)
509      REAL(wp) ::   zcumulneg, zcumulpos   ! temporary scalars     
510      REAL(wp) ::   zcoef                  ! temporary scalar
511      REAL(wp), DIMENSION(jpi,jpj) ::   ztx, zty   ! 2D workspace
512      !!----------------------------------------------------------------------
513
514      IF( kt == nit000 )   CALL sbc_cpl_init( k_ice )          ! initialisation
515
516      !                                                 ! Receive all the atmos. fields (including ice information)
517      isec = ( kt - nit000 ) * NINT( rdttra(1) )             ! date of exchanges
518      DO jn = 1, jprcv                                       ! received fields sent by the atmosphere
519         IF( srcv(jn)%laction )   CALL cpl_prism_rcv( jn, isec, frcv(:,:,jn), nrcvinfo(jn) )
520      END DO
521
522      !                                                      ! ========================= !
523      IF( srcv(jpr_otx1)%laction ) THEN                      !       ocean stress        !
524         !                                                   ! ========================= !
525         ! define frcv(:,:,jpr_otx1) and frcv(:,:,jpr_oty1): stress at U/V point along model grid
526         ! => need to be done only when we receive the field
527         IF(  nrcvinfo(jpr_otx1) == PRISM_Recvd   .OR. nrcvinfo(jpr_otx1) == PRISM_FromRest .OR.   &
528            & nrcvinfo(jpr_otx1) == PRISM_RecvOut .OR. nrcvinfo(jpr_otx1) == PRISM_FromRestOut ) THEN
529            !
530            IF( TRIM( cn_rcv_tau(2) ) == 'cartesian' ) THEN            ! 2 components on the sphere
531               !                                                       ! (cartesian to spherical -> 3 to 2 components)
532               !
533               CALL geo2oce( frcv(:,:,jpr_otx1), frcv(:,:,jpr_oty1), frcv(:,:,jpr_otz1),   &
534                  &          srcv(jpr_otx1)%clgrid, ztx, zty )
535               frcv(:,:,jpr_otx1) = ztx(:,:)   ! overwrite 1st comp. on the 1st grid
536               frcv(:,:,jpr_oty1) = zty(:,:)   ! overwrite 2nd comp. on the 1st grid
537               !
538               IF( srcv(jpr_otx2)%laction ) THEN
539                  CALL geo2oce( frcv(:,:,jpr_otx2), frcv(:,:,jpr_oty2), frcv(:,:,jpr_otz2),   &
540                     &          srcv(jpr_otx2)%clgrid, ztx, zty )
541                  frcv(:,:,jpr_otx2) = ztx(:,:)   ! overwrite 1st comp. on the 2nd grid
542                  frcv(:,:,jpr_oty2) = zty(:,:)   ! overwrite 2nd comp. on the 2nd grid
543               ENDIF
544               !
545            ENDIF
546            !
547            IF( TRIM( cn_rcv_tau(3) ) == 'eastward-northward' ) THEN   ! 2 components oriented along the local grid
548               !                                                       ! (geographical to local grid -> rotate the components)
549               CALL rot_rep( frcv(:,:,jpr_otx1), frcv(:,:,jpr_oty1), srcv(jpr_otx1)%clgrid, 'en->i', ztx )   
550               frcv(:,:,jpr_otx1) = ztx(:,:)      ! overwrite 1st component on the 1st grid
551               IF( srcv(jpr_otx2)%laction ) THEN
552                  CALL rot_rep( frcv(:,:,jpr_otx2), frcv(:,:,jpr_oty2), srcv(jpr_otx2)%clgrid, 'en->j', zty )   
553               ELSE
554                  CALL rot_rep( frcv(:,:,jpr_otx1), frcv(:,:,jpr_oty1), srcv(jpr_otx1)%clgrid, 'en->j', zty ) 
555               ENDIF
556               frcv(:,:,jpr_oty1) = zty(:,:)      ! overwrite 2nd component on the 2nd grid
557            ENDIF
558            !                             
559            IF( srcv(jpr_otx1)%clgrid == 'T' ) THEN
560               DO jj = 2, jpjm1                                          ! T ==> (U,V)
561                  DO ji = fs_2, fs_jpim1   ! vector opt.
562                     frcv(ji,jj,jpr_otx1) = 0.5 * ( frcv(ji+1,jj  ,jpr_otx1) + frcv(ji,jj,jpr_otx1) )
563                     frcv(ji,jj,jpr_oty1) = 0.5 * ( frcv(ji  ,jj+1,jpr_oty1) + frcv(ji,jj,jpr_oty1) )
564                  END DO
565               END DO
566               CALL lbc_lnk( frcv(:,:,jpr_otx1), 'U',  -1. )   ;   CALL lbc_lnk( frcv(:,:,jpr_oty1), 'V',  -1. )
567            ENDIF
568         ENDIF
569         !                                                   ! ========================= !
570      ELSE                                                   !   No dynamical coupling   !
571         !                                                   ! ========================= !
572         frcv(:,:,jpr_otx1) = 0.e0                               ! here simply set to zero
573         frcv(:,:,jpr_oty1) = 0.e0                               ! an external read in a file can be added instead
574         !
575      ENDIF
576     
577      ! u(v)tau will be modified by ice model -> need to be reset before each call of the ice/fsbc     
578      IF( MOD( kt-1, k_fsbc ) == 0 ) THEN
579         utau(:,:) = frcv(:,:,jpr_otx1)                   
580         vtau(:,:) = frcv(:,:,jpr_oty1)
581         IF( .NOT. srcv(jpr_w10m)%laction )   CALL sbc_tau2wnd
582      ENDIF
583      !                                                      ! ========================= !
584      IF( k_ice <= 1 ) THEN                                  !  heat & freshwater fluxes ! (Ocean only case)
585         !                                                   ! ========================= !
586         !
587         !                                                       ! non solar heat flux over the ocean (qns)
588         IF( srcv(jpr_qnsoce)%laction )   qns(:,:) = frcv(:,:,jpr_qnsoce)
589         IF( srcv(jpr_qnsmix)%laction )   qns(:,:) = frcv(:,:,jpr_qnsmix)       
590         !   energy for melting solid precipitation over free ocean
591         zcoef = xlsn / rhosn
592         qns(:,:) = qns(:,:) - frcv(:,:,jpr_snow) * zcoef
593         !                                                       ! solar flux over the ocean          (qsr)
594         IF( srcv(jpr_qsroce)%laction )   qsr(:,:) = frcv(:,:,jpr_qsroce) 
595         IF( srcv(jpr_qsrmix)%laction )   qsr(:,:) = frcv(:,:,jpr_qsrmix)
596         !
597         !                                                       ! total freshwater fluxes over the ocean (emp, emps)
598         SELECT CASE( TRIM( cn_rcv_emp ) )                                    ! evaporation - precipitation
599         CASE( 'conservative' )
600            emp(:,:) = frcv(:,:,jpr_tevp) - ( frcv(:,:,jpr_rain) + frcv(:,:,jpr_snow) )
601         CASE( 'oce only', 'oce and ice' )
602            emp(:,:) = frcv(:,:,jpr_oemp)
603         CASE default
604            CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of cn_rcv_emp' )
605         END SELECT
606         !
607         !                                                        ! runoffs and calving (added in emp)
608         IF( srcv(jpr_rnf)%laction )   emp(:,:) = emp(:,:) - frcv(:,:,jpr_rnf)       
609         IF( srcv(jpr_cal)%laction )   emp(:,:) = emp(:,:) - frcv(:,:,jpr_cal)
610         !
611!!gm :  this seems to be internal cooking, not sure to need that in a generic interface
612!!gm                                       at least should be optional...
613!!         IF( TRIM( cn_rcv_rnf ) == 'coupled' ) THEN     ! add to the total freshwater budget
614!!            ! remove negative runoff
615!!            zcumulpos = SUM( MAX( frcv(:,:,jpr_rnf), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )
616!!            zcumulneg = SUM( MIN( frcv(:,:,jpr_rnf), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )
617!!            IF( lk_mpp )   CALL mpp_sum( zcumulpos )   ! sum over the global domain
618!!            IF( lk_mpp )   CALL mpp_sum( zcumulneg )
619!!            IF( zcumulpos /= 0. ) THEN                 ! distribute negative runoff on positive runoff grid points
620!!               zcumulneg = 1.e0 + zcumulneg / zcumulpos
621!!               frcv(:,:,jpr_rnf) = MAX( frcv(:,:,jpr_rnf), 0.e0 ) * zcumulneg
622!!            ENDIF     
623!!            ! add runoff to e-p
624!!            emp(:,:) = emp(:,:) - frcv(:,:,jpr_rnf)
625!!         ENDIF
626!!gm  end of internal cooking
627         !
628         emps(:,:) = emp(:,:)                                        ! concentration/dilution = emp
629 
630         !                                                           ! 10 m wind speed
631         IF( srcv(jpr_w10m)%laction )   wndm(:,:) = frcv(:,:,jpr_w10m)
632         ! it not, we call sbc_tau2wnd in sbc_cpl_rcv (or later, after the ice???)
633         !
634      ENDIF
635      !
636   END SUBROUTINE sbc_cpl_rcv
637   
638
639   SUBROUTINE sbc_cpl_ice_tau( p_taui, p_tauj )     
640      !!----------------------------------------------------------------------
641      !!             ***  ROUTINE sbc_cpl_ice_tau  ***
642      !!
643      !! ** Purpose :   provide the stress over sea-ice in coupled mode
644      !!
645      !! ** Method  :   transform the received stress from the atmosphere into
646      !!             an atmosphere-ice stress in the (i,j) ocean referencial
647      !!             and at the velocity point of the sea-ice model (cigr_type):
648      !!                'C'-grid : i- (j-) components given at U- (V-) point
649      !!                'B'-grid : both components given at I-point
650      !!
651      !!                The received stress are :
652      !!                 - defined by 3 components (if cartesian coordinate)
653      !!                        or by 2 components (if spherical)
654      !!                 - oriented along geographical   coordinate (if eastward-northward)
655      !!                        or  along the local grid coordinate (if local grid)
656      !!                 - given at U- and V-point, resp.   if received on 2 grids
657      !!                        or at a same point (T or I) if received on 1 grid
658      !!                Therefore and if necessary, they are successively
659      !!             processed in order to obtain them
660      !!                 first  as  2 components on the sphere
661      !!                 second as  2 components oriented along the local grid
662      !!                 third  as  2 components on the cigr_type point
663      !!
664      !!                In 'oce and ice' case, only one vector stress field
665      !!             is received. It has already been processed in sbc_cpl_rcv
666      !!             so that it is now defined as (i,j) components given at U-
667      !!             and V-points, respectively. Therefore, here only the third
668      !!             transformation is done and only if the ice-grid is a 'B'-grid.
669      !!
670      !! ** Action  :   return ptau_i, ptau_j, the stress over the ice at cigr_type point
671      !!----------------------------------------------------------------------
672      REAL(wp), INTENT(out), DIMENSION(jpi,jpj) ::   p_taui   ! i- & j-components of atmos-ice stress [N/m2]
673      REAL(wp), INTENT(out), DIMENSION(jpi,jpj) ::   p_tauj   ! at I-point (B-grid) or U & V-point (C-grid)
674      !!
675      INTEGER ::   ji, jj                          ! dummy loop indices
676      INTEGER ::   itx                             ! index of taux over ice
677      REAL(wp), DIMENSION(jpi,jpj) ::   ztx, zty   ! 2D workspace
678      !!----------------------------------------------------------------------
679
680      IF( srcv(jpr_itx1)%laction ) THEN   ;   itx =  jpr_itx1   
681      ELSE                                ;   itx =  jpr_otx1
682      ENDIF
683
684      ! do something only if we just received the stress from atmosphere
685      IF(  nrcvinfo(itx) == PRISM_Recvd   .OR. nrcvinfo(itx) == PRISM_FromRest .OR.   &
686         & nrcvinfo(itx) == PRISM_RecvOut .OR. nrcvinfo(itx) == PRISM_FromRestOut ) THEN
687
688         !                                                      ! ======================= !
689         IF( srcv(jpr_itx1)%laction ) THEN                      !   ice stress received   !
690            !                                                   ! ======================= !
691           
692            IF( TRIM( cn_rcv_tau(2) ) == 'cartesian' ) THEN            ! 2 components on the sphere
693               !                                                       ! (cartesian to spherical -> 3 to 2 components)
694               CALL geo2oce( frcv(:,:,jpr_itx1), frcv(:,:,jpr_ity1), frcv(:,:,jpr_itz1),   &
695                  &          srcv(jpr_itx1)%clgrid, ztx, zty )
696               frcv(:,:,jpr_itx1) = ztx(:,:)   ! overwrite 1st comp. on the 1st grid
697               frcv(:,:,jpr_itx1) = zty(:,:)   ! overwrite 2nd comp. on the 1st grid
698               !
699               IF( srcv(jpr_itx2)%laction ) THEN
700                  CALL geo2oce( frcv(:,:,jpr_itx2), frcv(:,:,jpr_ity2), frcv(:,:,jpr_itz2),   &
701                     &          srcv(jpr_itx2)%clgrid, ztx, zty )
702                  frcv(:,:,jpr_itx2) = ztx(:,:)   ! overwrite 1st comp. on the 2nd grid
703                  frcv(:,:,jpr_ity2) = zty(:,:)   ! overwrite 2nd comp. on the 2nd grid
704               ENDIF
705               !
706            ENDIF
707            !
708            IF( TRIM( cn_rcv_tau(3) ) == 'eastward-northward' ) THEN   ! 2 components oriented along the local grid
709               !                                                       ! (geographical to local grid -> rotate the components)
710               CALL rot_rep( frcv(:,:,jpr_itx1), frcv(:,:,jpr_ity1), srcv(jpr_itx1)%clgrid, 'en->i', ztx )   
711               frcv(:,:,jpr_itx1) = ztx(:,:)      ! overwrite 1st component on the 1st grid
712               IF( srcv(jpr_itx2)%laction ) THEN
713                  CALL rot_rep( frcv(:,:,jpr_itx2), frcv(:,:,jpr_ity2), srcv(jpr_itx2)%clgrid, 'en->j', zty )   
714               ELSE
715                  CALL rot_rep( frcv(:,:,jpr_itx1), frcv(:,:,jpr_ity1), srcv(jpr_itx1)%clgrid, 'en->j', zty ) 
716               ENDIF
717               frcv(:,:,jpr_ity1) = zty(:,:)      ! overwrite 2nd component on the 1st grid
718            ENDIF
719            !                                                   ! ======================= !
720         ELSE                                                   !     use ocean stress    !
721            !                                                   ! ======================= !
722            frcv(:,:,jpr_itx1) = frcv(:,:,jpr_otx1)
723            frcv(:,:,jpr_ity1) = frcv(:,:,jpr_oty1)
724            !
725         ENDIF
726
727         !                                                      ! ======================= !
728         !                                                      !     put on ice grid     !
729         !                                                      ! ======================= !
730         !   
731         !                                                  j+1   j     -----V---F
732         ! ice stress on ice velocity point (cigr_type)                  !       |
733         ! (C-grid ==>(U,V) or B-grid ==> I or F)                 j      |   T   U
734         !                                                               |       |
735         !                                                   j    j-1   -I-------|
736         !                                               (for I)         |       |
737         !                                                              i-1  i   i
738         !                                                               i      i+1 (for I)
739         SELECT CASE ( cigr_type )
740            !
741         CASE( 'I' )                                         ! B-grid ==> I
742            SELECT CASE ( srcv(jpr_itx1)%clgrid )
743            CASE( 'U' )
744               DO jj = 2, jpjm1                                   ! (U,V) ==> I
745                  DO ji = fs_2, fs_jpim1   ! vector opt.
746                     p_taui(ji,jj) = 0.5 * ( frcv(ji-1,jj  ,jpr_itx1) + frcv(ji-1,jj-1,jpr_itx1) )
747                     p_tauj(ji,jj) = 0.5 * ( frcv(ji  ,jj-1,jpr_ity1) + frcv(ji-1,jj-1,jpr_ity1) )
748                  END DO
749               END DO
750            CASE( 'F' )
751               DO jj = 2, jpjm1                                   ! F ==> I
752                  DO ji = fs_2, fs_jpim1   ! vector opt.
753                     p_taui(ji,jj) = frcv(ji-1,jj-1,jpr_itx1) 
754                     p_tauj(ji,jj) = frcv(ji-1,jj-1,jpr_ity1) 
755                  END DO
756               END DO
757            CASE( 'T' )
758               DO jj = 2, jpjm1                                   ! T ==> I
759                  DO ji = fs_2, fs_jpim1   ! vector opt.
760                     p_taui(ji,jj) = 0.25 * ( frcv(ji,jj  ,jpr_itx1) + frcv(ji-1,jj  ,jpr_itx1)   &
761                        &                   + frcv(ji,jj-1,jpr_itx1) + frcv(ji-1,jj-1,jpr_itx1) ) 
762                     p_tauj(ji,jj) = 0.25 * ( frcv(ji,jj  ,jpr_ity1) + frcv(ji-1,jj  ,jpr_ity1)   &
763                        &                   + frcv(ji,jj-1,jpr_ity1) + frcv(ji-1,jj-1,jpr_ity1) )
764                  END DO
765               END DO
766            CASE( 'I' )
767               p_taui(:,:) = frcv(:,:,jpr_itx1)                   ! I ==> I
768               p_tauj(:,:) = frcv(:,:,jpr_ity1)
769            END SELECT
770            IF( srcv(jpr_itx1)%clgrid /= 'I' ) THEN
771               CALL lbc_lnk( p_taui, 'I',  -1. )   ;   CALL lbc_lnk( p_tauj, 'I',  -1. )
772            ENDIF
773            !
774         CASE( 'F' )                                         ! B-grid ==> F
775            SELECT CASE ( srcv(jpr_itx1)%clgrid )
776            CASE( 'U' )
777               DO jj = 2, jpjm1                                   ! (U,V) ==> F
778                  DO ji = fs_2, fs_jpim1   ! vector opt.
779                     p_taui(ji,jj) = 0.5 * ( frcv(ji,jj,jpr_itx1) + frcv(ji  ,jj+1,jpr_itx1) )
780                     p_tauj(ji,jj) = 0.5 * ( frcv(ji,jj,jpr_ity1) + frcv(ji+1,jj  ,jpr_ity1) )
781                  END DO
782               END DO
783            CASE( 'I' )
784               DO jj = 2, jpjm1                                   ! I ==> F
785                  DO ji = fs_2, fs_jpim1   ! vector opt.
786                     p_taui(ji,jj) = frcv(ji+1,jj+1,jpr_itx1) 
787                     p_tauj(ji,jj) = frcv(ji+1,jj+1,jpr_ity1) 
788                  END DO
789               END DO
790            CASE( 'T' )
791               DO jj = 2, jpjm1                                   ! T ==> F
792                  DO ji = fs_2, fs_jpim1   ! vector opt.
793                     p_taui(ji,jj) = 0.25 * ( frcv(ji,jj  ,jpr_itx1) + frcv(ji+1,jj  ,jpr_itx1)   &
794                        &                   + frcv(ji,jj+1,jpr_itx1) + frcv(ji+1,jj+1,jpr_itx1) ) 
795                     p_tauj(ji,jj) = 0.25 * ( frcv(ji,jj  ,jpr_ity1) + frcv(ji+1,jj  ,jpr_ity1)   &
796                        &                   + frcv(ji,jj+1,jpr_ity1) + frcv(ji+1,jj+1,jpr_ity1) )
797                  END DO
798               END DO
799            CASE( 'F' )
800               p_taui(:,:) = frcv(:,:,jpr_itx1)                   ! F ==> F
801               p_tauj(:,:) = frcv(:,:,jpr_ity1)
802            END SELECT
803            IF( srcv(jpr_itx1)%clgrid /= 'F' ) THEN
804               CALL lbc_lnk( p_taui, 'F',  -1. )   ;   CALL lbc_lnk( p_tauj, 'F',  -1. )
805            ENDIF
806            !
807         CASE( 'C' )                                         ! C-grid ==> U,V
808            SELECT CASE ( srcv(jpr_itx1)%clgrid )
809            CASE( 'U' )
810               p_taui(:,:) = frcv(:,:,jpr_itx1)                   ! (U,V) ==> (U,V)
811               p_tauj(:,:) = frcv(:,:,jpr_ity1)
812            CASE( 'F' )
813               DO jj = 2, jpjm1                                   ! F ==> (U,V)
814                  DO ji = fs_2, fs_jpim1   ! vector opt.
815                     p_taui(ji,jj) = 0.5 * ( frcv(ji,jj,jpr_itx1) + frcv(ji  ,jj-1,jpr_itx1) )
816                     p_tauj(ji,jj) = 0.5 * ( frcv(ji,jj,jpr_ity1) + frcv(ji-1,jj  ,jpr_ity1) )
817                  END DO
818               END DO
819            CASE( 'T' )
820               DO jj = 2, jpjm1                                   ! T ==> (U,V)
821                  DO ji = fs_2, fs_jpim1   ! vector opt.
822                     p_taui(ji,jj) = 0.5 * ( frcv(ji+1,jj  ,jpr_itx1) + frcv(ji,jj,jpr_itx1) )
823                     p_tauj(ji,jj) = 0.5 * ( frcv(ji  ,jj+1,jpr_ity1) + frcv(ji,jj,jpr_ity1) )
824                  END DO
825               END DO
826            CASE( 'I' )
827               DO jj = 2, jpjm1                                   ! I ==> (U,V)
828                  DO ji = fs_2, fs_jpim1   ! vector opt.
829                     p_taui(ji,jj) = 0.5 * ( frcv(ji+1,jj+1,jpr_itx1) + frcv(ji+1,jj  ,jpr_itx1) )
830                     p_tauj(ji,jj) = 0.5 * ( frcv(ji+1,jj+1,jpr_ity1) + frcv(ji  ,jj+1,jpr_ity1) )
831                  END DO
832               END DO
833            END SELECT
834            IF( srcv(jpr_itx1)%clgrid /= 'U' ) THEN
835               CALL lbc_lnk( p_taui, 'U',  -1. )   ;   CALL lbc_lnk( p_tauj, 'V',  -1. )
836            ENDIF
837         END SELECT
838
839         !!gm Should be useless as sbc_cpl_ice_tau only called at coupled frequency
840         ! The receive stress are transformed such that in all case frcv(:,:,jpr_itx1), frcv(:,:,jpr_ity1)
841         ! become the i-component and j-component of the stress at the right grid point
842         !!gm  frcv(:,:,jpr_itx1) = p_taui(:,:)
843         !!gm  frcv(:,:,jpr_ity1) = p_tauj(:,:)
844         !!gm
845      ENDIF
846      !   
847   END SUBROUTINE sbc_cpl_ice_tau
848   
849
850   SUBROUTINE sbc_cpl_ice_flx( p_frld  ,                                  &
851      &                        pqns_tot, pqns_ice, pqsr_tot , pqsr_ice,   &
852      &                        pemp_tot, pemp_ice, pdqns_ice, psprecip,   &
853      &                        palbi   , psst    , pist                 )
854      !!----------------------------------------------------------------------
855      !!             ***  ROUTINE sbc_cpl_ice_flx_rcv  ***
856      !!
857      !! ** Purpose :   provide the heat and freshwater fluxes of the
858      !!              ocean-ice system.
859      !!
860      !! ** Method  :   transform the fields received from the atmosphere into
861      !!             surface heat and fresh water boundary condition for the
862      !!             ice-ocean system. The following fields are provided:
863      !!              * total non solar, solar and freshwater fluxes (qns_tot,
864      !!             qsr_tot and emp_tot) (total means weighted ice-ocean flux)
865      !!             NB: emp_tot include runoffs and calving.
866      !!              * fluxes over ice (qns_ice, qsr_ice, emp_ice) where
867      !!             emp_ice = sublimation - solid precipitation as liquid
868      !!             precipitation are re-routed directly to the ocean and
869      !!             runoffs and calving directly enter the ocean.
870      !!              * solid precipitation (sprecip), used to add to qns_tot
871      !!             the heat lost associated to melting solid precipitation
872      !!             over the ocean fraction.
873      !!       ===>> CAUTION here this changes the net heat flux received from
874      !!             the atmosphere
875      !!              * 10m wind module (wndm)   
876      !!
877      !!             N.B. - fields over sea-ice are passed in argument so that
878      !!                 the module can be compile without sea-ice.
879      !!                  - the fluxes have been separated from the stress as
880      !!                 (a) they are updated at each ice time step compare to
881      !!                 an update at each coupled time step for the stress, and
882      !!                 (b) the conservative computation of the fluxes over the
883      !!                 sea-ice area requires the knowledge of the ice fraction
884      !!                 after the ice advection and before the ice thermodynamics,
885      !!                 so that the stress is updated before the ice dynamics
886      !!                 while the fluxes are updated after it.
887      !!
888      !! ** Action  :   update at each nf_ice time step:
889      !!                   pqns_tot, pqsr_tot  non-solar and solar total heat fluxes
890      !!                   pqns_ice, pqsr_ice  non-solar and solar heat fluxes over the ice
891      !!                   pemp_tot            total evaporation - precipitation(liquid and solid) (-runoff)(-calving)
892      !!                   pemp_ice            ice sublimation - solid precipitation over the ice
893      !!                   pdqns_ice           d(non-solar heat flux)/d(Temperature) over the ice
894      !!                   sprecip             solid precipitation over the ocean 
895      !!                   wndm                10m wind module
896      !!----------------------------------------------------------------------
897      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj,jpl) ::   p_frld     ! lead fraction                [0 to 1]
898      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj    ) ::   pqns_tot   ! total non solar heat flux    [W/m2]
899      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj,jpl) ::   pqns_ice   ! ice   non solar heat flux    [W/m2]
900      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj    ) ::   pqsr_tot   ! total     solar heat flux    [W/m2]
901      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj,jpl) ::   pqsr_ice   ! ice       solar heat flux    [W/m2]
902      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj    ) ::   pemp_tot   ! total     freshwater budget        [Kg/m2/s]
903      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj    ) ::   pemp_ice   ! solid freshwater budget over ice   [Kg/m2/s]
904      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj    ) ::   psprecip   ! Net solid precipitation (=emp_ice) [Kg/m2/s]
905      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj,jpl) ::   pdqns_ice  ! d(Q non solar)/d(Temperature) over ice
906      ! optional arguments, used only in 'mixed oce-ice' case
907      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj,jpl), OPTIONAL ::   palbi   ! ice albedo
908      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj    ), OPTIONAL ::   psst    ! sea surface temperature     [Celcius]
909      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj,jpl), OPTIONAL ::   pist    ! ice surface temperature     [Kelvin]
910     !!
911      INTEGER ::   ji, jj           ! dummy loop indices
912      INTEGER ::   isec, info       ! temporary integer
913      REAL(wp)::   zcoef, ztsurf    ! temporary scalar
914      REAL(wp), DIMENSION(jpi,jpj)::   zsnow    ! snow precipitation
915      !!----------------------------------------------------------------------
916      !
917      !                                                      ! ========================= !
918      !                                                      !    freshwater budget      !   (emp)
919      !                                                      ! ========================= !
920      !
921      !                                                           ! total Precipitations - total Evaporation (emp_tot)
922      !                                                           ! solid precipitation  - sublimation       (emp_ice)
923      !                                                           ! solid Precipitation                      (sprecip)
924      SELECT CASE( TRIM( cn_rcv_emp ) )
925      CASE( 'conservative'  )   ! received fields: jpr_rain, jpr_snow, jpr_ievp, jpr_tevp
926         pemp_tot(:,:) = frcv(:,:,jpr_tevp) - frcv(:,:,jpr_rain) - frcv(:,:,jpr_snow)
927         pemp_ice(:,:) = frcv(:,:,jpr_ievp) - frcv(:,:,jpr_snow)
928         zsnow   (:,:) = frcv(:,:,jpr_snow)
929      CASE( 'oce and ice'   )   ! received fields: jpr_sbpr, jpr_semp, jpr_oemp
930         pemp_tot(:,:) = p_frld(:,:,1) * frcv(:,:,jpr_oemp) + (1.- p_frld(:,:,1)) * frcv(:,:,jpr_sbpr) 
931         pemp_ice(:,:) = frcv(:,:,jpr_semp)
932         zsnow   (:,:) = - frcv(:,:,jpr_semp) + frcv(:,:,jpr_ievp)
933      END SELECT
934      psprecip(:,:) = - pemp_ice(:,:)
935      !   
936      !                                                           ! runoffs and calving (put in emp_tot)
937      IF( srcv(jpr_rnf)%laction )   pemp_tot(:,:) = pemp_tot(:,:) -      frcv(:,:,jpr_rnf)
938      IF( srcv(jpr_cal)%laction )   pemp_tot(:,:) = pemp_tot(:,:) - ABS( frcv(:,:,jpr_cal) )
939      !
940!!gm :  this seems to be internal cooking, not sure to need that in a generic interface
941!!gm                                       at least should be optional...
942!!       ! remove negative runoff                            ! sum over the global domain
943!!       zcumulpos = SUM( MAX( frcv(:,:,jpr_rnf), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )
944!!       zcumulneg = SUM( MIN( frcv(:,:,jpr_rnf), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )
945!!       IF( lk_mpp )   CALL mpp_sum( zcumulpos )
946!!       IF( lk_mpp )   CALL mpp_sum( zcumulneg )
947!!       IF( zcumulpos /= 0. ) THEN                          ! distribute negative runoff on positive runoff grid points
948!!          zcumulneg = 1.e0 + zcumulneg / zcumulpos
949!!          frcv(:,:,jpr_rnf) = MAX( frcv(:,:,jpr_rnf), 0.e0 ) * zcumulneg
950!!       ENDIF     
951!!       pemp_tot(:,:) = pemp_tot(:,:) - frcv(:,:,jpr_rnf)   ! add runoff to e-p
952!!
953!!gm  end of internal cooking
954
955
956      !                                                      ! ========================= !
957      SELECT CASE( TRIM( cn_rcv_qns ) )                      !   non solar heat fluxes   !   (qns)
958      !                                                      ! ========================= !
959      CASE( 'conservative' )                                      ! the required fields are directly provided
960         pqns_tot(:,:  ) = frcv(:,:,jpr_qnsmix)
961         pqns_ice(:,:,1) = frcv(:,:,jpr_qnsice)
962      CASE( 'oce and ice' )       ! the total flux is computed from ocean and ice fluxes
963         pqns_tot(:,:  ) =  p_frld(:,:,1) * frcv(:,:,jpr_qnsoce) + ( 1.- p_frld(:,:,1) ) * frcv(:,:,jpr_qnsice)
964         pqns_ice(:,:,1) =  frcv(:,:,jpr_qnsice)
965      CASE( 'mixed oce-ice' )     ! the ice flux is cumputed from the total flux, the SST and ice informations
966         pqns_tot(:,:  ) = frcv(:,:,jpr_qnsmix)
967         pqns_ice(:,:,1) = frcv(:,:,jpr_qnsmix)    &
968            &            + frcv(:,:,jpr_dqnsdt) * ( pist(:,:,1) - ( (rt0 + psst(:,:  ) ) *        p_frld(:,:,1)   &
969            &                                                   +          pist(:,:,1)   * ( 1. - p_frld(:,:,1) ) ) )
970      END SELECT
971      !                                                           ! snow melting heat flux ....
972      !   energy for melting solid precipitation over free ocean
973      zcoef = xlsn / rhosn
974      pqns_tot(:,:) = pqns_tot(:,:) - p_frld(:,:,1) * zsnow(:,:) * zcoef
975!!gm
976!!    currently it is taken into account in leads budget but not in the qns_tot, and thus not in
977!!    the flux that enter the ocean....
978!!    moreover 1 - it is not diagnose anywhere....
979!!             2 - it is unclear for me whether this heat lost is taken into account in the atmosphere or not...
980!!
981!! similar job should be done for snow and precipitation temperature
982
983      !                                                      ! ========================= !
984      SELECT CASE( TRIM( cn_rcv_qsr ) )                      !      solar heat fluxes    !   (qsr)
985      !                                                      ! ========================= !
986      CASE( 'conservative' )
987         pqsr_tot(:,:  ) = frcv(:,:,jpr_qsrmix)
988         pqsr_ice(:,:,1) = frcv(:,:,jpr_qsrice)
989      CASE( 'oce and ice' )
990         pqsr_tot(:,:  ) =  p_frld(:,:,1) * frcv(:,:,jpr_qsroce) + ( 1.- p_frld(:,:,1) ) * frcv(:,:,jpr_qsrice)
991         pqsr_ice(:,:,1) =  frcv(:,:,jpr_qsrice)
992      CASE( 'mixed oce-ice' )
993         pqsr_tot(:,:  ) = frcv(:,:,jpr_qsrmix)
994!       Create solar heat flux over ice using incoming solar heat flux and albedos
995!       ( see OASIS3 user guide, 5th edition, p39 )
996         pqsr_ice(:,:,1) = frcv(:,:,jpr_qsrmix) * ( 1.- palbi(:,:,1) )   &
997            &            / (  1.- ( albedo_oce_mix(:,:  ) * ( 1.- p_frld(:,:,1) )   &
998            &                     + palbi         (:,:,1) *       p_frld(:,:,1)   )  )
999      END SELECT
1000
1001
1002      SELECT CASE( TRIM( cn_rcv_dqnsdt ) )
1003      CASE ('coupled')
1004          pdqns_ice(:,:,1) = frcv(:,:,jpr_dqnsdt)
1005      END SELECT
1006
1007
1008      !                                                      ! ========================= !
1009      !                                                      !      10 m wind speed      !   (wndm)
1010      !                                                      ! ========================= !
1011      !
1012      IF( srcv(jpr_w10m  )%laction )   wndm(:,:) = frcv(:,:,jpr_w10m)
1013      !
1014   END SUBROUTINE sbc_cpl_ice_flx
1015   
1016   
1017   SUBROUTINE sbc_cpl_snd( kt )
1018      !!----------------------------------------------------------------------
1019      !!             ***  ROUTINE sbc_cpl_snd  ***
1020      !!
1021      !! ** Purpose :   provide the ocean-ice informations to the atmosphere
1022      !!
1023      !! ** Method  :   send to the atmosphere through a call to cpl_prism_snd
1024      !!              all the needed fields (as defined in sbc_cpl_init)
1025      !!----------------------------------------------------------------------
1026      INTEGER, INTENT(in) ::   kt
1027      !!
1028      INTEGER ::   ji, jj          ! dummy loop indices
1029      INTEGER ::   isec, info      ! temporary integer
1030      REAL(wp), DIMENSION(jpi,jpj) ::   zfr_l   ! 1. - fr_i(:,:)
1031      REAL(wp), DIMENSION(jpi,jpj) ::   ztmp1, ztmp2
1032      REAL(wp), DIMENSION(jpi,jpj) ::   zotx1 , zoty1 , zotz1, zitx1, zity1, zitz1
1033      !!----------------------------------------------------------------------
1034
1035      isec = ( kt - nit000 ) * NINT(rdttra(1))        ! date of exchanges
1036
1037      zfr_l(:,:) = 1.- fr_i(:,:)
1038
1039      !                                                      ! ------------------------- !
1040      !                                                      !    Surface temperature    !   in Kelvin
1041      !                                                      ! ------------------------- !
1042      SELECT CASE( cn_snd_temperature)
1043      CASE( 'oce only'             )   ;   ztmp1(:,:) =   tn(:,:,1) + rt0
1044      CASE( 'weighted oce and ice' )   ;   ztmp1(:,:) = ( tn(:,:,1) + rt0 ) * zfr_l(:,:)   
1045                                           ztmp2(:,:) =   tn_ice(:,:,1)     *  fr_i(:,:)
1046      CASE( 'mixed oce-ice'        )   ;   ztmp1(:,:) = ( tn(:,:,1) + rt0 ) * zfr_l(:,:) + tn_ice(:,:,1) * fr_i(:,:)
1047      CASE default                     ;   CALL ctl_stop( 'sbc_cpl_snd: wrong definition of cn_snd_temperature' )
1048      END SELECT
1049      IF( ssnd(jps_toce)%laction )   CALL cpl_prism_snd( jps_toce, isec, ztmp1, info )
1050      IF( ssnd(jps_tice)%laction )   CALL cpl_prism_snd( jps_tice, isec, ztmp2, info )
1051      IF( ssnd(jps_tmix)%laction )   CALL cpl_prism_snd( jps_tmix, isec, ztmp1, info )
1052      !
1053      !                                                      ! ------------------------- !
1054      !                                                      !           Albedo          !
1055      !                                                      ! ------------------------- !
1056      IF( ssnd(jps_albice)%laction ) THEN                         ! ice
1057         ztmp1(:,:) = alb_ice(:,:,1) * fr_i(:,:)
1058         CALL cpl_prism_snd( jps_albice, isec, ztmp1, info )
1059      ENDIF
1060      IF( ssnd(jps_albmix)%laction ) THEN                         ! mixed ice-ocean
1061         ztmp1(:,:) = albedo_oce_mix(:,:) * zfr_l(:,:) + alb_ice(:,:,1) * fr_i(:,:)
1062         CALL cpl_prism_snd( jps_albmix, isec, ztmp1, info )
1063      ENDIF
1064      !                                                      ! ------------------------- !
1065      !                                                      !  Ice fraction & Thickness !
1066      !                                                      ! ------------------------- !
1067      IF( ssnd(jps_fice)%laction )   CALL cpl_prism_snd( jps_fice, isec, fr_i                  , info )
1068      IF( ssnd(jps_hice)%laction )   CALL cpl_prism_snd( jps_hice, isec, hicif(:,:) * fr_i(:,:), info )
1069      IF( ssnd(jps_hsnw)%laction )   CALL cpl_prism_snd( jps_hsnw, isec, hsnif(:,:) * fr_i(:,:), info )
1070      !
1071      !                                                      ! ------------------------- !
1072      IF( ssnd(jps_ocx1)%laction ) THEN                      !      Surface current      !
1073         !                                                   ! ------------------------- !
1074         !   
1075         !                                                  j+1   j     -----V---F
1076         ! surfacz velocity always sent from T point                     !       |
1077         !                                                        j      |   T   U
1078         !                                                               |       |
1079         !                                                   j    j-1   -I-------|
1080         !                                               (for I)         |       |
1081         !                                                              i-1  i   i
1082         !                                                               i      i+1 (for I)
1083         SELECT CASE( TRIM( cn_snd_crt(1) ) )
1084         CASE( 'oce only'             )      ! C-grid ==> T
1085            DO jj = 2, jpjm1
1086               DO ji = fs_2, fs_jpim1   ! vector opt.
1087                  zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj  ,1) )
1088                  zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji  ,jj-1,1) ) 
1089               END DO
1090            END DO
1091         CASE( 'weighted oce and ice' )   
1092            SELECT CASE ( cigr_type )
1093            CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
1094               DO jj = 2, jpjm1
1095                  DO ji = fs_2, fs_jpim1   ! vector opt.
1096                     zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj) 
1097                     zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj)
1098                     zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj)
1099                     zity1(ji,jj) = 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj)
1100                  END DO
1101               END DO
1102            CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
1103               DO jj = 2, jpjm1
1104                  DO ji = fs_2, fs_jpim1   ! vector opt.
1105                     zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj-1,1) ) * zfr_l(ji,jj) 
1106                     zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj) 
1107                     zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     &
1108                        &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
1109                     zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     &
1110                        &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
1111                  END DO
1112               END DO
1113            CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
1114               DO jj = 2, jpjm1
1115                  DO ji = fs_2, fs_jpim1   ! vector opt.
1116                     zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj-1,1) ) * zfr_l(ji,jj) 
1117                     zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj) 
1118                     zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     &
1119                        &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
1120                     zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     &
1121                        &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
1122                  END DO
1123               END DO
1124            END SELECT
1125            CALL lbc_lnk( zitx1, 'T', -1. )   ;   CALL lbc_lnk( zity1, 'T', -1. )
1126         CASE( 'mixed oce-ice'        )
1127            SELECT CASE ( cigr_type )
1128            CASE( 'C' )                      ! Ocean and Ice on C-grid ==> T
1129               DO jj = 2, jpjm1
1130                  DO ji = fs_2, fs_jpim1   ! vector opt.
1131                     zotx1(ji,jj) = 0.5 * ( un   (ji,jj,1) + un   (ji-1,jj  ,1) ) * zfr_l(ji,jj)   &
1132                        &         + 0.5 * ( u_ice(ji,jj  ) + u_ice(ji-1,jj    ) ) *  fr_i(ji,jj)
1133                     zoty1(ji,jj) = 0.5 * ( vn   (ji,jj,1) + vn   (ji  ,jj-1,1) ) * zfr_l(ji,jj)   &
1134                        &         + 0.5 * ( v_ice(ji,jj  ) + v_ice(ji  ,jj-1  ) ) *  fr_i(ji,jj)
1135                  END DO
1136               END DO
1137            CASE( 'I' )                      ! Ocean on C grid, Ice on I-point (B-grid) ==> T
1138               DO jj = 2, jpjm1
1139                  DO ji = fs_2, fs_jpim1   ! vector opt.
1140                     zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj-1,1) ) * zfr_l(ji,jj)   &   
1141                        &         + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1)                     &
1142                        &                  + u_ice(ji+1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
1143                     zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
1144                        &         + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1)                     &
1145                        &                  + v_ice(ji+1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
1146                  END DO
1147               END DO
1148            CASE( 'F' )                      ! Ocean on C grid, Ice on F-point (B-grid) ==> T
1149               DO jj = 2, jpjm1
1150                  DO ji = fs_2, fs_jpim1   ! vector opt.
1151                     zotx1(ji,jj) = 0.5  * ( un(ji,jj,1)      + un(ji-1,jj-1,1) ) * zfr_l(ji,jj)   &   
1152                        &         + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1)                     &
1153                        &                  + u_ice(ji-1,jj  ) + u_ice(ji,jj  )  ) *  fr_i(ji,jj)
1154                     zoty1(ji,jj) = 0.5  * ( vn(ji,jj,1)      + vn(ji  ,jj-1,1) ) * zfr_l(ji,jj)   & 
1155                        &         + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1)                     &
1156                        &                  + v_ice(ji-1,jj  ) + v_ice(ji,jj  )  ) *  fr_i(ji,jj)
1157                  END DO
1158               END DO
1159            END SELECT
1160         END SELECT
1161         CALL lbc_lnk( zotx1, 'T', -1. )   ;   CALL lbc_lnk( zoty1, 'T', -1. )
1162         !
1163         !
1164         IF( TRIM( cn_snd_crt(3) ) == 'eastward-northward' ) THEN             ! Rotation of the components
1165            !                                                                     ! Ocean component
1166            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocx1)%clgrid, 'ij->e', ztmp1 )       ! 1st component
1167            CALL rot_rep( zotx1, zoty1, ssnd(jps_ocx1)%clgrid, 'ij->n', ztmp2 )       ! 2nd component
1168            zotx1(:,:) = ztmp1(:,:)                                                   ! overwrite the components
1169            zoty1(:,:) = ztmp2(:,:)
1170            IF( ssnd(jps_ivx1)%laction ) THEN                                     ! Ice component
1171               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->e', ztmp1 )    ! 1st component
1172               CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->n', ztmp2 )    ! 2nd component
1173               zitx1(:,:) = ztmp1(:,:)                                                ! overwrite the components
1174               zity1(:,:) = ztmp2(:,:)
1175            ENDIF
1176         ENDIF
1177         !
1178         ! spherical coordinates to cartesian -> 2 components to 3 components
1179         IF( TRIM( cn_snd_crt(2) ) == 'cartesian' ) THEN
1180            ztmp1(:,:) = zotx1(:,:)                     ! ocean currents
1181            ztmp2(:,:) = zoty1(:,:)
1182            CALL oce2geo ( ztmp1, ztmp2, 'T', zotx1, zoty1, zotz1 )
1183            !
1184            IF( ssnd(jps_ivx1)%laction ) THEN           ! ice velocities
1185               ztmp1(:,:) = zitx1(:,:)
1186               ztmp1(:,:) = zity1(:,:)
1187               CALL oce2geo ( ztmp1, ztmp2, 'T', zitx1, zity1, zitz1 )
1188            ENDIF
1189         ENDIF
1190         !
1191         IF( ssnd(jps_ocx1)%laction )   CALL cpl_prism_snd( jps_ocx1, isec, zotx1, info )   ! ocean x current 1st grid
1192         IF( ssnd(jps_ocy1)%laction )   CALL cpl_prism_snd( jps_ocy1, isec, zoty1, info )   ! ocean y current 1st grid
1193         IF( ssnd(jps_ocz1)%laction )   CALL cpl_prism_snd( jps_ocz1, isec, zotz1, info )   ! ocean z current 1st grid
1194         !
1195         IF( ssnd(jps_ivx1)%laction )   CALL cpl_prism_snd( jps_ivx1, isec, zitx1, info )   ! ice   x current 1st grid
1196         IF( ssnd(jps_ivy1)%laction )   CALL cpl_prism_snd( jps_ivy1, isec, zity1, info )   ! ice   y current 1st grid
1197         IF( ssnd(jps_ivz1)%laction )   CALL cpl_prism_snd( jps_ivz1, isec, zitz1, info )   ! ice   z current 1st grid
1198         !
1199      ENDIF
1200   !
1201   END SUBROUTINE sbc_cpl_snd
1202   
1203#else
1204   !!----------------------------------------------------------------------
1205   !!   Dummy module                                            NO coupling
1206   !!----------------------------------------------------------------------
1207   USE par_kind        ! kind definition
1208CONTAINS
1209   SUBROUTINE sbc_cpl_snd( kt )
1210      WRITE(*,*) 'sbc_cpl_snd: You should not have seen this print! error?', kt
1211   END SUBROUTINE sbc_cpl_snd
1212   !
1213   SUBROUTINE sbc_cpl_rcv( kt, k_fsbc, k_ice )     
1214      WRITE(*,*) 'sbc_cpl_snd: You should not have seen this print! error?', kt, k_fsbc, k_ice
1215   END SUBROUTINE sbc_cpl_rcv
1216   !
1217   SUBROUTINE sbc_cpl_ice_tau( p_taui, p_tauj )     
1218      REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_taui   ! i- & j-components of atmos-ice stress [N/m2]
1219      REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_tauj   ! at I-point (B-grid) or U & V-point (C-grid)
1220      p_taui(:,:) = 0.   ;   p_tauj(:,:) = 0. ! stupid definition to avoid warning message when compiling...
1221      WRITE(*,*) 'sbc_cpl_snd: You should not have seen this print! error?'
1222   END SUBROUTINE sbc_cpl_ice_tau
1223   !
1224   SUBROUTINE sbc_cpl_ice_flx( p_frld  ,                                  &
1225      &                        pqns_tot, pqns_ice, pqsr_tot , pqsr_ice,   &
1226      &                        pemp_tot, pemp_ice, pdqns_ice, psprecip,   &
1227      &                        palbi   , psst    , pist                )
1228      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   p_frld     ! lead fraction                [0 to 1]
1229      REAL(wp), INTENT(  out), DIMENSION(:,:  ) ::   pqns_tot   ! total non solar heat flux    [W/m2]
1230      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pqns_ice   ! ice   non solar heat flux    [W/m2]
1231      REAL(wp), INTENT(  out), DIMENSION(:,:  ) ::   pqsr_tot   ! total     solar heat flux    [W/m2]
1232      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pqsr_ice   ! ice       solar heat flux    [W/m2]
1233      REAL(wp), INTENT(  out), DIMENSION(:,:  ) ::   pemp_tot   ! total     freshwater budget  [Kg/m2/s]
1234      REAL(wp), INTENT(  out), DIMENSION(:,:  ) ::   pemp_ice   ! ice solid freshwater budget  [Kg/m2/s]
1235      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pdqns_ice  ! d(Q non solar)/d(Temperature) over ice
1236      REAL(wp), INTENT(  out), DIMENSION(:,:  ) ::   psprecip   ! solid precipitation          [Kg/m2/s]
1237      REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   palbi   ! ice albedo
1238      REAL(wp), INTENT(in   ), DIMENSION(:,:  ), OPTIONAL ::   psst    ! sea surface temperature      [Celcius]
1239      REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   pist    ! ice surface temperature      [Kelvin]
1240      WRITE(*,*) 'sbc_cpl_snd: You should not have seen this print! error?', p_frld(1,1,1)
1241      ! stupid definition to avoid warning message when compiling...
1242      pqns_tot(:,:) = 0. ; pqns_ice(:,:,:) = 0. ; pdqns_ice(:,:,:) = 0.
1243      pqsr_tot(:,:) = 0. ; pqsr_ice(:,:,:) = 0. 
1244      pemp_tot(:,:) = 0. ; pemp_ice(:,:)   = 0. ; psprecip(:,:) = 0.
1245   END SUBROUTINE sbc_cpl_ice_flx
1246   
1247#endif
1248
1249   !!======================================================================
1250END MODULE sbccpl
Note: See TracBrowser for help on using the repository browser.