NCEPLIBS-ip 4.0.0
ip_equid_cylind_grid_mod.f90
Go to the documentation of this file.
1
5
12 use ip_grid_mod
13 use earth_radius_mod
14 implicit none
15
16 private
17 public :: ip_equid_cylind_grid
18
19 type, extends(ip_grid) :: ip_equid_cylind_grid
20 real :: hi, rlat1, rlon1, rlat2, rlon2
21 real :: dlat, dlon
22 contains
23 procedure :: init_grib1
24 procedure :: init_grib2
25 procedure :: gdswzd => gdswzd_equid_cylind
27
28 REAL :: DLAT ! GRID RESOLUTION IN DEGREES N/S DIRECTION
29 REAL :: DLON ! GRID RESOLUTION IN DEGREES E/W DIRECTION
30 REAL :: RERTH
31
32contains
33
41 subroutine init_grib1(self, g1_desc)
42 class(ip_equid_cylind_grid), intent(inout) :: self
43 type(grib1_descriptor), intent(in) :: g1_desc
44
45 integer :: iscan
46
47 associate(kgds => g1_desc%gds)
48 self%IM=kgds(2)
49 self%JM=kgds(3)
50 self%RLAT1=kgds(4)*1.e-3
51 self%RLON1=kgds(5)*1.e-3
52 self%RLAT2=kgds(7)*1.e-3
53 self%RLON2=kgds(8)*1.e-3
54 iscan=mod(kgds(11)/128,2)
55 self%HI=(-1.)**iscan
56 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
57 self%DLAT=(self%RLAT2-self%RLAT1)/(self%JM-1)
58
59 ! defaults
60 self%iwrap = 0
61 self%jwrap1 = 0
62 self%jwrap2 = 0
63 self%nscan = mod(kgds(11) / 32, 2)
64 self%nscan_field_pos = self%nscan
65 self%kscan = 0
66
67 self%iwrap = nint(360/abs(self%dlon))
68
69 if(self%im < self%iwrap) self%iwrap=0
70 self%jwrap1 = 0
71 self%jwrap2 = 0
72 if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0) then
73 if(abs(self%rlat1) > 90-0.25*self%dlat) then
74 self%jwrap1 = 2
75 elseif(abs(self%rlat1) > 90-0.75*self%dlat) then
76 self%jwrap1 = 1
77 endif
78 if(abs(self%rlat2) > 90-0.25*self%dlat) then
79 self%jwrap2 = 2 * self%jm
80 elseif(abs(self%rlat2) > 90-0.75*self%dlat) then
81 self%jwrap2 = 2 * self%jm+1
82 endif
83 endif
84
85 self%rerth = 6.3712e6
86 self%eccen_squared = 0.0
87 end associate
88
89 end subroutine init_grib1
90
91
98 subroutine init_grib2(self, g2_desc)
99 class(ip_equid_cylind_grid), intent(inout) :: self
100 type(grib2_descriptor), intent(in) :: g2_desc
101
102 integer :: iscale, iscan
103
104 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
105 self%IM=igdtmpl(8)
106 self%JM=igdtmpl(9)
107 iscale=igdtmpl(10)*igdtmpl(11)
108 IF(iscale==0) iscale=10**6
109 self%RLAT1=float(igdtmpl(12))/float(iscale)
110 self%RLON1=float(igdtmpl(13))/float(iscale)
111 self%RLAT2=float(igdtmpl(15))/float(iscale)
112 self%RLON2=float(igdtmpl(16))/float(iscale)
113 iscan=mod(igdtmpl(19)/128,2)
114 self%HI=(-1.)**iscan
115 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
116 self%DLAT=(self%RLAT2-self%RLAT1)/(self%JM-1)
117
118 self%nscan = mod(igdtmpl(19)/32,2)
119 self%nscan_field_pos = self%nscan
120 self%kscan = 0
121 self%iwrap = nint(360/abs(self%DLON))
122
123 if(self%im.lt.self%iwrap) self%iwrap=0
124 self%jwrap1=0
125 self%jwrap2=0
126
127 if(self%im < self%iwrap) self%iwrap=0
128 self%jwrap1 = 0
129 self%jwrap2 = 0
130 if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0) then
131 if(abs(self%rlat1) > 90-0.25*self%dlat) then
132 self%jwrap1 = 2
133 elseif(abs(self%rlat1) > 90-0.75*self%dlat) then
134 self%jwrap1 = 1
135 endif
136 if(abs(self%rlat2) > 90-0.25*self%dlat) then
137 self%jwrap2 = 2 * self%jm
138 elseif(abs(self%rlat2) > 90-0.75*self%dlat) then
139 self%jwrap2 = 2 * self%jm+1
140 endif
141 endif
142
143 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
144
145 end associate
146 end subroutine init_grib2
147
148
193 SUBROUTINE gdswzd_equid_cylind(self,IOPT,NPTS,FILL, &
194 XPTS,YPTS,RLON,RLAT,NRET, &
195 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
196 IMPLICIT NONE
197 !
198 class(ip_equid_cylind_grid), intent(in) :: self
199 INTEGER, INTENT(IN ) :: IOPT, NPTS
200 INTEGER, INTENT( OUT) :: NRET
201 !
202 REAL, INTENT(IN ) :: FILL
203 REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
204 REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
205 REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
206 REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
207 REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
208 !
209 INTEGER :: IM, JM, N
210 !
211 LOGICAL :: LROT, LMAP, LAREA
212 !
213 REAL :: HI, RLAT1, RLON1, RLAT2, RLON2
214 REAL :: XMAX, XMIN, YMAX, YMIN
215 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216 IF(PRESENT(crot)) crot=fill
217 IF(PRESENT(srot)) srot=fill
218 IF(PRESENT(xlon)) xlon=fill
219 IF(PRESENT(xlat)) xlat=fill
220 IF(PRESENT(ylon)) ylon=fill
221 IF(PRESENT(ylat)) ylat=fill
222 IF(PRESENT(area)) area=fill
223
224 im=self%im
225 jm=self%jm
226
227 rlat1=self%rlat1
228 rlon1=self%rlon1
229 rlat2=self%rlat2
230 rlon2=self%rlon2
231
232 hi=self%hi
233
234 rerth = self%rerth
235 dlat = self%dlat
236 dlon = self%dlon
237
238 xmin=0
239 xmax=im+1
240 IF(im.EQ.nint(360/abs(dlon))) xmax=im+2
241 ymin=0
242 ymax=jm+1
243 nret=0
244 IF(PRESENT(crot).AND.PRESENT(srot))THEN
245 lrot=.true.
246 ELSE
247 lrot=.false.
248 ENDIF
249 IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
250 lmap=.true.
251 ELSE
252 lmap=.false.
253 ENDIF
254 IF(PRESENT(area))THEN
255 larea=.true.
256 ELSE
257 larea=.false.
258 ENDIF
259 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260 ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
261 IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
262 !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
263 DO n=1,npts
264 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
265 ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
266 rlon(n)=mod(rlon1+dlon*(xpts(n)-1)+3600,360.)
267 rlat(n)=min(max(rlat1+dlat*(ypts(n)-1),-90.),90.)
268 nret=nret+1
269 IF(lrot) CALL equid_cylind_vect_rot(crot(n),srot(n))
270 IF(lmap) CALL equid_cylind_map_jacob(xlon(n),xlat(n),ylon(n),ylat(n))
271 IF(larea) CALL equid_cylind_grid_area(rlat(n),area(n))
272 ELSE
273 rlon(n)=fill
274 rlat(n)=fill
275 ENDIF
276 ENDDO
277 !$OMP END PARALLEL DO
278 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279 ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
280 ELSEIF(iopt.EQ.-1) THEN
281 !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
282 DO n=1,npts
283 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90) THEN
284 xpts(n)=1+hi*mod(hi*(rlon(n)-rlon1)+3600,360.)/dlon
285 ypts(n)=1+(rlat(n)-rlat1)/dlat
286 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
287 ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
288 nret=nret+1
289 IF(lrot) CALL equid_cylind_vect_rot(crot(n),srot(n))
290 IF(lmap) CALL equid_cylind_map_jacob(xlon(n),xlat(n),ylon(n),ylat(n))
291 IF(larea) CALL equid_cylind_grid_area(rlat(n),area(n))
292 ELSE
293 xpts(n)=fill
294 ypts(n)=fill
295 ENDIF
296 ELSE
297 xpts(n)=fill
298 ypts(n)=fill
299 ENDIF
300 ENDDO
301 !$OMP END PARALLEL DO
302 ENDIF
303 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
304 END SUBROUTINE gdswzd_equid_cylind
305 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306 SUBROUTINE equid_cylind_error(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
307 !$$$ SUBPROGRAM DOCUMENTATION BLOCK
308 !
309 ! SUBPROGRAM: EQUID_CYLIND_ERROR ERROR HANDLER
310 ! PRGMMR: GAYNO ORG: W/NMC23 DATE: 2015-07-13
311 !
312 ! ABSTRACT: UPON AN ERROR, THIS SUBPROGRAM ASSIGNS
313 ! A "FILL" VALUE TO THE OUTPUT FIELDS.
314 !
315 ! PROGRAM HISTORY LOG:
316 ! 2015-07-13 GAYNO INITIAL VERSION
317 !
318 ! USAGE: CALL EQUID_CYLIND_ERROR(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
319 !
320 ! INPUT ARGUMENT LIST:
321 ! IOPT - INTEGER OPTION FLAG
322 ! (+1 TO COMPUTE EARTH COORDS OF SELECTED GRID COORDS)
323 ! (-1 TO COMPUTE GRID COORDS OF SELECTED EARTH COORDS)
324 ! NPTS - INTEGER MAXIMUM NUMBER OF COORDINATES
325 ! FILL - REAL FILL VALUE TO SET INVALID OUTPUT DATA
326 ! (MUST BE IMPOSSIBLE VALUE; SUGGESTED VALUE: -9999.)
327 ! OUTPUT ARGUMENT LIST:
328 ! RLON - REAL (NPTS) EARTH LONGITUDES IN DEGREES E IF IOPT<0
329 ! RLAT - REAL (NPTS) EARTH LATITUDES IN DEGREES N IF IOPT<0
330 ! XPTS - REAL (NPTS) GRID X POINT COORDINATES IF IOPT>0
331 ! YPTS - REAL (NPTS) GRID Y POINT COORDINATES IF IOPT>0
332 !
333 ! ATTRIBUTES:
334 ! LANGUAGE: FORTRAN 90
335 !
336 !$$$
337 IMPLICIT NONE
338 !
339 INTEGER, INTENT(IN ) :: IOPT, NPTS
340 !
341 REAL, INTENT(IN ) :: FILL
342 REAL, INTENT( OUT) :: RLAT(NPTS),RLON(NPTS)
343 REAL, INTENT( OUT) :: XPTS(NPTS),YPTS(NPTS)
344 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
345 IF(iopt>=0) THEN
346 rlon=fill
347 rlat=fill
348 ENDIF
349 IF(iopt<=0) THEN
350 xpts=fill
351 ypts=fill
352 ENDIF
353 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
354 END SUBROUTINE equid_cylind_error
355
368 SUBROUTINE equid_cylind_vect_rot(CROT,SROT)
369 IMPLICIT NONE
370
371 REAL, INTENT( OUT) :: CROT, SROT
372
373 crot=1.0
374 srot=0.0
375
376 END SUBROUTINE equid_cylind_vect_rot
377
388 SUBROUTINE equid_cylind_map_jacob(XLON,XLAT,YLON,YLAT)
389 REAL, INTENT( OUT) :: XLON,XLAT,YLON,YLAT
390
391 xlon=1.0/dlon
392 xlat=0.
393 ylon=0.
394 ylat=1.0/dlat
395
396 END SUBROUTINE equid_cylind_map_jacob
397
405 SUBROUTINE equid_cylind_grid_area(RLAT,AREA)
406 IMPLICIT NONE
407
408 REAL, INTENT(IN ) :: RLAT
409 REAL, INTENT( OUT) :: AREA
410
411 REAL, PARAMETER :: PI=3.14159265358979
412 REAL, PARAMETER :: DPR=180./pi
413
414 REAL :: DSLAT, RLATU, RLATD
415
416 rlatu=min(max(rlat+dlat/2,-90.),90.)
417 rlatd=min(max(rlat-dlat/2,-90.),90.)
418 dslat=sin(rlatu/dpr)-sin(rlatd/dpr)
419 area=rerth**2*abs(dslat*dlon)/dpr
420
421 END SUBROUTINE equid_cylind_grid_area
422
424
Equidistant cylindrical grib decoder and grid coordinate transformations.
subroutine equid_cylind_map_jacob(XLON, XLAT, YLON, YLAT)
Computes the map jacobians for a equidistant cylindrical grid.
subroutine init_grib2(self, g2_desc)
Initializes an equidistant cylindrical grid given a grib2_descriptor object.
subroutine init_grib1(self, g1_desc)
Initializes an equidistant cylindrical grid given a grib1_descriptor object.
subroutine equid_cylind_grid_area(RLAT, AREA)
Computes the grid box area for a equidistant cylindrical grid.
subroutine equid_cylind_vect_rot(CROT, SROT)
Computes the vector rotation sines and cosines for a equidistant cylindrical grid.
subroutine gdswzd_equid_cylind(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
Calculates Earth coordinates (iopt = 1) or grid coorindates (iopt = -1) for equidistant cylindrical g...
Uses derived type grid descriptor objects to abstract away the raw Grib-1 and Grib-2 grid definitions...
Abstract ip_grid type.
Definition: ip_grid_mod.f90:8
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Abstract grid that holds fields and methods common to all grids.
Definition: ip_grid_mod.f90:45