source: LMDZ6/branches/Ocean_skin/libf/phylmd/surf_land_orchidee_noopenmp_mod.F90 @ 3627

Last change on this file since 3627 was 2952, checked in by Laurent Fairhead, 7 years ago

Parametrization of drag by copses
Need version 4465 of ORCHIDEE at least

  1. Cheruy
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 27.5 KB
Line 
1!
2! $Header$
3!
4MODULE surf_land_orchidee_noopenmp_mod
5!
6! This module is compiled only if CPP key ORCHIDEE_NOOPENMP is defined.
7! This module should be used with ORCHIDEE sequentiel or parallele MPI version
8! (not MPI-OpenMP mixte) until revision 1077 in the ORCHIDEE trunk.
9
10#ifdef ORCHIDEE_NOOPENMP
11!
12! This module controles the interface towards the model ORCHIDEE
13!
14! Subroutines in this module : surf_land_orchidee
15!                              Init_orchidee_index
16!                              Get_orchidee_communicator
17!                              Init_neighbours
18  USE dimphy
19#ifdef CPP_VEGET
20  USE intersurf     ! module d'ORCHIDEE
21#endif
22  USE cpl_mod,      ONLY : cpl_send_land_fields
23  USE surface_data, ONLY : type_ocean
24  USE geometry_mod, ONLY : dx, dy
25  USE mod_grid_phy_lmdz
26  USE mod_phys_lmdz_para
27
28  IMPLICIT NONE
29
30  PRIVATE
31  PUBLIC  :: surf_land_orchidee
32
33CONTAINS
34!
35!****************************************************************************************
36
37  SUBROUTINE surf_land_orchidee(itime, dtime, date0, knon, &
38       knindex, rlon, rlat, yrmu0, pctsrf, &
39       debut, lafin, &
40       plev,  u1_lay, v1_lay, gustiness, temp_air, spechum, epot_air, ccanopy, &
41       tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
42       precip_rain, precip_snow, lwdown, swnet, swdown, &
43       ps, q2m, t2m, &
44       evap, fluxsens, fluxlat, &             
45       tsol_rad, tsurf_new, alb1_new, alb2_new, &
46       emis_new, z0_new, z0h_new, qsurf, &
47       veget, lai, height)
48!   
49! Cette routine sert d'interface entre le modele atmospherique et le
50! modele de sol continental. Appel a sechiba
51!
52! L. Fairhead 02/2000
53!
54! input:
55!   itime        numero du pas de temps
56!   dtime        pas de temps de la physique (en s)
57!   nisurf       index de la surface a traiter (1 = sol continental)
58!   knon         nombre de points de la surface a traiter
59!   knindex      index des points de la surface a traiter
60!   rlon         longitudes de la grille entiere
61!   rlat         latitudes de la grille entiere
62!   pctsrf       tableau des fractions de surface de chaque maille
63!   debut        logical: 1er appel a la physique (lire les restart)
64!   lafin        logical: dernier appel a la physique (ecrire les restart)
65!                     (si false calcul simplifie des fluxs sur les continents)
66!   plev         hauteur de la premiere couche (Pa)     
67!   u1_lay       vitesse u 1ere couche
68!   v1_lay       vitesse v 1ere couche
69!   temp_air     temperature de l'air 1ere couche
70!   spechum      humidite specifique 1ere couche
71!   epot_air     temp pot de l'air
72!   ccanopy      concentration CO2 canopee, correspond au co2_send de
73!                carbon_cycle_mod ou valeur constant co2_ppm
74!   tq_cdrag     cdrag
75!   petAcoef     coeff. A de la resolution de la CL pour t
76!   peqAcoef     coeff. A de la resolution de la CL pour q
77!   petBcoef     coeff. B de la resolution de la CL pour t
78!   peqBcoef     coeff. B de la resolution de la CL pour q
79!   precip_rain  precipitation liquide
80!   precip_snow  precipitation solide
81!   lwdown       flux IR descendant a la surface
82!   swnet        flux solaire net
83!   swdown       flux solaire entrant a la surface
84!   ps           pression au sol
85!   radsol       rayonnement net aus sol (LW + SW)
86!   
87!
88! output:
89!   evap         evaporation totale
90!   fluxsens     flux de chaleur sensible
91!   fluxlat      flux de chaleur latente
92!   tsol_rad     
93!   tsurf_new    temperature au sol
94!   alb1_new     albedo in visible SW interval
95!   alb2_new     albedo in near IR interval
96!   emis_new     emissivite
97!   z0_new       surface roughness
98!   z0h_new      surface roughness, it is the same as z0_new
99!   qsurf        air moisture at surface
100!
101    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl, fco2_land_inst, fco2_lu_inst
102    USE indice_sol_mod
103    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
104    USE print_control_mod, ONLY: lunout
105#ifdef CPP_VEGET
106    USE time_phylmdz_mod, ONLY: itau_phy
107#endif
108    IMPLICIT NONE
109
110    INCLUDE "YOMCST.h"
111    INCLUDE "dimpft.h" 
112!
113! Parametres d'entree
114!****************************************************************************************
115    INTEGER, INTENT(IN)                       :: itime
116    REAL, INTENT(IN)                          :: dtime
117    REAL, INTENT(IN)                          :: date0
118    INTEGER, INTENT(IN)                       :: knon
119    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
120    LOGICAL, INTENT(IN)                       :: debut, lafin
121    REAL, DIMENSION(klon,nbsrf), INTENT(IN)   :: pctsrf
122    REAL, DIMENSION(klon), INTENT(IN)         :: rlon, rlat
123    REAL, DIMENSION(klon), INTENT(IN)         :: yrmu0
124    REAL, DIMENSION(klon), INTENT(IN)         :: plev
125    REAL, DIMENSION(klon), INTENT(IN)         :: u1_lay, v1_lay, gustiness
126    REAL, DIMENSION(klon), INTENT(IN)         :: temp_air, spechum
127    REAL, DIMENSION(klon), INTENT(IN)         :: epot_air, ccanopy
128    REAL, DIMENSION(klon), INTENT(IN)         :: tq_cdrag
129    REAL, DIMENSION(klon), INTENT(IN)         :: petAcoef, peqAcoef
130    REAL, DIMENSION(klon), INTENT(IN)         :: petBcoef, peqBcoef
131    REAL, DIMENSION(klon), INTENT(IN)         :: precip_rain, precip_snow
132    REAL, DIMENSION(klon), INTENT(IN)         :: lwdown, swnet, swdown, ps
133    REAL, DIMENSION(klon), INTENT(IN)         :: q2m, t2m
134
135! Parametres de sortie
136!****************************************************************************************
137    REAL, DIMENSION(klon), INTENT(OUT)        :: evap, fluxsens, fluxlat, qsurf
138    REAL, DIMENSION(klon), INTENT(OUT)        :: tsol_rad, tsurf_new
139    REAL, DIMENSION(klon), INTENT(OUT)        :: alb1_new, alb2_new
140    REAL, DIMENSION(klon), INTENT(OUT)        :: emis_new, z0_new, z0h_new
141    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: veget ! dummy variables
142    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: lai   ! dummy variables
143    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: height! dummy variables
144
145! Local
146!****************************************************************************************
147    INTEGER                                   :: ij, jj, igrid, ireal, index
148    INTEGER                                   :: error
149    INTEGER, SAVE                             :: nb_fields_cpl ! number of fields for the climate-carbon coupling (between ATM and ORCHIDEE).
150    !$OMP THREADPRIVATE(nb_fields_cpl)
151    REAL, SAVE, ALLOCATABLE, DIMENSION(:,:)   :: fields_cpl    ! Fluxes for the climate-carbon coupling
152    !$OMP THREADPRIVATE(fields_cpl)
153    REAL, DIMENSION(klon)                     :: swdown_vrai
154    CHARACTER (len = 20)                      :: modname = 'surf_land_orchidee'
155    CHARACTER (len = 80)                      :: abort_message
156    LOGICAL,SAVE                              :: check = .FALSE.
157    !$OMP THREADPRIVATE(check)
158
159! type de couplage dans sechiba
160!  character (len=10)   :: coupling = 'implicit'
161! drapeaux controlant les appels dans SECHIBA
162!  type(control_type), save   :: control_in
163! Preserved albedo
164    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: albedo_keep, zlev
165    !$OMP THREADPRIVATE(albedo_keep,zlev)
166! coordonnees geographiques
167    REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: lalo
168    !$OMP THREADPRIVATE(lalo)
169! pts voisins
170    INTEGER,ALLOCATABLE, DIMENSION(:,:), SAVE :: neighbours
171    !$OMP THREADPRIVATE(neighbours)
172! fractions continents
173    REAL,ALLOCATABLE, DIMENSION(:), SAVE      :: contfrac
174    !$OMP THREADPRIVATE(contfrac)
175! resolution de la grille
176    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: resolution
177    !$OMP THREADPRIVATE(resolution)
178
179    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: lon_scat, lat_scat 
180    !$OMP THREADPRIVATE(lon_scat,lat_scat)
181
182    LOGICAL, SAVE                             :: lrestart_read = .TRUE.
183    !$OMP THREADPRIVATE(lrestart_read)
184    LOGICAL, SAVE                             :: lrestart_write = .FALSE.
185    !$OMP THREADPRIVATE(lrestart_write)
186
187    REAL, DIMENSION(knon,2)                   :: albedo_out
188    !$OMP THREADPRIVATE(albedo_out)
189
190! Pb de nomenclature
191    REAL, DIMENSION(klon)                     :: petA_orc, peqA_orc
192    REAL, DIMENSION(klon)                     :: petB_orc, peqB_orc
193! Pb de correspondances de grilles
194    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  :: ig, jg
195    !$OMP THREADPRIVATE(ig,jg)
196    INTEGER :: indi, indj
197    INTEGER, SAVE, ALLOCATABLE,DIMENSION(:)   :: ktindex
198    !$OMP THREADPRIVATE(ktindex)
199
200! Essai cdrag
201    REAL, DIMENSION(klon)                     :: cdrag
202    INTEGER,SAVE                              :: offset
203    !$OMP THREADPRIVATE(offset)
204
205    REAL, DIMENSION(klon_glo)                 :: rlon_g,rlat_g
206    INTEGER, SAVE                             :: orch_comm
207    !$OMP THREADPRIVATE(orch_comm)
208
209    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: coastalflow
210    !$OMP THREADPRIVATE(coastalflow)
211    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: riverflow
212    !$OMP THREADPRIVATE(riverflow)
213!
214! Fin definition
215!****************************************************************************************
216#ifdef CPP_VEGET
217
218    IF (check) WRITE(lunout,*)'Entree ', modname
219 
220! Initialisation
221 
222    IF (debut) THEN
223! Test de coherence
224#ifndef ORCH_NEW
225       ! Compilation avec orchidee nouvelle version necessaire avec carbon_cycle_cpl=y
226       IF (carbon_cycle_cpl) THEN
227          abort_message='You must define preprossing key ORCH_NEW when running carbon_cycle_cpl=y'
228          CALL abort_physic(modname,abort_message,1)
229       END IF
230#endif
231       ALLOCATE(ktindex(knon))
232       IF ( .NOT. ALLOCATED(albedo_keep)) THEN
233          ALLOCATE(albedo_keep(klon))
234          ALLOCATE(zlev(knon))
235       ENDIF
236! Pb de correspondances de grilles
237       ALLOCATE(ig(klon))
238       ALLOCATE(jg(klon))
239       ig(1) = 1
240       jg(1) = 1
241       indi = 0
242       indj = 2
243       DO igrid = 2, klon - 1
244          indi = indi + 1
245          IF ( indi > nbp_lon) THEN
246             indi = 1
247             indj = indj + 1
248          ENDIF
249          ig(igrid) = indi
250          jg(igrid) = indj
251       ENDDO
252       ig(klon) = 1
253       jg(klon) = nbp_lat
254
255       IF ((.NOT. ALLOCATED(lalo))) THEN
256          ALLOCATE(lalo(knon,2), stat = error)
257          IF (error /= 0) THEN
258             abort_message='Pb allocation lalo'
259             CALL abort_physic(modname,abort_message,1)
260          ENDIF
261       ENDIF
262       IF ((.NOT. ALLOCATED(lon_scat))) THEN
263          ALLOCATE(lon_scat(nbp_lon,nbp_lat), stat = error)
264          IF (error /= 0) THEN
265             abort_message='Pb allocation lon_scat'
266             CALL abort_physic(modname,abort_message,1)
267          ENDIF
268       ENDIF
269       IF ((.NOT. ALLOCATED(lat_scat))) THEN
270          ALLOCATE(lat_scat(nbp_lon,nbp_lat), stat = error)
271          IF (error /= 0) THEN
272             abort_message='Pb allocation lat_scat'
273             CALL abort_physic(modname,abort_message,1)
274          ENDIF
275       ENDIF
276       lon_scat = 0.
277       lat_scat = 0.
278       DO igrid = 1, knon
279          index = knindex(igrid)
280          lalo(igrid,2) = rlon(index)
281          lalo(igrid,1) = rlat(index)
282       ENDDO
283
284       
285       
286       CALL Gather(rlon,rlon_g)
287       CALL Gather(rlat,rlat_g)
288
289       IF (is_mpi_root) THEN
290          index = 1
291          DO jj = 2, nbp_lat-1
292             DO ij = 1, nbp_lon
293                index = index + 1
294                lon_scat(ij,jj) = rlon_g(index)
295                lat_scat(ij,jj) = rlat_g(index)
296             ENDDO
297          ENDDO
298          lon_scat(:,1) = lon_scat(:,2)
299          lat_scat(:,1) = rlat_g(1)
300          lon_scat(:,nbp_lat) = lon_scat(:,2)
301          lat_scat(:,nbp_lat) = rlat_g(klon_glo)
302       ENDIF
303
304       CALL bcast(lon_scat)
305       CALL bcast(lat_scat)
306
307!
308! Allouer et initialiser le tableau des voisins et des fraction de continents
309!
310       IF ( (.NOT.ALLOCATED(neighbours))) THEN
311          ALLOCATE(neighbours(knon,8), stat = error)
312          IF (error /= 0) THEN
313             abort_message='Pb allocation neighbours'
314             CALL abort_physic(modname,abort_message,1)
315          ENDIF
316       ENDIF
317       neighbours = -1.
318       IF (( .NOT. ALLOCATED(contfrac))) THEN
319          ALLOCATE(contfrac(knon), stat = error)
320          IF (error /= 0) THEN
321             abort_message='Pb allocation contfrac'
322             CALL abort_physic(modname,abort_message,1)
323          ENDIF
324       ENDIF
325
326       DO igrid = 1, knon
327          ireal = knindex(igrid)
328          contfrac(igrid) = pctsrf(ireal,is_ter)
329       ENDDO
330
331
332       CALL Init_neighbours(knon,neighbours,knindex,pctsrf(:,is_ter))
333
334!
335!  Allocation et calcul resolutions
336       IF ( (.NOT.ALLOCATED(resolution))) THEN
337          ALLOCATE(resolution(knon,2), stat = error)
338          IF (error /= 0) THEN
339             abort_message='Pb allocation resolution'
340             CALL abort_physic(modname,abort_message,1)
341          ENDIF
342       ENDIF
343       DO igrid = 1, knon
344          ij = knindex(igrid)
345          resolution(igrid,1) = dx(ij)
346          resolution(igrid,2) = dy(ij)
347       ENDDO
348     
349       ALLOCATE(coastalflow(klon), stat = error)
350       IF (error /= 0) THEN
351          abort_message='Pb allocation coastalflow'
352          CALL abort_physic(modname,abort_message,1)
353       ENDIF
354       
355       ALLOCATE(riverflow(klon), stat = error)
356       IF (error /= 0) THEN
357          abort_message='Pb allocation riverflow'
358          CALL abort_physic(modname,abort_message,1)
359       ENDIF
360
361!
362! Allocate variables needed for carbon_cycle_mod
363       IF ( carbon_cycle_cpl ) THEN
364          nb_fields_cpl=2
365       ELSE
366          nb_fields_cpl=1
367       END IF
368
369
370       IF (carbon_cycle_cpl) THEN
371          ALLOCATE(fco2_land_inst(klon),stat=error)
372          IF (error /= 0)  CALL abort_physic(modname,'Pb in allocation fco2_land_inst',1)
373         
374          ALLOCATE(fco2_lu_inst(klon),stat=error)
375          IF(error /=0) CALL abort_physic(modname,'Pb in allocation fco2_lu_inst',1)
376       END IF
377
378       ALLOCATE(fields_cpl(klon,nb_fields_cpl), stat = error)
379       IF (error /= 0) CALL abort_physic(modname,'Pb in allocation fields_cpl',1)
380
381    ENDIF                          ! (fin debut)
382
383!
384! Appel a la routine sols continentaux
385!
386    IF (lafin) lrestart_write = .TRUE.
387    IF (check) WRITE(lunout,*)'lafin ',lafin,lrestart_write
388   
389    petA_orc(1:knon) = petBcoef(1:knon) * dtime
390    petB_orc(1:knon) = petAcoef(1:knon)
391    peqA_orc(1:knon) = peqBcoef(1:knon) * dtime
392    peqB_orc(1:knon) = peqAcoef(1:knon)
393
394    cdrag = 0.
395    cdrag(1:knon) = tq_cdrag(1:knon)
396
397! zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/287.05*temp_air(1:knon))*9.80665)
398    zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/RD*temp_air(1:knon))*RG)
399
400
401! PF et PASB
402!   where(cdrag > 0.01)
403!     cdrag = 0.01
404!   endwhere
405!  write(*,*)'Cdrag = ',minval(cdrag),maxval(cdrag)
406
407!
408! Init Orchidee
409!
410!  if (pole_nord) then
411!    offset=0
412!    ktindex(:)=ktindex(:)+nbp_lon-1
413!  else
414!    offset = klon_mpi_begin-1+nbp_lon-1
415!    ktindex(:)=ktindex(:)+MOD(offset,nbp_lon)
416!    offset=offset-MOD(offset,nbp_lon)
417!  endif
418 
419    IF (debut) THEN
420       CALL Get_orchidee_communicator(knon,orch_comm)
421       IF (knon /=0) THEN
422          CALL Init_orchidee_index(knon,orch_comm,knindex,offset,ktindex)
423
424#ifndef CPP_MPI
425          ! Interface for ORCHIDEE compiled in sequential mode(without preprocessing flag CPP_MPI)
426          CALL intersurf_main (itime+itau_phy-1, nbp_lon, nbp_lat, knon, ktindex, dtime, &
427               lrestart_read, lrestart_write, lalo, &
428               contfrac, neighbours, resolution, date0, &
429               zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
430               cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
431               precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
432               evap, fluxsens, fluxlat, coastalflow, riverflow, &
433               tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
434               lon_scat, lat_scat, q2m, t2m &
435#ifdef ORCH_NEW
436               , nb_fields_cpl, fields_cpl)
437#else
438               )
439#endif
440
441#else         
442          ! Interface for ORCHIDEE version 1.9 or later(1.9.2, 1.9.3, 1.9.4, 1.9.5) compiled in parallel mode(with preprocessing flag CPP_MPI)
443          CALL intersurf_main (itime+itau_phy-1, nbp_lon, nbp_lat, offset, knon, ktindex, &
444               orch_comm, dtime, lrestart_read, lrestart_write, lalo, &
445               contfrac, neighbours, resolution, date0, &
446               zlev,  u1_lay(1:knon), v1_lay(1:knon), spechum(1:knon), temp_air(1:knon), epot_air(1:knon), ccanopy(1:knon), &
447               cdrag(1:knon), petA_orc(1:knon), peqA_orc(1:knon), petB_orc(1:knon), peqB_orc(1:knon), &
448               precip_rain(1:knon), precip_snow(1:knon), lwdown(1:knon), swnet(1:knon), swdown(1:knon), ps(1:knon), &
449               evap(1:knon), fluxsens(1:knon), fluxlat(1:knon), coastalflow(1:knon), riverflow(1:knon), &
450               tsol_rad(1:knon), tsurf_new(1:knon), qsurf(1:knon), albedo_out(1:knon,:), emis_new(1:knon), z0_new(1:knon), &
451               lon_scat, lat_scat, q2m, t2m &
452#ifdef ORCH_NEW
453               , nb_fields_cpl, fields_cpl(1:knon,:))
454#else
455               )
456#endif
457#endif
458         
459       ENDIF
460
461       albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
462
463    ENDIF
464
465!  swdown_vrai(1:knon) = swnet(1:knon)/(1. - albedo_keep(1:knon))
466    swdown_vrai(1:knon) = swdown(1:knon)
467
468    IF (knon /=0) THEN
469#ifndef CPP_MPI
470       ! Interface for ORCHIDEE compiled in sequential mode(without preprocessing flag CPP_MPI)
471       CALL intersurf_main (itime+itau_phy, nbp_lon, nbp_lat, knon, ktindex, dtime, &
472            lrestart_read, lrestart_write, lalo, &
473            contfrac, neighbours, resolution, date0, &
474            zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
475            cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
476            precip_rain, precip_snow, lwdown, swnet, swdown_vrai, ps, &
477            evap, fluxsens, fluxlat, coastalflow, riverflow, &
478            tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
479            lon_scat, lat_scat, q2m, t2m &
480#ifdef ORCH_NEW
481            , nb_fields_cpl, fields_cpl)
482#else
483            )
484#endif
485#else
486       ! Interface for ORCHIDEE version 1.9 or later compiled in parallel mode(with preprocessing flag CPP_MPI)
487       CALL intersurf_main (itime+itau_phy, nbp_lon, nbp_lat,offset, knon, ktindex, &
488            orch_comm,dtime, lrestart_read, lrestart_write, lalo, &
489            contfrac, neighbours, resolution, date0, &
490            zlev,  u1_lay(1:knon), v1_lay(1:knon), spechum(1:knon), temp_air(1:knon), epot_air(1:knon), ccanopy(1:knon), &
491            cdrag(1:knon), petA_orc(1:knon), peqA_orc(1:knon), petB_orc(1:knon), peqB_orc(1:knon), &
492            precip_rain(1:knon), precip_snow(1:knon), lwdown(1:knon), swnet(1:knon), swdown_vrai(1:knon), ps(1:knon), &
493            evap(1:knon), fluxsens(1:knon), fluxlat(1:knon), coastalflow(1:knon), riverflow(1:knon), &
494            tsol_rad(1:knon), tsurf_new(1:knon), qsurf(1:knon), albedo_out(1:knon,:), emis_new(1:knon), z0_new(1:knon), &
495            lon_scat, lat_scat, q2m, t2m &
496#ifdef ORCH_NEW
497            , nb_fields_cpl, fields_cpl(1:knon,:))
498#else
499            )
500#endif
501#endif
502    ENDIF
503
504    albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
505   
506    ! ORCHIDEE only gives one value for z0_new. Copy it into z0h_new.
507    z0h_new(:)=z0_new(:)
508
509!* Send to coupler
510!
511    IF (type_ocean=='couple') THEN
512       CALL cpl_send_land_fields(itime, knon, knindex, &
513            riverflow, coastalflow)
514    ENDIF
515
516    alb1_new(1:knon) = albedo_out(1:knon,1)
517    alb2_new(1:knon) = albedo_out(1:knon,2)
518
519! Convention orchidee: positif vers le haut
520    fluxsens(1:knon) = -1. * fluxsens(1:knon)
521    fluxlat(1:knon)  = -1. * fluxlat(1:knon)
522   
523!  evap     = -1. * evap
524
525    IF (debut) lrestart_read = .FALSE.
526
527! Decompress variables for the module carbon_cycle_mod
528    IF (carbon_cycle_cpl) THEN
529       fco2_land_inst(:)=0.
530       fco2_lu_inst(:)=0.
531       
532       DO igrid = 1, knon
533          ireal = knindex(igrid)
534          fco2_land_inst(ireal) = fields_cpl(igrid,1)
535          fco2_lu_inst(ireal)   = fields_cpl(igrid,2)
536       END DO
537    END IF
538
539#endif   
540  END SUBROUTINE surf_land_orchidee
541!
542!****************************************************************************************
543!
544  SUBROUTINE Init_orchidee_index(knon,orch_comm,knindex,offset,ktindex)
545   
546    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
547
548#ifdef CPP_MPI
549    INCLUDE 'mpif.h'
550#endif   
551
552
553! Input arguments
554!****************************************************************************************
555    INTEGER, INTENT(IN)                   :: knon
556    INTEGER, INTENT(IN)                   :: orch_comm
557    INTEGER, DIMENSION(klon), INTENT(IN)  :: knindex
558
559! Output arguments
560!****************************************************************************************
561    INTEGER, INTENT(OUT)                  :: offset
562    INTEGER, DIMENSION(knon), INTENT(OUT) :: ktindex
563
564! Local varables
565!****************************************************************************************
566#ifdef CPP_MPI
567    INTEGER, DIMENSION(MPI_STATUS_SIZE)   :: status
568#endif
569
570    INTEGER                               :: MyLastPoint
571    INTEGER                               :: LastPoint
572    INTEGER                               :: mpi_rank_orch
573    INTEGER                               :: mpi_size_orch
574    INTEGER                               :: ierr
575!
576! End definition
577!****************************************************************************************
578
579    MyLastPoint=klon_mpi_begin-1+knindex(knon)+nbp_lon-1
580   
581    IF (is_parallel) THEN
582#ifdef CPP_MPI   
583       CALL MPI_COMM_SIZE(orch_comm,mpi_size_orch,ierr)
584       CALL MPI_COMM_RANK(orch_comm,mpi_rank_orch,ierr)
585#endif
586    ELSE
587       mpi_rank_orch=0
588       mpi_size_orch=1
589    ENDIF
590
591    IF (is_parallel) THEN
592       IF (mpi_rank_orch /= 0) THEN
593#ifdef CPP_MPI
594          CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank_orch-1,1234,orch_comm,status,ierr)
595#endif
596       ENDIF
597       
598       IF (mpi_rank_orch /= mpi_size_orch-1) THEN
599#ifdef CPP_MPI
600          CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank_orch+1,1234,orch_comm,ierr) 
601#endif
602       ENDIF
603    ENDIF
604   
605    IF (mpi_rank_orch == 0) THEN
606       offset=0
607    ELSE
608       offset=LastPoint-MOD(LastPoint,nbp_lon)
609    ENDIF
610   
611    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin+nbp_lon-1)-offset-1
612   
613
614  END SUBROUTINE  Init_orchidee_index
615!
616!****************************************************************************************
617!
618  SUBROUTINE Get_orchidee_communicator(knon,orch_comm)
619   
620#ifdef CPP_MPI
621    INCLUDE 'mpif.h'
622#endif   
623
624
625    INTEGER,INTENT(IN)  :: knon
626    INTEGER,INTENT(OUT) :: orch_comm
627   
628    INTEGER             :: color
629    INTEGER             :: ierr
630!
631! End definition
632!****************************************************************************************
633
634    IF (knon==0) THEN
635       color = 0
636    ELSE
637       color = 1
638    ENDIF
639   
640#ifdef CPP_MPI   
641    CALL MPI_COMM_SPLIT(COMM_LMDZ_PHY,color,mpi_rank,orch_comm,ierr)
642#endif
643   
644  END SUBROUTINE Get_orchidee_communicator
645!
646!****************************************************************************************
647
648  SUBROUTINE Init_neighbours(knon,neighbours,ktindex,pctsrf)
649   
650    USE indice_sol_mod
651    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
652
653#ifdef CPP_MPI
654    INCLUDE 'mpif.h'
655#endif   
656
657! Input arguments
658!****************************************************************************************
659    INTEGER, INTENT(IN)                     :: knon
660    INTEGER, DIMENSION(klon), INTENT(IN)    :: ktindex
661    REAL, DIMENSION(klon), INTENT(IN)       :: pctsrf
662   
663! Output arguments
664!****************************************************************************************
665    INTEGER, DIMENSION(knon,8), INTENT(OUT) :: neighbours
666
667! Local variables
668!****************************************************************************************
669    INTEGER                              :: knon_g
670    INTEGER                              :: i, igrid, jj, ij, iglob
671    INTEGER                              :: ierr, ireal, index
672    INTEGER                              :: var_tmp
673    INTEGER, DIMENSION(0:mpi_size-1)     :: knon_nb
674    INTEGER, DIMENSION(0:mpi_size-1)     :: displs
675    INTEGER, DIMENSION(8,3)              :: off_ini
676    INTEGER, DIMENSION(8)                :: offset 
677    INTEGER, DIMENSION(knon)             :: ktindex_p
678    INTEGER, DIMENSION(nbp_lon,nbp_lat)        :: correspond
679    INTEGER, ALLOCATABLE, DIMENSION(:)   :: ktindex_g
680    INTEGER, ALLOCATABLE, DIMENSION(:,:) :: neighbours_g
681    REAL, DIMENSION(klon_glo)            :: pctsrf_g
682   
683!
684! End definition
685!****************************************************************************************
686
687    IF (is_sequential) THEN
688       knon_nb(:)=knon
689    ELSE 
690       
691#ifdef CPP_MPI 
692       CALL MPI_GATHER(knon,1,MPI_INTEGER,knon_nb,1,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
693#endif
694       
695    ENDIF
696   
697    IF (is_mpi_root) THEN
698       knon_g=SUM(knon_nb(:))
699       ALLOCATE(ktindex_g(knon_g))
700       ALLOCATE(neighbours_g(knon_g,8))
701       neighbours_g(:,:)=-1
702       displs(0)=0
703       DO i=1,mpi_size-1
704          displs(i)=displs(i-1)+knon_nb(i-1)
705       ENDDO
706   ELSE
707       ALLOCATE(ktindex_g(1))
708       ALLOCATE(neighbours_g(1,8))
709   ENDIF
710   
711    ktindex_p(1:knon)=ktindex(1:knon)+klon_mpi_begin-1+nbp_lon-1
712   
713    IF (is_sequential) THEN
714       ktindex_g(:)=ktindex_p(:)
715    ELSE
716       
717#ifdef CPP_MPI 
718       CALL MPI_GATHERV(ktindex_p,knon,MPI_INTEGER,ktindex_g,knon_nb,&
719            displs,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
720#endif
721       
722    ENDIF
723   
724    CALL Gather(pctsrf,pctsrf_g)
725   
726    IF (is_mpi_root) THEN
727!  Initialisation des offset   
728!
729! offset bord ouest
730       off_ini(1,1) = - nbp_lon  ; off_ini(2,1) = - nbp_lon + 1; off_ini(3,1) = 1
731       off_ini(4,1) = nbp_lon + 1; off_ini(5,1) = nbp_lon      ; off_ini(6,1) = 2 * nbp_lon - 1
732       off_ini(7,1) = nbp_lon -1 ; off_ini(8,1) = - 1
733! offset point normal
734       off_ini(1,2) = - nbp_lon  ; off_ini(2,2) = - nbp_lon + 1; off_ini(3,2) = 1
735       off_ini(4,2) = nbp_lon + 1; off_ini(5,2) = nbp_lon      ; off_ini(6,2) = nbp_lon - 1
736       off_ini(7,2) = -1     ; off_ini(8,2) = - nbp_lon - 1
737! offset bord   est
738       off_ini(1,3) = - nbp_lon; off_ini(2,3) = - 2 * nbp_lon + 1; off_ini(3,3) = - nbp_lon + 1
739       off_ini(4,3) =  1   ; off_ini(5,3) = nbp_lon          ; off_ini(6,3) = nbp_lon - 1
740       off_ini(7,3) = -1   ; off_ini(8,3) = - nbp_lon - 1
741!
742!
743! Attention aux poles
744!
745       DO igrid = 1, knon_g
746          index = ktindex_g(igrid)
747          jj = INT((index - 1)/nbp_lon) + 1
748          ij = index - (jj - 1) * nbp_lon
749          correspond(ij,jj) = igrid
750       ENDDO
751       
752       DO igrid = 1, knon_g
753          iglob = ktindex_g(igrid)
754          IF (MOD(iglob, nbp_lon) == 1) THEN
755             offset = off_ini(:,1)
756          ELSE IF(MOD(iglob, nbp_lon) == 0) THEN
757             offset = off_ini(:,3)
758          ELSE
759             offset = off_ini(:,2)
760          ENDIF
761          DO i = 1, 8
762             index = iglob + offset(i)
763             ireal = (MIN(MAX(1, index - nbp_lon + 1), klon_glo))
764             IF (pctsrf_g(ireal) > EPSFRA) THEN
765                jj = INT((index - 1)/nbp_lon) + 1
766                ij = index - (jj - 1) * nbp_lon
767                neighbours_g(igrid, i) = correspond(ij, jj)
768             ENDIF
769          ENDDO
770       ENDDO
771
772    ENDIF
773   
774    DO i=1,8
775       IF (is_sequential) THEN
776          neighbours(:,i)=neighbours_g(:,i)
777       ELSE
778#ifdef CPP_MPI
779          IF (knon > 0) THEN
780             ! knon>0, scattter global field neighbours_g from master process to local process
781             CALL MPI_SCATTERV(neighbours_g(:,i),knon_nb,displs,MPI_INTEGER,neighbours(:,i),knon,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
782          ELSE
783             ! knon=0, no need to save the field for this process
784             CALL MPI_SCATTERV(neighbours_g(:,i),knon_nb,displs,MPI_INTEGER,var_tmp,knon,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
785          END IF
786#endif
787       ENDIF
788    ENDDO
789   
790  END SUBROUTINE Init_neighbours
791!
792!****************************************************************************************
793!
794
795#endif
796END MODULE surf_land_orchidee_noopenmp_mod
Note: See TracBrowser for help on using the repository browser.