56 integer :: iscan, jscan
59 associate(kgds => g1_desc%gds)
61 self%eccen_squared = 0.0
66 self%RLAT1=kgds(4)*1.e-3
67 self%RLON1=kgds(5)*1.e-3
68 self%RLON2=kgds(8)*1.e-3
69 self%RLATI=kgds(9)*1.e-3
71 iscan=mod(kgds(11)/128,2)
72 jscan=mod(kgds(11)/64,2)
77 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
78 self%DPHI=hj*dy/(self%RERTH*cos(self%RLATI/
dpr))
84 self%nscan = mod(kgds(11) / 32, 2)
85 self%nscan_field_pos = self%nscan
88 self%iwrap = nint(360 / abs(self%dlon))
89 if (self%im < self%iwrap) self%iwrap = 0
102 type(grib2_descriptor),
intent(in) :: g2_desc
104 integer :: iscan, jscan
107 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
109 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
114 self%RLAT1=float(igdtmpl(10))*1.0e-6
115 self%RLON1=float(igdtmpl(11))*1.0e-6
116 self%RLON2=float(igdtmpl(15))*1.0e-6
117 self%RLATI=float(igdtmpl(13))*1.0e-6
119 iscan=mod(igdtmpl(16)/128,2)
120 jscan=mod(igdtmpl(16)/64,2)
122 dy=float(igdtmpl(19))*1.0e-3
125 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
126 self%DPHI=hj*dy/(self%RERTH*cos(self%RLATI/dpr))
131 self%nscan=mod(igdtmpl(16) / 32,2)
132 self%nscan_field_pos = self%nscan
134 self%iwrap = nint(360 / abs(self%dlon))
135 if(self%im < self%iwrap) self%iwrap = 0
199 XPTS,YPTS,RLON,RLAT,NRET, &
200 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
204 INTEGER,
INTENT(IN ) :: IOPT, NPTS
205 INTEGER,
INTENT( OUT) :: NRET
207 REAL,
INTENT(IN ) :: FILL
208 REAL,
INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
209 REAL,
INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
210 REAL,
OPTIONAL,
INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
211 REAL,
OPTIONAL,
INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
212 REAL,
OPTIONAL,
INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
216 LOGICAL :: LROT, LMAP, LAREA
219 REAL :: RLAT1, RLON1, RLON2, RLATI
220 REAL :: XMAX, XMIN, YMAX, YMIN
223 IF(
PRESENT(crot)) crot=fill
224 IF(
PRESENT(srot)) srot=fill
225 IF(
PRESENT(xlon)) xlon=fill
226 IF(
PRESENT(xlat)) xlat=fill
227 IF(
PRESENT(ylon)) ylon=fill
228 IF(
PRESENT(ylat)) ylat=fill
229 IF(
PRESENT(area)) area=fill
246 ye=1-log(tan((rlat1+90)/2/dpr))/
dphi
249 IF(im.EQ.nint(360/abs(
dlon))) xmax=im+2
253 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
258 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
263 IF(
PRESENT(area))
THEN
270 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
273 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
274 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
275 rlon(n)=mod(rlon1+
dlon*(xpts(n)-1)+3600,360.)
276 rlat(n)=2*atan(exp(
dphi*(ypts(n)-ye)))*dpr-90
289 ELSEIF(iopt.EQ.-1)
THEN
292 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LT.90)
THEN
293 xpts(n)=1+hi*mod(hi*(rlon(n)-rlon1)+3600,360.)/
dlon
294 ypts(n)=ye+log(tan((rlat(n)+90)/2/dpr))/
dphi
295 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
296 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
334 REAL,
INTENT( OUT) :: CROT, SROT
362 REAL,
INTENT(IN ) :: RLAT
363 REAL,
INTENT( OUT) :: XLON, XLAT, YLON, YLAT
368 ylat=1./
dphi/cos(rlat/dpr)/dpr
390 REAL,
INTENT(IN ) :: RLAT
391 REAL,
INTENT( OUT) :: AREA
void gdswzd(int igdtnum, int *igdtmpl, int igdtlen, int iopt, int npts, float fill, float *xpts, float *ypts, float *rlon, float *rlat, int *nret, float *crot, float *srot, float *xlon, float *xlat, float *ylon, float *ylat, float *area)
gdswzd() interface for C for _4 build of library.
Module containing common constants.
real, parameter dpr
Radians to degrees.
Determine earth radius and shape.
Users derived type grid descriptor objects to abstract away the raw GRIB1 and GRIB2 grid definitions.
GDS wizard for mercator cylindrical.
real dlon
Longitudinal direction grid length.
real dphi
Latitudinal direction grid length.
subroutine init_grib1(self, g1_desc)
Initializes a mercator grid given a grib1_descriptor object.
subroutine mercator_grid_area(RLAT, AREA)
Grid box area for mercator cylindrical grids.
subroutine init_grib2(self, g2_desc)
Init GRIB2.
subroutine mercator_map_jacob(RLAT, XLON, XLAT, YLON, YLAT)
Map jacobians for mercator cylindrical grids.
subroutine gdswzd_mercator(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for mercator cylindrical.
real rerth
Radius of the Earth.
subroutine mercator_vect_rot(CROT, SROT)
Vector rotation fields for mercator cylindrical grids.
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Abstract grid that holds fields and methods common to all grids.