20 real :: hi, rlat1, rlon1, rlat2, rlon2
23 procedure :: init_grib1
24 procedure :: init_grib2
47 associate(kgds => g1_desc%gds)
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)
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)
63 self%nscan = mod(kgds(11) / 32, 2)
64 self%nscan_field_pos = self%nscan
67 self%iwrap = nint(360/abs(self%dlon))
69 if(self%im < self%iwrap) self%iwrap=0
72 if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0)
then
73 if(abs(self%rlat1) > 90-0.25*self%dlat)
then
75 elseif(abs(self%rlat1) > 90-0.75*self%dlat)
then
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
86 self%eccen_squared = 0.0
100 type(grib2_descriptor),
intent(in) :: g2_desc
102 integer :: iscale, iscan
104 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
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)
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)
118 self%nscan = mod(igdtmpl(19)/32,2)
119 self%nscan_field_pos = self%nscan
121 self%iwrap = nint(360/abs(self%DLON))
123 if(self%im.lt.self%iwrap) self%iwrap=0
127 if(self%im < self%iwrap) self%iwrap=0
130 if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0)
then
131 if(abs(self%rlat1) > 90-0.25*self%dlat)
then
133 elseif(abs(self%rlat1) > 90-0.75*self%dlat)
then
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
143 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
194 XPTS,YPTS,RLON,RLAT,NRET, &
195 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
199 INTEGER,
INTENT(IN ) :: IOPT, NPTS
200 INTEGER,
INTENT( OUT) :: NRET
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)
211 LOGICAL :: LROT, LMAP, LAREA
213 REAL :: HI, RLAT1, RLON1, RLAT2, RLON2
214 REAL :: XMAX, XMIN, YMAX, YMIN
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
240 IF(im.EQ.nint(360/abs(dlon))) xmax=im+2
244 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
249 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
254 IF(
PRESENT(area))
THEN
261 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
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.)
280 ELSEIF(iopt.EQ.-1)
THEN
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
306 SUBROUTINE equid_cylind_error(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
339 INTEGER,
INTENT(IN ) :: IOPT, NPTS
341 REAL,
INTENT(IN ) :: FILL
342 REAL,
INTENT( OUT) :: RLAT(NPTS),RLON(NPTS)
343 REAL,
INTENT( OUT) :: XPTS(NPTS),YPTS(NPTS)
354 END SUBROUTINE equid_cylind_error
371 REAL,
INTENT( OUT) :: CROT, SROT
389 REAL,
INTENT( OUT) :: XLON,XLAT,YLON,YLAT
408 REAL,
INTENT(IN ) :: RLAT
409 REAL,
INTENT( OUT) :: AREA
411 REAL,
PARAMETER :: PI=3.14159265358979
412 REAL,
PARAMETER :: DPR=180./pi
414 REAL :: DSLAT, RLATU, RLATD
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
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...
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Abstract grid that holds fields and methods common to all grids.