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

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

coupling interface, bugfix regarding changeset:1248, see ticket:155

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