source: LMDZ6/branches/Ocean_skin/libf/phylmd/cpl_mod.F90 @ 3767

Last change on this file since 3767 was 3767, checked in by lguez, 4 years ago

Store delta_sal instead of s_int

Revision analoguous to revision [3744], for salinity. Store as a state
variable and send to the ocean delta_sal, the difference between
ocean-air interface salinity and bulk salinity, instead of s_int,
the interface salinity.

So replace dummy argument s_int of procedure cpl_send_ocean_fields
by dummy argument delta_sal. Replace dummy argument s_int of
procedures ocean_cpl_noice and surf_ocean by dummy argument
delta_sal. Replace variable s_int of module phys_state_var_mod
by variable delta_sal. Rename local variable ys_int of procedure
pbl_surface to ydelta_sal. Set variable delta_sal of module
phys_state_var_mod to 0 for an appearing ocean fraction and a
missing startup field. Replace variable o_s_int of module
phys_output_ctrlout_mod by variable o_delta_sal.

Rename variables cpl_s_int and cpl_s_int_2D of module cpl_mod to
cpl_delta_sal and cpl_delta_sal_2D. Rename variable ids_s_int of
module oasis to ids_delta_sal. Change infosend(ids_delta_sal)%name
to "CODELSSS".

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 64.7 KB
RevLine 
[782]1!
2MODULE cpl_mod
3!
4! This module excahanges and transforms all fields that should be recieved or sent to
5! coupler. The transformation of the fields are done from the grid 1D-array in phylmd
6! to the regular 2D grid accepted by the coupler. Cumulation of the fields for each
7! timestep is done in here.
8!
9! Each type of surface that recevie fields from the coupler have a subroutine named
10! cpl_receive_XXX_fields and each surface that have fields to be sent to the coupler
11! have a subroutine named cpl_send_XXX_fields.
12!
13!*************************************************************************************
14
15! Use statements
16!*************************************************************************************
[836]17  USE dimphy, ONLY : klon
[782]18  USE mod_phys_lmdz_para
19  USE ioipsl
20  USE iophy
21
22! The module oasis is always used. Without the cpp key CPP_COUPLE only the parameters
23! in the module are compiled and not the subroutines.
24  USE oasis
25  USE write_field_phy
[2344]26  USE time_phylmdz_mod, ONLY: day_step_phy
[782]27 
28! Global attributes
29!*************************************************************************************
30  IMPLICIT NONE
31  PRIVATE
32
[996]33  ! All subroutine are public except cpl_send_all
34  PUBLIC :: cpl_init, cpl_receive_frac, cpl_receive_ocean_fields, cpl_receive_seaice_fields, &
[782]35       cpl_send_ocean_fields, cpl_send_seaice_fields, cpl_send_land_fields, &
36       cpl_send_landice_fields, gath2cpl
37 
38
39! Declaration of module variables
40!*************************************************************************************
41! variable for coupling period
[1279]42  INTEGER, SAVE :: nexca
[782]43  !$OMP THREADPRIVATE(nexca)
44
45! variables for cumulating fields during a coupling periode :
46  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_sols, cpl_nsol, cpl_rain
47  !$OMP THREADPRIVATE(cpl_sols,cpl_nsol,cpl_rain)
48  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_snow, cpl_evap, cpl_tsol
49  !$OMP THREADPRIVATE(cpl_snow,cpl_evap,cpl_tsol)
[3627]50
[3767]51  REAL, ALLOCATABLE, SAVE:: cpl_delta_sst(:), cpl_delta_sal(:)
52  !$OMP THREADPRIVATE(cpl_delta_sst, cpl_delta_sal)
[3627]53 
[1279]54  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_fder, cpl_albe, cpl_taux, cpl_tauy
55  !$OMP THREADPRIVATE(cpl_fder,cpl_albe,cpl_taux,cpl_tauy)
[782]56  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_windsp
57  !$OMP THREADPRIVATE(cpl_windsp)
[2872]58  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_sens_rain, cpl_sens_snow
59  !$OMP THREADPRIVATE(cpl_sens_rain, cpl_sens_snow)
[1279]60  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_taumod
61  !$OMP THREADPRIVATE(cpl_taumod)
62  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_atm_co2
63  !$OMP THREADPRIVATE(cpl_atm_co2)
[782]64  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_rriv2D, cpl_rcoa2D, cpl_rlic2D
65  !$OMP THREADPRIVATE(cpl_rriv2D,cpl_rcoa2D,cpl_rlic2D)
66
67! variables read from coupler :
68  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sst     ! sea surface temperature
69  !$OMP THREADPRIVATE(read_sst)
70  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sit     ! sea ice temperature
71  !$OMP THREADPRIVATE(read_sit)
[3627]72
73  REAL, ALLOCATABLE, SAVE:: read_sss(:, :)
74  ! bulk salinity of the surface layer of the ocean, in ppt
75  !$OMP THREADPRIVATE(read_sss)
76
[782]77  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sic     ! sea ice fraction
78  !$OMP THREADPRIVATE(read_sic)
79  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_alb_sic ! albedo at sea ice
80  !$OMP THREADPRIVATE(read_alb_sic)
[1067]81  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_u0, read_v0 ! ocean surface current
82  !$OMP THREADPRIVATE(read_u0,read_v0)
[1279]83  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_co2     ! ocean co2 flux
84  !$OMP THREADPRIVATE(read_co2)
[782]85  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: unity
86  !$OMP THREADPRIVATE(unity)
87  INTEGER, SAVE                             :: nidct, nidcs
88  !$OMP THREADPRIVATE(nidct,nidcs)
89
90! variables to be sent to the coupler
91  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_sols2D, cpl_nsol2D, cpl_rain2D
92  !$OMP THREADPRIVATE(cpl_sols2D, cpl_nsol2D, cpl_rain2D)
93  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_snow2D, cpl_evap2D, cpl_tsol2D
94  !$OMP THREADPRIVATE(cpl_snow2D, cpl_evap2D, cpl_tsol2D)
[3627]95
[3767]96  REAL, ALLOCATABLE, SAVE:: cpl_delta_sst_2D(:,:), cpl_delta_sal_2D(:,:)
97  !$OMP THREADPRIVATE(cpl_delta_sst_2D, cpl_delta_sal_2D)
[3627]98
[782]99  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_fder2D, cpl_albe2D
100  !$OMP THREADPRIVATE(cpl_fder2D, cpl_albe2D)
101  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_taux2D, cpl_tauy2D
102  !$OMP THREADPRIVATE(cpl_taux2D, cpl_tauy2D)
[1279]103  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_taumod2D
104  !$OMP THREADPRIVATE(cpl_taumod2D)
[782]105  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_windsp2D
106  !$OMP THREADPRIVATE(cpl_windsp2D)
[2872]107  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE   :: cpl_sens_rain2D, cpl_sens_snow2D
108  !$OMP THREADPRIVATE(cpl_sens_rain2D, cpl_sens_snow2D)
[1279]109  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_atm_co22D
110  !$OMP THREADPRIVATE(cpl_atm_co22D)
[1001]111
[3605]112!!!!!!!!!! variable for calving
113  INTEGER, PARAMETER :: nb_zone_calving = 3
114  REAL,ALLOCATABLE, DIMENSION(:,:,:),SAVE :: area_calving
115  !$OMP THREADPRIVATE(area_calving)
116  REAL,ALLOCATABLE, DIMENSION(:,:),SAVE :: cell_area2D
117  !$OMP THREADPRIVATE(cell_area2D)
118  INTEGER, SAVE :: ind_calving(nb_zone_calving)
119  !$OMP THREADPRIVATE(ind_calving)
120
121  LOGICAL,SAVE :: cpl_old_calving
122  !$OMP THREADPRIVATE(cpl_old_calving)
123 
[782]124CONTAINS
125!
126!************************************************************************************
127!
128  SUBROUTINE cpl_init(dtime, rlon, rlat)
[1279]129    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl, fco2_ocn_day
[1454]130    USE surface_data
[1785]131    USE indice_sol_mod
[3605]132    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, grid1dTo2d_glo, klon_glo, grid_type, unstructured, regular_lonlat
[2344]133    USE time_phylmdz_mod, ONLY: annee_ref, day_ini, itau_phy, itaufin_phy
[2311]134    USE print_control_mod, ONLY: lunout
[3605]135    USE geometry_mod, ONLY : longitude_deg, latitude_deg, ind_cell_glo, cell_area
136    USE ioipsl_getin_p_mod, ONLY: getin_p
[3627]137    use config_ocean_skin_m, only: activate_ocean_skin
[782]138
139! Input arguments
140!*************************************************************************************
141    REAL, INTENT(IN)                  :: dtime
142    REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat
143
144! Local variables
145!*************************************************************************************
146    INTEGER                           :: error, sum_error, ig, i
147    INTEGER                           :: jf, nhoridct
148    INTEGER                           :: nhoridcs
149    INTEGER                           :: idtime
150    INTEGER                           :: idayref
151    INTEGER                           :: npas ! only for OASIS2
152    REAL                              :: zjulian
[2344]153    REAL, DIMENSION(nbp_lon,nbp_lat)  :: zx_lon, zx_lat
[782]154    CHARACTER(len = 20)               :: modname = 'cpl_init'
155    CHARACTER(len = 80)               :: abort_message
156    CHARACTER(len=80)                 :: clintocplnam, clfromcplnam
[3605]157    REAL, DIMENSION(klon_mpi)         :: rlon_mpi, rlat_mpi, cell_area_mpi
158    INTEGER, DIMENSION(klon_mpi)           :: ind_cell_glo_mpi
159    REAL, DIMENSION(nbp_lon,jj_nb)         :: lon2D, lat2D
160    INTEGER :: mask_calving(nbp_lon,jj_nb,nb_zone_calving)
161    REAL :: pos
[782]162
[3605]163!***************************************
164! Use old calving or not (default new calving method)
165! New calving method should be used with DYNAMICO and when using new coupling
166! weights.
167    cpl_old_calving=.FALSE.
168    CALL getin_p("cpl_old_calving",cpl_old_calving)
169
170
[782]171!*************************************************************************************
172! Calculate coupling period
173!
174!*************************************************************************************
175     
[2344]176    npas = itaufin_phy
[2075]177!    nexca = 86400 / dtime
178    nexca = t_coupl / dtime
[782]179    WRITE(lunout,*)' ##### Ocean couple #####'
180    WRITE(lunout,*)' Valeurs des pas de temps'
181    WRITE(lunout,*)' npas = ', npas
182    WRITE(lunout,*)' nexca = ', nexca
183   
184!*************************************************************************************
185! Allocate variables
186!
187!*************************************************************************************
188    error = 0
189    sum_error = 0
190
191    ALLOCATE(unity(klon), stat = error)
192    sum_error = sum_error + error
193    ALLOCATE(cpl_sols(klon,2), stat = error)
194    sum_error = sum_error + error
195    ALLOCATE(cpl_nsol(klon,2), stat = error)
196    sum_error = sum_error + error
197    ALLOCATE(cpl_rain(klon,2), stat = error)
198    sum_error = sum_error + error
199    ALLOCATE(cpl_snow(klon,2), stat = error)
200    sum_error = sum_error + error
201    ALLOCATE(cpl_evap(klon,2), stat = error)
202    sum_error = sum_error + error
203    ALLOCATE(cpl_tsol(klon,2), stat = error)
204    sum_error = sum_error + error
205    ALLOCATE(cpl_fder(klon,2), stat = error)
206    sum_error = sum_error + error
207    ALLOCATE(cpl_albe(klon,2), stat = error)
208    sum_error = sum_error + error
209    ALLOCATE(cpl_taux(klon,2), stat = error)
210    sum_error = sum_error + error
[1279]211    ALLOCATE(cpl_tauy(klon,2), stat = error)
212    sum_error = sum_error + error
[782]213    ALLOCATE(cpl_windsp(klon,2), stat = error)
214    sum_error = sum_error + error
[2545]215    ALLOCATE(cpl_taumod(klon,2), stat = error)
[782]216    sum_error = sum_error + error
[2872]217    ALLOCATE(cpl_sens_rain(klon,2), stat = error)
218    sum_error = sum_error + error
219    ALLOCATE(cpl_sens_snow(klon,2), stat = error)
220    sum_error = sum_error + error
[2344]221    ALLOCATE(cpl_rriv2D(nbp_lon,jj_nb), stat=error)
[782]222    sum_error = sum_error + error
[2344]223    ALLOCATE(cpl_rcoa2D(nbp_lon,jj_nb), stat=error)
[782]224    sum_error = sum_error + error
[2344]225    ALLOCATE(cpl_rlic2D(nbp_lon,jj_nb), stat=error)
[782]226    sum_error = sum_error + error
[2344]227    ALLOCATE(read_sst(nbp_lon, jj_nb), stat = error)
[782]228    sum_error = sum_error + error
[2344]229    ALLOCATE(read_sic(nbp_lon, jj_nb), stat = error)
[782]230    sum_error = sum_error + error
[2344]231    ALLOCATE(read_sit(nbp_lon, jj_nb), stat = error)
[782]232    sum_error = sum_error + error
[3627]233
234    if (activate_ocean_skin >= 1) then
235       ALLOCATE(read_sss(nbp_lon, jj_nb), stat = error)
236       sum_error = sum_error + error
237   
238       if (activate_ocean_skin == 2) then
[3767]239          ALLOCATE(cpl_delta_sst(klon), cpl_delta_sal(klon), stat = error)
[3627]240          sum_error = sum_error + error
241       end if
242    end if
243
[2344]244    ALLOCATE(read_alb_sic(nbp_lon, jj_nb), stat = error)
[782]245    sum_error = sum_error + error
[2344]246    ALLOCATE(read_u0(nbp_lon, jj_nb), stat = error)
[1067]247    sum_error = sum_error + error
[2344]248    ALLOCATE(read_v0(nbp_lon, jj_nb), stat = error)
[1067]249    sum_error = sum_error + error
250
[1279]251    IF (carbon_cycle_cpl) THEN
[2344]252       ALLOCATE(read_co2(nbp_lon, jj_nb), stat = error)
[1279]253       sum_error = sum_error + error
254       ALLOCATE(cpl_atm_co2(klon,2), stat = error)
255       sum_error = sum_error + error
256
257! Allocate variable in carbon_cycle_mod
[3605]258       IF (.NOT.ALLOCATED(fco2_ocn_day)) ALLOCATE(fco2_ocn_day(klon), stat = error)
[1279]259       sum_error = sum_error + error
[3605]260    ENDIF
[1279]261
[3605]262! calving initialization
263    ALLOCATE(area_calving(nbp_lon, jj_nb, nb_zone_calving), stat = error)
264    sum_error = sum_error + error
265    ALLOCATE(cell_area2D(nbp_lon, jj_nb), stat = error)   
266    sum_error = sum_error + error
267
268
269    CALL gather_omp(longitude_deg,rlon_mpi)
270    CALL gather_omp(latitude_deg,rlat_mpi)
271    CALL gather_omp(ind_cell_glo,ind_cell_glo_mpi)
272    CALL gather_omp(cell_area,cell_area_mpi)
273     
274    IF (is_omp_master) THEN
275      CALL Grid1DTo2D_mpi(rlon_mpi,lon2D)
276      CALL Grid1DTo2D_mpi(rlat_mpi,lat2D)
277      CALL Grid1DTo2D_mpi(cell_area_mpi,cell_area2D)
278      mask_calving(:,:,:) = 0
279      WHERE ( lat2D >= 40) mask_calving(:,:,1) = 1
280      WHERE ( lat2D < 40 .AND. lat2D > -50) mask_calving(:,:,2) = 1
281      WHERE ( lat2D <= -50) mask_calving(:,:,3) = 1
282   
283   
284      DO i=1,nb_zone_calving
285        area_calving(:,:,i)=mask_calving(:,:,i)*cell_area2D(:,:)
286        pos=1
287        IF (i>1) pos = 1 + ((nbp_lon*nbp_lat-1)*(i-1))/(nb_zone_calving-1)
288     
289        ind_calving(i)=0
290        IF (grid_type==unstructured) THEN
291
292          DO ig=1,klon_mpi
293            IF (ind_cell_glo_mpi(ig)==pos) ind_calving(i)=ig
294          ENDDO
295
296        ELSE IF (grid_type==regular_lonlat) THEN
297          IF ((ij_begin<=pos .AND. ij_end>=pos) .OR. (ij_begin<=pos .AND. is_south_pole_dyn )) THEN
298            ind_calving(i)=pos-(jj_begin-1)*nbp_lon
299          ENDIF
300        ENDIF
301     
302      ENDDO
303    ENDIF
304   
305           
[782]306    IF (sum_error /= 0) THEN
307       abort_message='Pb allocation variables couplees'
[2311]308       CALL abort_physic(modname,abort_message,1)
[782]309    ENDIF
310!*************************************************************************************
311! Initialize the allocated varaibles
312!
313!*************************************************************************************
314    DO ig = 1, klon
315       unity(ig) = ig
316    ENDDO
317
318!*************************************************************************************
319! Initialize coupling
320!
321!*************************************************************************************
322    idtime = INT(dtime)
323#ifdef CPP_COUPLE
324    CALL inicma
325#endif
326
327!*************************************************************************************
328! initialize NetCDF output
329!
330!*************************************************************************************
331    IF (is_sequential) THEN
332       idayref = day_ini
333       CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
[3605]334       CALL grid1dTo2d_glo(rlon,zx_lon)
[2344]335       DO i = 1, nbp_lon
[782]336          zx_lon(i,1) = rlon(i+1)
[2344]337          zx_lon(i,nbp_lat) = rlon(i+1)
[782]338       ENDDO
[3605]339       CALL grid1dTo2d_glo(rlat,zx_lat)
[782]340       clintocplnam="cpl_atm_tauflx"
[2344]341       CALL histbeg(clintocplnam,nbp_lon,zx_lon(:,1),nbp_lat,zx_lat(1,:),&
342            1,nbp_lon,1,nbp_lat, itau_phy,zjulian,dtime,nhoridct,nidct)
[782]343! no vertical axis
344       CALL histdef(nidct, 'tauxe','tauxe', &
[2344]345            "-",nbp_lon,nbp_lat, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
[782]346       CALL histdef(nidct, 'tauyn','tauyn', &
[2344]347            "-",nbp_lon,nbp_lat, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
[782]348       CALL histdef(nidct, 'tmp_lon','tmp_lon', &
[2344]349            "-",nbp_lon,nbp_lat, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
[782]350       CALL histdef(nidct, 'tmp_lat','tmp_lat', &
[2344]351            "-",nbp_lon,nbp_lat, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
[1279]352       DO jf=1,maxsend
353         IF (infosend(i)%action) THEN
354             CALL histdef(nidct, infosend(i)%name ,infosend(i)%name , &
[2344]355                "-",nbp_lon,nbp_lat,nhoridct,1,1,1,-99,32,"inst",dtime,dtime)
[1279]356         ENDIF
[3605]357       ENDDO
[782]358       CALL histend(nidct)
359       CALL histsync(nidct)
360       
361       clfromcplnam="cpl_atm_sst"
[2344]362       CALL histbeg(clfromcplnam,nbp_lon,zx_lon(:,1),nbp_lat,zx_lat(1,:),1,nbp_lon,1,nbp_lat, &
[782]363            0,zjulian,dtime,nhoridcs,nidcs)
364! no vertical axis
[1279]365       DO jf=1,maxrecv
366         IF (inforecv(i)%action) THEN
367             CALL histdef(nidcs,inforecv(i)%name ,inforecv(i)%name , &
[2344]368                "-",nbp_lon,nbp_lat,nhoridcs,1,1,1,-99,32,"inst",dtime,dtime)
[1279]369         ENDIF
[3605]370       ENDDO
[782]371       CALL histend(nidcs)
372       CALL histsync(nidcs)
373
374    ENDIF    ! is_sequential
375   
[1454]376
377!*************************************************************************************
378! compatibility test
379!
380!*************************************************************************************
381    IF (carbon_cycle_cpl .AND. version_ocean=='opa8') THEN
382       abort_message='carbon_cycle_cpl does not work with opa8'
[2311]383       CALL abort_physic(modname,abort_message,1)
[3605]384    ENDIF
[1454]385
[782]386  END SUBROUTINE cpl_init
387 
388!
389!*************************************************************************************
390!
[996]391 
392  SUBROUTINE cpl_receive_frac(itime, dtime, pctsrf, is_modified)
393! This subroutine receives from coupler for both ocean and seaice
[782]394! 4 fields : read_sst, read_sic, read_sit and read_alb_sic.
[996]395! The new sea-ice-land-landice fraction is returned. The others fields
396! are stored in this module.
397    USE surface_data
[2399]398    USE geometry_mod, ONLY : longitude_deg, latitude_deg
[1279]399    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl
[1785]400    USE indice_sol_mod
[2344]401    USE time_phylmdz_mod, ONLY: start_time, itau_phy
402    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[3627]403    use config_ocean_skin_m, only: activate_ocean_skin
[1785]404
[793]405    INCLUDE "YOMCST.h"
[782]406
[996]407! Arguments
[782]408!************************************************************************************
[996]409    INTEGER, INTENT(IN)                        :: itime
410    REAL, INTENT(IN)                           :: dtime
411    REAL, DIMENSION(klon,nbsrf), INTENT(INOUT) :: pctsrf
412    LOGICAL, INTENT(OUT)                       :: is_modified
[782]413
414! Local variables
415!************************************************************************************
[996]416    INTEGER                                 :: j, i, time_sec
[782]417    INTEGER                                 :: itau_w
[2344]418    INTEGER, DIMENSION(nbp_lon*nbp_lat)     :: ndexcs
[996]419    CHARACTER(len = 20)                     :: modname = 'cpl_receive_frac'
[782]420    CHARACTER(len = 80)                     :: abort_message
421    REAL, DIMENSION(klon)                   :: read_sic1D
[2344]422    REAL, DIMENSION(nbp_lon,jj_nb,maxrecv)      :: tab_read_flds
[996]423    REAL, DIMENSION(klon,nbsrf)             :: pctsrf_old
[1067]424    REAL, DIMENSION(klon_mpi)               :: rlon_mpi, rlat_mpi
[2344]425    REAL, DIMENSION(nbp_lon,jj_nb)             :: tmp_lon, tmp_lat
426    REAL, DIMENSION(nbp_lon,jj_nb)             :: tmp_r0
[782]427
428!*************************************************************************************
429! Start calculation
430! Get fields from coupler
431!
432!*************************************************************************************
433
[996]434    is_modified=.FALSE.
435
[1279]436! Check if right moment to receive from coupler
[996]437    IF (MOD(itime, nexca) == 1) THEN
438       is_modified=.TRUE.
439 
440       time_sec=(itime-1)*dtime
[782]441#ifdef CPP_COUPLE
[987]442!$OMP MASTER
[1010]443    CALL fromcpl(time_sec, tab_read_flds)
[987]444!$OMP END MASTER
[782]445#endif
446   
447! NetCDF output of received fields
[996]448       IF (is_sequential) THEN
449          ndexcs(:) = 0
[2344]450          itau_w = itau_phy + itime + start_time * day_step_phy
[1279]451          DO i = 1, maxrecv
452            IF (inforecv(i)%action) THEN
[2344]453                CALL histwrite(nidcs,inforecv(i)%name,itau_w,tab_read_flds(:,:,i),nbp_lon*(nbp_lat),ndexcs)
[1279]454            ENDIF
[3605]455          ENDDO
[996]456       ENDIF
[782]457
[1001]458
[996]459! Save each field in a 2D array.
[987]460!$OMP MASTER
[1279]461       read_sst(:,:)     = tab_read_flds(:,:,idr_sisutw)  ! Sea surface temperature
462       read_sic(:,:)     = tab_read_flds(:,:,idr_icecov)  ! Sea ice concentration
463       read_alb_sic(:,:) = tab_read_flds(:,:,idr_icealw)  ! Albedo at sea ice
464       read_sit(:,:)     = tab_read_flds(:,:,idr_icetem)  ! Sea ice temperature
[3627]465       if (activate_ocean_skin >= 1) read_sss(:,:) = tab_read_flds(:,:,idr_sss)
[996]466!$OMP END MASTER
[987]467
[1067]468       IF (cpl_current) THEN
469
470! Transform the longitudes and latitudes on 2D arrays
[2399]471          CALL gather_omp(longitude_deg,rlon_mpi)
472          CALL gather_omp(latitude_deg,rlat_mpi)
[1067]473!$OMP MASTER
474          CALL Grid1DTo2D_mpi(rlon_mpi,tmp_lon)
475          CALL Grid1DTo2D_mpi(rlat_mpi,tmp_lat)
476
477! Transform the currents from cartesian to spheric coordinates
478! tmp_r0 should be zero
[2344]479          CALL geo2atm(nbp_lon, jj_nb, tab_read_flds(:,:,idr_curenx), &
[1279]480             tab_read_flds(:,:,idr_cureny), tab_read_flds(:,:,idr_curenz), &
[1067]481               tmp_lon, tmp_lat, &
482               read_u0(:,:), read_v0(:,:), tmp_r0(:,:))
483!$OMP END MASTER
[1146]484
[1279]485      ELSE
[1067]486          read_u0(:,:) = 0.
487          read_v0(:,:) = 0.
[1279]488      ENDIF
489
490       IF (carbon_cycle_cpl) THEN
491!$OMP MASTER
492           read_co2(:,:) = tab_read_flds(:,:,idr_oceco2) ! CO2 flux
493!$OMP END MASTER
[1067]494       ENDIF
495
[782]496!*************************************************************************************
[996]497!  Transform seaice fraction (read_sic : ocean-seaice mask) into global
498!  fraction (pctsrf : ocean-seaice-land-landice mask)
[782]499!
500!*************************************************************************************
[996]501       CALL cpl2gath(read_sic, read_sic1D, klon, unity)
502
503       pctsrf_old(:,:) = pctsrf(:,:)
504       DO i = 1, klon
505          ! treatment only of points with ocean and/or seaice
[1279]506          ! old land-ocean mask can not be changed
[996]507          IF (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic) > 0.) THEN
508             pctsrf(i,is_sic) = (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic)) &
509                  * read_sic1D(i)
510             pctsrf(i,is_oce) = (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic)) &
511                  - pctsrf(i,is_sic)
[782]512          ENDIF
513       ENDDO
[987]514
[3605]515    ENDIF ! if time to receive
[782]516
[996]517  END SUBROUTINE cpl_receive_frac
[782]518
519!
520!*************************************************************************************
521!
[996]522
[3627]523  SUBROUTINE cpl_receive_ocean_fields(knon, knindex, tsurf_new, u0_new, &
524       v0_new, sss)
[782]525!
[996]526! This routine returns the field for the ocean that has been read from the coupler
527! (done earlier with cpl_receive_frac). The field is the temperature.
528! The temperature is transformed into 1D array with valid points from index 1 to knon.
[782]529!
[1279]530    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl, fco2_ocn_day
[1785]531    USE indice_sol_mod
[3627]532    use config_ocean_skin_m, only: activate_ocean_skin
[782]533
534! Input arguments
535!*************************************************************************************
536    INTEGER, INTENT(IN)                     :: knon
537    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
538
539! Output arguments
540!*************************************************************************************
541    REAL, DIMENSION(klon), INTENT(OUT)      :: tsurf_new
[3627]542
543    REAL, INTENT(OUT):: sss(:) ! (klon)
544    ! bulk salinity of the surface layer of the ocean, in ppt
545
[1067]546    REAL, DIMENSION(klon), INTENT(OUT)      :: u0_new
547    REAL, DIMENSION(klon), INTENT(OUT)      :: v0_new
[782]548
[996]549! Local variables
[782]550!*************************************************************************************
[1279]551    INTEGER                  :: i
552    INTEGER, DIMENSION(klon) :: index
553    REAL, DIMENSION(klon)    :: sic_new
[782]554
555!*************************************************************************************
556! Transform read_sst into compressed 1D variable tsurf_new
557!
558!*************************************************************************************
559    CALL cpl2gath(read_sst, tsurf_new, knon, knindex)
[3627]560    if (activate_ocean_skin >= 1) CALL cpl2gath(read_sss, sss, knon, knindex)
[996]561    CALL cpl2gath(read_sic, sic_new, knon, knindex)
[1067]562    CALL cpl2gath(read_u0, u0_new, knon, knindex)
563    CALL cpl2gath(read_v0, v0_new, knon, knindex)
[782]564
[996]565!*************************************************************************************
[1279]566! Transform read_co2 into uncompressed 1D variable fco2_ocn_day added directly in
567! the module carbon_cycle_mod
568!
569!*************************************************************************************
570    IF (carbon_cycle_cpl) THEN
571       DO i=1,klon
572          index(i)=i
[3605]573       ENDDO
[1279]574       CALL cpl2gath(read_co2, fco2_ocn_day, klon, index)
[3605]575    ENDIF
[1279]576
577!*************************************************************************************
[996]578! The fields received from the coupler have to be weighted with the fraction of ocean
579! in relation to the total sea-ice+ocean
580!
581!*************************************************************************************
582    DO i=1, knon
583       tsurf_new(i) = tsurf_new(i)/(1. - sic_new(i))
[3605]584    ENDDO
[782]585
586  END SUBROUTINE cpl_receive_ocean_fields
[996]587
[782]588!
589!*************************************************************************************
590!
[996]591
[782]592  SUBROUTINE cpl_receive_seaice_fields(knon, knindex, &
[1146]593       tsurf_new, alb_new, u0_new, v0_new)
[782]594!
595! This routine returns the fields for the seaice that have been read from the coupler
[996]596! (done earlier with cpl_receive_frac). These fields are the temperature and
[782]597! albedo at sea ice surface and fraction of sea ice.
[996]598! The fields are transformed into 1D arrays with valid points from index 1 to knon.
[782]599!
600
601! Input arguments
602!*************************************************************************************
603    INTEGER, INTENT(IN)                     :: knon
604    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
605
606! Output arguments
607!*************************************************************************************
608    REAL, DIMENSION(klon), INTENT(OUT)      :: tsurf_new
609    REAL, DIMENSION(klon), INTENT(OUT)      :: alb_new
[1146]610    REAL, DIMENSION(klon), INTENT(OUT)      :: u0_new
611    REAL, DIMENSION(klon), INTENT(OUT)      :: v0_new
[782]612
[996]613! Local variables
614!*************************************************************************************
615    INTEGER               :: i
616    REAL, DIMENSION(klon) :: sic_new
[782]617
618!*************************************************************************************
619! Transform fields read from coupler from 2D into compressed 1D variables
620!
621!*************************************************************************************
622    CALL cpl2gath(read_sit, tsurf_new, knon, knindex)
623    CALL cpl2gath(read_alb_sic, alb_new, knon, knindex)
[996]624    CALL cpl2gath(read_sic, sic_new, knon, knindex)
[1146]625    CALL cpl2gath(read_u0, u0_new, knon, knindex)
626    CALL cpl2gath(read_v0, v0_new, knon, knindex)
[782]627
[996]628!*************************************************************************************
629! The fields received from the coupler have to be weighted with the sea-ice
630! concentration (in relation to the total sea-ice + ocean).
631!
632!*************************************************************************************
633    DO i= 1, knon
634       tsurf_new(i) = tsurf_new(i) / sic_new(i)
635       alb_new(i)   = alb_new(i)   / sic_new(i)
[3605]636    ENDDO
[996]637
[782]638  END SUBROUTINE cpl_receive_seaice_fields
639
640!
641!*************************************************************************************
642!
643
644  SUBROUTINE cpl_send_ocean_fields(itime, knon, knindex, &
645       swdown, lwdown, fluxlat, fluxsens, &
[2872]646       precip_rain, precip_snow, evap, tsurf, fder, albsol, taux, tauy, windsp,&
[3744]647       sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol, delta_sst, &
[3767]648       delta_sal)
[3627]649
650    ! This subroutine cumulates some fields for each time-step during
651    ! a coupling period. At last time-step in a coupling period the
652    ! fields are transformed to the grid accepted by the coupler. No
653    ! sending to the coupler will be done from here (it is done in
654    ! cpl_send_seaice_fields). Crucial hypothesis is that the surface
655    ! fractions do not change between coupling time-steps.
656
[1279]657    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl, co2_send
[1785]658    USE indice_sol_mod
[2344]659    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[3627]660    use config_ocean_skin_m, only: activate_ocean_skin
[782]661
662! Input arguments
663!*************************************************************************************
664    INTEGER, INTENT(IN)                     :: itime
665    INTEGER, INTENT(IN)                     :: knon
666    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
667    REAL, DIMENSION(klon), INTENT(IN)       :: swdown, lwdown
668    REAL, DIMENSION(klon), INTENT(IN)       :: fluxlat, fluxsens
669    REAL, DIMENSION(klon), INTENT(IN)       :: precip_rain, precip_snow
670    REAL, DIMENSION(klon), INTENT(IN)       :: evap, tsurf, fder, albsol
671    REAL, DIMENSION(klon), INTENT(IN)       :: taux, tauy, windsp
[3687]672    REAL, INTENT(IN):: sens_prec_liq(:), sens_prec_sol(:) ! (knon)
[2872]673    REAL, DIMENSION(klon), INTENT(IN)       :: lat_prec_liq, lat_prec_sol
[3740]674   
[3744]675    REAL, intent(in):: delta_sst(:) ! (knon)
676    ! Ocean-air interface temperature minus bulk SST, in
677    ! K. Defined only if activate_ocean_skin >= 1.
[3740]678
[3767]679    real, intent(in):: delta_sal(:) ! (knon)
680    ! Ocean-air interface salinity minus bulk salinity, in ppt.
[782]681
682! Local variables
683!*************************************************************************************
684    INTEGER                                 :: cpl_index, ig
685    INTEGER                                 :: error, sum_error
686    CHARACTER(len = 25)                     :: modname = 'cpl_send_ocean_fields'
687    CHARACTER(len = 80)                     :: abort_message
688
689!*************************************************************************************
690! Start calculation
691! The ocean points are saved with second array index=1
692!
693!*************************************************************************************
694    cpl_index = 1
695
696!*************************************************************************************
697! Reset fields to zero in the beginning of a new coupling period
698!
699!*************************************************************************************
700    IF (MOD(itime, nexca) == 1) THEN
[996]701       cpl_sols(1:knon,cpl_index) = 0.0
702       cpl_nsol(1:knon,cpl_index) = 0.0
703       cpl_rain(1:knon,cpl_index) = 0.0
704       cpl_snow(1:knon,cpl_index) = 0.0
705       cpl_evap(1:knon,cpl_index) = 0.0
706       cpl_tsol(1:knon,cpl_index) = 0.0
707       cpl_fder(1:knon,cpl_index) = 0.0
708       cpl_albe(1:knon,cpl_index) = 0.0
709       cpl_taux(1:knon,cpl_index) = 0.0
710       cpl_tauy(1:knon,cpl_index) = 0.0
711       cpl_windsp(1:knon,cpl_index) = 0.0
[2872]712       cpl_sens_rain(1:knon,cpl_index) = 0.0
713       cpl_sens_snow(1:knon,cpl_index) = 0.0
[1279]714       cpl_taumod(1:knon,cpl_index) = 0.0
715       IF (carbon_cycle_cpl) cpl_atm_co2(1:knon,cpl_index) = 0.0
[3628]716
717       if (activate_ocean_skin == 2) then
[3767]718          cpl_delta_sst = 0.
719          cpl_delta_sal = 0.
[3628]720       end if
[782]721    ENDIF
722       
723!*************************************************************************************
724! Cumulate at each time-step
725!
726!*************************************************************************************   
727    DO ig = 1, knon
728       cpl_sols(ig,cpl_index) = cpl_sols(ig,cpl_index) + &
[1403]729            swdown(ig)      / REAL(nexca)
[782]730       cpl_nsol(ig,cpl_index) = cpl_nsol(ig,cpl_index) + &
[1403]731            (lwdown(ig) + fluxlat(ig) +fluxsens(ig)) / REAL(nexca)
[782]732       cpl_rain(ig,cpl_index) = cpl_rain(ig,cpl_index) + &
[1403]733            precip_rain(ig) / REAL(nexca)
[782]734       cpl_snow(ig,cpl_index) = cpl_snow(ig,cpl_index) + &
[1403]735            precip_snow(ig) / REAL(nexca)
[782]736       cpl_evap(ig,cpl_index) = cpl_evap(ig,cpl_index) + &
[1403]737            evap(ig)        / REAL(nexca)
[782]738       cpl_tsol(ig,cpl_index) = cpl_tsol(ig,cpl_index) + &
[1403]739            tsurf(ig)       / REAL(nexca)
[782]740       cpl_fder(ig,cpl_index) = cpl_fder(ig,cpl_index) + &
[1403]741            fder(ig)        / REAL(nexca)
[782]742       cpl_albe(ig,cpl_index) = cpl_albe(ig,cpl_index) + &
[1403]743            albsol(ig)      / REAL(nexca)
[782]744       cpl_taux(ig,cpl_index) = cpl_taux(ig,cpl_index) + &
[1403]745            taux(ig)        / REAL(nexca)
[782]746       cpl_tauy(ig,cpl_index) = cpl_tauy(ig,cpl_index) + &
[1403]747            tauy(ig)        / REAL(nexca)     
[782]748       cpl_windsp(ig,cpl_index) = cpl_windsp(ig,cpl_index) + &
[1403]749            windsp(ig)      / REAL(nexca)
[2872]750       cpl_sens_rain(ig,cpl_index) = cpl_sens_rain(ig,cpl_index) + &
751            sens_prec_liq(ig)      / REAL(nexca)
752       cpl_sens_snow(ig,cpl_index) = cpl_sens_snow(ig,cpl_index) + &
753            sens_prec_sol(ig)      / REAL(nexca)
[1279]754       cpl_taumod(ig,cpl_index) =   cpl_taumod(ig,cpl_index) + &
[1403]755          SQRT ( taux(ig)*taux(ig)+tauy(ig)*tauy(ig) ) / REAL (nexca)
[782]756
[1279]757       IF (carbon_cycle_cpl) THEN
758          cpl_atm_co2(ig,cpl_index) = cpl_atm_co2(ig,cpl_index) + &
[1403]759               co2_send(knindex(ig))/ REAL(nexca)
[3605]760!!---OB: this is correct but why knindex ??
761       ENDIF
[3627]762
[3628]763       if (activate_ocean_skin == 2) then
[3744]764          cpl_delta_sst(ig) = cpl_delta_sst(ig) + delta_sst(ig) / REAL(nexca)
[3767]765          cpl_delta_sal(ig) = cpl_delta_sal(ig) + delta_sal(ig) / REAL(nexca)
[3628]766       end if
[1279]767     ENDDO
768
[782]769!*************************************************************************************
770! If the time-step corresponds to the end of coupling period the
771! fields are transformed to the 2D grid.
772! No sending to the coupler (it is done from cpl_send_seaice_fields).
773!
774!*************************************************************************************
775    IF (MOD(itime, nexca) == 0) THEN
776
777       IF (.NOT. ALLOCATED(cpl_sols2D)) THEN
778          sum_error = 0
[2344]779          ALLOCATE(cpl_sols2D(nbp_lon,jj_nb,2), stat=error)
[782]780          sum_error = sum_error + error
[2344]781          ALLOCATE(cpl_nsol2D(nbp_lon,jj_nb,2), stat=error)
[782]782          sum_error = sum_error + error
[2344]783          ALLOCATE(cpl_rain2D(nbp_lon,jj_nb,2), stat=error)
[782]784          sum_error = sum_error + error
[2344]785          ALLOCATE(cpl_snow2D(nbp_lon,jj_nb,2), stat=error)
[782]786          sum_error = sum_error + error
[2344]787          ALLOCATE(cpl_evap2D(nbp_lon,jj_nb,2), stat=error)
[782]788          sum_error = sum_error + error
[2344]789          ALLOCATE(cpl_tsol2D(nbp_lon,jj_nb,2), stat=error)
[782]790          sum_error = sum_error + error
[2344]791          ALLOCATE(cpl_fder2D(nbp_lon,jj_nb,2), stat=error)
[782]792          sum_error = sum_error + error
[2344]793          ALLOCATE(cpl_albe2D(nbp_lon,jj_nb,2), stat=error)
[782]794          sum_error = sum_error + error
[2344]795          ALLOCATE(cpl_taux2D(nbp_lon,jj_nb,2), stat=error)
[782]796          sum_error = sum_error + error
[2344]797          ALLOCATE(cpl_tauy2D(nbp_lon,jj_nb,2), stat=error)
[782]798          sum_error = sum_error + error
[2344]799          ALLOCATE(cpl_windsp2D(nbp_lon,jj_nb), stat=error)
[782]800          sum_error = sum_error + error
[2872]801          ALLOCATE(cpl_sens_rain2D(nbp_lon,jj_nb,2), stat=error)
802          sum_error = sum_error + error
803          ALLOCATE(cpl_sens_snow2D(nbp_lon,jj_nb,2), stat=error)
804          sum_error = sum_error + error
[2344]805          ALLOCATE(cpl_taumod2D(nbp_lon,jj_nb,2), stat=error)
[1279]806          sum_error = sum_error + error
[782]807         
[1279]808          IF (carbon_cycle_cpl) THEN
[2344]809             ALLOCATE(cpl_atm_co22D(nbp_lon,jj_nb), stat=error)
[1279]810             sum_error = sum_error + error
[3605]811          ENDIF
[1279]812
[3627]813          if (activate_ocean_skin == 2) then
[3744]814             ALLOCATE(cpl_delta_sst_2D(nbp_lon, jj_nb), &
[3767]815                  cpl_delta_sal_2D(nbp_lon, jj_nb), stat = error)
[3627]816             sum_error = sum_error + error
817          end if
818
[782]819          IF (sum_error /= 0) THEN
820             abort_message='Pb allocation variables couplees pour l''ecriture'
[2311]821             CALL abort_physic(modname,abort_message,1)
[782]822          ENDIF
823       ENDIF
824       
825
[1146]826       CALL gath2cpl(cpl_sols(:,cpl_index), cpl_sols2D(:,:,cpl_index), &
[782]827            knon, knindex)
828
[1146]829       CALL gath2cpl(cpl_nsol(:,cpl_index), cpl_nsol2D(:,:,cpl_index), &
[782]830            knon, knindex)
831
[1146]832       CALL gath2cpl(cpl_rain(:,cpl_index), cpl_rain2D(:,:,cpl_index), &
[782]833            knon, knindex)
834
[1146]835       CALL gath2cpl(cpl_snow(:,cpl_index), cpl_snow2D(:,:,cpl_index), &
[782]836            knon, knindex)
837
[1146]838       CALL gath2cpl(cpl_evap(:,cpl_index), cpl_evap2D(:,:,cpl_index), &
[782]839            knon, knindex)
840
841! cpl_tsol2D(:,:,:) not used!
[1146]842       CALL gath2cpl(cpl_tsol(:,cpl_index), cpl_tsol2D(:,:, cpl_index), &
[782]843            knon, knindex)
844
845! cpl_fder2D(:,:,1) not used, only cpl_fder(:,:,2)!
[1146]846       CALL gath2cpl(cpl_fder(:,cpl_index), cpl_fder2D(:,:,cpl_index), &
[782]847            knon, knindex)
848
849! cpl_albe2D(:,:,:) not used!
[1146]850       CALL gath2cpl(cpl_albe(:,cpl_index), cpl_albe2D(:,:,cpl_index), &
[782]851            knon, knindex)
852
[1146]853       CALL gath2cpl(cpl_taux(:,cpl_index), cpl_taux2D(:,:,cpl_index), &
[782]854            knon, knindex)
855
[1146]856       CALL gath2cpl(cpl_tauy(:,cpl_index), cpl_tauy2D(:,:,cpl_index), &
[782]857            knon, knindex)
858
[1146]859       CALL gath2cpl(cpl_windsp(:,cpl_index), cpl_windsp2D(:,:), &
[782]860            knon, knindex)
861
[2872]862       CALL gath2cpl(cpl_sens_rain(:,cpl_index), cpl_sens_rain2D(:,:,cpl_index), &
863            knon, knindex)
864
865       CALL gath2cpl(cpl_sens_snow(:,cpl_index), cpl_sens_snow2D(:,:,cpl_index), &
866            knon, knindex)
867
[1279]868       CALL gath2cpl(cpl_taumod(:,cpl_index), cpl_taumod2D(:,:,cpl_index), &
869            knon, knindex)
[782]870
[1279]871       IF (carbon_cycle_cpl) &
872            CALL gath2cpl(cpl_atm_co2(:,cpl_index), cpl_atm_co22D(:,:), knon, knindex)
[3628]873       if (activate_ocean_skin == 2) then
[3744]874          CALL gath2cpl(cpl_delta_sst, cpl_delta_sst_2D, knon, knindex)
[3767]875          CALL gath2cpl(cpl_delta_sal, cpl_delta_sal_2D, knon, knindex)
[3628]876       end if
[3627]877    ENDIF
[1279]878
[782]879  END SUBROUTINE cpl_send_ocean_fields
880
881!
882!*************************************************************************************
883!
884
885  SUBROUTINE cpl_send_seaice_fields(itime, dtime, knon, knindex, &
886       pctsrf, lafin, rlon, rlat, &
887       swdown, lwdown, fluxlat, fluxsens, &
[2872]888       precip_rain, precip_snow, evap, tsurf, fder, albsol, taux, tauy,&
889       sens_prec_liq, sens_prec_sol, lat_prec_liq, lat_prec_sol)
[782]890!
891! This subroutine cumulates some fields for each time-step during a coupling
892! period. At last time-step in a coupling period the fields are transformed to the
893! grid accepted by the coupler. All fields for all types of surfaces are sent to
894! the coupler.
895!
[1279]896    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl
[1785]897    USE indice_sol_mod
[2344]898    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[782]899
900! Input arguments
901!*************************************************************************************
902    INTEGER, INTENT(IN)                     :: itime
903    INTEGER, INTENT(IN)                     :: knon
904    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
905    REAL, INTENT(IN)                        :: dtime
906    REAL, DIMENSION(klon), INTENT(IN)       :: rlon, rlat
907    REAL, DIMENSION(klon), INTENT(IN)       :: swdown, lwdown
908    REAL, DIMENSION(klon), INTENT(IN)       :: fluxlat, fluxsens
909    REAL, DIMENSION(klon), INTENT(IN)       :: precip_rain, precip_snow
910    REAL, DIMENSION(klon), INTENT(IN)       :: evap, tsurf, fder
911    REAL, DIMENSION(klon), INTENT(IN)       :: albsol, taux, tauy
912    REAL, DIMENSION(klon,nbsrf), INTENT(IN) :: pctsrf
[3687]913    REAL, INTENT(IN):: sens_prec_liq(:), sens_prec_sol(:) ! (knon)
[2872]914    REAL, DIMENSION(klon), INTENT(IN)       :: lat_prec_liq, lat_prec_sol
[782]915    LOGICAL, INTENT(IN)                     :: lafin
916
917! Local variables
918!*************************************************************************************
919    INTEGER                                 :: cpl_index, ig
920    INTEGER                                 :: error, sum_error
921    CHARACTER(len = 25)                     :: modname = 'cpl_send_seaice_fields'
922    CHARACTER(len = 80)                     :: abort_message
[1146]923    REAL, DIMENSION(klon)                   :: cpl_fder_tmp
[782]924
925!*************************************************************************************
926! Start calulation
927! The sea-ice points are saved with second array index=2
928!
929!*************************************************************************************
930    cpl_index = 2
931
932!*************************************************************************************
933! Reset fields to zero in the beginning of a new coupling period
934!
935!*************************************************************************************
936    IF (MOD(itime, nexca) == 1) THEN
[996]937       cpl_sols(1:knon,cpl_index) = 0.0
938       cpl_nsol(1:knon,cpl_index) = 0.0
939       cpl_rain(1:knon,cpl_index) = 0.0
940       cpl_snow(1:knon,cpl_index) = 0.0
941       cpl_evap(1:knon,cpl_index) = 0.0
942       cpl_tsol(1:knon,cpl_index) = 0.0
943       cpl_fder(1:knon,cpl_index) = 0.0
944       cpl_albe(1:knon,cpl_index) = 0.0
945       cpl_taux(1:knon,cpl_index) = 0.0
946       cpl_tauy(1:knon,cpl_index) = 0.0
[2872]947       cpl_sens_rain(1:knon,cpl_index) = 0.0
948       cpl_sens_snow(1:knon,cpl_index) = 0.0
[1279]949       cpl_taumod(1:knon,cpl_index) = 0.0
[782]950    ENDIF
951       
952!*************************************************************************************
953! Cumulate at each time-step
954!
955!*************************************************************************************   
956    DO ig = 1, knon
957       cpl_sols(ig,cpl_index) = cpl_sols(ig,cpl_index) + &
[1403]958            swdown(ig)      / REAL(nexca)
[782]959       cpl_nsol(ig,cpl_index) = cpl_nsol(ig,cpl_index) + &
[1403]960            (lwdown(ig) + fluxlat(ig) +fluxsens(ig)) / REAL(nexca)
[782]961       cpl_rain(ig,cpl_index) = cpl_rain(ig,cpl_index) + &
[1403]962            precip_rain(ig) / REAL(nexca)
[782]963       cpl_snow(ig,cpl_index) = cpl_snow(ig,cpl_index) + &
[1403]964            precip_snow(ig) / REAL(nexca)
[782]965       cpl_evap(ig,cpl_index) = cpl_evap(ig,cpl_index) + &
[1403]966            evap(ig)        / REAL(nexca)
[782]967       cpl_tsol(ig,cpl_index) = cpl_tsol(ig,cpl_index) + &
[1403]968            tsurf(ig)       / REAL(nexca)
[782]969       cpl_fder(ig,cpl_index) = cpl_fder(ig,cpl_index) + &
[1403]970            fder(ig)        / REAL(nexca)
[782]971       cpl_albe(ig,cpl_index) = cpl_albe(ig,cpl_index) + &
[1403]972            albsol(ig)      / REAL(nexca)
[782]973       cpl_taux(ig,cpl_index) = cpl_taux(ig,cpl_index) + &
[1403]974            taux(ig)        / REAL(nexca)
[782]975       cpl_tauy(ig,cpl_index) = cpl_tauy(ig,cpl_index) + &
[1403]976            tauy(ig)        / REAL(nexca)     
[2872]977       cpl_sens_rain(ig,cpl_index) = cpl_sens_rain(ig,cpl_index) + &
978            sens_prec_liq(ig)      / REAL(nexca)
979       cpl_sens_snow(ig,cpl_index) = cpl_sens_snow(ig,cpl_index) + &
980            sens_prec_sol(ig)      / REAL(nexca)
[1279]981       cpl_taumod(ig,cpl_index) = cpl_taumod(ig,cpl_index) + &
[1403]982            SQRT ( taux(ig)*taux(ig)+tauy(ig)*tauy(ig) ) / REAL(nexca)
[782]983    ENDDO
984
985!*************************************************************************************
986! If the time-step corresponds to the end of coupling period the
987! fields are transformed to the 2D grid and all fields are sent to coupler.
988!
989!*************************************************************************************
990    IF (MOD(itime, nexca) == 0) THEN
991       IF (.NOT. ALLOCATED(cpl_sols2D)) THEN
992          sum_error = 0
[2344]993          ALLOCATE(cpl_sols2D(nbp_lon,jj_nb,2), stat=error)
[782]994          sum_error = sum_error + error
[2344]995          ALLOCATE(cpl_nsol2D(nbp_lon,jj_nb,2), stat=error)
[782]996          sum_error = sum_error + error
[2344]997          ALLOCATE(cpl_rain2D(nbp_lon,jj_nb,2), stat=error)
[782]998          sum_error = sum_error + error
[2344]999          ALLOCATE(cpl_snow2D(nbp_lon,jj_nb,2), stat=error)
[782]1000          sum_error = sum_error + error
[2344]1001          ALLOCATE(cpl_evap2D(nbp_lon,jj_nb,2), stat=error)
[782]1002          sum_error = sum_error + error
[2344]1003          ALLOCATE(cpl_tsol2D(nbp_lon,jj_nb,2), stat=error)
[782]1004          sum_error = sum_error + error
[2344]1005          ALLOCATE(cpl_fder2D(nbp_lon,jj_nb,2), stat=error)
[782]1006          sum_error = sum_error + error
[2344]1007          ALLOCATE(cpl_albe2D(nbp_lon,jj_nb,2), stat=error)
[782]1008          sum_error = sum_error + error
[2344]1009          ALLOCATE(cpl_taux2D(nbp_lon,jj_nb,2), stat=error)
[782]1010          sum_error = sum_error + error
[2344]1011          ALLOCATE(cpl_tauy2D(nbp_lon,jj_nb,2), stat=error)
[782]1012          sum_error = sum_error + error
[2344]1013          ALLOCATE(cpl_windsp2D(nbp_lon,jj_nb), stat=error)
[782]1014          sum_error = sum_error + error
[2872]1015          ALLOCATE(cpl_sens_rain2D(nbp_lon,jj_nb,2), stat=error)
1016          sum_error = sum_error + error
1017          ALLOCATE(cpl_sens_snow2D(nbp_lon,jj_nb,2), stat=error)
1018          sum_error = sum_error + error
[2344]1019          ALLOCATE(cpl_taumod2D(nbp_lon,jj_nb,2), stat=error)
[1279]1020          sum_error = sum_error + error
1021
1022          IF (carbon_cycle_cpl) THEN
[2344]1023             ALLOCATE(cpl_atm_co22D(nbp_lon,jj_nb), stat=error)
[1279]1024             sum_error = sum_error + error
[3605]1025          ENDIF
[1279]1026
[782]1027          IF (sum_error /= 0) THEN
1028             abort_message='Pb allocation variables couplees pour l''ecriture'
[2311]1029             CALL abort_physic(modname,abort_message,1)
[782]1030          ENDIF
1031       ENDIF
1032
[1146]1033       CALL gath2cpl(cpl_sols(:,cpl_index), cpl_sols2D(:,:,cpl_index), &
[782]1034            knon, knindex)
1035
[1146]1036       CALL gath2cpl(cpl_nsol(:,cpl_index), cpl_nsol2D(:,:,cpl_index), &
[782]1037            knon, knindex)
1038
[1146]1039       CALL gath2cpl(cpl_rain(:,cpl_index), cpl_rain2D(:,:,cpl_index), &
[782]1040            knon, knindex)
1041
[1146]1042       CALL gath2cpl(cpl_snow(:,cpl_index), cpl_snow2D(:,:,cpl_index), &
[782]1043            knon, knindex)
1044
[1146]1045       CALL gath2cpl(cpl_evap(:,cpl_index), cpl_evap2D(:,:,cpl_index), &
[782]1046            knon, knindex)
1047
1048! cpl_tsol2D(:,:,:) not used!
[1146]1049       CALL gath2cpl(cpl_tsol(:,cpl_index), cpl_tsol2D(:,:, cpl_index), &
[782]1050            knon, knindex)
1051
[1146]1052       ! Set default value and decompress before gath2cpl
1053       cpl_fder_tmp(:) = -20.
1054       DO ig = 1, knon
1055          cpl_fder_tmp(knindex(ig))=cpl_fder(ig,cpl_index)
[3605]1056       ENDDO
[1146]1057       CALL gath2cpl(cpl_fder_tmp(:), cpl_fder2D(:,:,cpl_index), &
1058            klon, unity)
[782]1059
1060! cpl_albe2D(:,:,:) not used!
[1146]1061       CALL gath2cpl(cpl_albe(:,cpl_index), cpl_albe2D(:,:,cpl_index), &
[782]1062            knon, knindex)
1063
[1146]1064       CALL gath2cpl(cpl_taux(:,cpl_index), cpl_taux2D(:,:,cpl_index), &
[782]1065            knon, knindex)
1066
[1146]1067       CALL gath2cpl(cpl_tauy(:,cpl_index), cpl_tauy2D(:,:,cpl_index), &
[782]1068            knon, knindex)
1069
[2872]1070       CALL gath2cpl(cpl_sens_rain(:,cpl_index), cpl_sens_rain2D(:,:,cpl_index), &
1071            knon, knindex)
1072
1073       CALL gath2cpl(cpl_sens_snow(:,cpl_index), cpl_sens_snow2D(:,:,cpl_index), &
1074            knon, knindex)
1075
[1279]1076       CALL gath2cpl(cpl_taumod(:,cpl_index), cpl_taumod2D(:,:,cpl_index), &
1077            knon, knindex)
1078
[782]1079       ! Send all fields
1080       CALL cpl_send_all(itime, dtime, pctsrf, lafin, rlon, rlat)
1081    ENDIF
1082
1083  END SUBROUTINE cpl_send_seaice_fields
1084
1085!
1086!*************************************************************************************
1087!
1088
1089  SUBROUTINE cpl_send_land_fields(itime, knon, knindex, rriv_in, rcoa_in)
1090!
1091! This subroutine cumulates some fields for each time-step during a coupling
1092! period. At last time-step in a coupling period the fields are transformed to the
1093! grid accepted by the coupler. No sending to the coupler will be done from here
1094! (it is done in cpl_send_seaice_fields).
1095!
[2344]1096    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[782]1097
1098! Input arguments
1099!*************************************************************************************
1100    INTEGER, INTENT(IN)                       :: itime
1101    INTEGER, INTENT(IN)                       :: knon
1102    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
1103    REAL, DIMENSION(klon), INTENT(IN)         :: rriv_in
1104    REAL, DIMENSION(klon), INTENT(IN)         :: rcoa_in
1105
1106! Local variables
1107!*************************************************************************************
[2344]1108    REAL, DIMENSION(nbp_lon,jj_nb)             :: rriv2D
1109    REAL, DIMENSION(nbp_lon,jj_nb)             :: rcoa2D
[782]1110
1111!*************************************************************************************
1112! Rearrange fields in 2D variables
1113! First initialize to zero to avoid unvalid points causing problems
1114!
1115!*************************************************************************************
[987]1116!$OMP MASTER
[782]1117    rriv2D(:,:) = 0.0
1118    rcoa2D(:,:) = 0.0
[987]1119!$OMP END MASTER
[782]1120    CALL gath2cpl(rriv_in, rriv2D, knon, knindex)
1121    CALL gath2cpl(rcoa_in, rcoa2D, knon, knindex)
1122
1123!*************************************************************************************
1124! Reset cumulated fields to zero in the beginning of a new coupling period
1125!
1126!*************************************************************************************
1127    IF (MOD(itime, nexca) == 1) THEN
[987]1128!$OMP MASTER
[782]1129       cpl_rriv2D(:,:) = 0.0
1130       cpl_rcoa2D(:,:) = 0.0
[987]1131!$OMP END MASTER
[782]1132    ENDIF
1133
1134!*************************************************************************************
1135! Cumulate : Following fields should be cumulated at each time-step
1136!
1137!*************************************************************************************   
[987]1138!$OMP MASTER
[1403]1139    cpl_rriv2D(:,:) = cpl_rriv2D(:,:) + rriv2D(:,:) / REAL(nexca)
1140    cpl_rcoa2D(:,:) = cpl_rcoa2D(:,:) + rcoa2D(:,:) / REAL(nexca)
[987]1141!$OMP END MASTER
[782]1142
1143  END SUBROUTINE cpl_send_land_fields
1144
1145!
1146!*************************************************************************************
1147!
1148
1149  SUBROUTINE cpl_send_landice_fields(itime, knon, knindex, rlic_in)
1150! This subroutine cumulates the field for melting ice for each time-step
1151! during a coupling period. This routine will not send to coupler. Sending
1152! will be done in cpl_send_seaice_fields.
1153!
[1279]1154
[2344]1155    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[996]1156
[782]1157! Input varibales
1158!*************************************************************************************
1159    INTEGER, INTENT(IN)                       :: itime
1160    INTEGER, INTENT(IN)                       :: knon
1161    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
1162    REAL, DIMENSION(klon), INTENT(IN)         :: rlic_in
1163
1164! Local varibales
1165!*************************************************************************************
[2344]1166    REAL, DIMENSION(nbp_lon,jj_nb)             :: rlic2D
[782]1167
1168!*************************************************************************************
1169! Rearrange field in a 2D variable
1170! First initialize to zero to avoid unvalid points causing problems
1171!
1172!*************************************************************************************
[987]1173!$OMP MASTER
[782]1174    rlic2D(:,:) = 0.0
[987]1175!$OMP END MASTER
[782]1176    CALL gath2cpl(rlic_in, rlic2D, knon, knindex)
1177
1178!*************************************************************************************
1179! Reset field to zero in the beginning of a new coupling period
1180!
1181!*************************************************************************************
1182    IF (MOD(itime, nexca) == 1) THEN
[987]1183!$OMP MASTER
[782]1184       cpl_rlic2D(:,:) = 0.0
[987]1185!$OMP END MASTER
[782]1186    ENDIF
1187
1188!*************************************************************************************
1189! Cumulate : Melting ice should be cumulated at each time-step
1190!
1191!*************************************************************************************   
[987]1192!$OMP MASTER
[1403]1193    cpl_rlic2D(:,:) = cpl_rlic2D(:,:) + rlic2D(:,:) / REAL(nexca)
[987]1194!$OMP END MASTER
[782]1195
1196  END SUBROUTINE cpl_send_landice_fields
1197
1198!
1199!*************************************************************************************
1200!
1201
1202  SUBROUTINE cpl_send_all(itime, dtime, pctsrf, lafin, rlon, rlat)
1203! This routine will send fields for all different surfaces to the coupler.
1204! This subroutine should be executed after calculations by the last surface(sea-ice),
1205! all calculations at the different surfaces have to be done before.
1206!   
[996]1207    USE surface_data
[1279]1208    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl
[1785]1209    USE indice_sol_mod
[2344]1210    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
1211    USE time_phylmdz_mod, ONLY: start_time, itau_phy
[3627]1212    use config_ocean_skin_m, only: activate_ocean_skin
[782]1213! Some includes
[2344]1214!   
[782]1215! Input arguments
1216!*************************************************************************************
1217    INTEGER, INTENT(IN)                                  :: itime
1218    REAL, INTENT(IN)                                     :: dtime
1219    REAL, DIMENSION(klon), INTENT(IN)                    :: rlon, rlat
1220    REAL, DIMENSION(klon,nbsrf), INTENT(IN)              :: pctsrf
1221    LOGICAL, INTENT(IN)                                  :: lafin
1222   
1223! Local variables
1224!*************************************************************************************
[3605]1225    INTEGER                                              :: error, sum_error, i,j,k
[782]1226    INTEGER                                              :: itau_w
[996]1227    INTEGER                                              :: time_sec
[2344]1228    INTEGER, DIMENSION(nbp_lon*(nbp_lat))                      :: ndexct
[782]1229    REAL                                                 :: Up, Down
[2344]1230    REAL, DIMENSION(nbp_lon, jj_nb)                          :: tmp_lon, tmp_lat
1231    REAL, DIMENSION(nbp_lon, jj_nb, 4)                       :: pctsrf2D
1232    REAL, DIMENSION(nbp_lon, jj_nb)                          :: deno
[782]1233    CHARACTER(len = 20)                                  :: modname = 'cpl_send_all'
1234    CHARACTER(len = 80)                                  :: abort_message
1235   
1236! Variables with fields to coupler
[2344]1237    REAL, DIMENSION(nbp_lon, jj_nb)                          :: tmp_taux
1238    REAL, DIMENSION(nbp_lon, jj_nb)                          :: tmp_tauy
1239    REAL, DIMENSION(nbp_lon, jj_nb)                          :: tmp_calv
[782]1240! Table with all fields to send to coupler
[2344]1241    REAL, DIMENSION(nbp_lon, jj_nb, maxsend)                 :: tab_flds
[3605]1242    REAL, DIMENSION(klon_mpi)                                :: rlon_mpi, rlat_mpi
1243    REAL  :: calving(nb_zone_calving)
1244    REAL  :: calving_glo(nb_zone_calving)
1245   
[1001]1246#ifdef CPP_MPI
[782]1247    INCLUDE 'mpif.h'
1248    INTEGER, DIMENSION(MPI_STATUS_SIZE)                  :: status
1249#endif
1250
1251! End definitions
1252!*************************************************************************************
1253   
1254
1255
1256!*************************************************************************************
1257! All fields are stored in a table tab_flds(:,:,:)
[1146]1258! First store the fields which are already on the right format
[782]1259!
1260!*************************************************************************************
[987]1261!$OMP MASTER
[1279]1262    tab_flds(:,:,ids_windsp) = cpl_windsp2D(:,:)
1263    tab_flds(:,:,ids_shfice) = cpl_sols2D(:,:,2)
1264    tab_flds(:,:,ids_nsfice) = cpl_nsol2D(:,:,2)
1265    tab_flds(:,:,ids_dflxdt) = cpl_fder2D(:,:,2)
[2872]1266    tab_flds(:,:,ids_qraioc) = cpl_sens_rain2D(:,:,1)
1267    tab_flds(:,:,ids_qsnooc) = cpl_sens_snow2D(:,:,1)
1268    tab_flds(:,:,ids_qraiic) = cpl_sens_rain2D(:,:,2)
1269    tab_flds(:,:,ids_qsnoic) = cpl_sens_snow2D(:,:,2)
[3628]1270
1271    if (activate_ocean_skin == 2) then
[3744]1272       tab_flds(:, :, ids_delta_sst) = cpl_delta_sst_2D
[3767]1273       tab_flds(:, :, ids_delta_sal) = cpl_delta_sal_2D
[3628]1274    end if
[1146]1275   
[996]1276    IF (version_ocean=='nemo') THEN
[3605]1277       tab_flds(:,:,ids_liqrun) = (cpl_rriv2D(:,:) + cpl_rcoa2D(:,:))
[1279]1278       IF (carbon_cycle_cpl) tab_flds(:,:,ids_atmco2)=cpl_atm_co22D(:,:)
[996]1279    ELSE IF (version_ocean=='opa8') THEN
[1279]1280       tab_flds(:,:,ids_shfoce) = cpl_sols2D(:,:,1)
1281       tab_flds(:,:,ids_nsfoce) = cpl_nsol2D(:,:,1)
1282       tab_flds(:,:,ids_icevap) = cpl_evap2D(:,:,2)
1283       tab_flds(:,:,ids_ocevap) = cpl_evap2D(:,:,1)
1284       tab_flds(:,:,ids_runcoa) = cpl_rcoa2D(:,:)
1285       tab_flds(:,:,ids_rivflu) = cpl_rriv2D(:,:)
[3605]1286    ENDIF
[1146]1287
[782]1288!*************************************************************************************
1289! Transform the fraction of sub-surfaces from 1D to 2D array
1290!
1291!*************************************************************************************
1292    pctsrf2D(:,:,:) = 0.
[987]1293!$OMP END MASTER
[782]1294    CALL gath2cpl(pctsrf(:,is_oce), pctsrf2D(:,:,is_oce), klon, unity)
1295    CALL gath2cpl(pctsrf(:,is_sic), pctsrf2D(:,:,is_sic), klon, unity)
1296    CALL gath2cpl(pctsrf(:,is_lic), pctsrf2D(:,:,is_lic), klon, unity)
1297
1298!*************************************************************************************
1299! Calculate the average calving per latitude
1300! Store calving in tab_flds(:,:,19)
1301!
1302!*************************************************************************************     
[987]1303    IF (is_omp_root) THEN
1304
[3605]1305      IF (cpl_old_calving) THEN   ! use old calving
1306
1307        DO j = 1, jj_nb
1308           tmp_calv(:,j) = DOT_PRODUCT (cpl_rlic2D(1:nbp_lon,j), &
1309                pctsrf2D(1:nbp_lon,j,is_lic)) / REAL(nbp_lon)
1310        ENDDO
[782]1311   
[1001]1312   
[3605]1313        IF (is_parallel) THEN
1314           IF (.NOT. is_north_pole_dyn) THEN
[1001]1315#ifdef CPP_MPI
[3605]1316              CALL MPI_RECV(Up,1,MPI_REAL_LMDZ,mpi_rank-1,1234,COMM_LMDZ_PHY,status,error)
1317              CALL MPI_SEND(tmp_calv(1,1),1,MPI_REAL_LMDZ,mpi_rank-1,1234,COMM_LMDZ_PHY,error)
[782]1318#endif
[3605]1319           ENDIF
[1001]1320       
[3605]1321           IF (.NOT. is_south_pole_dyn) THEN
[1001]1322#ifdef CPP_MPI
[3605]1323              CALL MPI_SEND(tmp_calv(1,jj_nb),1,MPI_REAL_LMDZ,mpi_rank+1,1234,COMM_LMDZ_PHY,error)
1324              CALL MPI_RECV(down,1,MPI_REAL_LMDZ,mpi_rank+1,1234,COMM_LMDZ_PHY,status,error)
[782]1325#endif
[3605]1326           ENDIF
[996]1327         
[3605]1328           IF (.NOT. is_north_pole_dyn .AND. ii_begin /=1) THEN
1329              Up=Up+tmp_calv(nbp_lon,1)
1330              tmp_calv(:,1)=Up
1331           ENDIF
1332           
1333           IF (.NOT. is_south_pole_dyn .AND. ii_end /= nbp_lon) THEN
1334              Down=Down+tmp_calv(1,jj_nb)
1335              tmp_calv(:,jj_nb)=Down
1336           ENDIF
1337        ENDIF
1338        tab_flds(:,:,ids_calvin) = tmp_calv(:,:)
1339
1340      ELSE
1341         ! cpl_old_calving=FALSE
1342         ! To be used with new method for calculation of coupling weights
1343         DO k=1,nb_zone_calving
1344            calving(k)=0
1345            DO j = 1, jj_nb
1346               calving(k)= calving(k)+DOT_PRODUCT(cpl_rlic2D(:,j)*area_calving(:,j,k),pctsrf2D(:,j,is_lic))
1347            ENDDO
1348         ENDDO
[996]1349         
[3605]1350#ifdef CPP_MPI
1351         CALL MPI_ALLREDUCE(calving, calving_glo, nb_zone_calving, MPI_REAL_LMDZ, MPI_SUM, COMM_LMDZ_PHY, error)
1352#endif
1353         
1354         tab_flds(:,:,ids_calvin) = 0
1355         DO k=1,nb_zone_calving
1356            IF (ind_calving(k)>0 ) THEN
1357               j=(ind_calving(k)-1)/nbp_lon + 1
1358               i=MOD(ind_calving(k)-1,nbp_lon)+1
1359               tab_flds(i,j,ids_calvin) = calving_glo(k)
1360            ENDIF
1361         ENDDO
1362         
[987]1363      ENDIF
[996]1364     
[782]1365!*************************************************************************************
1366! Calculate total flux for snow, rain and wind with weighted addition using the
1367! fractions of ocean and seaice.
1368!
1369!*************************************************************************************   
[996]1370       ! fraction oce+seaice
1371       deno =  pctsrf2D(:,:,is_oce) + pctsrf2D(:,:,is_sic)
[836]1372
[996]1373       IF (version_ocean=='nemo') THEN
[1279]1374          tab_flds(:,:,ids_shftot)  = 0.0
1375          tab_flds(:,:,ids_nsftot) = 0.0
1376          tab_flds(:,:,ids_totrai) = 0.0
1377          tab_flds(:,:,ids_totsno) = 0.0
1378          tab_flds(:,:,ids_toteva) = 0.0
1379          tab_flds(:,:,ids_taumod) = 0.0
[1146]1380 
[996]1381          tmp_taux(:,:)    = 0.0
1382          tmp_tauy(:,:)    = 0.0
1383          ! For all valid grid cells containing some fraction of ocean or sea-ice
1384          WHERE ( deno(:,:) /= 0 )
1385             tmp_taux = cpl_taux2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1386                  cpl_taux2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1387             tmp_tauy = cpl_tauy2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1388                  cpl_tauy2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1389
1390             tab_flds(:,:,ids_shftot) = cpl_sols2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[1146]1391                  cpl_sols2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1392             tab_flds(:,:,ids_nsftot) = cpl_nsol2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[1146]1393                  cpl_nsol2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1394             tab_flds(:,:,ids_totrai) = cpl_rain2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[1146]1395                  cpl_rain2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1396             tab_flds(:,:,ids_totsno) = cpl_snow2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[1146]1397                  cpl_snow2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1398             tab_flds(:,:,ids_toteva) = cpl_evap2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[1146]1399                  cpl_evap2D(:,:,2)  * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1400             tab_flds(:,:,ids_taumod) = cpl_taumod2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1401                  cpl_taumod2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1402             
[1146]1403         ENDWHERE
1404
[1279]1405          tab_flds(:,:,ids_icevap) = cpl_evap2D(:,:,2)
[996]1406         
1407       ELSE IF (version_ocean=='opa8') THEN
[1146]1408          ! Store fields for rain and snow in tab_flds(:,:,15) and tab_flds(:,:,16)
[1279]1409          tab_flds(:,:,ids_totrai) = 0.0
1410          tab_flds(:,:,ids_totsno) = 0.0
[996]1411          tmp_taux(:,:)    = 0.0
1412          tmp_tauy(:,:)    = 0.0
1413          ! For all valid grid cells containing some fraction of ocean or sea-ice
1414          WHERE ( deno(:,:) /= 0 )
[1279]1415             tab_flds(:,:,ids_totrai) = cpl_rain2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[996]1416                  cpl_rain2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
[1279]1417             tab_flds(:,:,ids_totsno) = cpl_snow2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
[996]1418                  cpl_snow2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1419             
1420             tmp_taux = cpl_taux2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1421                  cpl_taux2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1422             tmp_tauy = cpl_tauy2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1423                  cpl_tauy2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1424          ENDWHERE
[3605]1425       ENDIF
[987]1426
[996]1427    ENDIF ! is_omp_root
1428 
[782]1429!*************************************************************************************
1430! Transform the wind components from local atmospheric 2D coordinates to geocentric
1431! 3D coordinates.
1432! Store the resulting wind components in tab_flds(:,:,1:6)
1433!*************************************************************************************
1434
1435! Transform the longitudes and latitudes on 2D arrays
[1001]1436   
[987]1437    CALL gather_omp(rlon,rlon_mpi)
1438    CALL gather_omp(rlat,rlat_mpi)
1439!$OMP MASTER
1440    CALL Grid1DTo2D_mpi(rlon_mpi,tmp_lon)
1441    CALL Grid1DTo2D_mpi(rlat_mpi,tmp_lat)
1442!$OMP END MASTER   
1443
[782]1444    IF (is_sequential) THEN
[2429]1445       IF (is_north_pole_dyn) tmp_lon(:,1)     = tmp_lon(:,2)
1446       IF (is_south_pole_dyn) tmp_lon(:,nbp_lat) = tmp_lon(:,nbp_lat-1)
[782]1447    ENDIF
1448     
1449! NetCDF output of the wind before transformation of coordinate system
1450    IF (is_sequential) THEN
1451       ndexct(:) = 0
[2344]1452       itau_w = itau_phy + itime + start_time * day_step_phy
1453       CALL histwrite(nidct,'tauxe',itau_w,tmp_taux,nbp_lon*(nbp_lat),ndexct)
1454       CALL histwrite(nidct,'tauyn',itau_w,tmp_tauy,nbp_lon*(nbp_lat),ndexct)
1455       CALL histwrite(nidct,'tmp_lon',itau_w,tmp_lon,nbp_lon*(nbp_lat),ndexct)
1456       CALL histwrite(nidct,'tmp_lat',itau_w,tmp_lat,nbp_lon*(nbp_lat),ndexct)
[782]1457    ENDIF
1458
[1067]1459! Transform the wind from spherical atmospheric 2D coordinates to geocentric
1460! cartesian 3D coordinates
[987]1461!$OMP MASTER
[2344]1462    CALL atm2geo (nbp_lon, jj_nb, tmp_taux, tmp_tauy, tmp_lon, tmp_lat, &
[1279]1463         tab_flds(:,:,ids_tauxxu), tab_flds(:,:,ids_tauyyu), tab_flds(:,:,ids_tauzzu) )
[782]1464   
[1279]1465    tab_flds(:,:,ids_tauxxv)  = tab_flds(:,:,ids_tauxxu)
1466    tab_flds(:,:,ids_tauyyv)  = tab_flds(:,:,ids_tauyyu)
1467    tab_flds(:,:,ids_tauzzv)  = tab_flds(:,:,ids_tauzzu)
[987]1468!$OMP END MASTER
[782]1469
1470!*************************************************************************************
1471! NetCDF output of all fields just before sending to coupler.
1472!
1473!*************************************************************************************
1474    IF (is_sequential) THEN
[1279]1475        DO j=1,maxsend
1476          IF (infosend(j)%action) CALL histwrite(nidct,infosend(j)%name, itau_w, &
[2344]1477             tab_flds(:,:,j),nbp_lon*(nbp_lat),ndexct)
[1279]1478        ENDDO
[782]1479    ENDIF
1480!*************************************************************************************
1481! Send the table of all fields
1482!
1483!*************************************************************************************
[996]1484    time_sec=(itime-1)*dtime
[782]1485#ifdef CPP_COUPLE
[987]1486!$OMP MASTER
[1010]1487    CALL intocpl(time_sec, lafin, tab_flds(:,:,:))
[987]1488!$OMP END MASTER
[782]1489#endif
1490
1491!*************************************************************************************
1492! Finish with some dellocate
1493!
1494!************************************************************************************* 
1495    sum_error=0
1496    DEALLOCATE(cpl_sols2D, cpl_nsol2D, cpl_rain2D, cpl_snow2D, stat=error )
1497    sum_error = sum_error + error
1498    DEALLOCATE(cpl_evap2D, cpl_tsol2D, cpl_fder2D, cpl_albe2D, stat=error )
1499    sum_error = sum_error + error
[1279]1500    DEALLOCATE(cpl_taux2D, cpl_tauy2D, cpl_windsp2D, cpl_taumod2D, stat=error )
[782]1501    sum_error = sum_error + error
[2872]1502    DEALLOCATE(cpl_sens_rain2D, cpl_sens_snow2D, stat=error)
1503    sum_error = sum_error + error
1504
[1279]1505   
1506    IF (carbon_cycle_cpl) THEN
1507       DEALLOCATE(cpl_atm_co22D, stat=error )
1508       sum_error = sum_error + error
[3605]1509    ENDIF
[1279]1510
[3767]1511    if (activate_ocean_skin == 2) deallocate(cpl_delta_sst_2d, cpl_delta_sal_2d)
[3627]1512
[782]1513    IF (sum_error /= 0) THEN
1514       abort_message='Pb in deallocation of cpl_xxxx2D coupling variables'
[2311]1515       CALL abort_physic(modname,abort_message,1)
[782]1516    ENDIF
1517   
1518  END SUBROUTINE cpl_send_all
1519!
1520!*************************************************************************************
1521!
1522  SUBROUTINE cpl2gath(champ_in, champ_out, knon, knindex)
[1001]1523  USE mod_phys_lmdz_para
[1279]1524! Cette routine transforme un champs de la grille 2D recu du coupleur sur la grille
1525! 'gathered' (la grille physiq comprime).
[782]1526!
1527!
1528! input:         
[1279]1529!   champ_in     champ sur la grille 2D
[782]1530!   knon         nombre de points dans le domaine a traiter
1531!   knindex      index des points de la surface a traiter
1532!
1533! output:
[1279]1534!   champ_out    champ sur la grille 'gatherd'
[782]1535!
[2344]1536    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[782]1537
1538! Input
1539    INTEGER, INTENT(IN)                       :: knon
[2344]1540    REAL, DIMENSION(nbp_lon,jj_nb), INTENT(IN)    :: champ_in
[782]1541    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
1542
1543! Output
[987]1544    REAL, DIMENSION(klon_mpi), INTENT(OUT)        :: champ_out
[782]1545
1546! Local
1547    INTEGER                                   :: i, ig
[987]1548    REAL, DIMENSION(klon_mpi)                 :: temp_mpi
1549    REAL, DIMENSION(klon)                     :: temp_omp
[782]1550
1551!*************************************************************************************
1552!
[1001]1553   
1554
[2344]1555! Transform from 2 dimensions (nbp_lon,jj_nb) to 1 dimension (klon)
[987]1556!$OMP MASTER
1557    CALL Grid2Dto1D_mpi(champ_in,temp_mpi)
1558!$OMP END MASTER
[782]1559
[987]1560    CALL scatter_omp(temp_mpi,temp_omp)
1561   
[782]1562! Compress from klon to knon
1563    DO i = 1, knon
1564       ig = knindex(i)
[987]1565       champ_out(i) = temp_omp(ig)
[782]1566    ENDDO
[1001]1567
[782]1568  END SUBROUTINE cpl2gath
1569!
1570!*************************************************************************************
1571!
1572  SUBROUTINE gath2cpl(champ_in, champ_out, knon, knindex)
[987]1573  USE mod_phys_lmdz_para
[782]1574! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
1575! au coupleur.
1576!
1577! input:         
1578!   champ_in     champ sur la grille gathere       
1579!   knon         nombre de points dans le domaine a traiter
1580!   knindex      index des points de la surface a traiter
1581!
1582! output:
1583!   champ_out    champ sur la grille 2D
1584!
[2344]1585    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
[782]1586   
1587! Input arguments
1588!*************************************************************************************
1589    INTEGER, INTENT(IN)                    :: knon
1590    REAL, DIMENSION(klon), INTENT(IN)      :: champ_in
1591    INTEGER, DIMENSION(klon), INTENT(IN)   :: knindex
1592
1593! Output arguments
1594!*************************************************************************************
[2344]1595    REAL, DIMENSION(nbp_lon,jj_nb), INTENT(OUT) :: champ_out
[782]1596
1597! Local variables
1598!*************************************************************************************
1599    INTEGER                                :: i, ig
[987]1600    REAL, DIMENSION(klon)                  :: temp_omp
1601    REAL, DIMENSION(klon_mpi)              :: temp_mpi
[782]1602!*************************************************************************************
1603
1604! Decompress from knon to klon
[987]1605    temp_omp = 0.
[782]1606    DO i = 1, knon
1607       ig = knindex(i)
[987]1608       temp_omp(ig) = champ_in(i)
[782]1609    ENDDO
1610
[2344]1611! Transform from 1 dimension (klon) to 2 dimensions (nbp_lon,jj_nb)
[987]1612    CALL gather_omp(temp_omp,temp_mpi)
1613
1614!$OMP MASTER   
1615    CALL Grid1Dto2D_mpi(temp_mpi,champ_out)
[782]1616   
[2429]1617    IF (is_north_pole_dyn) champ_out(:,1)=temp_mpi(1)
1618    IF (is_south_pole_dyn) champ_out(:,jj_nb)=temp_mpi(klon)
[987]1619!$OMP END MASTER
[782]1620   
1621  END SUBROUTINE gath2cpl
1622!
1623!*************************************************************************************
1624!
1625END MODULE cpl_mod
1626
Note: See TracBrowser for help on using the repository browser.