50 REAL :: tinyreal=tiny(1.0)
64 real :: dx, dy, hi, hj
65 integer :: iproj, iscan, jscan
67 associate(kgds => g1_desc%gds)
69 self%eccen_squared = 0.0
74 self%RLAT1=kgds(4)*1.e-3
75 self%RLON1=kgds(5)*1.e-3
77 self%IROT=mod(kgds(6)/8,2)
78 self%ORIENT=kgds(7)*1.e-3
83 iproj=mod(kgds(10)/128,2)
84 iscan=mod(kgds(11)/128,2)
85 jscan=mod(kgds(11)/64,2)
87 self%RLATI1=kgds(12)*1.e-3
88 self%RLATI2=kgds(13)*1.e-3
99 self%nscan = mod(kgds(11) / 32, 2)
100 self%nscan_field_pos = self%nscan
114 type(grib2_descriptor),
intent(in) :: g2_desc
116 real :: dx, dy, hi, hj
117 integer :: iproj, iscan, jscan
120 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
121 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
126 self%RLAT1=float(igdtmpl(10))*1.0e-6
127 self%RLON1=float(igdtmpl(11))*1.0e-6
129 self%IROT=mod(igdtmpl(12)/8,2)
130 self%ORIENT=float(igdtmpl(14))*1.0e-6
132 dx=float(igdtmpl(15))*1.0e-3
133 dy=float(igdtmpl(16))*1.0e-3
135 iproj=mod(igdtmpl(17)/128,2)
136 iscan=mod(igdtmpl(18)/128,2)
137 jscan=mod(igdtmpl(18)/64,2)
139 self%RLATI1=float(igdtmpl(19))*1.0e-6
140 self%RLATI2=float(igdtmpl(20))*1.0e-6
148 self%nscan = mod(igdtmpl(18) / 32, 2)
149 self%nscan_field_pos = self%nscan
220 XPTS,YPTS,RLON,RLAT,NRET, &
221 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
225 INTEGER,
INTENT(IN ) :: IOPT, NPTS
226 INTEGER,
INTENT( OUT) :: NRET
228 REAL,
INTENT(IN ) :: FILL
229 REAL,
INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
230 REAL,
INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
231 REAL,
OPTIONAL,
INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
232 REAL,
OPTIONAL,
INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
233 REAL,
OPTIONAL,
INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
237 LOGICAL :: LROT, LMAP, LAREA
242 REAL :: ORIENT, RLAT1, RLON1
243 REAL :: RLATI1, RLATI2
244 REAL :: XMAX, XMIN, YMAX, YMIN, XP, YP
247 IF(
PRESENT(crot)) crot=fill
248 IF(
PRESENT(srot)) srot=fill
249 IF(
PRESENT(xlon)) xlon=fill
250 IF(
PRESENT(xlat)) xlat=fill
251 IF(
PRESENT(ylon)) ylon=fill
252 IF(
PRESENT(ylat)) ylat=fill
253 IF(
PRESENT(area)) area=fill
273 IF(abs(rlati1-rlati2).LT.tinyreal)
THEN
276 an=log(cos(rlati1/dpr)/cos(rlati2/dpr))/ &
277 log(tan((90-rlati1)/2/dpr)/tan((90-rlati2)/2/dpr))
279 de=
rerth*cos(rlati1/dpr)*tan((rlati1+90)/2/dpr)**
an/
an
280 IF(abs(
h*rlat1-90).LT.tinyreal)
THEN
284 dr=de/tan((rlat1+90)/2/dpr)**
an
285 dlon1=mod(rlon1-orient+180+3600,360.)-180
286 xp=1-sin(
an*dlon1/dpr)*dr/
dxs
287 yp=1+cos(
an*dlon1/dpr)*dr/
dys
296 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
301 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
306 IF(
PRESENT(area))
THEN
313 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
316 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
317 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
318 di=
h*(xpts(n)-xp)*
dxs
319 dj=
h*(ypts(n)-yp)*
dys
322 IF(dr2.LT.de2*1.e-6)
THEN
326 rlon(n)=mod(orient+1./
an*dpr*atan2(di,-dj)+3600,360.)
327 rlat(n)=(2*dpr*atan((de2/dr2)**antr)-90)
330 dlon=mod(rlon(n)-orient+180+3600,360.)-180
333 xlon(n),xlat(n),ylon(n),ylat(n))
342 ELSEIF(iopt.EQ.-1)
THEN
345 IF(abs(rlon(n)).LT.(360.+tinyreal).AND.abs(rlat(n)).LT.(90.+tinyreal).AND. &
346 abs(
h*rlat(n)+90).GT.tinyreal)
THEN
347 dr=
h*de*tan((90-rlat(n))/2/dpr)**
an
348 dlon=mod(rlon(n)-orient+180+3600,360.)-180
349 xpts(n)=xp+
h*sin(
an*dlon/dpr)*dr/
dxs
350 ypts(n)=yp-
h*cos(
an*dlon/dpr)*dr/
dys
351 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
352 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
356 xlon(n),xlat(n),ylon(n),ylat(n))
392 REAL,
INTENT( IN) :: DLON
393 REAL,
INTENT( OUT) :: CROT, SROT
396 crot=cos(
an*dlon/dpr)
397 srot=sin(
an*dlon/dpr)
430 REAL,
INTENT(IN ) :: RLAT, FILL, DLON, DR
431 REAL,
INTENT( OUT) :: XLON, XLAT, YLON, YLAT
436 IF(clat.LE.0.OR.dr.LE.0)
THEN
442 xlon=
h*cos(
an*dlon/dpr)*
an/dpr*dr/
dxs
443 xlat=-
h*sin(
an*dlon/dpr)*
an/dpr*dr/
dxs/clat
444 ylon=
h*sin(
an*dlon/dpr)*
an/dpr*dr/
dys
445 ylat=
h*cos(
an*dlon/dpr)*
an/dpr*dr/
dys/clat
471 REAL,
INTENT(IN ) :: RLAT
472 REAL,
INTENT(IN ) :: FILL
473 REAL,
INTENT(IN ) :: DR
474 REAL,
INTENT( OUT) :: AREA
479 IF(clat.LE.0.OR.dr.LE.0)
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.
Determine earth radius and shape.
Users derived type grid descriptor objects to abstract away the raw GRIB1 and GRIB2 grid definitions.
Lambert conformal grib decoder and grid coordinate transformations.
subroutine lambert_conf_grid_area(RLAT, FILL, DR, AREA)
Grid box area for lambert conformal conical.
real dxs
x-direction grid length adjusted for scan mode.
subroutine init_grib1(self, g1_desc)
Initializes a Lambert Conformal grid given a grib1_descriptor object.
subroutine lambert_conf_vect_rot(DLON, CROT, SROT)
Vector rotation fields for lambert conformal conical.
subroutine lambert_conf_map_jacob(RLAT, FILL, DLON, DR, XLON, XLAT, YLON, YLAT)
Map jacobians for lambert conformal conical.
subroutine gdswzd_lambert_conf(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for lambert conformal conical.
integer irot
vector rotation flag.
real rerth
Radius of the earth.
real dys
y-direction grid length adjusted for scan model.
subroutine init_grib2(self, g2_desc)
Initializes a Lambert Conformal grid given a grib2_descriptor object.
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Abstract grid that holds fields and methods common to all grids.