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

Last change on this file since 1134 was 1067, checked in by Laurent Fairhead, 16 years ago
  • Modifications lie au premier niveau du modele pour la diffusion turbulent

du vent.

  • Preparation pour un couplage des courrant oceaniques.

JG

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