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

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

Merge revisions 3427:3600 of trunk into branch Ocean_skin

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