source: LMDZ6/branches/Ocean_skin/libf/phylmd/surf_land_orchidee_nounstruct_mod.F90 @ 3605

Last change on this file since 3605 was 3605, checked in by lguez, 4 years ago

Merge revisions 3427:3600 of trunk into branch Ocean_skin

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