source: LMDZ6/trunk/libf/phylmd/cpl_mod.F90 @ 3428

Last change on this file since 3428 was 3102, checked in by oboucher, 6 years ago

Removing x permission from these files

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