63 real :: dx, dy, hi, hj
64 integer :: iproj, iscan, jscan
66 associate(kgds => g1_desc%gds)
68 self%eccen_squared = 0.0
73 self%RLAT1=kgds(4)*1.e-3
74 self%RLON1=kgds(5)*1.e-3
76 self%IROT=mod(kgds(6)/8,2)
77 self%ORIENT=kgds(7)*1.e-3
82 iproj=mod(kgds(10)/128,2)
83 iscan=mod(kgds(11)/128,2)
84 jscan=mod(kgds(11)/64,2)
86 self%RLATI1=kgds(12)*1.e-3
87 self%RLATI2=kgds(13)*1.e-3
98 self%nscan = mod(kgds(11) / 32, 2)
99 self%nscan_field_pos = self%nscan
113 type(grib2_descriptor),
intent(in) :: g2_desc
115 real :: dx, dy, hi, hj
116 integer :: iproj, iscan, jscan
119 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
120 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
125 self%RLAT1=float(igdtmpl(10))*1.0e-6
126 self%RLON1=float(igdtmpl(11))*1.0e-6
128 self%IROT=mod(igdtmpl(12)/8,2)
129 self%ORIENT=float(igdtmpl(14))*1.0e-6
131 dx=float(igdtmpl(15))*1.0e-3
132 dy=float(igdtmpl(16))*1.0e-3
134 iproj=mod(igdtmpl(17)/128,2)
135 iscan=mod(igdtmpl(18)/128,2)
136 jscan=mod(igdtmpl(18)/64,2)
138 self%RLATI1=float(igdtmpl(19))*1.0e-6
139 self%RLATI2=float(igdtmpl(20))*1.0e-6
147 self%nscan = mod(igdtmpl(18) / 32, 2)
148 self%nscan_field_pos = self%nscan
219 XPTS,YPTS,RLON,RLAT,NRET, &
220 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
224 INTEGER,
INTENT(IN ) :: IOPT, NPTS
225 INTEGER,
INTENT( OUT) :: NRET
227 REAL,
INTENT(IN ) :: FILL
228 REAL,
INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
229 REAL,
INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
230 REAL,
OPTIONAL,
INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
231 REAL,
OPTIONAL,
INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
232 REAL,
OPTIONAL,
INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
236 LOGICAL :: LROT, LMAP, LAREA
241 REAL :: ORIENT, RLAT1, RLON1
242 REAL :: RLATI1, RLATI2
243 REAL :: XMAX, XMIN, YMAX, YMIN, XP, YP
246 IF(
PRESENT(crot)) crot=fill
247 IF(
PRESENT(srot)) srot=fill
248 IF(
PRESENT(xlon)) xlon=fill
249 IF(
PRESENT(xlat)) xlat=fill
250 IF(
PRESENT(ylon)) ylon=fill
251 IF(
PRESENT(ylat)) ylat=fill
252 IF(
PRESENT(area)) area=fill
272 IF(rlati1.EQ.rlati2)
THEN
275 an=log(cos(rlati1/dpr)/cos(rlati2/dpr))/ &
276 log(tan((90-rlati1)/2/dpr)/tan((90-rlati2)/2/dpr))
278 de=
rerth*cos(rlati1/dpr)*tan((rlati1+90)/2/dpr)**
an/
an
279 IF(
h*rlat1.EQ.90)
THEN
283 dr=de/tan((rlat1+90)/2/dpr)**
an
284 dlon1=mod(rlon1-orient+180+3600,360.)-180
285 xp=1-sin(
an*dlon1/dpr)*dr/
dxs
286 yp=1+cos(
an*dlon1/dpr)*dr/
dys
295 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
300 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
305 IF(
PRESENT(area))
THEN
312 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
315 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
316 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
317 di=
h*(xpts(n)-xp)*
dxs
318 dj=
h*(ypts(n)-yp)*
dys
321 IF(dr2.LT.de2*1.e-6)
THEN
325 rlon(n)=mod(orient+1./
an*dpr*atan2(di,-dj)+3600,360.)
326 rlat(n)=(2*dpr*atan((de2/dr2)**antr)-90)
329 dlon=mod(rlon(n)-orient+180+3600,360.)-180
332 xlon(n),xlat(n),ylon(n),ylat(n))
341 ELSEIF(iopt.EQ.-1)
THEN
344 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90.AND. &
345 h*rlat(n).NE.-90)
THEN
346 dr=
h*de*tan((90-rlat(n))/2/dpr)**
an
347 dlon=mod(rlon(n)-orient+180+3600,360.)-180
348 xpts(n)=xp+
h*sin(
an*dlon/dpr)*dr/
dxs
349 ypts(n)=yp-
h*cos(
an*dlon/dpr)*dr/
dys
350 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
351 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
355 xlon(n),xlat(n),ylon(n),ylat(n))
391 REAL,
INTENT( IN) :: DLON
392 REAL,
INTENT( OUT) :: CROT, SROT
395 crot=cos(
an*dlon/dpr)
396 srot=sin(
an*dlon/dpr)
429 REAL,
INTENT(IN ) :: RLAT, FILL, DLON, DR
430 REAL,
INTENT( OUT) :: XLON, XLAT, YLON, YLAT
435 IF(clat.LE.0.OR.dr.LE.0)
THEN
441 xlon=
h*cos(
an*dlon/dpr)*
an/dpr*dr/
dxs
442 xlat=-
h*sin(
an*dlon/dpr)*
an/dpr*dr/
dxs/clat
443 ylon=
h*sin(
an*dlon/dpr)*
an/dpr*dr/
dys
444 ylat=
h*cos(
an*dlon/dpr)*
an/dpr*dr/
dys/clat
470 REAL,
INTENT(IN ) :: RLAT
471 REAL,
INTENT(IN ) :: FILL
472 REAL,
INTENT(IN ) :: DR
473 REAL,
INTENT( OUT) :: AREA
478 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.