source: LMDZ6/trunk/libf/phylmd/surf_land_orchidee_noopenmp_mod.F90 @ 5319

Last change on this file since 5319 was 5305, checked in by abarral, 8 days ago

Turn YOMCST2.h.h into module

  • 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.4 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
108USE dimpft_mod_h
109        USE yomcst_mod_h
110IMPLICIT NONE
111
112
113!
114! Parametres d'entree
115!****************************************************************************************
116    INTEGER, INTENT(IN)                       :: itime
117    REAL, INTENT(IN)                          :: dtime
118    REAL, INTENT(IN)                          :: date0
119    INTEGER, INTENT(IN)                       :: knon
120    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
121    LOGICAL, INTENT(IN)                       :: debut, lafin
122    REAL, DIMENSION(klon,nbsrf), INTENT(IN)   :: pctsrf
123    REAL, DIMENSION(klon), INTENT(IN)         :: rlon, rlat
124    REAL, DIMENSION(klon), INTENT(IN)         :: yrmu0
125    REAL, DIMENSION(klon), INTENT(IN)         :: plev
126    REAL, DIMENSION(klon), INTENT(IN)         :: u1_lay, v1_lay, gustiness
127    REAL, DIMENSION(klon), INTENT(IN)         :: temp_air, spechum
128    REAL, DIMENSION(klon), INTENT(IN)         :: epot_air, ccanopy
129    REAL, DIMENSION(klon), INTENT(IN)         :: tq_cdrag
130    REAL, DIMENSION(klon), INTENT(IN)         :: petAcoef, peqAcoef
131    REAL, DIMENSION(klon), INTENT(IN)         :: petBcoef, peqBcoef
132    REAL, DIMENSION(klon), INTENT(IN)         :: precip_rain, precip_snow
133    REAL, DIMENSION(klon), INTENT(IN)         :: lwdown, swnet, swdown, ps
134    REAL, DIMENSION(klon), INTENT(IN)         :: q2m, t2m
135
136! Parametres de sortie
137!****************************************************************************************
138    REAL, DIMENSION(klon), INTENT(OUT)        :: evap, fluxsens, fluxlat, qsurf
139    REAL, DIMENSION(klon), INTENT(OUT)        :: tsol_rad, tsurf_new
140    REAL, DIMENSION(klon), INTENT(OUT)        :: alb1_new, alb2_new
141    REAL, DIMENSION(klon), INTENT(OUT)        :: emis_new, z0_new, z0h_new
142    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: veget ! dummy variables
143    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: lai   ! dummy variables
144    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: height! dummy variables
145
146! Local
147!****************************************************************************************
148    INTEGER                                   :: ij, jj, igrid, ireal, index
149    INTEGER                                   :: error
150    INTEGER, SAVE                             :: nb_fields_cpl ! number of fields for the climate-carbon coupling (between ATM and ORCHIDEE).
151    !$OMP THREADPRIVATE(nb_fields_cpl)
152    REAL, SAVE, ALLOCATABLE, DIMENSION(:,:)   :: fields_cpl    ! Fluxes for the climate-carbon coupling
153    !$OMP THREADPRIVATE(fields_cpl)
154    REAL, DIMENSION(klon)                     :: swdown_vrai
155    CHARACTER (len = 20)                      :: modname = 'surf_land_orchidee'
156    CHARACTER (len = 80)                      :: abort_message
157    LOGICAL,SAVE                              :: check = .FALSE.
158    !$OMP THREADPRIVATE(check)
159
160! type de couplage dans sechiba
161!  character (len=10)   :: coupling = 'implicit'
162! drapeaux controlant les appels dans SECHIBA
163!  type(control_type), save   :: control_in
164! Preserved albedo
165    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: albedo_keep, zlev
166    !$OMP THREADPRIVATE(albedo_keep,zlev)
167! coordonnees geographiques
168    REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: lalo
169    !$OMP THREADPRIVATE(lalo)
170! pts voisins
171    INTEGER,ALLOCATABLE, DIMENSION(:,:), SAVE :: neighbours
172    !$OMP THREADPRIVATE(neighbours)
173! fractions continents
174    REAL,ALLOCATABLE, DIMENSION(:), SAVE      :: contfrac
175    !$OMP THREADPRIVATE(contfrac)
176! resolution de la grille
177    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: resolution
178    !$OMP THREADPRIVATE(resolution)
179
180    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: lon_scat, lat_scat 
181    !$OMP THREADPRIVATE(lon_scat,lat_scat)
182
183    LOGICAL, SAVE                             :: lrestart_read = .TRUE.
184    !$OMP THREADPRIVATE(lrestart_read)
185    LOGICAL, SAVE                             :: lrestart_write = .FALSE.
186    !$OMP THREADPRIVATE(lrestart_write)
187
188    REAL, DIMENSION(knon,2)                   :: albedo_out
189    !$OMP THREADPRIVATE(albedo_out)
190
191! Pb de nomenclature
192    REAL, DIMENSION(klon)                     :: petA_orc, peqA_orc
193    REAL, DIMENSION(klon)                     :: petB_orc, peqB_orc
194! Pb de correspondances de grilles
195    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  :: ig, jg
196    !$OMP THREADPRIVATE(ig,jg)
197    INTEGER :: indi, indj
198    INTEGER, SAVE, ALLOCATABLE,DIMENSION(:)   :: ktindex
199    !$OMP THREADPRIVATE(ktindex)
200
201! Essai cdrag
202    REAL, DIMENSION(klon)                     :: cdrag
203    INTEGER,SAVE                              :: offset
204    !$OMP THREADPRIVATE(offset)
205
206    REAL, DIMENSION(klon_glo)                 :: rlon_g,rlat_g
207    INTEGER, SAVE                             :: orch_comm
208    !$OMP THREADPRIVATE(orch_comm)
209
210    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: coastalflow
211    !$OMP THREADPRIVATE(coastalflow)
212    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: riverflow
213    !$OMP THREADPRIVATE(riverflow)
214!
215! Fin definition
216!****************************************************************************************
217#ifdef CPP_VEGET
218
219    IF (check) WRITE(lunout,*)'Entree ', modname
220 
221! Initialisation
222 
223    IF (debut) THEN
224! Test de coherence
225#ifndef ORCH_NEW
226       ! Compilation avec orchidee nouvelle version necessaire avec carbon_cycle_cpl=y
227       IF (carbon_cycle_cpl) THEN
228          abort_message='You must define preprossing key ORCH_NEW when running carbon_cycle_cpl=y'
229          CALL abort_physic(modname,abort_message,1)
230       END IF
231#endif
232       ALLOCATE(ktindex(knon))
233       IF ( .NOT. ALLOCATED(albedo_keep)) THEN
234          ALLOCATE(albedo_keep(klon))
235          ALLOCATE(zlev(knon))
236       ENDIF
237! Pb de correspondances de grilles
238       ALLOCATE(ig(klon))
239       ALLOCATE(jg(klon))
240       ig(1) = 1
241       jg(1) = 1
242       indi = 0
243       indj = 2
244       DO igrid = 2, klon - 1
245          indi = indi + 1
246          IF ( indi > nbp_lon) THEN
247             indi = 1
248             indj = indj + 1
249          ENDIF
250          ig(igrid) = indi
251          jg(igrid) = indj
252       ENDDO
253       ig(klon) = 1
254       jg(klon) = nbp_lat
255
256       IF ((.NOT. ALLOCATED(lalo))) THEN
257          ALLOCATE(lalo(knon,2), stat = error)
258          IF (error /= 0) THEN
259             abort_message='Pb allocation lalo'
260             CALL abort_physic(modname,abort_message,1)
261          ENDIF
262       ENDIF
263       IF ((.NOT. ALLOCATED(lon_scat))) THEN
264          ALLOCATE(lon_scat(nbp_lon,nbp_lat), stat = error)
265          IF (error /= 0) THEN
266             abort_message='Pb allocation lon_scat'
267             CALL abort_physic(modname,abort_message,1)
268          ENDIF
269       ENDIF
270       IF ((.NOT. ALLOCATED(lat_scat))) THEN
271          ALLOCATE(lat_scat(nbp_lon,nbp_lat), stat = error)
272          IF (error /= 0) THEN
273             abort_message='Pb allocation lat_scat'
274             CALL abort_physic(modname,abort_message,1)
275          ENDIF
276       ENDIF
277       lon_scat = 0.
278       lat_scat = 0.
279       DO igrid = 1, knon
280          index = knindex(igrid)
281          lalo(igrid,2) = rlon(index)
282          lalo(igrid,1) = rlat(index)
283       ENDDO
284
285       
286       
287       CALL Gather(rlon,rlon_g)
288       CALL Gather(rlat,rlat_g)
289
290       IF (is_mpi_root) THEN
291          index = 1
292          DO jj = 2, nbp_lat-1
293             DO ij = 1, nbp_lon
294                index = index + 1
295                lon_scat(ij,jj) = rlon_g(index)
296                lat_scat(ij,jj) = rlat_g(index)
297             ENDDO
298          ENDDO
299          lon_scat(:,1) = lon_scat(:,2)
300          lat_scat(:,1) = rlat_g(1)
301          lon_scat(:,nbp_lat) = lon_scat(:,2)
302          lat_scat(:,nbp_lat) = rlat_g(klon_glo)
303       ENDIF
304
305       CALL bcast(lon_scat)
306       CALL bcast(lat_scat)
307
308!
309! Allouer et initialiser le tableau des voisins et des fraction de continents
310!
311       IF ( (.NOT.ALLOCATED(neighbours))) THEN
312          ALLOCATE(neighbours(knon,8), stat = error)
313          IF (error /= 0) THEN
314             abort_message='Pb allocation neighbours'
315             CALL abort_physic(modname,abort_message,1)
316          ENDIF
317       ENDIF
318       neighbours = -1.
319       IF (( .NOT. ALLOCATED(contfrac))) THEN
320          ALLOCATE(contfrac(knon), stat = error)
321          IF (error /= 0) THEN
322             abort_message='Pb allocation contfrac'
323             CALL abort_physic(modname,abort_message,1)
324          ENDIF
325       ENDIF
326
327       DO igrid = 1, knon
328          ireal = knindex(igrid)
329          contfrac(igrid) = pctsrf(ireal,is_ter)
330       ENDDO
331
332
333       CALL Init_neighbours(knon,neighbours,knindex,pctsrf(:,is_ter))
334
335!
336!  Allocation et calcul resolutions
337       IF ( (.NOT.ALLOCATED(resolution))) THEN
338          ALLOCATE(resolution(knon,2), stat = error)
339          IF (error /= 0) THEN
340             abort_message='Pb allocation resolution'
341             CALL abort_physic(modname,abort_message,1)
342          ENDIF
343       ENDIF
344       DO igrid = 1, knon
345          ij = knindex(igrid)
346          resolution(igrid,1) = dx(ij)
347          resolution(igrid,2) = dy(ij)
348       ENDDO
349     
350       ALLOCATE(coastalflow(klon), stat = error)
351       IF (error /= 0) THEN
352          abort_message='Pb allocation coastalflow'
353          CALL abort_physic(modname,abort_message,1)
354       ENDIF
355       
356       ALLOCATE(riverflow(klon), stat = error)
357       IF (error /= 0) THEN
358          abort_message='Pb allocation riverflow'
359          CALL abort_physic(modname,abort_message,1)
360       ENDIF
361
362!
363! Allocate variables needed for carbon_cycle_mod
364       IF ( carbon_cycle_cpl ) THEN
365          nb_fields_cpl=2
366       ELSE
367          nb_fields_cpl=1
368       END IF
369
370
371       IF (carbon_cycle_cpl) THEN
372          ALLOCATE(fco2_land_inst(klon),stat=error)
373          IF (error /= 0)  CALL abort_physic(modname,'Pb in allocation fco2_land_inst',1)
374         
375          ALLOCATE(fco2_lu_inst(klon),stat=error)
376          IF(error /=0) CALL abort_physic(modname,'Pb in allocation fco2_lu_inst',1)
377       END IF
378
379       ALLOCATE(fields_cpl(klon,nb_fields_cpl), stat = error)
380       IF (error /= 0) CALL abort_physic(modname,'Pb in allocation fields_cpl',1)
381
382    ENDIF                          ! (fin debut)
383
384!
385! Appel a la routine sols continentaux
386!
387    IF (lafin) lrestart_write = .TRUE.
388    IF (check) WRITE(lunout,*)'lafin ',lafin,lrestart_write
389   
390    petA_orc(1:knon) = petBcoef(1:knon) * dtime
391    petB_orc(1:knon) = petAcoef(1:knon)
392    peqA_orc(1:knon) = peqBcoef(1:knon) * dtime
393    peqB_orc(1:knon) = peqAcoef(1:knon)
394
395    cdrag = 0.
396    cdrag(1:knon) = tq_cdrag(1:knon)
397
398! zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/287.05*temp_air(1:knon))*9.80665)
399    zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/RD*temp_air(1:knon))*RG)
400
401
402! PF et PASB
403!   where(cdrag > 0.01)
404!     cdrag = 0.01
405!   endwhere
406!  write(*,*)'Cdrag = ',minval(cdrag),maxval(cdrag)
407
408!
409! Init Orchidee
410!
411!  if (pole_nord) then
412!    offset=0
413!    ktindex(:)=ktindex(:)+nbp_lon-1
414!  else
415!    offset = klon_mpi_begin-1+nbp_lon-1
416!    ktindex(:)=ktindex(:)+MOD(offset,nbp_lon)
417!    offset=offset-MOD(offset,nbp_lon)
418!  endif
419 
420    IF (debut) THEN
421       CALL Get_orchidee_communicator(knon,orch_comm)
422       IF (knon /=0) THEN
423          CALL Init_orchidee_index(knon,orch_comm,knindex,offset,ktindex)
424         
425          IF (.NOT. using_mpi) THEN
426            ! Interface for ORCHIDEE compiled in sequential mode(without preprocessing flag CPP_MPI)
427            CALL intersurf_main (itime+itau_phy-1, nbp_lon, nbp_lat, knon, ktindex, dtime, &
428                 lrestart_read, lrestart_write, lalo, &
429                 contfrac, neighbours, resolution, date0, &
430                 zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
431                 cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
432                 precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
433                 evap, fluxsens, fluxlat, coastalflow, riverflow, &
434                 tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
435                 lon_scat, lat_scat, q2m, t2m &
436#ifdef ORCH_NEW
437                 , nb_fields_cpl, fields_cpl)
438#else
439                 )
440#endif
441
442          ELSE         
443            ! 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)
444            CALL intersurf_main (itime+itau_phy-1, nbp_lon, nbp_lat, offset, knon, ktindex, &
445                 orch_comm, dtime, lrestart_read, lrestart_write, lalo, &
446                 contfrac, neighbours, resolution, date0, &
447                 zlev,  u1_lay(1:knon), v1_lay(1:knon), spechum(1:knon), temp_air(1:knon), epot_air(1:knon), ccanopy(1:knon), &
448                 cdrag(1:knon), petA_orc(1:knon), peqA_orc(1:knon), petB_orc(1:knon), peqB_orc(1:knon), &
449                 precip_rain(1:knon), precip_snow(1:knon), lwdown(1:knon), swnet(1:knon), swdown(1:knon), ps(1:knon), &
450                 evap(1:knon), fluxsens(1:knon), fluxlat(1:knon), coastalflow(1:knon), riverflow(1:knon), &
451                 tsol_rad(1:knon), tsurf_new(1:knon), qsurf(1:knon), albedo_out(1:knon,:), emis_new(1:knon), z0_new(1:knon), &
452                 lon_scat, lat_scat, q2m, t2m &
453#ifdef ORCH_NEW
454                 , nb_fields_cpl, fields_cpl(1:knon,:))
455#else
456                 )
457#endif
458          ENDIF
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       IF (.NOT. using_mpi) THEN
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    USE lmdz_mpi
548
549! Input arguments
550!****************************************************************************************
551    INTEGER, INTENT(IN)                   :: knon
552    INTEGER, INTENT(IN)                   :: orch_comm
553    INTEGER, DIMENSION(klon), INTENT(IN)  :: knindex
554
555! Output arguments
556!****************************************************************************************
557    INTEGER, INTENT(OUT)                  :: offset
558    INTEGER, DIMENSION(knon), INTENT(OUT) :: ktindex
559
560! Local varables
561!****************************************************************************************
562    INTEGER, DIMENSION(MPI_STATUS_SIZE)   :: status
563
564    INTEGER                               :: MyLastPoint
565    INTEGER                               :: LastPoint
566    INTEGER                               :: mpi_rank_orch
567    INTEGER                               :: mpi_size_orch
568    INTEGER                               :: ierr
569!
570! End definition
571!****************************************************************************************
572
573    MyLastPoint=klon_mpi_begin-1+knindex(knon)+nbp_lon-1
574   
575    IF (is_parallel) THEN
576       CALL MPI_COMM_SIZE(orch_comm,mpi_size_orch,ierr)
577       CALL MPI_COMM_RANK(orch_comm,mpi_rank_orch,ierr)
578    ELSE
579       mpi_rank_orch=0
580       mpi_size_orch=1
581    ENDIF
582
583    IF (is_parallel) THEN
584       IF (mpi_rank_orch /= 0) THEN
585          CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank_orch-1,1234,orch_comm,status,ierr)
586       ENDIF
587       
588       IF (mpi_rank_orch /= mpi_size_orch-1) THEN
589          CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank_orch+1,1234,orch_comm,ierr) 
590       ENDIF
591    ENDIF
592   
593    IF (mpi_rank_orch == 0) THEN
594       offset=0
595    ELSE
596       offset=LastPoint-MOD(LastPoint,nbp_lon)
597    ENDIF
598   
599    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin+nbp_lon-1)-offset-1
600   
601
602  END SUBROUTINE  Init_orchidee_index
603!
604!****************************************************************************************
605!
606  SUBROUTINE Get_orchidee_communicator(knon,orch_comm)
607  USE lmdz_mpi
608
609    INTEGER,INTENT(IN)  :: knon
610    INTEGER,INTENT(OUT) :: orch_comm
611   
612    INTEGER             :: color
613    INTEGER             :: ierr
614!
615! End definition
616!****************************************************************************************
617
618    IF (knon==0) THEN
619       color = 0
620    ELSE
621       color = 1
622    ENDIF
623   
624    IF (using_mpi) THEN
625      CALL MPI_COMM_SPLIT(COMM_LMDZ_PHY,color,mpi_rank,orch_comm,ierr)
626    ENDIF
627   
628  END SUBROUTINE Get_orchidee_communicator
629!
630!****************************************************************************************
631
632  SUBROUTINE Init_neighbours(knon,neighbours,ktindex,pctsrf)
633   
634    USE indice_sol_mod
635    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
636    USE lmdz_mpi
637
638! Input arguments
639!****************************************************************************************
640    INTEGER, INTENT(IN)                     :: knon
641    INTEGER, DIMENSION(klon), INTENT(IN)    :: ktindex
642    REAL, DIMENSION(klon), INTENT(IN)       :: pctsrf
643   
644! Output arguments
645!****************************************************************************************
646    INTEGER, DIMENSION(knon,8), INTENT(OUT) :: neighbours
647
648! Local variables
649!****************************************************************************************
650    INTEGER                              :: knon_g
651    INTEGER                              :: i, igrid, jj, ij, iglob
652    INTEGER                              :: ierr, ireal, index
653    INTEGER                              :: var_tmp
654    INTEGER, DIMENSION(0:mpi_size-1)     :: knon_nb
655    INTEGER, DIMENSION(0:mpi_size-1)     :: displs
656    INTEGER, DIMENSION(8,3)              :: off_ini
657    INTEGER, DIMENSION(8)                :: offset 
658    INTEGER, DIMENSION(knon)             :: ktindex_p
659    INTEGER, DIMENSION(nbp_lon,nbp_lat)        :: correspond
660    INTEGER, ALLOCATABLE, DIMENSION(:)   :: ktindex_g
661    INTEGER, ALLOCATABLE, DIMENSION(:,:) :: neighbours_g
662    REAL, DIMENSION(klon_glo)            :: pctsrf_g
663   
664!
665! End definition
666!****************************************************************************************
667
668    IF (is_sequential) THEN
669       knon_nb(:)=knon
670    ELSE 
671        CALL MPI_GATHER(knon,1,MPI_INTEGER,knon_nb,1,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
672    ENDIF
673   
674    IF (is_mpi_root) THEN
675       knon_g=SUM(knon_nb(:))
676       ALLOCATE(ktindex_g(knon_g))
677       ALLOCATE(neighbours_g(knon_g,8))
678       neighbours_g(:,:)=-1
679       displs(0)=0
680       DO i=1,mpi_size-1
681          displs(i)=displs(i-1)+knon_nb(i-1)
682       ENDDO
683   ELSE
684       ALLOCATE(ktindex_g(1))
685       ALLOCATE(neighbours_g(1,8))
686   ENDIF
687   
688    ktindex_p(1:knon)=ktindex(1:knon)+klon_mpi_begin-1+nbp_lon-1
689   
690    IF (is_sequential) THEN
691       ktindex_g(:)=ktindex_p(:)
692    ELSE
693       CALL MPI_GATHERV(ktindex_p,knon,MPI_INTEGER,ktindex_g,knon_nb,&
694            displs,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
695    ENDIF
696   
697    CALL Gather(pctsrf,pctsrf_g)
698   
699    IF (is_mpi_root) THEN
700!  Initialisation des offset   
701!
702! offset bord ouest
703       off_ini(1,1) = - nbp_lon  ; off_ini(2,1) = - nbp_lon + 1; off_ini(3,1) = 1
704       off_ini(4,1) = nbp_lon + 1; off_ini(5,1) = nbp_lon      ; off_ini(6,1) = 2 * nbp_lon - 1
705       off_ini(7,1) = nbp_lon -1 ; off_ini(8,1) = - 1
706! offset point normal
707       off_ini(1,2) = - nbp_lon  ; off_ini(2,2) = - nbp_lon + 1; off_ini(3,2) = 1
708       off_ini(4,2) = nbp_lon + 1; off_ini(5,2) = nbp_lon      ; off_ini(6,2) = nbp_lon - 1
709       off_ini(7,2) = -1     ; off_ini(8,2) = - nbp_lon - 1
710! offset bord   est
711       off_ini(1,3) = - nbp_lon; off_ini(2,3) = - 2 * nbp_lon + 1; off_ini(3,3) = - nbp_lon + 1
712       off_ini(4,3) =  1   ; off_ini(5,3) = nbp_lon          ; off_ini(6,3) = nbp_lon - 1
713       off_ini(7,3) = -1   ; off_ini(8,3) = - nbp_lon - 1
714!
715!
716! Attention aux poles
717!
718       DO igrid = 1, knon_g
719          index = ktindex_g(igrid)
720          jj = INT((index - 1)/nbp_lon) + 1
721          ij = index - (jj - 1) * nbp_lon
722          correspond(ij,jj) = igrid
723       ENDDO
724       
725       DO igrid = 1, knon_g
726          iglob = ktindex_g(igrid)
727          IF (MOD(iglob, nbp_lon) == 1) THEN
728             offset = off_ini(:,1)
729          ELSE IF(MOD(iglob, nbp_lon) == 0) THEN
730             offset = off_ini(:,3)
731          ELSE
732             offset = off_ini(:,2)
733          ENDIF
734          DO i = 1, 8
735             index = iglob + offset(i)
736             ireal = (MIN(MAX(1, index - nbp_lon + 1), klon_glo))
737             IF (pctsrf_g(ireal) > EPSFRA) THEN
738                jj = INT((index - 1)/nbp_lon) + 1
739                ij = index - (jj - 1) * nbp_lon
740                neighbours_g(igrid, i) = correspond(ij, jj)
741             ENDIF
742          ENDDO
743       ENDDO
744
745    ENDIF
746   
747    DO i=1,8
748       IF (is_sequential) THEN
749          neighbours(:,i)=neighbours_g(:,i)
750       ELSE
751          IF (knon > 0) THEN
752             ! knon>0, scattter global field neighbours_g from master process to local process
753             CALL MPI_SCATTERV(neighbours_g(:,i),knon_nb,displs,MPI_INTEGER,neighbours(:,i),knon,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
754          ELSE
755             ! knon=0, no need to save the field for this process
756             CALL MPI_SCATTERV(neighbours_g(:,i),knon_nb,displs,MPI_INTEGER,var_tmp,knon,MPI_INTEGER,0,COMM_LMDZ_PHY,ierr)
757          END IF
758       ENDIF
759    ENDDO
760   
761  END SUBROUTINE Init_neighbours
762!
763!****************************************************************************************
764!
765
766#endif
767END MODULE surf_land_orchidee_noopenmp_mod
Note: See TracBrowser for help on using the repository browser.