SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
get_adj_mes_gauss.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
5 ! #########
6  SUBROUTINE get_adj_mes_gauss(KGRID_PAR,KL,PGRID_PAR,KLEFT,KRIGHT,KTOP,KBOTTOM)
7 ! ##############################################################
8 !
9 !!**** *GET_ADJACENT_MESHES_GAUSS* get the near grid mesh indices
10 !!
11 !! PURPOSE
12 !! -------
13 !!
14 !! METHOD
15 !! ------
16 !!
17 !! REFERENCE
18 !! ---------
19 !!
20 !! AUTHOR
21 !! ------
22 !!
23 !! V. Masson Meteo-France
24 !!
25 !! MODIFICATION
26 !! ------------
27 !!
28 !! Original 03/2004
29 !!
30 !----------------------------------------------------------------------------
31 !
32 !* 0. DECLARATION
33 ! -----------
34 !
36 !
37 !
38 USE yomhook ,ONLY : lhook, dr_hook
39 USE parkind1 ,ONLY : jprb
40 !
41 IMPLICIT NONE
42 !
43 !* 0.1 Declaration of arguments
44 ! ------------------------
45 !
46 INTEGER, INTENT(IN) :: kgrid_par ! size of PGRID_PAR
47 INTEGER, INTENT(IN) :: kl ! number of points
48 REAL, DIMENSION(KGRID_PAR), INTENT(IN) :: pgrid_par ! grid parameters
49 INTEGER, DIMENSION(KL), INTENT(OUT) :: kleft ! left mesh index
50 INTEGER, DIMENSION(KL), INTENT(OUT) :: kright ! right mesh index
51 INTEGER, DIMENSION(KL), INTENT(OUT) :: ktop ! top mesh index
52 INTEGER, DIMENSION(KL), INTENT(OUT) :: kbottom ! bottom mesh index
53 !
54 !* 0.2 Declaration of other local variables
55 ! ------------------------------------
56 !
57 INTEGER :: inlati ! number of pseudo-latitudes
58 REAL :: zlapo ! latitude of the rotated pole (deg)
59 REAL :: zlopo ! logitude of the rotated pole (deg)
60 REAL :: zcodil ! stretching factor
61 INTEGER, DIMENSION(:),ALLOCATABLE :: inlopa ! number of pseudo-longitudes
62 ! ! on each pseudo-latitude circle
63 ! ! on pseudo-northern hemisphere
64 ! ! (starting from the rotated pole)
65 REAL, DIMENSION(:),ALLOCATABLE :: zxcen
66 !
67 INTEGER :: jlat, jlon, il, jl, ilgrid, jlon2, id, jl0
68 !
69 REAL :: zdis, zinter
70 !
71 REAL(KIND=JPRB) :: zhook_handle
72 !
73 !----------------------------------------------------------------------------
74 !
75 IF (lhook) CALL dr_hook('GET_ADJ_MES_GAUSS',0,zhook_handle)
76  CALL get_gridtype_gauss(pgrid_par,inlati,zlapo,zlopo,zcodil)
77 !
78 ALLOCATE(inlopa(0:inlati))
79 ALLOCATE(zxcen(kl))
80 !
81  CALL get_gridtype_gauss(pgrid_par,inlati,zlapo,zlopo,zcodil,inlopa(1:inlati), &
82  plon_xy=zxcen )
83 !
84 inlopa(0) = 0
85 !
86 kleft(:) = 0
87 kright(:) = 0
88 ktop(:) = 0
89 kbottom(:) = 0
90 !
91 il=0
92 DO jlat=1,inlati
93  DO jlon=1,inlopa(jlat)
94  il=il+1
95  ENDDO
96 ENDDO
97 !
98 jl = 0.0
99 IF (il==kl) THEN
100  DO jlat=1,inlati
101  !
102  jl0 = jl
103  !
104  DO jlon=1,inlopa(jlat)
105  !
106  jl = jl + 1
107  !
108  IF (jlon>1 ) kleft(jl) = jl-1
109  IF (jlon<inlopa(jlat) ) kright(jl) = jl+1
110  !
111  IF (jlon==1 ) kleft(jl) = jl+inlopa(jlat)-1
112  IF (jlon==inlopa(jlat)) kright(jl) = jl-inlopa(jlat)+1
113  !
114  IF (jlat>1 ) THEN
115  zdis = abs(zxcen(jl) - zxcen(jl0 - inlopa(jlat-1) + 1))
116  id = 1
117  DO jlon2 = 1,inlopa(jlat-1)
118  zinter = abs(zxcen(jl) - zxcen(jl0 - inlopa(jlat-1) + jlon2))
119  IF (zinter<zdis) THEN
120  zdis = zinter
121  id = jlon2
122  ENDIF
123  ENDDO
124  ktop(jl) = jl0 - inlopa(jlat-1) + id
125  ENDIF
126  !
127  IF (jlat<inlati ) THEN
128  zdis = abs(zxcen(jl) - zxcen(jl0 + inlopa(jlat) + 1))
129  id = 1
130  DO jlon2 = 1,inlopa(jlat+1)
131  zinter = abs(zxcen(jl) - zxcen(jl0 + inlopa(jlat) + jlon2))
132  IF (zinter<zdis) THEN
133  zdis = zinter
134  id = jlon2
135  ENDIF
136  ENDDO
137  kbottom(jl) = jl0 + inlopa(jlat) + id
138  ENDIF
139  !
140  END DO
141  END DO
142 END IF
143 !
144 DEALLOCATE(inlopa)
145 DEALLOCATE(zxcen)
146 !
147 IF (lhook) CALL dr_hook('GET_ADJ_MES_GAUSS',1,zhook_handle)
148 !
149 !-------------------------------------------------------------------------------
150 !
151 END SUBROUTINE get_adj_mes_gauss
subroutine get_gridtype_gauss(PGRID_PAR, KNLATI, PLAPO, PLOPO, PCODIL, KNLOPA, KL, PLAT, PLON, PLAT_XY, PLON_XY, PMESH_SIZE, PLONINF, PLATINF, PLONSUP, PLATSUP)
subroutine get_adj_mes_gauss(KGRID_PAR, KL, PGRID_PAR, KLEFT, KRIGHT, KTOP, KBOTTOM)