source: LMDZ6/trunk/libf/phylmd/surf_land_orchidee_nofrein_mod.F90 @ 5301

Last change on this file since 5301 was 5296, checked in by abarral, 46 hours ago

Turn compbl.h into a module
Move calcul_REGDYN.h to obsolete
Create phys_constants_mod.f90

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