1module ip_mercator_grid_mod
12 real :: rlat1, rlon1, rlon2, rlati, hi, dlon, dphi
14 procedure :: init_grib1
15 procedure :: init_grib2
16 procedure :: gdswzd => gdswzd_mercator
25 subroutine init_grib1(self, g1_desc)
26 class(ip_mercator_grid),
intent(inout) :: self
27 type(grib1_descriptor),
intent(in) :: g1_desc
29 integer :: iscan, jscan
32 associate(kgds => g1_desc%gds)
34 self%eccen_squared = 0.0
39 self%RLAT1=kgds(4)*1.e-3
40 self%RLON1=kgds(5)*1.e-3
41 self%RLON2=kgds(8)*1.e-3
42 self%RLATI=kgds(9)*1.e-3
44 iscan=mod(kgds(11)/128,2)
45 jscan=mod(kgds(11)/64,2)
50 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
51 self%DPHI=hj*dy/(self%RERTH*cos(self%RLATI/
dpr))
57 self%nscan = mod(kgds(11) / 32, 2)
58 self%nscan_field_pos = self%nscan
61 self%iwrap = nint(360 / abs(self%dlon))
62 if (self%im < self%iwrap) self%iwrap = 0
65 end subroutine init_grib1
67 subroutine init_grib2(self, g2_desc)
68 class(ip_mercator_grid),
intent(inout) :: self
69 type(grib2_descriptor),
intent(in) :: g2_desc
71 integer :: iscan, jscan
74 associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
76 call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
81 self%RLAT1=float(igdtmpl(10))*1.0e-6
82 self%RLON1=float(igdtmpl(11))*1.0e-6
83 self%RLON2=float(igdtmpl(15))*1.0e-6
84 self%RLATI=float(igdtmpl(13))*1.0e-6
86 iscan=mod(igdtmpl(16)/128,2)
87 jscan=mod(igdtmpl(16)/64,2)
89 dy=float(igdtmpl(19))*1.0e-3
92 self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
93 self%DPHI=hj*dy/(self%RERTH*cos(self%RLATI/dpr))
98 self%nscan=mod(igdtmpl(16) / 32,2)
99 self%nscan_field_pos = self%nscan
101 self%iwrap = nint(360 / abs(self%dlon))
102 if(self%im < self%iwrap) self%iwrap = 0
105 end subroutine init_grib2
107 SUBROUTINE gdswzd_mercator(self,IOPT,NPTS,FILL, &
108 XPTS,YPTS,RLON,RLAT,NRET, &
109 CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
224 class(ip_mercator_grid),
intent(in) :: self
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
240 REAL :: RLAT1, RLON1, RLON2, RLATI
241 REAL :: XMAX, XMIN, YMAX, YMIN
244 IF(
PRESENT(crot)) crot=fill
245 IF(
PRESENT(srot)) srot=fill
246 IF(
PRESENT(xlon)) xlon=fill
247 IF(
PRESENT(xlat)) xlat=fill
248 IF(
PRESENT(ylon)) ylon=fill
249 IF(
PRESENT(ylat)) ylat=fill
250 IF(
PRESENT(area)) area=fill
267 ye=1-log(tan((rlat1+90)/2/dpr))/dphi
270 IF(im.EQ.nint(360/abs(dlon))) xmax=im+2
274 IF(
PRESENT(crot).AND.
PRESENT(srot))
THEN
279 IF(
PRESENT(xlon).AND.
PRESENT(xlat).AND.
PRESENT(ylon).AND.
PRESENT(ylat))
THEN
284 IF(
PRESENT(area))
THEN
291 IF(iopt.EQ.0.OR.iopt.EQ.1)
THEN
294 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
295 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
296 rlon(n)=mod(rlon1+dlon*(xpts(n)-1)+3600,360.)
297 rlat(n)=2*atan(exp(dphi*(ypts(n)-ye)))*dpr-90
299 IF(lrot)
CALL mercator_vect_rot(crot(n),srot(n))
300 IF(lmap)
CALL mercator_map_jacob(rlat(n),xlon(n),xlat(n),ylon(n),ylat(n))
301 IF(larea)
CALL mercator_grid_area(rlat(n),area(n))
310 ELSEIF(iopt.EQ.-1)
THEN
313 IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LT.90)
THEN
314 xpts(n)=1+hi*mod(hi*(rlon(n)-rlon1)+3600,360.)/dlon
315 ypts(n)=ye+log(tan((rlat(n)+90)/2/dpr))/dphi
316 IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
317 ypts(n).GE.ymin.AND.ypts(n).LE.ymax)
THEN
319 IF(lrot)
CALL mercator_vect_rot(crot(n),srot(n))
320 IF(lmap)
CALL mercator_map_jacob(rlat(n),xlon(n),xlat(n),ylon(n),ylat(n))
321 IF(larea)
CALL mercator_grid_area(rlat(n),area(n))
334 END SUBROUTINE gdswzd_mercator
337 SUBROUTINE mercator_error(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
370 INTEGER,
INTENT(IN ) :: IOPT, NPTS
372 REAL,
INTENT(IN ) :: FILL
373 REAL,
INTENT( OUT) :: RLAT(NPTS),RLON(NPTS)
374 REAL,
INTENT( OUT) :: XPTS(NPTS),YPTS(NPTS)
385 END SUBROUTINE mercator_error
387 SUBROUTINE mercator_vect_rot(CROT,SROT)
420 REAL,
INTENT( OUT) :: CROT, SROT
425 END SUBROUTINE mercator_vect_rot
427 SUBROUTINE mercator_map_jacob(RLAT,XLON,XLAT,YLON,YLAT)
460 REAL,
INTENT(IN ) :: RLAT
461 REAL,
INTENT( OUT) :: XLON, XLAT, YLON, YLAT
466 ylat=1./dphi/cos(rlat/dpr)/dpr
468 END SUBROUTINE mercator_map_jacob
470 SUBROUTINE mercator_grid_area(RLAT,AREA)
500 REAL,
INTENT(IN ) :: RLAT
501 REAL,
INTENT( OUT) :: AREA
503 area=rerth**2*cos(rlat/dpr)**2*dphi*dlon/dpr
505 END SUBROUTINE mercator_grid_area
507end module ip_mercator_grid_mod
Module containing common constants.
real(real64), parameter pi
PI.
real(real64), parameter dpr
Radians to degrees.
Uses derived type grid descriptor objects to abstract away the raw Grib-1 and Grib-2 grid definitions...
Abstract grid that holds fields and methods common to all grids.