source: LMDZ6/trunk/libf/phylmd/surf_land_orchidee_nounstruct_mod.F90 @ 5324

Last change on this file since 5324 was 5305, checked in by abarral, 2 weeks ago

Turn YOMCST2.h.h into module

  • 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    USE yomcst_mod_h
59    USE dimpft_mod_h
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!   z0m_new      surface roughness for momentum
110!   z0h_new      surface roughness for heat
111!   qsurf        air moisture at surface
112!
113
114!
115! Parametres d'entree
116!****************************************************************************************
117    INTEGER, INTENT(IN)                       :: itime
118    REAL, INTENT(IN)                          :: dtime
119    REAL, INTENT(IN)                          :: date0
120    INTEGER, INTENT(IN)                       :: knon
121    INTEGER, DIMENSION(klon), INTENT(IN)      :: knindex
122    LOGICAL, INTENT(IN)                       :: debut, lafin
123    REAL, DIMENSION(klon,nbsrf), INTENT(IN)   :: pctsrf
124    REAL, DIMENSION(klon), INTENT(IN)         :: rlon, rlat
125    REAL, DIMENSION(klon), INTENT(IN)         :: yrmu0 ! cosine of solar zenith angle
126    REAL, DIMENSION(klon), INTENT(IN)         :: plev
127    REAL, DIMENSION(klon), INTENT(IN)         :: u1_lay, v1_lay, gustiness
128    REAL, DIMENSION(klon), INTENT(IN)         :: temp_air, spechum
129    REAL, DIMENSION(klon), INTENT(IN)         :: epot_air, ccanopy
130    REAL, DIMENSION(klon), INTENT(IN)         :: tq_cdrag
131    REAL, DIMENSION(klon), INTENT(IN)         :: petAcoef, peqAcoef
132    REAL, DIMENSION(klon), INTENT(IN)         :: petBcoef, peqBcoef
133    REAL, DIMENSION(klon), INTENT(IN)         :: precip_rain, precip_snow
134    REAL, DIMENSION(klon), INTENT(IN)         :: lwdown, swnet, swdown, ps
135    REAL, DIMENSION(klon), INTENT(IN)         :: q2m, t2m
136
137! Parametres de sortie
138!****************************************************************************************
139    REAL, DIMENSION(klon), INTENT(OUT)        :: evap, fluxsens, fluxlat, qsurf
140    REAL, DIMENSION(klon), INTENT(OUT)        :: tsol_rad, tsurf_new
141    REAL, DIMENSION(klon), INTENT(OUT)        :: alb1_new, alb2_new
142    REAL, DIMENSION(klon), INTENT(OUT)        :: emis_new, z0m_new, z0h_new
143    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: veget
144    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: lai
145    REAL, DIMENSION(klon,nvm_lmdz), INTENT(OUT) :: height
146
147
148! Local
149!****************************************************************************************
150    INTEGER                                   :: ij, jj, igrid, ireal, index
151    INTEGER                                   :: error
152    REAL, DIMENSION(klon)                     :: swdown_vrai
153    CHARACTER (len = 20)                      :: modname = 'surf_land_orchidee'
154    CHARACTER (len = 80)                      :: abort_message
155    LOGICAL,SAVE                              :: check = .FALSE.
156    !$OMP THREADPRIVATE(check)
157
158! type de couplage dans sechiba
159!  character (len=10)   :: coupling = 'implicit'
160! drapeaux controlant les appels dans SECHIBA
161!  type(control_type), save   :: control_in
162! Preserved albedo
163    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: albedo_keep, zlev
164    !$OMP THREADPRIVATE(albedo_keep,zlev)
165! coordonnees geographiques
166    REAL, ALLOCATABLE, DIMENSION(:,:), SAVE   :: lalo
167    !$OMP THREADPRIVATE(lalo)
168! pts voisins
169    INTEGER,ALLOCATABLE, DIMENSION(:,:), SAVE :: neighbours
170    !$OMP THREADPRIVATE(neighbours)
171! fractions continents
172    REAL,ALLOCATABLE, DIMENSION(:), SAVE      :: contfrac
173    !$OMP THREADPRIVATE(contfrac)
174! resolution de la grille
175    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: resolution
176    !$OMP THREADPRIVATE(resolution)
177
178    REAL, ALLOCATABLE, DIMENSION (:,:), SAVE  :: lon_scat, lat_scat 
179    !$OMP THREADPRIVATE(lon_scat,lat_scat)
180
181    LOGICAL, SAVE                             :: lrestart_read = .TRUE.
182    !$OMP THREADPRIVATE(lrestart_read)
183    LOGICAL, SAVE                             :: lrestart_write = .FALSE.
184    !$OMP THREADPRIVATE(lrestart_write)
185
186    REAL, DIMENSION(knon,2)                   :: albedo_out
187
188! Pb de nomenclature
189    REAL, DIMENSION(klon)                     :: petA_orc, peqA_orc
190    REAL, DIMENSION(klon)                     :: petB_orc, peqB_orc
191! Pb de correspondances de grilles
192    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  :: ig, jg
193    !$OMP THREADPRIVATE(ig,jg)
194    INTEGER :: indi, indj
195    INTEGER, SAVE, ALLOCATABLE,DIMENSION(:)   :: ktindex
196    !$OMP THREADPRIVATE(ktindex)
197
198! Essai cdrag
199    REAL, DIMENSION(klon)                     :: cdrag
200    INTEGER,SAVE                              :: offset
201    !$OMP THREADPRIVATE(offset)
202
203    REAL, DIMENSION(klon_glo)                 :: rlon_g,rlat_g
204    INTEGER, SAVE                             :: orch_comm
205    !$OMP THREADPRIVATE(orch_comm)
206
207    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: coastalflow
208    !$OMP THREADPRIVATE(coastalflow)
209    REAL, ALLOCATABLE, DIMENSION(:), SAVE     :: riverflow
210    !$OMP THREADPRIVATE(riverflow)
211   
212    INTEGER :: orch_omp_rank
213    INTEGER :: orch_omp_size
214!
215! Fin definition
216!****************************************************************************************
217
218    IF (check) WRITE(lunout,*)'Entree ', modname
219 
220! Initialisation
221 
222    IF (debut) THEN
223! Test of coherence between variable ok_veget and cpp key CPP_VEGET
224#ifndef CPP_VEGET
225       abort_message='Pb de coherence: ok_veget = .true. mais CPP_VEGET = .false.'
226       CALL abort_physic(modname,abort_message,1)
227#endif
228
229       CALL Init_surf_para(knon)
230       ALLOCATE(ktindex(knon))
231       IF ( .NOT. ALLOCATED(albedo_keep)) THEN
232!ym          ALLOCATE(albedo_keep(klon))
233!ym bizarre que non allou� en knon precedement
234          ALLOCATE(albedo_keep(knon))
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! 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! carbon_cycle_cpl not possible with this interface and version of ORHCHIDEE
362!
363       IF (carbon_cycle_cpl) THEN
364          abort_message='carbon_cycle_cpl not yet possible with this interface of ORCHIDEE'
365          CALL abort_physic(modname,abort_message,1)
366       END IF
367       
368    ENDIF                          ! (fin debut)
369 
370
371!
372! Appel a la routine sols continentaux
373!
374    IF (lafin) lrestart_write = .TRUE.
375    IF (check) WRITE(lunout,*)'lafin ',lafin,lrestart_write
376     
377    petA_orc(1:knon) = petBcoef(1:knon) * dtime
378    petB_orc(1:knon) = petAcoef(1:knon)
379    peqA_orc(1:knon) = peqBcoef(1:knon) * dtime
380    peqB_orc(1:knon) = peqAcoef(1:knon)
381
382    cdrag = 0.
383    cdrag(1:knon) = tq_cdrag(1:knon)
384
385! zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/287.05*temp_air(1:knon))*9.80665)
386!    zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/RD*temp_air(1:knon))*RG)
387     zlev(1:knon) = plev(1:knon)*RD*temp_air(1:knon)/((ps(1:knon)*100.0)*RG)
388
389
390! PF et PASB
391!   where(cdrag > 0.01)
392!     cdrag = 0.01
393!   endwhere
394!  write(*,*)'Cdrag = ',minval(cdrag),maxval(cdrag)
395
396 
397    IF (debut) THEN
398       CALL Init_orchidee_index(knon,knindex,offset,ktindex)
399       CALL Get_orchidee_communicator(orch_comm,orch_omp_size,orch_omp_rank)
400       CALL Init_synchro_omp
401       
402       IF (knon > 0) THEN
403#ifdef CPP_VEGET
404         CALL Init_intersurf(nbp_lon,nbp_lat,knon,ktindex,offset,orch_omp_size,orch_omp_rank,orch_comm)
405#endif
406       ENDIF
407
408       
409       IF (knon > 0) THEN
410
411#ifdef CPP_VEGET
412         CALL intersurf_initialize_gathered (itime+itau_phy-1, nbp_lon, nbp_lat, knon, ktindex, dtime, &
413               lrestart_read, lrestart_write, lalo, contfrac, neighbours, resolution, date0, &
414               zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, &
415               cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
416               precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
417               evap, fluxsens, fluxlat, coastalflow, riverflow, &
418               tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0m_new, &   
419               lon_scat, lat_scat, q2m, t2m, z0h_new, nvm_orch)
420#endif         
421       ENDIF
422
423       CALL Synchro_omp
424
425       albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
426
427    ENDIF
428
429   
430!  swdown_vrai(1:knon) = swnet(1:knon)/(1. - albedo_keep(1:knon))
431    swdown_vrai(1:knon) = swdown(1:knon)
432
433    IF (knon > 0) THEN
434#ifdef CPP_VEGET   
435       IF (nvm_orch .NE. nvm_lmdz ) THEN
436          abort_message='Pb de dimensiosn PFT: nvm_orch et nvm_lmdz differents.'
437          CALL abort_physic(modname,abort_message,1)
438       ENDIF
439
440       CALL intersurf_main_gathered (itime+itau_phy, nbp_lon, nbp_lat, knon, ktindex, dtime,  &
441            lrestart_read, lrestart_write, lalo, &
442            contfrac, neighbours, resolution, date0, &
443            zlev,  u1_lay(1:knon), v1_lay(1:knon), spechum(1:knon), temp_air(1:knon), epot_air(1:knon), ccanopy(1:knon), &
444            cdrag(1:knon), petA_orc(1:knon), peqA_orc(1:knon), petB_orc(1:knon), peqB_orc(1:knon), &
445            precip_rain(1:knon), precip_snow(1:knon), lwdown(1:knon), swnet(1:knon), swdown_vrai(1:knon), ps(1:knon), &
446            evap(1:knon), fluxsens(1:knon), fluxlat(1:knon), coastalflow(1:knon), riverflow(1:knon), &
447            tsol_rad(1:knon), tsurf_new(1:knon), qsurf(1:knon), albedo_out(1:knon,:), emis_new(1:knon), z0m_new(1:knon), &
448            lon_scat, lat_scat, q2m, t2m, z0h_new(1:knon),&
449            veget(1:knon,:),lai(1:knon,:),height(1:knon,:),&
450            coszang=yrmu0(1:knon))
451#endif       
452    ENDIF
453
454    CALL Synchro_omp
455   
456    albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
457
458!* Send to coupler
459!
460    IF (type_ocean=='couple') THEN
461       CALL cpl_send_land_fields(itime, knon, knindex, &
462            riverflow, coastalflow)
463    ENDIF
464
465    alb1_new(1:knon) = albedo_out(1:knon,1)
466    alb2_new(1:knon) = albedo_out(1:knon,2)
467
468! Convention orchidee: positif vers le haut
469    fluxsens(1:knon) = -1. * fluxsens(1:knon)
470    fluxlat(1:knon)  = -1. * fluxlat(1:knon)
471   
472!  evap     = -1. * evap
473
474    IF (debut) lrestart_read = .FALSE.
475   
476    IF (debut) CALL Finalize_surf_para
477
478   
479  END SUBROUTINE surf_land_orchidee
480!
481!****************************************************************************************
482!
483  SUBROUTINE Init_orchidee_index(knon,knindex,offset,ktindex)
484  USE mod_surf_para
485  USE mod_grid_phy_lmdz
486 
487    INTEGER,INTENT(IN)    :: knon
488    INTEGER,INTENT(IN)    :: knindex(klon)   
489    INTEGER,INTENT(OUT)   :: offset
490    INTEGER,INTENT(OUT)   :: ktindex(klon)
491   
492    INTEGER               :: ktindex_glo(knon_glo)
493    INTEGER               :: offset_para(0:omp_size*mpi_size-1)
494    INTEGER               :: LastPoint
495    INTEGER               :: task
496   
497    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin-1)+(klon_omp_begin-1)+nbp_lon-1
498   
499    CALL gather_surf(ktindex(1:knon),ktindex_glo)
500   
501    IF (is_mpi_root .AND. is_omp_root) THEN
502      LastPoint=0
503      DO Task=0,mpi_size*omp_size-1
504        IF (knon_glo_para(Task)>0) THEN
505           offset_para(task)= LastPoint-MOD(LastPoint,nbp_lon)
506           LastPoint=ktindex_glo(knon_glo_end_para(task))
507        ENDIF
508      ENDDO
509    ENDIF
510   
511    CALL bcast(offset_para)
512   
513    offset=offset_para(omp_size*mpi_rank+omp_rank)
514   
515    ktindex(1:knon)=ktindex(1:knon)-offset
516
517  END SUBROUTINE Init_orchidee_index
518
519!
520!************************* ***************************************************************
521!
522
523  SUBROUTINE Get_orchidee_communicator(orch_comm,orch_omp_size,orch_omp_rank)
524  USE  mod_surf_para
525  USE lmdz_mpi
526     
527    INTEGER,INTENT(OUT) :: orch_comm
528    INTEGER,INTENT(OUT) :: orch_omp_size
529    INTEGER,INTENT(OUT) :: orch_omp_rank
530    INTEGER             :: color
531    INTEGER             :: i,ierr
532!
533! End definition
534!****************************************************************************************
535   
536   
537    IF (is_omp_root) THEN         
538     
539      IF (knon_mpi==0) THEN
540         color = 0
541      ELSE
542         color = 1
543      ENDIF
544   
545      IF (using_mpi) THEN
546        CALL MPI_COMM_SPLIT(COMM_LMDZ_PHY,color,mpi_rank,orch_comm,ierr)
547      ENDIF
548   
549    ENDIF
550    CALL bcast_omp(orch_comm)
551   
552    IF (knon_mpi /= 0) THEN
553      orch_omp_size=0
554      DO i=0,omp_size-1
555        IF (knon_omp_para(i) /=0) THEN
556          orch_omp_size=orch_omp_size+1
557          IF (i==omp_rank) orch_omp_rank=orch_omp_size-1
558        ENDIF
559      ENDDO
560    ENDIF
561   
562   
563  END SUBROUTINE Get_orchidee_communicator
564!
565!****************************************************************************************
566
567
568  SUBROUTINE Init_neighbours(knon,neighbours,knindex,pctsrf)
569    USE mod_grid_phy_lmdz
570    USE mod_surf_para   
571    USE indice_sol_mod
572    USE lmdz_mpi
573
574! Input arguments
575!****************************************************************************************
576    INTEGER, INTENT(IN)                     :: knon
577    INTEGER, DIMENSION(klon), INTENT(IN)    :: knindex
578    REAL, DIMENSION(klon), INTENT(IN)       :: pctsrf
579   
580! Output arguments
581!****************************************************************************************
582    INTEGER, DIMENSION(knon,8), INTENT(OUT) :: neighbours
583
584! Local variables
585!****************************************************************************************
586    INTEGER                              :: i, igrid, jj, ij, iglob
587    INTEGER                              :: ierr, ireal, index
588    INTEGER, DIMENSION(8,3)              :: off_ini
589    INTEGER, DIMENSION(8)                :: offset 
590    INTEGER, DIMENSION(nbp_lon,nbp_lat)  :: correspond
591    INTEGER, DIMENSION(knon_glo)         :: ktindex_glo
592    INTEGER, DIMENSION(knon_glo,8)       :: neighbours_glo
593    REAL, DIMENSION(klon_glo)            :: pctsrf_glo
594    INTEGER                              :: ktindex(klon)
595!
596! End definition
597!****************************************************************************************
598
599    ktindex(1:knon)=knindex(1:knon)+(klon_mpi_begin-1)+(klon_omp_begin-1)+nbp_lon-1
600   
601    CALL gather_surf(ktindex(1:knon),ktindex_glo)
602    CALL gather(pctsrf,pctsrf_glo)
603   
604    IF (is_mpi_root .AND. is_omp_root) THEN
605      neighbours_glo(:,:)=-1
606!  Initialisation des offset   
607!
608! offset bord ouest
609       off_ini(1,1) = - nbp_lon   ; off_ini(2,1) = - nbp_lon + 1     ; off_ini(3,1) = 1
610       off_ini(4,1) = nbp_lon + 1 ; off_ini(5,1) = nbp_lon           ; off_ini(6,1) = 2 * nbp_lon - 1
611       off_ini(7,1) = nbp_lon -1  ; off_ini(8,1) = - 1
612! offset point normal
613       off_ini(1,2) = - nbp_lon   ; off_ini(2,2) = - nbp_lon + 1     ; off_ini(3,2) = 1
614       off_ini(4,2) = nbp_lon + 1 ; off_ini(5,2) = nbp_lon           ; off_ini(6,2) = nbp_lon - 1
615       off_ini(7,2) = -1          ; off_ini(8,2) = - nbp_lon - 1
616! offset bord   est
617       off_ini(1,3) = - nbp_lon   ; off_ini(2,3) = - 2 * nbp_lon + 1 ; off_ini(3,3) = - nbp_lon + 1
618       off_ini(4,3) =  1          ; off_ini(5,3) = nbp_lon           ; off_ini(6,3) = nbp_lon - 1
619       off_ini(7,3) = -1          ; off_ini(8,3) = - nbp_lon - 1
620!
621!
622! Attention aux poles
623!
624       DO igrid = 1, knon_glo
625          index = ktindex_glo(igrid)
626          jj = INT((index - 1)/nbp_lon) + 1
627          ij = index - (jj - 1) * nbp_lon
628          correspond(ij,jj) = igrid
629       ENDDO
630!sonia : Les mailles des voisines doivent etre toutes egales (pour couplage orchidee)
631       IF (knon_glo == 1) THEN
632         igrid = 1
633         DO i = 1,8
634           neighbours_glo(igrid, i) = igrid
635         ENDDO
636       ELSE
637       print*,'sonia : knon_glo,ij,jj', knon_glo, ij,jj
638       
639       DO igrid = 1, knon_glo
640          iglob = ktindex_glo(igrid)
641         
642          IF (MOD(iglob, nbp_lon) == 1) THEN
643             offset = off_ini(:,1)
644          ELSE IF(MOD(iglob, nbp_lon) == 0) THEN
645             offset = off_ini(:,3)
646          ELSE
647             offset = off_ini(:,2)
648          ENDIF
649         
650          DO i = 1, 8
651             index = iglob + offset(i)
652             ireal = (MIN(MAX(1, index - nbp_lon + 1), klon_glo))
653             IF (pctsrf_glo(ireal) > EPSFRA) THEN
654                jj = INT((index - 1)/nbp_lon) + 1
655                ij = index - (jj - 1) * nbp_lon
656                neighbours_glo(igrid, i) = correspond(ij, jj)
657             ENDIF
658          ENDDO
659       ENDDO
660       ENDIF !fin knon_glo == 1
661
662    ENDIF
663   
664    DO i = 1, 8
665      CALL scatter_surf(neighbours_glo(:,i),neighbours(1:knon,i))
666    ENDDO
667  END SUBROUTINE Init_neighbours
668
669!
670!****************************************************************************************
671!
672#endif
673#endif
674END MODULE surf_land_orchidee_nounstruct_mod
Note: See TracBrowser for help on using the repository browser.