source: LMDZ4/trunk/libf/phylmd/cpl_mod.F90 @ 1054

Last change on this file since 1054 was 1010, checked in by lsce, 16 years ago

Corrected bug due to merge.
JG

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Author Date Id Revision
File size: 51.3 KB
Line 
1!
2! $Header$
3!
4MODULE cpl_mod
5!
6! This module excahanges and transforms all fields that should be recieved or sent to
7! coupler. The transformation of the fields are done from the grid 1D-array in phylmd
8! to the regular 2D grid accepted by the coupler. Cumulation of the fields for each
9! timestep is done in here.
10!
11! Each type of surface that recevie fields from the coupler have a subroutine named
12! cpl_receive_XXX_fields and each surface that have fields to be sent to the coupler
13! have a subroutine named cpl_send_XXX_fields.
14!
15!*************************************************************************************
16
17! Use statements
18!*************************************************************************************
19  USE dimphy, ONLY : klon
20  USE mod_phys_lmdz_para
21  USE ioipsl
22  USE iophy
23
24! The module oasis is always used. Without the cpp key CPP_COUPLE only the parameters
25! in the module are compiled and not the subroutines.
26  USE oasis
27  USE write_field_phy
28 
29! Global attributes
30!*************************************************************************************
31  IMPLICIT NONE
32  PRIVATE
33
34  ! All subroutine are public except cpl_send_all
35  PUBLIC :: cpl_init, cpl_receive_frac, cpl_receive_ocean_fields, cpl_receive_seaice_fields, &
36       cpl_send_ocean_fields, cpl_send_seaice_fields, cpl_send_land_fields, &
37       cpl_send_landice_fields, gath2cpl
38 
39
40! Declaration of module variables
41!*************************************************************************************
42! variable for coupling period
43  INTEGER, SAVE                             :: nexca
44  !$OMP THREADPRIVATE(nexca)
45
46! variables for cumulating fields during a coupling periode :
47  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_sols, cpl_nsol, cpl_rain
48  !$OMP THREADPRIVATE(cpl_sols,cpl_nsol,cpl_rain)
49  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_snow, cpl_evap, cpl_tsol
50  !$OMP THREADPRIVATE(cpl_snow,cpl_evap,cpl_tsol)
51  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_fder, cpl_albe, cpl_taux
52  !$OMP THREADPRIVATE(cpl_fder,cpl_albe,cpl_taux)
53  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_windsp
54  !$OMP THREADPRIVATE(cpl_windsp)
55  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_tauy
56  !$OMP THREADPRIVATE(cpl_tauy)
57  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_rriv2D, cpl_rcoa2D, cpl_rlic2D
58  !$OMP THREADPRIVATE(cpl_rriv2D,cpl_rcoa2D,cpl_rlic2D)
59
60! variables read from coupler :
61  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sst     ! sea surface temperature
62  !$OMP THREADPRIVATE(read_sst)
63  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sit     ! sea ice temperature
64  !$OMP THREADPRIVATE(read_sit)
65  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_sic     ! sea ice fraction
66  !$OMP THREADPRIVATE(read_sic)
67  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: read_alb_sic ! albedo at sea ice
68  !$OMP THREADPRIVATE(read_alb_sic)
69 
70  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE  :: unity
71  !$OMP THREADPRIVATE(unity)
72  INTEGER, SAVE                             :: nidct, nidcs
73  !$OMP THREADPRIVATE(nidct,nidcs)
74
75! variables to be sent to the coupler
76  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_sols2D, cpl_nsol2D, cpl_rain2D
77  !$OMP THREADPRIVATE(cpl_sols2D, cpl_nsol2D, cpl_rain2D)
78  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_snow2D, cpl_evap2D, cpl_tsol2D
79  !$OMP THREADPRIVATE(cpl_snow2D, cpl_evap2D, cpl_tsol2D)
80  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_fder2D, cpl_albe2D
81  !$OMP THREADPRIVATE(cpl_fder2D, cpl_albe2D)
82  REAL, ALLOCATABLE, DIMENSION(:,:,:), SAVE :: cpl_taux2D, cpl_tauy2D
83  !$OMP THREADPRIVATE(cpl_taux2D, cpl_tauy2D)
84  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: cpl_windsp2D
85  !$OMP THREADPRIVATE(cpl_windsp2D)
86 
87! variable for OPENMP parallelisation
88
89  INTEGER,ALLOCATABLE,DIMENSION(:),SAVE :: knon_omp
90  REAL,ALLOCATABLE,DIMENSION(:,:),SAVE ::  buffer_omp
91 
92CONTAINS
93!
94!************************************************************************************
95!
96  SUBROUTINE cpl_init(dtime, rlon, rlat)
97
98    INCLUDE "dimensions.h"
99    INCLUDE "indicesol.h"
100    INCLUDE "control.h"
101    INCLUDE "temps.h"
102    INCLUDE "iniprint.h"
103
104! Input arguments
105!*************************************************************************************
106    REAL, INTENT(IN)                  :: dtime
107    REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat
108
109! Local variables
110!*************************************************************************************
111    INTEGER                           :: error, sum_error, ig, i
112    INTEGER                           :: jf, nhoridct
113    INTEGER                           :: nhoridcs
114    INTEGER                           :: idtime
115    INTEGER                           :: idayref
116    INTEGER                           :: npas ! only for OASIS2
117    REAL                              :: zjulian
118    REAL, DIMENSION(iim,jjm+1)        :: zx_lon, zx_lat
119    CHARACTER(len = 20)               :: modname = 'cpl_init'
120    CHARACTER(len = 80)               :: abort_message
121    CHARACTER(len=80)                 :: clintocplnam, clfromcplnam
122
123!*************************************************************************************
124! Calculate coupling period
125!
126!*************************************************************************************
127     
128    npas = itaufin/ iphysiq
129    nexca = 86400 / dtime
130    WRITE(lunout,*)' ##### Ocean couple #####'
131    WRITE(lunout,*)' Valeurs des pas de temps'
132    WRITE(lunout,*)' npas = ', npas
133    WRITE(lunout,*)' nexca = ', nexca
134   
135!*************************************************************************************
136! Allocate variables
137!
138!*************************************************************************************
139    error = 0
140    sum_error = 0
141
142    ALLOCATE(unity(klon), stat = error)
143    sum_error = sum_error + error
144    ALLOCATE(cpl_sols(klon,2), stat = error)
145    sum_error = sum_error + error
146    ALLOCATE(cpl_nsol(klon,2), stat = error)
147    sum_error = sum_error + error
148    ALLOCATE(cpl_rain(klon,2), stat = error)
149    sum_error = sum_error + error
150    ALLOCATE(cpl_snow(klon,2), stat = error)
151    sum_error = sum_error + error
152    ALLOCATE(cpl_evap(klon,2), stat = error)
153    sum_error = sum_error + error
154    ALLOCATE(cpl_tsol(klon,2), stat = error)
155    sum_error = sum_error + error
156    ALLOCATE(cpl_fder(klon,2), stat = error)
157    sum_error = sum_error + error
158    ALLOCATE(cpl_albe(klon,2), stat = error)
159    sum_error = sum_error + error
160    ALLOCATE(cpl_taux(klon,2), stat = error)
161    sum_error = sum_error + error
162    ALLOCATE(cpl_windsp(klon,2), stat = error)
163    sum_error = sum_error + error
164    ALLOCATE(cpl_tauy(klon,2), stat = error)
165    sum_error = sum_error + error
166    ALLOCATE(cpl_rriv2D(iim,jj_nb), stat=error)
167    sum_error = sum_error + error
168    ALLOCATE(cpl_rcoa2D(iim,jj_nb), stat=error)
169    sum_error = sum_error + error
170    ALLOCATE(cpl_rlic2D(iim,jj_nb), stat=error)
171    sum_error = sum_error + error
172    ALLOCATE(read_sst(iim, jj_nb), stat = error)
173    sum_error = sum_error + error
174    ALLOCATE(read_sic(iim, jj_nb), stat = error)
175    sum_error = sum_error + error
176    ALLOCATE(read_sit(iim, jj_nb), stat = error)
177    sum_error = sum_error + error
178    ALLOCATE(read_alb_sic(iim, jj_nb), stat = error)
179    sum_error = sum_error + error
180
181    IF (sum_error /= 0) THEN
182       abort_message='Pb allocation variables couplees'
183       CALL abort_gcm(modname,abort_message,1)
184    ENDIF
185!*************************************************************************************
186! Initialize the allocated varaibles
187!
188!*************************************************************************************
189    DO ig = 1, klon
190       unity(ig) = ig
191    ENDDO
192
193!    cpl_sols = 0.   ; cpl_nsol = 0.  ; cpl_rain = 0.   ; cpl_snow = 0.
194!    cpl_evap = 0.   ; cpl_tsol = 0.  ; cpl_fder = 0.   ; cpl_albe = 0.
195!    cpl_taux = 0.   ; cpl_tauy = 0.  ; cpl_rriv2D = 0. ; cpl_rcoa2D = 0.
196!    cpl_rlic2D = 0. ; cpl_windsp = 0.
197
198!*************************************************************************************
199! Initialize coupling
200!
201!*************************************************************************************
202    idtime = INT(dtime)
203#ifdef CPP_COUPLE
204!$OMP MASTER   
205    CALL inicma
206!$OMP END MASTER
207#endif
208
209!*************************************************************************************
210! initialize NetCDF output
211!
212!*************************************************************************************
213    IF (is_sequential) THEN
214       idayref = day_ini
215       CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
216       CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlon,zx_lon)
217       DO i = 1, iim
218          zx_lon(i,1) = rlon(i+1)
219          zx_lon(i,jjm+1) = rlon(i+1)
220       ENDDO
221       CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlat,zx_lat)
222       clintocplnam="cpl_atm_tauflx"
223       CALL histbeg(clintocplnam, iim,zx_lon(:,1),jjm+1,zx_lat(1,:),&
224            1,iim,1,jjm+1, itau_phy,zjulian,dtime,nhoridct,nidct)
225! no vertical axis
226       CALL histdef(nidct, 'tauxe','tauxe', &
227            "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
228       CALL histdef(nidct, 'tauyn','tauyn', &
229            "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
230       CALL histdef(nidct, 'tmp_lon','tmp_lon', &
231            "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
232       CALL histdef(nidct, 'tmp_lat','tmp_lat', &
233            "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
234       DO jf=1,jpflda2o1 + jpflda2o2
235          CALL histdef(nidct, cl_writ(jf),cl_writ(jf), &
236               "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
237       END DO
238       CALL histend(nidct)
239       CALL histsync(nidct)
240       
241       clfromcplnam="cpl_atm_sst"
242       CALL histbeg(clfromcplnam, iim,zx_lon(:,1),jjm+1,zx_lat(1,:),1,iim,1,jjm+1, &
243            0,zjulian,dtime,nhoridcs,nidcs)
244! no vertical axis
245       DO jf=1,jpfldo2a
246          CALL histdef(nidcs, cl_read(jf),cl_read(jf), &
247               "-",iim, jjm+1, nhoridcs, 1, 1, 1, -99, 32, "inst", dtime,dtime)
248       END DO
249       CALL histend(nidcs)
250       CALL histsync(nidcs)
251
252    ENDIF    ! is_sequential
253   
254! OPENMP Initialization
255
256!$OMP MASTER
257  ALLOCATE(knon_omp(0:omp_size-1))
258  ALLOCATE(buffer_omp(klon_mpi,0:omp_size-1))       
259!$OMP END MASTER
260!$OMP BARRIER
261   
262  END SUBROUTINE cpl_init
263 
264!
265!*************************************************************************************
266!
267 
268  SUBROUTINE cpl_receive_frac(itime, dtime, pctsrf, is_modified)
269! This subroutine receives from coupler for both ocean and seaice
270! 4 fields : read_sst, read_sic, read_sit and read_alb_sic.
271! The new sea-ice-land-landice fraction is returned. The others fields
272! are stored in this module.
273    USE surface_data
274
275    INCLUDE "indicesol.h"
276    INCLUDE "temps.h"
277    INCLUDE "iniprint.h"
278    INCLUDE "YOMCST.h"
279    INCLUDE "dimensions.h"
280
281! Arguments
282!************************************************************************************
283    INTEGER, INTENT(IN)                        :: itime
284    REAL, INTENT(IN)                           :: dtime
285    REAL, DIMENSION(klon,nbsrf), INTENT(INOUT) :: pctsrf
286    LOGICAL, INTENT(OUT)                       :: is_modified
287
288! Local variables
289!************************************************************************************
290    INTEGER                                 :: j, i, time_sec
291    INTEGER                                 :: itau_w
292    INTEGER, DIMENSION(iim*(jjm+1))         :: ndexcs
293    CHARACTER(len = 20)                     :: modname = 'cpl_receive_frac'
294    CHARACTER(len = 80)                     :: abort_message
295    REAL, DIMENSION(klon)                   :: read_sic1D
296    REAL, DIMENSION(iim,jj_nb,jpfldo2a)     :: tab_read_flds
297    REAL, DIMENSION(klon,nbsrf)             :: pctsrf_old
298
299!*************************************************************************************
300! Start calculation
301! Get fields from coupler
302!
303!*************************************************************************************
304
305    is_modified=.FALSE.
306
307! Check if right moment to recevie from coupler
308    IF (MOD(itime, nexca) == 1) THEN
309       is_modified=.TRUE.
310 
311       time_sec=(itime-1)*dtime
312#ifdef CPP_COUPLE
313    time_sec=(itime-1)*dtime
314!$OMP MASTER
315    CALL fromcpl(time_sec, tab_read_flds)
316!$OMP END MASTER
317#endif
318   
319! NetCDF output of received fields
320       IF (is_sequential) THEN
321          ndexcs(:) = 0
322          itau_w = itau_phy + itime
323          DO i = 1, jpfldo2a
324             CALL histwrite(nidcs,cl_read(i),itau_w,tab_read_flds(:,:,i),iim*(jjm+1),ndexcs)
325          END DO
326       ENDIF
327
328
329! Save each field in a 2D array.
330!$OMP MASTER
331       IF (version_ocean=='nemo') THEN
332          read_sst(:,:)     = tab_read_flds(:,:,1)  ! Sea surface temperature
333          read_sic(:,:)     = tab_read_flds(:,:,2)  ! Sea ice concentration
334          read_sit(:,:)     = tab_read_flds(:,:,3)  ! Sea ice temperature
335          read_alb_sic(:,:) = tab_read_flds(:,:,4)  ! Albedo at sea ice
336       ELSE IF (version_ocean=='opa8') THEN
337          read_sst(:,:)     = tab_read_flds(:,:,1)  ! Sea surface temperature (multiplicated by fraction)
338          read_sic(:,:)     = tab_read_flds(:,:,2)  ! Sea ice concentration
339          read_alb_sic(:,:) = tab_read_flds(:,:,3)  ! Albedo at sea ice (multiplicated by fraction)
340          read_sit(:,:)     = tab_read_flds(:,:,4)  ! Sea ice temperature (multiplicated by fraction)
341       END IF
342!$OMP END MASTER
343
344!*************************************************************************************
345!  Transform seaice fraction (read_sic : ocean-seaice mask) into global
346!  fraction (pctsrf : ocean-seaice-land-landice mask)
347!
348!*************************************************************************************
349       CALL cpl2gath(read_sic, read_sic1D, klon, unity)
350
351       pctsrf_old(:,:) = pctsrf(:,:)
352       DO i = 1, klon
353          ! treatment only of points with ocean and/or seaice
354          IF (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic) > 0.) THEN
355             pctsrf(i,is_sic) = (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic)) &
356                  * read_sic1D(i)
357             pctsrf(i,is_oce) = (pctsrf_old(i,is_oce) + pctsrf_old(i,is_sic)) &
358                  - pctsrf(i,is_sic)
359          ENDIF
360       ENDDO
361
362    END IF ! if time to receive
363
364  END SUBROUTINE cpl_receive_frac
365
366!
367!*************************************************************************************
368!
369
370  SUBROUTINE cpl_receive_ocean_fields(knon, knindex, tsurf_new)
371!
372! This routine returns the field for the ocean that has been read from the coupler
373! (done earlier with cpl_receive_frac). The field is the temperature.
374! The temperature is transformed into 1D array with valid points from index 1 to knon.
375!
376    INCLUDE "indicesol.h"
377
378! Input arguments
379!*************************************************************************************
380    INTEGER, INTENT(IN)                     :: knon
381    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
382
383! Output arguments
384!*************************************************************************************
385    REAL, DIMENSION(klon), INTENT(OUT)      :: tsurf_new
386
387! Local variables
388!*************************************************************************************
389    INTEGER               :: i
390    REAL, DIMENSION(klon) :: sic_new
391
392!*************************************************************************************
393! Transform read_sst into compressed 1D variable tsurf_new
394!
395!*************************************************************************************
396    CALL cpl2gath(read_sst, tsurf_new, knon, knindex)
397    CALL cpl2gath(read_sic, sic_new, knon, knindex)
398
399!*************************************************************************************
400! The fields received from the coupler have to be weighted with the fraction of ocean
401! in relation to the total sea-ice+ocean
402!
403!*************************************************************************************
404    DO i=1, knon
405       tsurf_new(i) = tsurf_new(i)/(1. - sic_new(i))
406    END DO
407
408  END SUBROUTINE cpl_receive_ocean_fields
409
410!
411!*************************************************************************************
412!
413
414  SUBROUTINE cpl_receive_seaice_fields(knon, knindex, &
415       tsurf_new, alb_new)
416!
417! This routine returns the fields for the seaice that have been read from the coupler
418! (done earlier with cpl_receive_frac). These fields are the temperature and
419! albedo at sea ice surface and fraction of sea ice.
420! The fields are transformed into 1D arrays with valid points from index 1 to knon.
421!
422
423! Input arguments
424!*************************************************************************************
425    INTEGER, INTENT(IN)                     :: knon
426    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
427
428! Output arguments
429!*************************************************************************************
430    REAL, DIMENSION(klon), INTENT(OUT)      :: tsurf_new
431    REAL, DIMENSION(klon), INTENT(OUT)      :: alb_new
432
433! Local variables
434!*************************************************************************************
435    INTEGER               :: i
436    REAL, DIMENSION(klon) :: sic_new
437
438!*************************************************************************************
439! Transform fields read from coupler from 2D into compressed 1D variables
440!
441!*************************************************************************************
442    CALL cpl2gath(read_sit, tsurf_new, knon, knindex)
443    CALL cpl2gath(read_alb_sic, alb_new, knon, knindex)
444    CALL cpl2gath(read_sic, sic_new, knon, knindex)
445
446!*************************************************************************************
447! The fields received from the coupler have to be weighted with the sea-ice
448! concentration (in relation to the total sea-ice + ocean).
449!
450!*************************************************************************************
451    DO i= 1, knon
452       tsurf_new(i) = tsurf_new(i) / sic_new(i)
453       alb_new(i)   = alb_new(i)   / sic_new(i)
454    END DO
455
456  END SUBROUTINE cpl_receive_seaice_fields
457
458!
459!*************************************************************************************
460!
461
462  SUBROUTINE cpl_send_ocean_fields(itime, knon, knindex, &
463       swdown, lwdown, fluxlat, fluxsens, &
464       precip_rain, precip_snow, evap, tsurf, fder, albsol, taux, tauy, windsp)
465!
466! This subroutine cumulates some fields for each time-step during a coupling
467! period. At last time-step in a coupling period the fields are transformed to the
468! grid accepted by the coupler. No sending to the coupler will be done from here
469! (it is done in cpl_send_seaice_fields).
470!
471    INCLUDE "indicesol.h"
472    INCLUDE "dimensions.h"
473
474! Input arguments
475!*************************************************************************************
476    INTEGER, INTENT(IN)                     :: itime
477    INTEGER, INTENT(IN)                     :: knon
478    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
479    REAL, DIMENSION(klon), INTENT(IN)       :: swdown, lwdown
480    REAL, DIMENSION(klon), INTENT(IN)       :: fluxlat, fluxsens
481    REAL, DIMENSION(klon), INTENT(IN)       :: precip_rain, precip_snow
482    REAL, DIMENSION(klon), INTENT(IN)       :: evap, tsurf, fder, albsol
483    REAL, DIMENSION(klon), INTENT(IN)       :: taux, tauy, windsp
484
485! Local variables
486!*************************************************************************************
487    INTEGER                                 :: cpl_index, ig
488    INTEGER                                 :: error, sum_error
489    CHARACTER(len = 25)                     :: modname = 'cpl_send_ocean_fields'
490    CHARACTER(len = 80)                     :: abort_message
491
492!*************************************************************************************
493! Start calculation
494! The ocean points are saved with second array index=1
495!
496!*************************************************************************************
497    cpl_index = 1
498
499!*************************************************************************************
500! Reset fields to zero in the beginning of a new coupling period
501!
502!*************************************************************************************
503    IF (MOD(itime, nexca) == 1) THEN
504       cpl_sols(1:knon,cpl_index) = 0.0
505       cpl_nsol(1:knon,cpl_index) = 0.0
506       cpl_rain(1:knon,cpl_index) = 0.0
507       cpl_snow(1:knon,cpl_index) = 0.0
508       cpl_evap(1:knon,cpl_index) = 0.0
509       cpl_tsol(1:knon,cpl_index) = 0.0
510       cpl_fder(1:knon,cpl_index) = 0.0
511       cpl_albe(1:knon,cpl_index) = 0.0
512       cpl_taux(1:knon,cpl_index) = 0.0
513       cpl_tauy(1:knon,cpl_index) = 0.0
514       cpl_windsp(1:knon,cpl_index) = 0.0
515    ENDIF
516       
517!*************************************************************************************
518! Cumulate at each time-step
519!
520!*************************************************************************************   
521    DO ig = 1, knon
522       cpl_sols(ig,cpl_index) = cpl_sols(ig,cpl_index) + &
523            swdown(ig)      / FLOAT(nexca)
524       cpl_nsol(ig,cpl_index) = cpl_nsol(ig,cpl_index) + &
525            (lwdown(ig) + fluxlat(ig) +fluxsens(ig)) / FLOAT(nexca)
526       cpl_rain(ig,cpl_index) = cpl_rain(ig,cpl_index) + &
527            precip_rain(ig) / FLOAT(nexca)
528       cpl_snow(ig,cpl_index) = cpl_snow(ig,cpl_index) + &
529            precip_snow(ig) / FLOAT(nexca)
530       cpl_evap(ig,cpl_index) = cpl_evap(ig,cpl_index) + &
531            evap(ig)        / FLOAT(nexca)
532       cpl_tsol(ig,cpl_index) = cpl_tsol(ig,cpl_index) + &
533            tsurf(ig)       / FLOAT(nexca)
534       cpl_fder(ig,cpl_index) = cpl_fder(ig,cpl_index) + &
535            fder(ig)        / FLOAT(nexca)
536       cpl_albe(ig,cpl_index) = cpl_albe(ig,cpl_index) + &
537            albsol(ig)      / FLOAT(nexca)
538       cpl_taux(ig,cpl_index) = cpl_taux(ig,cpl_index) + &
539            taux(ig)        / FLOAT(nexca)
540       cpl_tauy(ig,cpl_index) = cpl_tauy(ig,cpl_index) + &
541            tauy(ig)        / FLOAT(nexca)     
542       cpl_windsp(ig,cpl_index) = cpl_windsp(ig,cpl_index) + &
543            windsp(ig)      / FLOAT(nexca)
544    ENDDO
545
546!*************************************************************************************
547! If the time-step corresponds to the end of coupling period the
548! fields are transformed to the 2D grid.
549! No sending to the coupler (it is done from cpl_send_seaice_fields).
550!
551!*************************************************************************************
552    IF (MOD(itime, nexca) == 0) THEN
553
554       IF (.NOT. ALLOCATED(cpl_sols2D)) THEN
555          sum_error = 0
556          ALLOCATE(cpl_sols2D(iim,jj_nb,2), stat=error)
557          sum_error = sum_error + error
558          ALLOCATE(cpl_nsol2D(iim,jj_nb,2), stat=error)
559          sum_error = sum_error + error
560          ALLOCATE(cpl_rain2D(iim,jj_nb,2), stat=error)
561          sum_error = sum_error + error
562          ALLOCATE(cpl_snow2D(iim,jj_nb,2), stat=error)
563          sum_error = sum_error + error
564          ALLOCATE(cpl_evap2D(iim,jj_nb,2), stat=error)
565          sum_error = sum_error + error
566          ALLOCATE(cpl_tsol2D(iim,jj_nb,2), stat=error)
567          sum_error = sum_error + error
568          ALLOCATE(cpl_fder2D(iim,jj_nb,2), stat=error)
569          sum_error = sum_error + error
570          ALLOCATE(cpl_albe2D(iim,jj_nb,2), stat=error)
571          sum_error = sum_error + error
572          ALLOCATE(cpl_taux2D(iim,jj_nb,2), stat=error)
573          sum_error = sum_error + error
574          ALLOCATE(cpl_tauy2D(iim,jj_nb,2), stat=error)
575          sum_error = sum_error + error
576          ALLOCATE(cpl_windsp2D(iim,jj_nb), stat=error)
577          sum_error = sum_error + error
578         
579          IF (sum_error /= 0) THEN
580             abort_message='Pb allocation variables couplees pour l''ecriture'
581             CALL abort_gcm(modname,abort_message,1)
582          ENDIF
583       ENDIF
584       
585
586       CALL gath2cpl(cpl_sols(1,cpl_index), cpl_sols2D(1,1,cpl_index), &
587            knon, knindex)
588
589       CALL gath2cpl(cpl_nsol(1,cpl_index), cpl_nsol2D(1,1,cpl_index), &
590            knon, knindex)
591
592       CALL gath2cpl(cpl_rain(1,cpl_index), cpl_rain2D(1,1,cpl_index), &
593            knon, knindex)
594
595       CALL gath2cpl(cpl_snow(1,cpl_index), cpl_snow2D(1,1,cpl_index), &
596            knon, knindex)
597
598       CALL gath2cpl(cpl_evap(1,cpl_index), cpl_evap2D(1,1,cpl_index), &
599            knon, knindex)
600
601! cpl_tsol2D(:,:,:) not used!
602       CALL gath2cpl(cpl_tsol(1,cpl_index), cpl_tsol2D(1,1, cpl_index), &
603            knon, knindex)
604
605! cpl_fder2D(:,:,1) not used, only cpl_fder(:,:,2)!
606       CALL gath2cpl(cpl_fder(1,cpl_index), cpl_fder2D(1,1,cpl_index), &
607            knon, knindex)
608
609! cpl_albe2D(:,:,:) not used!
610       CALL gath2cpl(cpl_albe(1,cpl_index), cpl_albe2D(1,1,cpl_index), &
611            knon, knindex)
612
613       CALL gath2cpl(cpl_taux(1,cpl_index), cpl_taux2D(1,1,cpl_index), &
614            knon, knindex)
615
616       CALL gath2cpl(cpl_tauy(1,cpl_index), cpl_tauy2D(1,1,cpl_index), &
617            knon, knindex)
618
619       CALL gath2cpl(cpl_windsp(1,cpl_index), cpl_windsp2D(1,1), &
620            knon, knindex)
621
622    ENDIF
623
624  END SUBROUTINE cpl_send_ocean_fields
625
626!
627!*************************************************************************************
628!
629
630  SUBROUTINE cpl_send_seaice_fields(itime, dtime, knon, knindex, &
631       pctsrf, lafin, rlon, rlat, &
632       swdown, lwdown, fluxlat, fluxsens, &
633       precip_rain, precip_snow, evap, tsurf, fder, albsol, taux, tauy)
634!
635! This subroutine cumulates some fields for each time-step during a coupling
636! period. At last time-step in a coupling period the fields are transformed to the
637! grid accepted by the coupler. All fields for all types of surfaces are sent to
638! the coupler.
639!
640    INCLUDE "indicesol.h"
641    INCLUDE "dimensions.h"
642
643! Input arguments
644!*************************************************************************************
645    INTEGER, INTENT(IN)                     :: itime
646    INTEGER, INTENT(IN)                     :: knon
647    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
648    REAL, INTENT(IN)                        :: dtime
649    REAL, DIMENSION(klon), INTENT(IN)       :: rlon, rlat
650    REAL, DIMENSION(klon), INTENT(IN)       :: swdown, lwdown
651    REAL, DIMENSION(klon), INTENT(IN)       :: fluxlat, fluxsens
652    REAL, DIMENSION(klon), INTENT(IN)       :: precip_rain, precip_snow
653    REAL, DIMENSION(klon), INTENT(IN)       :: evap, tsurf, fder
654    REAL, DIMENSION(klon), INTENT(IN)       :: albsol, taux, tauy
655    REAL, DIMENSION(klon,nbsrf), INTENT(IN) :: pctsrf
656    LOGICAL, INTENT(IN)                     :: lafin
657
658! Local variables
659!*************************************************************************************
660    INTEGER                                 :: cpl_index, ig
661    INTEGER                                 :: error, sum_error
662    CHARACTER(len = 25)                     :: modname = 'cpl_send_seaice_fields'
663    CHARACTER(len = 80)                     :: abort_message
664
665
666!*************************************************************************************
667! Start calulation
668! The sea-ice points are saved with second array index=2
669!
670!*************************************************************************************
671    cpl_index = 2
672
673!*************************************************************************************
674! Reset fields to zero in the beginning of a new coupling period
675!
676!*************************************************************************************
677    IF (MOD(itime, nexca) == 1) THEN
678       cpl_sols(1:knon,cpl_index) = 0.0
679       cpl_nsol(1:knon,cpl_index) = 0.0
680       cpl_rain(1:knon,cpl_index) = 0.0
681       cpl_snow(1:knon,cpl_index) = 0.0
682       cpl_evap(1:knon,cpl_index) = 0.0
683       cpl_tsol(1:knon,cpl_index) = 0.0
684       cpl_fder(1:knon,cpl_index) = 0.0
685       cpl_albe(1:knon,cpl_index) = 0.0
686       cpl_taux(1:knon,cpl_index) = 0.0
687       cpl_tauy(1:knon,cpl_index) = 0.0
688    ENDIF
689       
690!*************************************************************************************
691! Cumulate at each time-step
692!
693!*************************************************************************************   
694    DO ig = 1, knon
695       cpl_sols(ig,cpl_index) = cpl_sols(ig,cpl_index) + &
696            swdown(ig)      / FLOAT(nexca)
697       cpl_nsol(ig,cpl_index) = cpl_nsol(ig,cpl_index) + &
698            (lwdown(ig) + fluxlat(ig) +fluxsens(ig)) / FLOAT(nexca)
699       cpl_rain(ig,cpl_index) = cpl_rain(ig,cpl_index) + &
700            precip_rain(ig) / FLOAT(nexca)
701       cpl_snow(ig,cpl_index) = cpl_snow(ig,cpl_index) + &
702            precip_snow(ig) / FLOAT(nexca)
703       cpl_evap(ig,cpl_index) = cpl_evap(ig,cpl_index) + &
704            evap(ig)        / FLOAT(nexca)
705       cpl_tsol(ig,cpl_index) = cpl_tsol(ig,cpl_index) + &
706            tsurf(ig)       / FLOAT(nexca)
707       cpl_fder(ig,cpl_index) = cpl_fder(ig,cpl_index) + &
708            fder(ig)        / FLOAT(nexca)
709       cpl_albe(ig,cpl_index) = cpl_albe(ig,cpl_index) + &
710            albsol(ig)      / FLOAT(nexca)
711       cpl_taux(ig,cpl_index) = cpl_taux(ig,cpl_index) + &
712            taux(ig)        / FLOAT(nexca)
713       cpl_tauy(ig,cpl_index) = cpl_tauy(ig,cpl_index) + &
714            tauy(ig)        / FLOAT(nexca)     
715    ENDDO
716
717!*************************************************************************************
718! If the time-step corresponds to the end of coupling period the
719! fields are transformed to the 2D grid and all fields are sent to coupler.
720!
721!*************************************************************************************
722    IF (MOD(itime, nexca) == 0) THEN
723       IF (.NOT. ALLOCATED(cpl_sols2D)) THEN
724          sum_error = 0
725          ALLOCATE(cpl_sols2D(iim,jj_nb,2), stat=error)
726          sum_error = sum_error + error
727          ALLOCATE(cpl_nsol2D(iim,jj_nb,2), stat=error)
728          sum_error = sum_error + error
729          ALLOCATE(cpl_rain2D(iim,jj_nb,2), stat=error)
730          sum_error = sum_error + error
731          ALLOCATE(cpl_snow2D(iim,jj_nb,2), stat=error)
732          sum_error = sum_error + error
733          ALLOCATE(cpl_evap2D(iim,jj_nb,2), stat=error)
734          sum_error = sum_error + error
735          ALLOCATE(cpl_tsol2D(iim,jj_nb,2), stat=error)
736          sum_error = sum_error + error
737          ALLOCATE(cpl_fder2D(iim,jj_nb,2), stat=error)
738          sum_error = sum_error + error
739          ALLOCATE(cpl_albe2D(iim,jj_nb,2), stat=error)
740          sum_error = sum_error + error
741          ALLOCATE(cpl_taux2D(iim,jj_nb,2), stat=error)
742          sum_error = sum_error + error
743          ALLOCATE(cpl_tauy2D(iim,jj_nb,2), stat=error)
744          sum_error = sum_error + error
745          ALLOCATE(cpl_windsp2D(iim,jj_nb), stat=error)
746          sum_error = sum_error + error
747         
748          IF (sum_error /= 0) THEN
749             abort_message='Pb allocation variables couplees pour l''ecriture'
750             CALL abort_gcm(modname,abort_message,1)
751          ENDIF
752       ENDIF
753
754       CALL gath2cpl(cpl_sols(1,cpl_index), cpl_sols2D(1,1,cpl_index), &
755            knon, knindex)
756
757       CALL gath2cpl(cpl_nsol(1,cpl_index), cpl_nsol2D(1,1,cpl_index), &
758            knon, knindex)
759
760       CALL gath2cpl(cpl_rain(1,cpl_index), cpl_rain2D(1,1,cpl_index), &
761            knon, knindex)
762
763       CALL gath2cpl(cpl_snow(1,cpl_index), cpl_snow2D(1,1,cpl_index), &
764            knon, knindex)
765
766       CALL gath2cpl(cpl_evap(1,cpl_index), cpl_evap2D(1,1,cpl_index), &
767            knon, knindex)
768
769! cpl_tsol2D(:,:,:) not used!
770       CALL gath2cpl(cpl_tsol(1,cpl_index), cpl_tsol2D(1,1, cpl_index), &
771            knon, knindex)
772
773       CALL gath2cpl(cpl_fder(1,cpl_index), cpl_fder2D(1,1,cpl_index), &
774            knon, knindex)
775
776! cpl_albe2D(:,:,:) not used!
777       CALL gath2cpl(cpl_albe(1,cpl_index), cpl_albe2D(1,1,cpl_index), &
778            knon, knindex)
779
780       CALL gath2cpl(cpl_taux(1,cpl_index), cpl_taux2D(1,1,cpl_index), &
781            knon, knindex)
782
783       CALL gath2cpl(cpl_tauy(1,cpl_index), cpl_tauy2D(1,1,cpl_index), &
784            knon, knindex)
785
786       ! Send all fields
787       CALL cpl_send_all(itime, dtime, pctsrf, lafin, rlon, rlat)
788    ENDIF
789
790  END SUBROUTINE cpl_send_seaice_fields
791
792!
793!*************************************************************************************
794!
795
796  SUBROUTINE cpl_send_land_fields(itime, knon, knindex, rriv_in, rcoa_in)
797!
798! This subroutine cumulates some fields for each time-step during a coupling
799! period. At last time-step in a coupling period the fields are transformed to the
800! grid accepted by the coupler. No sending to the coupler will be done from here
801! (it is done in cpl_send_seaice_fields).
802!
803    INCLUDE "dimensions.h"
804
805! Input arguments
806!*************************************************************************************
807    INTEGER, INTENT(IN)                       :: itime
808    INTEGER, INTENT(IN)                       :: knon
809    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
810    REAL, DIMENSION(klon), INTENT(IN)         :: rriv_in
811    REAL, DIMENSION(klon), INTENT(IN)         :: rcoa_in
812
813! Local variables
814!*************************************************************************************
815    REAL, DIMENSION(iim,jj_nb)             :: rriv2D
816    REAL, DIMENSION(iim,jj_nb)             :: rcoa2D
817
818!*************************************************************************************
819! Rearrange fields in 2D variables
820! First initialize to zero to avoid unvalid points causing problems
821!
822!*************************************************************************************
823!$OMP MASTER
824    rriv2D(:,:) = 0.0
825    rcoa2D(:,:) = 0.0
826!$OMP END MASTER
827    CALL gath2cpl(rriv_in, rriv2D, knon, knindex)
828    CALL gath2cpl(rcoa_in, rcoa2D, knon, knindex)
829
830!*************************************************************************************
831! Reset cumulated fields to zero in the beginning of a new coupling period
832!
833!*************************************************************************************
834    IF (MOD(itime, nexca) == 1) THEN
835!$OMP MASTER
836       cpl_rriv2D(:,:) = 0.0
837       cpl_rcoa2D(:,:) = 0.0
838!$OMP END MASTER
839    ENDIF
840
841!*************************************************************************************
842! Cumulate : Following fields should be cumulated at each time-step
843!
844!*************************************************************************************   
845!$OMP MASTER
846    cpl_rriv2D(:,:) = cpl_rriv2D(:,:) + rriv2D(:,:) / FLOAT(nexca)
847    cpl_rcoa2D(:,:) = cpl_rcoa2D(:,:) + rcoa2D(:,:) / FLOAT(nexca)
848!$OMP END MASTER
849
850  END SUBROUTINE cpl_send_land_fields
851
852!
853!*************************************************************************************
854!
855
856  SUBROUTINE cpl_send_landice_fields(itime, knon, knindex, rlic_in)
857! This subroutine cumulates the field for melting ice for each time-step
858! during a coupling period. This routine will not send to coupler. Sending
859! will be done in cpl_send_seaice_fields.
860!
861    INCLUDE "dimensions.h"
862
863! Input varibales
864!*************************************************************************************
865    INTEGER, INTENT(IN)                       :: itime
866    INTEGER, INTENT(IN)                       :: knon
867    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
868    REAL, DIMENSION(klon), INTENT(IN)         :: rlic_in
869
870! Local varibales
871!*************************************************************************************
872    REAL, DIMENSION(iim,jj_nb)             :: rlic2D
873
874!*************************************************************************************
875! Rearrange field in a 2D variable
876! First initialize to zero to avoid unvalid points causing problems
877!
878!*************************************************************************************
879!$OMP MASTER
880    rlic2D(:,:) = 0.0
881!$OMP END MASTER
882    CALL gath2cpl(rlic_in, rlic2D, knon, knindex)
883
884!*************************************************************************************
885! Reset field to zero in the beginning of a new coupling period
886!
887!*************************************************************************************
888    IF (MOD(itime, nexca) == 1) THEN
889!$OMP MASTER
890       cpl_rlic2D(:,:) = 0.0
891!$OMP END MASTER
892    ENDIF
893
894!*************************************************************************************
895! Cumulate : Melting ice should be cumulated at each time-step
896!
897!*************************************************************************************   
898!$OMP MASTER
899    cpl_rlic2D(:,:) = cpl_rlic2D(:,:) + rlic2D(:,:) / FLOAT(nexca)
900!$OMP END MASTER
901
902  END SUBROUTINE cpl_send_landice_fields
903
904!
905!*************************************************************************************
906!
907
908  SUBROUTINE cpl_send_all(itime, dtime, pctsrf, lafin, rlon, rlat)
909! This routine will send fields for all different surfaces to the coupler.
910! This subroutine should be executed after calculations by the last surface(sea-ice),
911! all calculations at the different surfaces have to be done before.
912!   
913    USE surface_data
914! Some includes
915!*************************************************************************************
916    INCLUDE "indicesol.h"
917    INCLUDE "temps.h"
918    INCLUDE "dimensions.h"
919   
920! Input arguments
921!*************************************************************************************
922    INTEGER, INTENT(IN)                                  :: itime
923    REAL, INTENT(IN)                                     :: dtime
924    REAL, DIMENSION(klon), INTENT(IN)                    :: rlon, rlat
925    REAL, DIMENSION(klon,nbsrf), INTENT(IN)              :: pctsrf
926    LOGICAL, INTENT(IN)                                  :: lafin
927   
928! Local variables
929!*************************************************************************************
930    INTEGER                                              :: error, sum_error, j
931    INTEGER                                              :: itau_w
932    INTEGER                                              :: time_sec
933    INTEGER, DIMENSION(iim*(jjm+1))                      :: ndexct
934    REAL                                                 :: Up, Down
935    REAL, DIMENSION(iim, jj_nb)                          :: tmp_lon, tmp_lat
936    REAL, DIMENSION(iim, jj_nb, 4)                       :: pctsrf2D
937    REAL, DIMENSION(iim, jj_nb)                          :: deno
938    CHARACTER(len = 20)                                  :: modname = 'cpl_send_all'
939    CHARACTER(len = 80)                                  :: abort_message
940   
941! Variables with fields to coupler
942    REAL, DIMENSION(iim, jj_nb)                          :: tmp_taux
943    REAL, DIMENSION(iim, jj_nb)                          :: tmp_tauy
944    REAL, DIMENSION(iim, jj_nb)                          :: tmp_calv
945! Table with all fields to send to coupler
946    REAL, DIMENSION(iim, jj_nb, jpflda2o1+jpflda2o2)     :: tab_flds
947    REAL, DIMENSION(klon_mpi)                            :: rlon_mpi, rlat_mpi
948
949#ifdef CPP_MPI
950    INCLUDE 'mpif.h'
951    INTEGER, DIMENSION(MPI_STATUS_SIZE)                  :: status
952#endif
953
954! End definitions
955!*************************************************************************************
956   
957
958
959!*************************************************************************************
960! All fields are stored in a table tab_flds(:,:,:)
961! First store the fields 7 to 18 which are already on the right format
962!
963!*************************************************************************************
964!$OMP MASTER
965    IF (version_ocean=='nemo') THEN
966       tab_flds(:,:,7)  = cpl_windsp2D(:,:)
967       tab_flds(:,:,14) = cpl_sols2D(:,:,2)
968       tab_flds(:,:,12) = cpl_sols2D(:,:,1)
969       tab_flds(:,:,15) = cpl_nsol2D(:,:,2)
970       tab_flds(:,:,13) = cpl_nsol2D(:,:,1)
971       tab_flds(:,:,16) = cpl_fder2D(:,:,2)
972       tab_flds(:,:,11) = cpl_evap2D(:,:,2)
973       tab_flds(:,:,18) = cpl_rriv2D(:,:)
974       tab_flds(:,:,19) = cpl_rcoa2D(:,:)
975    ELSE IF (version_ocean=='opa8') THEN
976       tab_flds(:,:,7)  = cpl_windsp2D(:,:)
977       tab_flds(:,:,8)  = cpl_sols2D(:,:,2)
978       tab_flds(:,:,9)  = cpl_sols2D(:,:,1)
979       tab_flds(:,:,10) = cpl_nsol2D(:,:,2)
980       tab_flds(:,:,11) = cpl_nsol2D(:,:,1)
981       tab_flds(:,:,12) = cpl_fder2D(:,:,2)
982       tab_flds(:,:,13) = cpl_evap2D(:,:,2)
983       tab_flds(:,:,14) = cpl_evap2D(:,:,1)
984       tab_flds(:,:,17) = cpl_rcoa2D(:,:)
985       tab_flds(:,:,18) = cpl_rriv2D(:,:)
986    END IF
987   
988!*************************************************************************************
989! Transform the fraction of sub-surfaces from 1D to 2D array
990!
991!*************************************************************************************
992    pctsrf2D(:,:,:) = 0.
993!$OMP END MASTER
994    CALL gath2cpl(pctsrf(:,is_oce), pctsrf2D(:,:,is_oce), klon, unity)
995    CALL gath2cpl(pctsrf(:,is_sic), pctsrf2D(:,:,is_sic), klon, unity)
996    CALL gath2cpl(pctsrf(:,is_lic), pctsrf2D(:,:,is_lic), klon, unity)
997
998!*************************************************************************************
999! Calculate the average calving per latitude
1000! Store calving in tab_flds(:,:,19)
1001!
1002!*************************************************************************************     
1003    IF (is_omp_root) THEN
1004
1005      DO j = 1, jj_nb
1006         tmp_calv(:,j) = DOT_PRODUCT (cpl_rlic2D(1:iim,j), &
1007              pctsrf2D(1:iim,j,is_lic)) / REAL(iim)
1008      ENDDO
1009   
1010   
1011      IF (is_parallel) THEN
1012         IF (.NOT. is_north_pole) THEN
1013#ifdef CPP_MPI
1014            CALL MPI_RECV(Up,1,MPI_REAL_LMDZ,mpi_rank-1,1234,COMM_LMDZ_PHY,status,error)
1015            CALL MPI_SEND(tmp_calv(1,1),1,MPI_REAL_LMDZ,mpi_rank-1,1234,COMM_LMDZ_PHY,error)
1016#endif
1017         ENDIF
1018       
1019         IF (.NOT. is_south_pole) THEN
1020#ifdef CPP_MPI
1021            CALL MPI_SEND(tmp_calv(1,jj_nb),1,MPI_REAL_LMDZ,mpi_rank+1,1234,COMM_LMDZ_PHY,error)
1022            CALL MPI_RECV(down,1,MPI_REAL_LMDZ,mpi_rank+1,1234,COMM_LMDZ_PHY,status,error)
1023#endif
1024         ENDIF
1025         
1026         IF (.NOT. is_north_pole .AND. ii_begin /=1) THEN
1027            Up=Up+tmp_calv(iim,1)
1028            tmp_calv(:,1)=Up
1029         ENDIF
1030         
1031         IF (.NOT. is_south_pole .AND. ii_end /= iim) THEN
1032            Down=Down+tmp_calv(1,jj_nb)
1033            tmp_calv(:,jj_nb)=Down       
1034         ENDIF
1035      ENDIF
1036     
1037      IF (version_ocean=='nemo') THEN
1038         tab_flds(:,:,17) = tmp_calv(:,:)
1039      ELSE IF (version_ocean=='opa8') THEN
1040         tab_flds(:,:,19) = tmp_calv(:,:)
1041      END IF
1042
1043!*************************************************************************************
1044! Calculate total flux for snow, rain and wind with weighted addition using the
1045! fractions of ocean and seaice.
1046!
1047! Store the fields for rain and snow directly in tab_flds(:,:,15) and
1048! tab_flds(:,:,16) respectively.
1049!
1050!*************************************************************************************   
1051       ! fraction oce+seaice
1052       deno =  pctsrf2D(:,:,is_oce) + pctsrf2D(:,:,is_sic)
1053
1054       IF (version_ocean=='nemo') THEN
1055          tab_flds(:,:,10) = 0.0
1056          tmp_taux(:,:)    = 0.0
1057          tmp_tauy(:,:)    = 0.0
1058          ! For all valid grid cells containing some fraction of ocean or sea-ice
1059          WHERE ( deno(:,:) /= 0 )
1060             tab_flds(:,:,10) = cpl_snow2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1061                  cpl_snow2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1062             
1063             tmp_taux = cpl_taux2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1064                  cpl_taux2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1065             tmp_tauy = cpl_tauy2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1066                  cpl_tauy2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1067          ENDWHERE
1068          tab_flds(:,:,8) = (cpl_evap2D(:,:,1) - ( cpl_rain2D(:,:,1) + cpl_snow2D(:,:,1)))
1069          tab_flds(:,:,9) = (cpl_evap2D(:,:,2) - ( cpl_rain2D(:,:,2) + cpl_snow2D(:,:,2)))
1070         
1071       ELSE IF (version_ocean=='opa8') THEN
1072          tab_flds(:,:,15) = 0.0
1073          tab_flds(:,:,16) = 0.0
1074          tmp_taux(:,:)    = 0.0
1075          tmp_tauy(:,:)    = 0.0
1076          ! For all valid grid cells containing some fraction of ocean or sea-ice
1077          WHERE ( deno(:,:) /= 0 )
1078             tab_flds(:,:,15) = cpl_rain2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1079                  cpl_rain2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1080             tab_flds(:,:,16) = cpl_snow2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1081                  cpl_snow2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1082             
1083             tmp_taux = cpl_taux2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1084                  cpl_taux2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1085             tmp_tauy = cpl_tauy2D(:,:,1) * pctsrf2D(:,:,is_oce) / deno(:,:) +    &
1086                  cpl_tauy2D(:,:,2) * pctsrf2D(:,:,is_sic) / deno(:,:)
1087          ENDWHERE
1088       END IF
1089
1090    ENDIF ! is_omp_root
1091 
1092!*************************************************************************************
1093! Transform the wind components from local atmospheric 2D coordinates to geocentric
1094! 3D coordinates.
1095! Store the resulting wind components in tab_flds(:,:,1:6)
1096!*************************************************************************************
1097
1098! Transform the longitudes and latitudes on 2D arrays
1099   
1100    CALL gather_omp(rlon,rlon_mpi)
1101    CALL gather_omp(rlat,rlat_mpi)
1102!$OMP MASTER
1103    CALL Grid1DTo2D_mpi(rlon_mpi,tmp_lon)
1104    CALL Grid1DTo2D_mpi(rlat_mpi,tmp_lat)
1105!$OMP END MASTER   
1106
1107    IF (is_sequential) THEN
1108       IF (is_north_pole) tmp_lon(:,1)     = tmp_lon(:,2)
1109       IF (is_south_pole) tmp_lon(:,jjm+1) = tmp_lon(:,jjm)
1110    ENDIF
1111     
1112! NetCDF output of the wind before transformation of coordinate system
1113    IF (is_sequential) THEN
1114       ndexct(:) = 0
1115       itau_w = itau_phy + itime
1116       CALL histwrite(nidct,'tauxe',itau_w,tmp_taux,iim*(jjm+1),ndexct)
1117       CALL histwrite(nidct,'tauyn',itau_w,tmp_tauy,iim*(jjm+1),ndexct)
1118       CALL histwrite(nidct,'tmp_lon',itau_w,tmp_lon,iim*(jjm+1),ndexct)
1119       CALL histwrite(nidct,'tmp_lat',itau_w,tmp_lat,iim*(jjm+1),ndexct)
1120    ENDIF
1121
1122! Transform the wind from local atmospheric 2D coordinates to geocentric
1123! 3D coordinates
1124!$OMP MASTER
1125    CALL atm2geo (iim, jj_nb, tmp_taux, tmp_tauy, tmp_lon, tmp_lat, &
1126         tab_flds(:,:,1), tab_flds(:,:,2), tab_flds(:,:,3) )
1127   
1128    tab_flds(:,:,4)  = tab_flds(:,:,1)
1129    tab_flds(:,:,5)  = tab_flds(:,:,2)
1130    tab_flds(:,:,6)  = tab_flds(:,:,3)
1131!$OMP END MASTER
1132
1133!*************************************************************************************
1134! NetCDF output of all fields just before sending to coupler.
1135!
1136!*************************************************************************************
1137    IF (is_sequential) THEN
1138       CALL histwrite(nidct,cl_writ(8), itau_w,tab_flds(:,:,8), iim*(jjm+1),ndexct)
1139       CALL histwrite(nidct,cl_writ(9), itau_w,tab_flds(:,:,9), iim*(jjm+1),ndexct)
1140       CALL histwrite(nidct,cl_writ(10),itau_w,tab_flds(:,:,10),iim*(jjm+1),ndexct)
1141       CALL histwrite(nidct,cl_writ(11),itau_w,tab_flds(:,:,11),iim*(jjm+1),ndexct)
1142       CALL histwrite(nidct,cl_writ(12),itau_w,tab_flds(:,:,12),iim*(jjm+1),ndexct)
1143       CALL histwrite(nidct,cl_writ(13),itau_w,tab_flds(:,:,13),iim*(jjm+1),ndexct)
1144       CALL histwrite(nidct,cl_writ(14),itau_w,tab_flds(:,:,14),iim*(jjm+1),ndexct)
1145       CALL histwrite(nidct,cl_writ(15),itau_w,tab_flds(:,:,15),iim*(jjm+1),ndexct)
1146       CALL histwrite(nidct,cl_writ(16),itau_w,tab_flds(:,:,16),iim*(jjm+1),ndexct)
1147       CALL histwrite(nidct,cl_writ(17),itau_w,tab_flds(:,:,17),iim*(jjm+1),ndexct)
1148       CALL histwrite(nidct,cl_writ(18),itau_w,tab_flds(:,:,18),iim*(jjm+1),ndexct)
1149       CALL histwrite(nidct,cl_writ(19),itau_w,tab_flds(:,:,19),iim*(jjm+1),ndexct)
1150       CALL histwrite(nidct,cl_writ(1), itau_w,tab_flds(:,:,1), iim*(jjm+1),ndexct)
1151       CALL histwrite(nidct,cl_writ(2), itau_w,tab_flds(:,:,2), iim*(jjm+1),ndexct)
1152       CALL histwrite(nidct,cl_writ(3), itau_w,tab_flds(:,:,3), iim*(jjm+1),ndexct)
1153       CALL histwrite(nidct,cl_writ(4), itau_w,tab_flds(:,:,4), iim*(jjm+1),ndexct)
1154       CALL histwrite(nidct,cl_writ(5), itau_w,tab_flds(:,:,5), iim*(jjm+1),ndexct)
1155       CALL histwrite(nidct,cl_writ(6), itau_w,tab_flds(:,:,6), iim*(jjm+1),ndexct)
1156       CALL histwrite(nidct,cl_writ(7), itau_w,tab_flds(:,:,7), iim*(jjm+1),ndexct)
1157       CALL histsync(nidct)
1158    ENDIF
1159
1160
1161!*************************************************************************************
1162! Send the table of all fields
1163!
1164!*************************************************************************************
1165    time_sec=(itime-1)*dtime
1166#ifdef CPP_COUPLE
1167    time_sec=(itime-1)*dtime
1168!$OMP MASTER
1169    CALL intocpl(time_sec, lafin, tab_flds(:,:,:))
1170!$OMP END MASTER
1171#endif
1172
1173!*************************************************************************************
1174! Finish with some dellocate
1175!
1176!************************************************************************************* 
1177    sum_error=0
1178    DEALLOCATE(cpl_sols2D, cpl_nsol2D, cpl_rain2D, cpl_snow2D, stat=error )
1179    sum_error = sum_error + error
1180    DEALLOCATE(cpl_evap2D, cpl_tsol2D, cpl_fder2D, cpl_albe2D, stat=error )
1181    sum_error = sum_error + error
1182    DEALLOCATE(cpl_taux2D, cpl_tauy2D, cpl_windsp2D, stat=error )
1183    sum_error = sum_error + error
1184    IF (sum_error /= 0) THEN
1185       abort_message='Pb in deallocation of cpl_xxxx2D coupling variables'
1186       CALL abort_gcm(modname,abort_message,1)
1187    ENDIF
1188   
1189  END SUBROUTINE cpl_send_all
1190!
1191!*************************************************************************************
1192!
1193  SUBROUTINE cpl2gath(champ_in, champ_out, knon, knindex)
1194  USE mod_phys_lmdz_para
1195! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
1196! au coupleur.
1197!
1198!
1199! input:         
1200!   champ_in     champ sur la grille gathere       
1201!   knon         nombre de points dans le domaine a traiter
1202!   knindex      index des points de la surface a traiter
1203!
1204! output:
1205!   champ_out    champ sur la grille 2D
1206!
1207    INCLUDE "dimensions.h"
1208
1209! Input
1210    INTEGER, INTENT(IN)                       :: knon
1211    REAL, DIMENSION(iim,jj_nb), INTENT(IN)    :: champ_in
1212    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
1213
1214! Output
1215    REAL, DIMENSION(klon_mpi), INTENT(OUT)        :: champ_out
1216
1217! Local
1218    INTEGER                                   :: i, ig
1219    REAL, DIMENSION(klon_mpi)                 :: temp_mpi
1220    REAL, DIMENSION(klon)                     :: temp_omp
1221
1222!*************************************************************************************
1223!
1224   
1225
1226! Transform from 2 dimensions (iim,jj_nb) to 1 dimension (klon)
1227!$OMP MASTER
1228    CALL Grid2Dto1D_mpi(champ_in,temp_mpi)
1229!$OMP END MASTER
1230
1231    CALL scatter_omp(temp_mpi,temp_omp)
1232   
1233! Compress from klon to knon
1234    DO i = 1, knon
1235       ig = knindex(i)
1236       champ_out(i) = temp_omp(ig)
1237    ENDDO
1238
1239  END SUBROUTINE cpl2gath
1240!
1241!*************************************************************************************
1242!
1243  SUBROUTINE gath2cpl(champ_in, champ_out, knon, knindex)
1244  USE mod_phys_lmdz_para
1245! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
1246! au coupleur.
1247!
1248! input:         
1249!   champ_in     champ sur la grille gathere       
1250!   knon         nombre de points dans le domaine a traiter
1251!   knindex      index des points de la surface a traiter
1252!
1253! output:
1254!   champ_out    champ sur la grille 2D
1255!
1256    INCLUDE "dimensions.h"
1257   
1258! Input arguments
1259!*************************************************************************************
1260    INTEGER, INTENT(IN)                    :: knon
1261    REAL, DIMENSION(klon), INTENT(IN)      :: champ_in
1262    INTEGER, DIMENSION(klon), INTENT(IN)   :: knindex
1263
1264! Output arguments
1265!*************************************************************************************
1266    REAL, DIMENSION(iim,jj_nb), INTENT(OUT) :: champ_out
1267
1268! Local variables
1269!*************************************************************************************
1270    INTEGER                                :: i, ig
1271    REAL, DIMENSION(klon)                  :: temp_omp
1272    REAL, DIMENSION(klon_mpi)              :: temp_mpi
1273!*************************************************************************************
1274
1275! Decompress from knon to klon
1276    temp_omp = 0.
1277    DO i = 1, knon
1278       ig = knindex(i)
1279       temp_omp(ig) = champ_in(i)
1280    ENDDO
1281
1282! Transform from 1 dimension (klon) to 2 dimensions (iim,jj_nb)
1283    CALL gather_omp(temp_omp,temp_mpi)
1284
1285!$OMP MASTER   
1286    CALL Grid1Dto2D_mpi(temp_mpi,champ_out)
1287   
1288    IF (is_north_pole) champ_out(:,1)=temp_mpi(1)
1289    IF (is_south_pole) champ_out(:,jj_nb)=temp_mpi(klon)
1290!$OMP END MASTER
1291   
1292  END SUBROUTINE gath2cpl
1293!
1294!*************************************************************************************
1295!
1296END MODULE cpl_mod
1297
Note: See TracBrowser for help on using the repository browser.