66 REAL,
PARAMETER :: SLAT=60.0
68 real :: dx, dy, hi, hj
69 integer :: iproj, iscan, jscan
71 associate(kgds => g1_desc%gds)
72 self%ELLIPTICAL=mod(kgds(6)/64,2).EQ.1
74 if (.not. self%elliptical)
then
76 self%eccen_squared = 0d0
85 self%RLAT1=kgds(4)*1.e-3
86 self%RLON1=kgds(5)*1.e-3
88 self%IROT=mod(kgds(6)/8,2)
92 self%ORIENT=kgds(7)*1.e-3
97 iproj=mod(kgds(10)/128,2)
98 iscan=mod(kgds(11)/128,2)
99 jscan=mod(kgds(11)/64,2)
105 IF(self%H.EQ.-1)self%ORIENT=self%ORIENT+180.
113 self%nscan = mod(kgds(11) / 32, 2)
114 self%nscan_field_pos = self%nscan
129 type(grib2_descriptor),
intent(in) :: g2_desc
131 real :: slat, dx, dy, hi, hj
132 integer :: iproj, iscan, jscan
134 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
135 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
137 self%ELLIPTICAL = self%eccen_squared > 0.0
142 self%RLAT1=float(igdtmpl(10))*1.e-6
143 self%RLON1=float(igdtmpl(11))*1.e-6
145 self%IROT=mod(igdtmpl(12)/8,2)
147 slat=float(abs(igdtmpl(13)))*1.e-6
150 self%ORIENT=float(igdtmpl(14))*1.e-6
152 dx=float(igdtmpl(15))*1.e-3
153 dy=float(igdtmpl(16))*1.e-3
155 iproj=mod(igdtmpl(17)/128,2)
156 iscan=mod(igdtmpl(18)/128,2)
157 jscan=mod(igdtmpl(18)/64,2)
166 self%nscan = mod(igdtmpl(18) / 32, 2)
167 self%nscan_field_pos = self%nscan
239 FILL,XPTS,YPTS,RLON,RLAT,NRET, &
240 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
245 INTEGER,
INTENT(IN ) :: IOPT, NPTS
246 INTEGER,
INTENT( OUT) :: NRET
248 REAL,
INTENT(IN ) :: FILL
249 REAL,
INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
250 REAL,
INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
251 REAL,
OPTIONAL,
INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
252 REAL,
OPTIONAL,
INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
253 REAL,
OPTIONAL,
INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
258 LOGICAL :: ELLIPTICAL, LROT, LMAP, LAREA
260 REAL :: ALAT, ALAT1, ALONG, DIFF
262 REAL :: DR, E, E_OVER_2
264 REAL :: RLAT1, RLON1, RHO, T, TC
265 REAL :: XMAX, XMIN, YMAX, YMIN
268 IF(
PRESENT(crot)) crot=fill
269 IF(
PRESENT(srot)) srot=fill
270 IF(
PRESENT(xlon)) xlon=fill
271 IF(
PRESENT(xlat)) xlat=fill
272 IF(
PRESENT(ylon)) ylon=fill
273 IF(
PRESENT(ylat)) ylat=fill
274 IF(
PRESENT(area)) area=fill
276 elliptical = self%elliptical
292 e2 = self%eccen_squared
295 IF (.NOT.elliptical)
THEN
296 de=(1.+sin(slatr))*
rerth
297 dr=de*cos(rlat1/dpr)/(1+
h*sin(rlat1/dpr))
305 along = (rlon1-
orient)/dpr
306 t=tan(pi4-alat/2.)/((1.-e*sin(alat))/ &
307 (1.+e*sin(alat)))**(e_over_2)
308 tc=tan(pi4-slatr/2.)/((1.-e*sin(slatr))/ &
309 (1.+e*sin(slatr)))**(e_over_2)
310 mc=cos(slatr)/sqrt(1.0-
e2*(sin(slatr)**2))
312 yp = 1.0 + rho*cos(
h*along)/
dys
313 xp = 1.0 - rho*sin(
h*along)/
dxs
320 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
325 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
330 IF(
PRESENT(area))
THEN
337 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
338 IF(.NOT.elliptical)
THEN
341 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
342 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
346 IF(dr2.LT.
de2*1.e-6)
THEN
350 rlon(n)=mod(
orient+
h*dpr*atan2(di,-dj)+3600,360.)
351 rlat(n)=
h*dpr*asin((
de2-dr2)/(
de2+dr2))
356 xlon(n),xlat(n),ylon(n),ylat(n))
368 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
369 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
372 rho=sqrt(di*di+dj*dj)
373 t=(rho*tc)/(
rerth*mc)
374 IF(abs(ypts(n)-yp)<0.01)
THEN
376 IF(di<=0.0) along=
orient-
h*90.0
378 along=
orient+
h*atan(di/(-dj))*dpr
379 IF(dj>0) along=along+180.
381 alat1=pi2-2.0*atan(t)
383 alat = pi2 - 2.0*atan(t*(((1.0-e*sin(alat1))/ &
384 (1.0+e*sin(alat1)))**(e_over_2)))
385 diff = abs(alat-alat1)*dpr
386 IF (diff < 0.000001)
EXIT
391 IF(rlon(n)<0.0) rlon(n)=rlon(n)+360.
392 IF(rlon(n)>360.0) rlon(n)=rlon(n)-360.0
404 ELSEIF(iopt.EQ.-1)
THEN
405 IF(.NOT.elliptical)
THEN
408 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90.AND. &
409 h*rlat(n).NE.-90)
THEN
410 dr=de*tan((90-
h*rlat(n))/2/dpr)
413 ypts(n)=yp-cos((rlon(n)-
orient)/dpr)*dr/
dys
414 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
415 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
419 xlon(n),xlat(n),ylon(n),ylat(n))
434 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90.AND. &
435 h*rlat(n).NE.-90)
THEN
437 along = (rlon(n)-
orient)/dpr
438 t=tan(pi4-alat*0.5)/((1.-e*sin(alat))/ &
439 (1.+e*sin(alat)))**(e_over_2)
441 xpts(n)= xp + rho*sin(
h*along) /
dxs
442 ypts(n)= yp - rho*cos(
h*along) /
dys
443 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
444 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
484 REAL,
INTENT(IN ) :: RLON
485 REAL,
INTENT( OUT) :: CROT, SROT
489 srot=sin((rlon-
orient)/dpr)
522 REAL,
INTENT(IN ) :: RLON, RLAT, DR2
523 REAL,
INTENT( OUT) :: XLON, XLAT, YLON, YLAT
527 IF(dr2.LT.
de2*1.e-6)
THEN
530 xlat=-sin((rlon-
orient)/dpr)/dpr*de/
dxs/2
537 xlat=-sin((rlon-
orient)/dpr)/dpr*dr/
dxs/clat
565 REAL,
INTENT(IN ) :: RLAT, DR2
566 REAL,
INTENT( OUT) :: AREA
570 IF(dr2.LT.
de2*1.e-6)
THEN
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 pi4
PI / 4.0.
real, parameter pi2
PI / 2.0.
real, parameter e2_wgs84
Eccentricity squared of Earth defined by WGS-84.
real, parameter rerth_wgs84
Radius of the Earth defined by WGS-84.
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 polar stereographic azimuthal.
integer irot
Local copy of irot.
subroutine polar_stereo_grid_area(RLAT, DR2, AREA)
Grid box area for polar stereographic grids.
real rerth
Radius of the Earth.
subroutine polar_stereo_map_jacob(RLON, RLAT, DR2, XLON, XLAT, YLON, YLAT)
Map jacobians for polar stereographic grids.
subroutine init_grib2(self, g2_desc)
Initializes a polar stereographic grid given a grib2_descriptor object.
subroutine init_grib1(self, g1_desc)
Initializes a polar stereographic grid given a grib1_descriptor object.
real orient
Local copy of orient.
subroutine polar_stereo_vect_rot(RLON, CROT, SROT)
Vector rotation fields for polar stereographic grids.
real e2
Eccentricity squared.
real dxs
Local copy of dxs.
subroutine gdswzd_polar_stereo(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for polar stereographic azimuthal.
real dys
Local copy of dys.
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Abstract grid that holds fields and methods common to all grids.