source: LMDZ5/branches/LMDZ_tree_FC/libf/phylmd/surf_land_orchidee_noz0h_mod.F90 @ 2924

Last change on this file since 2924 was 2571, checked in by jghattas, 8 years ago

Interface LMDZ/ORCHIDEE :

  • copied previous default module surf_land_orchidee_mod.f90 into surf_land_orchidee_noz0h.f90. This interface can still be compiled if adding cpp key ORCHIDEE_NOZ0H
  • modified default interface by adding z0h as output from ORCHIDEE.
  • added comments in each module surf_land_orchidee_xxx of compatiblity with ORCHIDEE.
  • all modules surf_land_orchidee_xxx now send back z0h and z0m to surf_land_mod. But note that z0m and zOh are different only in the new default version surf_land_orchidee_mod.f90. In the old interfaces, z0h is a copy of z0m.


cosp : some small changes to be able to compile with gfortran

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Author Date Id Revision
File size: 22.9 KB
Line 
1!
2MODULE surf_land_orchidee_noz0h_mod
3#ifdef ORCHIDEE_NOZ0H
4!
5! This module controles the interface towards the model ORCHIDEE.
6!
7! Compatibility with ORCHIDIEE :
8! This module is compiled only if CPP key ORCHIDEE_NOZ0H is defined.
9! This version can be used with ORCHIDEE/trunk revision 2961-3525.
10! It is possible to use this interface with later versions of ORCHIDEE/trunk than 3525 but
11! in that case the option ROUGH_DYN=y in ORCHIDEE should not be set.
12! This interface can also be used with ORCHIDEE/trunk revision 1078-2960 if changing
13! coszang=yrmu0 into sinang=yrmu0 at 2 places later below in this module.
14!
15! Subroutines in this module : surf_land_orchidee
16!                              Init_orchidee_index
17!                              Get_orchidee_communicator
18!                              Init_neighbours
19
20  USE dimphy
21#ifdef CPP_VEGET
22  USE intersurf     ! module d'ORCHIDEE
23#endif
24  USE cpl_mod,      ONLY : cpl_send_land_fields
25  USE surface_data, ONLY : type_ocean
26  USE geometry_mod, ONLY : dx, dy
27  USE mod_grid_phy_lmdz
28  USE mod_phys_lmdz_para, mpi_root_rank=>mpi_master
29
30  IMPLICIT NONE
31
32  PRIVATE
33  PUBLIC  :: surf_land_orchidee
34
35  LOGICAL, ALLOCATABLE, SAVE :: flag_omp(:)
36CONTAINS
37!
38!****************************************************************************************
39
40  SUBROUTINE surf_land_orchidee(itime, dtime, date0, knon, &
41       knindex, rlon, rlat, yrmu0, pctsrf, &
42       debut, lafin, &
43       plev,  u1_lay, v1_lay, gustiness, temp_air, spechum, epot_air, ccanopy, &
44       tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
45       precip_rain, precip_snow, lwdown, swnet, swdown, &
46       ps, q2m, t2m, &
47       evap, fluxsens, fluxlat, &             
48       tsol_rad, tsurf_new, alb1_new, alb2_new, &
49       emis_new, z0_new, z0h_new, qsurf)
50
51    USE mod_surf_para
52    USE mod_synchro_omp
53    USE carbon_cycle_mod, ONLY : carbon_cycle_cpl
54    USE indice_sol_mod
55    USE print_control_mod, ONLY: lunout
56    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
57#ifdef CPP_VEGET
58    USE time_phylmdz_mod, ONLY: itau_phy
59#endif
60!   
61! Cette routine sert d'interface entre le modele atmospherique et le
62! modele de sol continental. Appel a sechiba
63!
64! L. Fairhead 02/2000
65!
66! input:
67!   itime        numero du pas de temps
68!   dtime        pas de temps de la physique (en s)
69!   nisurf       index de la surface a traiter (1 = sol continental)
70!   knon         nombre de points de la surface a traiter
71!   knindex      index des points de la surface a traiter
72!   rlon         longitudes de la grille entiere
73!   rlat         latitudes de la grille entiere
74!   pctsrf       tableau des fractions de surface de chaque maille
75!   debut        logical: 1er appel a la physique (lire les restart)
76!   lafin        logical: dernier appel a la physique (ecrire les restart)
77!                     (si false calcul simplifie des fluxs sur les continents)
78!   plev         hauteur de la premiere couche (Pa)     
79!   u1_lay       vitesse u 1ere couche
80!   v1_lay       vitesse v 1ere couche
81!   temp_air     temperature de l'air 1ere couche
82!   spechum      humidite specifique 1ere couche
83!   epot_air     temp pot de l'air
84!   ccanopy      concentration CO2 canopee, correspond au co2_send de
85!                carbon_cycle_mod ou valeur constant co2_ppm
86!   tq_cdrag     cdrag
87!   petAcoef     coeff. A de la resolution de la CL pour t
88!   peqAcoef     coeff. A de la resolution de la CL pour q
89!   petBcoef     coeff. B de la resolution de la CL pour t
90!   peqBcoef     coeff. B de la resolution de la CL pour q
91!   precip_rain  precipitation liquide
92!   precip_snow  precipitation solide
93!   lwdown       flux IR descendant a la surface
94!   swnet        flux solaire net
95!   swdown       flux solaire entrant a la surface
96!   ps           pression au sol
97!   radsol       rayonnement net aus sol (LW + SW)
98!   
99!
100! output:
101!   evap         evaporation totale
102!   fluxsens     flux de chaleur sensible
103!   fluxlat      flux de chaleur latente
104!   tsol_rad     
105!   tsurf_new    temperature au sol
106!   alb1_new     albedo in visible SW interval
107!   alb2_new     albedo in near IR interval
108!   emis_new     emissivite
109!   z0_new       surface roughness
110!   z0h_new      surface roughness, it is a copy from z0_new
111!   qsurf        air moisture at surface
112!
113    INCLUDE "YOMCST.h"
114 
115!
116! Parametres d'entree
117!****************************************************************************************
118    INTEGER, INTENT(IN)                       :: itime
119    REAL, INTENT(IN)                          :: dtime
120    REAL, INTENT(IN)                          :: date0
121    INTEGER, INTENT(IN)                       :: knon
122    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
123    LOGICAL, INTENT(IN)                       :: debut, lafin
124    REAL, DIMENSION(klon,nbsrf), INTENT(IN)   :: pctsrf
125    REAL, DIMENSION(klon), INTENT(IN)         :: rlon, rlat
126    REAL, DIMENSION(klon), INTENT(IN)         :: yrmu0 ! cosine of solar zenith angle
127    REAL, DIMENSION(klon), INTENT(IN)         :: plev
128    REAL, DIMENSION(klon), INTENT(IN)         :: u1_lay, v1_lay, gustiness
129    REAL, DIMENSION(klon), INTENT(IN)         :: temp_air, spechum
130    REAL, DIMENSION(klon), INTENT(IN)         :: epot_air, ccanopy
131    REAL, DIMENSION(klon), INTENT(IN)         :: tq_cdrag
132    REAL, DIMENSION(klon), INTENT(IN)         :: petAcoef, peqAcoef
133    REAL, DIMENSION(klon), INTENT(IN)         :: petBcoef, peqBcoef
134    REAL, DIMENSION(klon), INTENT(IN)         :: precip_rain, precip_snow
135    REAL, DIMENSION(klon), INTENT(IN)         :: lwdown, swnet, swdown, ps
136    REAL, DIMENSION(klon), INTENT(IN)         :: q2m, t2m
137
138! Parametres de sortie
139!****************************************************************************************
140    REAL, DIMENSION(klon), INTENT(OUT)        :: evap, fluxsens, fluxlat, qsurf
141    REAL, DIMENSION(klon), INTENT(OUT)        :: tsol_rad, tsurf_new
142    REAL, DIMENSION(klon), INTENT(OUT)        :: alb1_new, alb2_new
143    REAL, DIMENSION(klon), INTENT(OUT)        :: emis_new, z0_new, z0h_new
144
145! Local
146!****************************************************************************************
147    INTEGER                                   :: ij, jj, igrid, ireal, index
148    INTEGER                                   :: error
149    REAL, DIMENSION(klon)                     :: swdown_vrai
150    CHARACTER (len = 20)                      :: modname = 'surf_land_orchidee'
151    CHARACTER (len = 80)                      :: abort_message
152    LOGICAL,SAVE                              :: check = .FALSE.
153    !$OMP THREADPRIVATE(check)
154
155! type de couplage dans sechiba
156!  character (len=10)   :: coupling = 'implicit'
157! drapeaux controlant les appels dans SECHIBA
158!  type(control_type), save   :: control_in
159! Preserved albedo
160    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: albedo_keep, zlev
161    !$OMP THREADPRIVATE(albedo_keep,zlev)
162! coordonnees geographiques
163    REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: lalo
164    !$OMP THREADPRIVATE(lalo)
165! pts voisins
166    INTEGER,ALLOCATABLE, DIMENSION(:,:), SAVE :: neighbours
167    !$OMP THREADPRIVATE(neighbours)
168! fractions continents
169    REAL,ALLOCATABLE, DIMENSION(:), SAVE      :: contfrac
170    !$OMP THREADPRIVATE(contfrac)
171! resolution de la grille
172    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: resolution
173    !$OMP THREADPRIVATE(resolution)
174
175    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: lon_scat, lat_scat 
176    !$OMP THREADPRIVATE(lon_scat,lat_scat)
177
178    LOGICAL, SAVE                             :: lrestart_read = .TRUE.
179    !$OMP THREADPRIVATE(lrestart_read)
180    LOGICAL, SAVE                             :: lrestart_write = .FALSE.
181    !$OMP THREADPRIVATE(lrestart_write)
182
183    REAL, DIMENSION(knon,2)                   :: albedo_out
184
185! Pb de nomenclature
186    REAL, DIMENSION(klon)                     :: petA_orc, peqA_orc
187    REAL, DIMENSION(klon)                     :: petB_orc, peqB_orc
188! Pb de correspondances de grilles
189    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  :: ig, jg
190    !$OMP THREADPRIVATE(ig,jg)
191    INTEGER :: indi, indj
192    INTEGER, SAVE, ALLOCATABLE,DIMENSION(:)   :: ktindex
193    !$OMP THREADPRIVATE(ktindex)
194
195! Essai cdrag
196    REAL, DIMENSION(klon)                     :: cdrag
197    INTEGER,SAVE                              :: offset
198    !$OMP THREADPRIVATE(offset)
199
200    REAL, DIMENSION(klon_glo)                 :: rlon_g,rlat_g
201    INTEGER, SAVE                             :: orch_comm
202    !$OMP THREADPRIVATE(orch_comm)
203
204    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: coastalflow
205    !$OMP THREADPRIVATE(coastalflow)
206    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: riverflow
207    !$OMP THREADPRIVATE(riverflow)
208   
209    INTEGER :: orch_omp_rank
210    INTEGER :: orch_omp_size
211!
212! Fin definition
213!****************************************************************************************
214
215    IF (check) WRITE(lunout,*)'Entree ', modname
216 
217! Initialisation
218 
219    IF (debut) THEN
220! Test of coherence between variable ok_veget and cpp key CPP_VEGET
221#ifndef CPP_VEGET
222       abort_message='Pb de coherence: ok_veget = .true. mais CPP_VEGET = .false.'
223       CALL abort_physic(modname,abort_message,1)
224#endif
225
226       CALL Init_surf_para(knon)
227       ALLOCATE(ktindex(knon))
228       IF ( .NOT. ALLOCATED(albedo_keep)) THEN
229!ym          ALLOCATE(albedo_keep(klon))
230!ym bizarre que non alloué en knon precedement
231          ALLOCATE(albedo_keep(knon))
232          ALLOCATE(zlev(knon))
233       ENDIF
234! Pb de correspondances de grilles
235       ALLOCATE(ig(klon))
236       ALLOCATE(jg(klon))
237       ig(1) = 1
238       jg(1) = 1
239       indi = 0
240       indj = 2
241       DO igrid = 2, klon - 1
242          indi = indi + 1
243          IF ( indi > nbp_lon) THEN
244             indi = 1
245             indj = indj + 1
246          ENDIF
247          ig(igrid) = indi
248          jg(igrid) = indj
249       ENDDO
250       ig(klon) = 1
251       jg(klon) = nbp_lat
252
253       IF ((.NOT. ALLOCATED(lalo))) THEN
254          ALLOCATE(lalo(knon,2), stat = error)
255          IF (error /= 0) THEN
256             abort_message='Pb allocation lalo'
257             CALL abort_physic(modname,abort_message,1)
258          ENDIF
259       ENDIF
260       IF ((.NOT. ALLOCATED(lon_scat))) THEN
261          ALLOCATE(lon_scat(nbp_lon,nbp_lat), stat = error)
262          IF (error /= 0) THEN
263             abort_message='Pb allocation lon_scat'
264             CALL abort_physic(modname,abort_message,1)
265          ENDIF
266       ENDIF
267       IF ((.NOT. ALLOCATED(lat_scat))) THEN
268          ALLOCATE(lat_scat(nbp_lon,nbp_lat), stat = error)
269          IF (error /= 0) THEN
270             abort_message='Pb allocation lat_scat'
271             CALL abort_physic(modname,abort_message,1)
272          ENDIF
273       ENDIF
274       lon_scat = 0.
275       lat_scat = 0.
276       DO igrid = 1, knon
277          index = knindex(igrid)
278          lalo(igrid,2) = rlon(index)
279          lalo(igrid,1) = rlat(index)
280       ENDDO
281
282       
283       
284       CALL Gather(rlon,rlon_g)
285       CALL Gather(rlat,rlat_g)
286
287       IF (is_mpi_root) THEN
288          index = 1
289          DO jj = 2, nbp_lat-1
290             DO ij = 1, nbp_lon
291                index = index + 1
292                lon_scat(ij,jj) = rlon_g(index)
293                lat_scat(ij,jj) = rlat_g(index)
294             ENDDO
295          ENDDO
296          lon_scat(:,1) = lon_scat(:,2)
297          lat_scat(:,1) = rlat_g(1)
298          lon_scat(:,nbp_lat) = lon_scat(:,2)
299          lat_scat(:,nbp_lat) = rlat_g(klon_glo)
300       ENDIF
301   
302       CALL bcast(lon_scat)
303       CALL bcast(lat_scat)
304!
305! Allouer et initialiser le tableau des voisins et des fraction de continents
306!
307       IF ( (.NOT.ALLOCATED(neighbours))) THEN
308          ALLOCATE(neighbours(knon,8), stat = error)
309          IF (error /= 0) THEN
310             abort_message='Pb allocation neighbours'
311             CALL abort_physic(modname,abort_message,1)
312          ENDIF
313       ENDIF
314       neighbours = -1.
315       IF (( .NOT. ALLOCATED(contfrac))) THEN
316          ALLOCATE(contfrac(knon), stat = error)
317          IF (error /= 0) THEN
318             abort_message='Pb allocation contfrac'
319             CALL abort_physic(modname,abort_message,1)
320          ENDIF
321       ENDIF
322
323       DO igrid = 1, knon
324          ireal = knindex(igrid)
325          contfrac(igrid) = pctsrf(ireal,is_ter)
326       ENDDO
327
328
329       CALL Init_neighbours(knon,neighbours,knindex,pctsrf(:,is_ter))
330
331!
332!  Allocation et calcul resolutions
333       IF ( (.NOT.ALLOCATED(resolution))) THEN
334          ALLOCATE(resolution(knon,2), stat = error)
335          IF (error /= 0) THEN
336             abort_message='Pb allocation resolution'
337             CALL abort_physic(modname,abort_message,1)
338          ENDIF
339       ENDIF
340       DO igrid = 1, knon
341          ij = knindex(igrid)
342          resolution(igrid,1) = dx(ij)
343          resolution(igrid,2) = dy(ij)
344       ENDDO
345     
346       ALLOCATE(coastalflow(klon), stat = error)
347       IF (error /= 0) THEN
348          abort_message='Pb allocation coastalflow'
349          CALL abort_physic(modname,abort_message,1)
350       ENDIF
351       
352       ALLOCATE(riverflow(klon), stat = error)
353       IF (error /= 0) THEN
354          abort_message='Pb allocation riverflow'
355          CALL abort_physic(modname,abort_message,1)
356       ENDIF
357!
358! carbon_cycle_cpl not possible with this interface and version of ORHCHIDEE
359!
360       IF (carbon_cycle_cpl) THEN
361          abort_message='carbon_cycle_cpl not yet possible with this interface of ORCHIDEE'
362          CALL abort_physic(modname,abort_message,1)
363       END IF
364       
365    ENDIF                          ! (fin debut)
366 
367
368!
369! Appel a la routine sols continentaux
370!
371    IF (lafin) lrestart_write = .TRUE.
372    IF (check) WRITE(lunout,*)'lafin ',lafin,lrestart_write
373     
374    petA_orc(1:knon) = petBcoef(1:knon) * dtime
375    petB_orc(1:knon) = petAcoef(1:knon)
376    peqA_orc(1:knon) = peqBcoef(1:knon) * dtime
377    peqB_orc(1:knon) = peqAcoef(1:knon)
378
379    cdrag = 0.
380    cdrag(1:knon) = tq_cdrag(1:knon)
381
382! zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/287.05*temp_air(1:knon))*9.80665)
383!    zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/RD*temp_air(1:knon))*RG)
384     zlev(1:knon) = plev(1:knon)*RD*temp_air(1:knon)/((ps(1:knon)*100.0)*RG)
385
386
387! PF et PASB
388!   where(cdrag > 0.01)
389!     cdrag = 0.01
390!   endwhere
391!  write(*,*)'Cdrag = ',minval(cdrag),maxval(cdrag)
392
393 
394    IF (debut) THEN
395       CALL Init_orchidee_index(knon,knindex,offset,ktindex)
396       CALL Get_orchidee_communicator(orch_comm,orch_omp_size,orch_omp_rank)
397       CALL Init_synchro_omp
398       
399       IF (knon > 0) THEN
400#ifdef CPP_VEGET
401         CALL Init_intersurf(nbp_lon,nbp_lat,knon,ktindex,offset,orch_omp_size,orch_omp_rank,orch_comm)
402#endif
403       ENDIF
404
405       
406       IF (knon > 0) THEN
407
408#ifdef CPP_VEGET
409          CALL intersurf_main (itime+itau_phy-1, nbp_lon, nbp_lat, knon, ktindex, dtime, &
410               lrestart_read, lrestart_write, lalo, &
411               contfrac, neighbours, resolution, date0, &
412               zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
413               cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
414               precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
415               evap, fluxsens, fluxlat, coastalflow, riverflow, &
416               tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
417               lon_scat, lat_scat, q2m, t2m, coszang=yrmu0)
418#endif         
419       ENDIF
420
421       CALL Synchro_omp
422
423       albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
424
425    ENDIF
426
427   
428!  swdown_vrai(1:knon) = swnet(1:knon)/(1. - albedo_keep(1:knon))
429    swdown_vrai(1:knon) = swdown(1:knon)
430
431    IF (knon > 0) THEN
432#ifdef CPP_VEGET   
433       CALL intersurf_main (itime+itau_phy, nbp_lon, nbp_lat, knon, ktindex, dtime,  &
434            lrestart_read, lrestart_write, lalo, &
435            contfrac, neighbours, resolution, date0, &
436            zlev,  u1_lay(1:knon), v1_lay(1:knon), spechum(1:knon), temp_air(1:knon), epot_air(1:knon), ccanopy(1:knon), &
437            cdrag(1:knon), petA_orc(1:knon), peqA_orc(1:knon), petB_orc(1:knon), peqB_orc(1:knon), &
438            precip_rain(1:knon), precip_snow(1:knon), lwdown(1:knon), swnet(1:knon), swdown_vrai(1:knon), ps(1:knon), &
439            evap(1:knon), fluxsens(1:knon), fluxlat(1:knon), coastalflow(1:knon), riverflow(1:knon), &
440            tsol_rad(1:knon), tsurf_new(1:knon), qsurf(1:knon), albedo_out(1:knon,:), emis_new(1:knon), z0_new(1:knon), &
441            lon_scat, lat_scat, q2m, t2m, coszang=yrmu0(1:knon))
442#endif       
443    ENDIF
444    z0h_new(:)=z0_new(:)
445
446    CALL Synchro_omp
447   
448    albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
449
450!* Send to coupler
451!
452    IF (type_ocean=='couple') THEN
453       CALL cpl_send_land_fields(itime, knon, knindex, &
454            riverflow, coastalflow)
455    ENDIF
456
457    alb1_new(1:knon) = albedo_out(1:knon,1)
458    alb2_new(1:knon) = albedo_out(1:knon,2)
459
460! Convention orchidee: positif vers le haut
461    fluxsens(1:knon) = -1. * fluxsens(1:knon)
462    fluxlat(1:knon)  = -1. * fluxlat(1:knon)
463   
464!  evap     = -1. * evap
465
466    IF (debut) lrestart_read = .FALSE.
467   
468    IF (debut) CALL Finalize_surf_para
469
470   
471  END SUBROUTINE surf_land_orchidee
472!
473!****************************************************************************************
474!
475  SUBROUTINE Init_orchidee_index(knon,knindex,offset,ktindex)
476  USE mod_surf_para
477  USE mod_grid_phy_lmdz
478 
479    INTEGER,INTENT(IN)    :: knon
480    INTEGER,INTENT(IN)    :: knindex(klon)   
481    INTEGER,INTENT(OUT)   :: offset
482    INTEGER,INTENT(OUT)   :: ktindex(klon)
483   
484    INTEGER               :: ktindex_glo(knon_glo)
485    INTEGER               :: offset_para(0:omp_size*mpi_size-1)
486    INTEGER               :: LastPoint
487    INTEGER               :: task
488   
489    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin-1)+(klon_omp_begin-1)+nbp_lon-1
490   
491    CALL gather_surf(ktindex(1:knon),ktindex_glo)
492   
493    IF (is_mpi_root .AND. is_omp_root) THEN
494      LastPoint=0
495      DO Task=0,mpi_size*omp_size-1
496        IF (knon_glo_para(Task)>0) THEN
497           offset_para(task)= LastPoint-MOD(LastPoint,nbp_lon)
498           LastPoint=ktindex_glo(knon_glo_end_para(task))
499        ENDIF
500      ENDDO
501    ENDIF
502   
503    CALL bcast(offset_para)
504   
505    offset=offset_para(omp_size*mpi_rank+omp_rank)
506   
507    ktindex(1:knon)=ktindex(1:knon)-offset
508
509  END SUBROUTINE Init_orchidee_index
510
511!
512!************************* ***************************************************************
513!
514
515  SUBROUTINE Get_orchidee_communicator(orch_comm,orch_omp_size,orch_omp_rank)
516  USE  mod_surf_para
517     
518#ifdef CPP_MPI
519    INCLUDE 'mpif.h'
520#endif   
521
522    INTEGER,INTENT(OUT) :: orch_comm
523    INTEGER,INTENT(OUT) :: orch_omp_size
524    INTEGER,INTENT(OUT) :: orch_omp_rank
525    INTEGER             :: color
526    INTEGER             :: i,ierr
527!
528! End definition
529!****************************************************************************************
530   
531   
532    IF (is_omp_root) THEN         
533     
534      IF (knon_mpi==0) THEN
535         color = 0
536      ELSE
537         color = 1
538      ENDIF
539   
540#ifdef CPP_MPI   
541      CALL MPI_COMM_SPLIT(COMM_LMDZ_PHY,color,mpi_rank,orch_comm,ierr)
542#endif
543   
544    ENDIF
545    CALL bcast_omp(orch_comm)
546   
547    IF (knon_mpi /= 0) THEN
548      orch_omp_size=0
549      DO i=0,omp_size-1
550        IF (knon_omp_para(i) /=0) THEN
551          orch_omp_size=orch_omp_size+1
552          IF (i==omp_rank) orch_omp_rank=orch_omp_size-1
553        ENDIF
554      ENDDO
555    ENDIF
556   
557   
558  END SUBROUTINE Get_orchidee_communicator
559!
560!****************************************************************************************
561
562
563  SUBROUTINE Init_neighbours(knon,neighbours,knindex,pctsrf)
564    USE mod_grid_phy_lmdz
565    USE mod_surf_para   
566    USE indice_sol_mod
567
568#ifdef CPP_MPI
569    INCLUDE 'mpif.h'
570#endif   
571
572! Input arguments
573!****************************************************************************************
574    INTEGER, INTENT(IN)                     :: knon
575    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
576    REAL, DIMENSION(klon), INTENT(IN)       :: pctsrf
577   
578! Output arguments
579!****************************************************************************************
580    INTEGER, DIMENSION(knon,8), INTENT(OUT) :: neighbours
581
582! Local variables
583!****************************************************************************************
584    INTEGER                              :: i, igrid, jj, ij, iglob
585    INTEGER                              :: ierr, ireal, index
586    INTEGER, DIMENSION(8,3)              :: off_ini
587    INTEGER, DIMENSION(8)                :: offset 
588    INTEGER, DIMENSION(nbp_lon,nbp_lat)  :: correspond
589    INTEGER, DIMENSION(knon_glo)         :: ktindex_glo
590    INTEGER, DIMENSION(knon_glo,8)       :: neighbours_glo
591    REAL, DIMENSION(klon_glo)            :: pctsrf_glo
592    INTEGER                              :: ktindex(klon)
593!
594! End definition
595!****************************************************************************************
596
597    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin-1)+(klon_omp_begin-1)+nbp_lon-1
598   
599    CALL gather_surf(ktindex(1:knon),ktindex_glo)
600    CALL gather(pctsrf,pctsrf_glo)
601   
602    IF (is_mpi_root .AND. is_omp_root) THEN
603      neighbours_glo(:,:)=-1
604!  Initialisation des offset   
605!
606! offset bord ouest
607       off_ini(1,1) = - nbp_lon   ; off_ini(2,1) = - nbp_lon + 1     ; off_ini(3,1) = 1
608       off_ini(4,1) = nbp_lon + 1 ; off_ini(5,1) = nbp_lon           ; off_ini(6,1) = 2 * nbp_lon - 1
609       off_ini(7,1) = nbp_lon -1  ; off_ini(8,1) = - 1
610! offset point normal
611       off_ini(1,2) = - nbp_lon   ; off_ini(2,2) = - nbp_lon + 1     ; off_ini(3,2) = 1
612       off_ini(4,2) = nbp_lon + 1 ; off_ini(5,2) = nbp_lon           ; off_ini(6,2) = nbp_lon - 1
613       off_ini(7,2) = -1          ; off_ini(8,2) = - nbp_lon - 1
614! offset bord   est
615       off_ini(1,3) = - nbp_lon   ; off_ini(2,3) = - 2 * nbp_lon + 1 ; off_ini(3,3) = - nbp_lon + 1
616       off_ini(4,3) =  1          ; off_ini(5,3) = nbp_lon           ; off_ini(6,3) = nbp_lon - 1
617       off_ini(7,3) = -1          ; off_ini(8,3) = - nbp_lon - 1
618!
619!
620! Attention aux poles
621!
622       DO igrid = 1, knon_glo
623          index = ktindex_glo(igrid)
624          jj = INT((index - 1)/nbp_lon) + 1
625          ij = index - (jj - 1) * nbp_lon
626          correspond(ij,jj) = igrid
627       ENDDO
628!sonia : Les mailles des voisines doivent etre toutes egales (pour couplage orchidee)
629       IF (knon_glo == 1) THEN
630         igrid = 1
631         DO i = 1,8
632           neighbours_glo(igrid, i) = igrid
633         ENDDO
634       ELSE
635       print*,'sonia : knon_glo,ij,jj', knon_glo, ij,jj
636       
637       DO igrid = 1, knon_glo
638          iglob = ktindex_glo(igrid)
639         
640          IF (MOD(iglob, nbp_lon) == 1) THEN
641             offset = off_ini(:,1)
642          ELSE IF(MOD(iglob, nbp_lon) == 0) THEN
643             offset = off_ini(:,3)
644          ELSE
645             offset = off_ini(:,2)
646          ENDIF
647         
648          DO i = 1, 8
649             index = iglob + offset(i)
650             ireal = (MIN(MAX(1, index - nbp_lon + 1), klon_glo))
651             IF (pctsrf_glo(ireal) > EPSFRA) THEN
652                jj = INT((index - 1)/nbp_lon) + 1
653                ij = index - (jj - 1) * nbp_lon
654                neighbours_glo(igrid, i) = correspond(ij, jj)
655             ENDIF
656          ENDDO
657       ENDDO
658       ENDIF !fin knon_glo == 1
659
660    ENDIF
661   
662    DO i = 1, 8
663      CALL scatter_surf(neighbours_glo(:,i),neighbours(1:knon,i))
664    ENDDO
665  END SUBROUTINE Init_neighbours
666
667!
668!****************************************************************************************
669!
670#endif
671END MODULE surf_land_orchidee_noz0h_mod
Note: See TracBrowser for help on using the repository browser.