SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
mode_gridtype_lonlat_rot.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 ! ##############################
7 ! ##############################
8 !
9 !############################################################################
10 !############################################################################
11 !############################################################################
12 !
13 USE yomhook ,ONLY : lhook, dr_hook
14 USE parkind1 ,ONLY : jprb
15 !
16  CONTAINS
17 !############################################################################
18 !############################################################################
19 !############################################################################
20 ! ####################################################################
21  SUBROUTINE put_gridtype_lonlat_rot(PGRID_PAR, &
22  pwest,psouth,pdlon,pdlat,ppolon,ppolat, &
23  klon,klat,kl,plon,plat )
24 ! ####################################################################
25 !
26 !!**** *PUT_GRIDTYPE_LONLAT_ROT* - routine to store in PGRID_PAR the horizontal grid
27 !!
28 !! AUTHOR
29 !! ------
30 !! P. Samuelsson SMHI
31 !!
32 !! MODIFICATIONS
33 !! -------------
34 !! Original 12/2012
35 !-------------------------------------------------------------------------------
36 !
37 !* 0. DECLARATIONS
38 ! ------------
39 !
40 IMPLICIT NONE
41 !
42 !
43 !* 0.1 Declarations of arguments
44 ! -------------------------
45 !
46 REAL, INTENT(IN) :: pwest ! West longitude in rotated grid (degrees)
47 REAL, INTENT(IN) :: psouth ! South latitude in rotated grid (degrees)
48 REAL, INTENT(IN) :: pdlon ! Longitudal grid spacing (degrees)
49 REAL, INTENT(IN) :: pdlat ! Latitudal grid spacing (degrees)
50 REAL, INTENT(IN) :: ppolon ! Longitude of rotated pole (degrees)
51 REAL, INTENT(IN) :: ppolat ! Latitude of rotated pole (degrees)
52 INTEGER, INTENT(IN) :: klon ! number of points in longitude
53 INTEGER, INTENT(IN) :: klat ! number of points in latitude
54 INTEGER, INTENT(IN) :: kl ! number of points used
55 REAL, DIMENSION(:), INTENT(IN) :: plon ! regular longitudes of all points
56 REAL, DIMENSION(:), INTENT(IN) :: plat ! regular latitudes of all points
57 REAL, DIMENSION(:),POINTER :: pgrid_par! parameters defining this grid
58 REAL(KIND=JPRB) :: zhook_handle
59 !
60 !
61 !* 0.2 Declarations of local variables
62 ! -------------------------------
63 !
64 !-------------------------------------------------------------------------------
65 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:PUT_GRIDTYPE_LONLAT_ROT',0,zhook_handle)
66 ALLOCATE(pgrid_par(9+2*kl))
67 pgrid_par(1) = pwest
68 pgrid_par(2) = psouth
69 pgrid_par(3) = pdlon
70 pgrid_par(4) = pdlat
71 pgrid_par(5) = ppolon
72 pgrid_par(6) = ppolat
73 pgrid_par(7) = float(klon)
74 pgrid_par(8) = float(klat)
75 pgrid_par(9) = float(kl)
76 pgrid_par(10:9+kl) = plon(:)
77 pgrid_par(10+kl:9+2*kl) = plat(:)
78 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:PUT_GRIDTYPE_LONLAT_ROT',1,zhook_handle)
79 !-------------------------------------------------------------------------------
80 END SUBROUTINE put_gridtype_lonlat_rot
81 !############################################################################
82 !############################################################################
83 !############################################################################
84 ! ####################################################################
85  SUBROUTINE get_gridtype_lonlat_rot(PGRID_PAR, &
86  pwest,psouth,pdlon,pdlat,ppolon,ppolat, &
87  klon,klat,kl,plon,plat )
88 ! ####################################################################
89 !
90 !!**** *GET_GRIDTYPE_LONLAT_ROT* - routine to get from PGRID_PAR the horizontal grid
91 !!
92 !! AUTHOR
93 !! ------
94 !! P. Samuelsson SMHI
95 !!
96 !! MODIFICATIONS
97 !! -------------
98 !! Original 12/2012
99 !-------------------------------------------------------------------------------
100 !
101 !* 0. DECLARATIONS
102 ! ------------
103 !
104 IMPLICIT NONE
105 !
106 !
107 !* 0.1 Declarations of arguments
108 ! -------------------------
109 !
110 REAL, DIMENSION(:), INTENT(IN) :: pgrid_par! parameters defining this grid
111 REAL, INTENT(OUT), OPTIONAL :: pwest ! West longitude in rotated grid (degrees)
112 REAL, INTENT(OUT), OPTIONAL :: psouth ! South latitude in rotated grid (degrees)
113 REAL, INTENT(OUT), OPTIONAL :: pdlon ! Longitudal grid spacing (degrees)
114 REAL, INTENT(OUT), OPTIONAL :: pdlat ! Latitudal grid spacing (degrees)
115 REAL, INTENT(OUT), OPTIONAL :: ppolon ! Longitude of rotated pole (degrees)
116 REAL, INTENT(OUT), OPTIONAL :: ppolat ! Latitude of rotated pole (degrees)
117 INTEGER, INTENT(OUT), OPTIONAL :: klon ! number of points in longitude
118 INTEGER, INTENT(OUT), OPTIONAL :: klat ! number of points in latitude
119 INTEGER, INTENT(OUT), OPTIONAL :: kl ! number of points used
120 REAL, DIMENSION(:), INTENT(OUT), OPTIONAL :: plon ! longitudes of all points
121 REAL, DIMENSION(:), INTENT(OUT), OPTIONAL :: plat ! latitudes of all points
122 !
123 !
124 !* 0.2 Declarations of local variables
125 ! -------------------------------
126 !
127 INTEGER :: ilon, ilat
128 INTEGER :: il
129 REAL(KIND=JPRB) :: zhook_handle
130 !-------------------------------------------------------------------------------
131 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:GET_GRIDTYPE_LONLAT_ROT',0,zhook_handle)
132 ilon = nint(pgrid_par(7))
133 ilat = nint(pgrid_par(8))
134 il = nint(pgrid_par(9))
135 !
136 IF (present(pwest )) pwest = pgrid_par(1)
137 IF (present(psouth )) psouth = pgrid_par(2)
138 IF (present(pdlon )) pdlon = pgrid_par(3)
139 IF (present(pdlat )) pdlat = pgrid_par(4)
140 IF (present(ppolon )) ppolon = pgrid_par(5)
141 IF (present(ppolat )) ppolat = pgrid_par(6)
142 IF (present(klon )) klon = ilon
143 IF (present(klat )) klat = ilat
144 IF (present(kl )) kl = il
145 IF (present(plon )) plon(:) = pgrid_par(10:9+il)
146 IF (present(plat )) plat(:) = pgrid_par(10+il:9+2*il)
147 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:GET_GRIDTYPE_LONLAT_ROT',1,zhook_handle)
148 !
149 !-------------------------------------------------------------------------------
150 END SUBROUTINE get_gridtype_lonlat_rot
151 !---------------------------------------------------------------------------------
152 !
153 ! ###################################################
154  SUBROUTINE latlon_lonlat_rot(PWEST,PSOUTH,PDLON,PDLAT,PPOLON,PPOLAT, &
155  klon,klat,plon,plat )
156 ! ###################################################
157 !
158 !!**** *LATLON_LONLAT_ROT * - Routine to compute regular geographical coordinates
159 !!
160 !! PURPOSE
161 !! -------
162 !!
163 !! AUTHOR
164 !! ------
165 !! P. Samuelsson SMHI
166 !!
167 !! MODIFICATIONS
168 !! -------------
169 !! Original 12/2012
170 !-------------------------------------------------------------------------------
171 !
172 !* 0. DECLARATIONS
173 ! ------------
174 !
175 !RJ: missing modi
176 USE modi_regrot_lonlat_rot
177 !
178 IMPLICIT NONE
179 !
180 !* 0.1 Declarations of arguments and results
181 !
182 REAL, INTENT(IN) :: pwest ! West longitude in rotated grid (degrees)
183 REAL, INTENT(IN) :: psouth ! South latitude in rotated grid (degrees)
184 REAL, INTENT(IN) :: pdlon ! Longitudal grid spacing (degrees)
185 REAL, INTENT(IN) :: pdlat ! Latitudal grid spacing (degrees)
186 REAL, INTENT(IN) :: ppolon ! Longitude of rotated south pole (degrees)
187 REAL, INTENT(IN) :: ppolat ! Latitude of rotated south pole (degrees)
188 INTEGER, INTENT(IN) :: klon ! number of points in longitude
189 INTEGER, INTENT(IN) :: klat ! number of points in latitude
190 REAL, DIMENSION(:), INTENT(OUT) :: plon,plat
191  ! returned geographic latitudes and
192  ! longitudes of the processed points
193  ! (degrees).
194 !
195 !* 0.2 Declarations of local variables
196 !
197 INTEGER :: jlon, jlat
198 INTEGER :: jl, kl
199 INTEGER :: jc
200 
201 REAL, DIMENSION(:), ALLOCATABLE :: zlatrot ! rotated latitude of all points
202 REAL, DIMENSION(:), ALLOCATABLE :: zlonrot ! rotated longitude of all points
203 
204 
205 REAL(KIND=JPRB) :: zhook_handle
206 !--------------------------------------------------------------------------------
207 !
208 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:LATLON_LONLAT_ROT',0,zhook_handle)
209 
210 kl = klon * klat
211 
212 ALLOCATE(zlatrot(kl))
213 ALLOCATE(zlonrot(kl))
214 
215 !
216 DO jlat=1,klat
217  DO jlon=1,klon
218  jl = jlon + klon * (jlat-1)
219  zlatrot(jl) = psouth + pdlat * (jlat-1)
220  zlonrot(jl) = pwest + pdlon * (jlon-1)
221  END DO
222 END DO
223 
224  CALL regrot_lonlat_rot(plon,plat,zlonrot,zlatrot, &
225  kl,1,kl,1, &
226  ppolon,ppolat,-1 )
227 
228 WHERE (plon(:)>180.) plon(:)=plon(:)-360.
229 WHERE (plon(:)<-180.) plon(:)=plon(:)+360.
230 
231 IF (lhook) CALL dr_hook('MODE_GRIDTYPE_LONLAT_ROT:LATLON_LONLAT_ROT',1,zhook_handle)
232 !---------------------------------------------------------------------------------
233 END SUBROUTINE latlon_lonlat_rot
234 !---------------------------------------------------------------------------------
235 !
236 !############################################################################
237 !############################################################################
238 !############################################################################
239 !############################################################################
240 !############################################################################
241 !############################################################################
242 !
243 END MODULE mode_gridtype_lonlat_rot
subroutine get_gridtype_lonlat_rot(PGRID_PAR, PWEST, PSOUTH, PDLON, PDLAT, PPOLON, PPOLAT, KLON, KLAT, KL, PLON, PLAT)
subroutine latlon_lonlat_rot(PWEST, PSOUTH, PDLON, PDLAT, PPOLON, PPOLAT, KLON, KLAT, PLON, PLAT)
subroutine regrot_lonlat_rot(PXREG, PYREG, PXROT, PYROT, KXDIM, KYDIM, KX, KY, PXCEN, PYCEN, KCALL)
subroutine put_gridtype_lonlat_rot(PGRID_PAR, PWEST, PSOUTH, PDLON, PDLAT, PPOLON, PPOLAT, KLON, KLAT, KL, PLON, PLAT)