NCEPLIBS-ip  5.0.0
ip_mercator_grid_mod.F90
Go to the documentation of this file.
1 
5 
14  use ip_grid_mod
15  use ip_constants_mod, only: dpr, pi
17  implicit none
18 
19  private
20  public :: ip_mercator_grid
21 
22  type, extends(ip_grid) :: ip_mercator_grid
23  real :: rlat1
24  real :: rlon1
25  real :: rlon2
26  real :: rlati
27  real :: hi
28  real :: dlon
29  real :: dphi
30  contains
32  procedure :: init_grib1
34  procedure :: init_grib2
37  procedure :: gdswzd => gdswzd_mercator
38  end type ip_mercator_grid
39 
40  REAL :: dlon
41  REAL :: dphi
42  REAL :: rerth
43 
44 CONTAINS
45 
52  subroutine init_grib1(self, g1_desc)
53  class(ip_mercator_grid), intent(inout) :: self
54  type(grib1_descriptor), intent(in) :: g1_desc
55 
56  integer :: iscan, jscan
57  real :: dy, hj
58 
59  associate(kgds => g1_desc%gds)
60  self%rerth = 6.3712e6
61  self%eccen_squared = 0.0
62 
63  self%IM=kgds(2)
64  self%JM=kgds(3)
65 
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
70 
71  iscan=mod(kgds(11)/128,2)
72  jscan=mod(kgds(11)/64,2)
73 
74  dy=kgds(13)
75  self%HI=(-1.)**iscan
76  hj=(-1.)**(1-jscan)
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))
79 
80  ! defaults
81  self%iwrap = 0
82  self%jwrap1 = 0
83  self%jwrap2 = 0
84  self%nscan = mod(kgds(11) / 32, 2)
85  self%nscan_field_pos = self%nscan
86  self%kscan = 0
87 
88  self%iwrap = nint(360 / abs(self%dlon))
89  if (self%im < self%iwrap) self%iwrap = 0
90  end associate
91 
92  end subroutine init_grib1
93 
100  subroutine init_grib2(self, g2_desc)
101  class(ip_mercator_grid), intent(inout) :: self
102  type(grib2_descriptor), intent(in) :: g2_desc
103 
104  integer :: iscan, jscan
105  real :: hj, dy
106 
107  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
108 
109  call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
110 
111  self%IM=igdtmpl(8)
112  self%JM=igdtmpl(9)
113 
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
118 
119  iscan=mod(igdtmpl(16)/128,2)
120  jscan=mod(igdtmpl(16)/64,2)
121 
122  dy=float(igdtmpl(19))*1.0e-3
123  self%HI=(-1.)**iscan
124  hj=(-1.)**(1-jscan)
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))
127 
128  self%jwrap1 = 0
129  self%jwrap2 = 0
130  self%kscan = 0
131  self%nscan=mod(igdtmpl(16) / 32,2)
132  self%nscan_field_pos = self%nscan
133 
134  self%iwrap = nint(360 / abs(self%dlon))
135  if(self%im < self%iwrap) self%iwrap = 0
136 
137  end associate
138  end subroutine init_grib2
139 
198  SUBROUTINE gdswzd_mercator(self,IOPT,NPTS,FILL, &
199  XPTS,YPTS,RLON,RLAT,NRET, &
200  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
201  IMPLICIT NONE
202  !
203  class(ip_mercator_grid), intent(in) :: self
204  INTEGER, INTENT(IN ) :: IOPT, NPTS
205  INTEGER, INTENT( OUT) :: NRET
206  !
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)
213  !
214  INTEGER :: IM, JM, N
215  !
216  LOGICAL :: LROT, LMAP, LAREA
217  !
218  REAL :: HI
219  REAL :: RLAT1, RLON1, RLON2, RLATI
220  REAL :: XMAX, XMIN, YMAX, YMIN
221  REAL :: YE
222  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
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
230  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231 
232  im=self%im
233  jm=self%jm
234 
235  rlat1=self%rlat1
236  rlon1=self%rlon1
237  rlon2=self%rlon2
238  rlati=self%rlati
239 
240  hi=self%hi
241 
242  dlon=self%dlon
243  dphi=self%dphi
244  rerth = self%rerth
245 
246  ye=1-log(tan((rlat1+90)/2/dpr))/dphi
247  xmin=0
248  xmax=im+1
249  IF(im.EQ.nint(360/abs(dlon))) xmax=im+2
250  ymin=0
251  ymax=jm+1
252  nret=0
253  IF(PRESENT(crot).AND.PRESENT(srot))THEN
254  lrot=.true.
255  ELSE
256  lrot=.false.
257  ENDIF
258  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
259  lmap=.true.
260  ELSE
261  lmap=.false.
262  ENDIF
263  IF(PRESENT(area))THEN
264  larea=.true.
265  ELSE
266  larea=.false.
267  ENDIF
268  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
270  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
271  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
272  DO n=1,npts
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
277  nret=nret+1
278  IF(lrot) CALL mercator_vect_rot(crot(n),srot(n))
279  IF(lmap) CALL mercator_map_jacob(rlat(n),xlon(n),xlat(n),ylon(n),ylat(n))
280  IF(larea) CALL mercator_grid_area(rlat(n),area(n))
281  ELSE
282  rlon(n)=fill
283  rlat(n)=fill
284  ENDIF
285  ENDDO
286  !$OMP END PARALLEL DO
287  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
289  ELSEIF(iopt.EQ.-1) THEN
290  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
291  DO n=1,npts
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
297  nret=nret+1
298  IF(lrot) CALL mercator_vect_rot(crot(n),srot(n))
299  IF(lmap) CALL mercator_map_jacob(rlat(n),xlon(n),xlat(n),ylon(n),ylat(n))
300  IF(larea) CALL mercator_grid_area(rlat(n),area(n))
301  ELSE
302  xpts(n)=fill
303  ypts(n)=fill
304  ENDIF
305  ELSE
306  xpts(n)=fill
307  ypts(n)=fill
308  ENDIF
309  ENDDO
310  !$OMP END PARALLEL DO
311  ENDIF
312  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
313  END SUBROUTINE gdswzd_mercator
314 
331  SUBROUTINE mercator_vect_rot(CROT,SROT)
332  IMPLICIT NONE
333 
334  REAL, INTENT( OUT) :: CROT, SROT
335 
336  crot=1.0
337  srot=0.0
338 
339  END SUBROUTINE mercator_vect_rot
340 
359  SUBROUTINE mercator_map_jacob(RLAT,XLON,XLAT,YLON,YLAT)
360  IMPLICIT NONE
361 
362  REAL, INTENT(IN ) :: RLAT
363  REAL, INTENT( OUT) :: XLON, XLAT, YLON, YLAT
364 
365  xlon=1./dlon
366  xlat=0.
367  ylon=0.
368  ylat=1./dphi/cos(rlat/dpr)/dpr
369 
370  END SUBROUTINE mercator_map_jacob
371 
387  SUBROUTINE mercator_grid_area(RLAT,AREA)
388  IMPLICIT NONE
389 
390  REAL, INTENT(IN ) :: RLAT
391  REAL, INTENT( OUT) :: AREA
392 
393  area=rerth**2*cos(rlat/dpr)**2*dphi*dlon/dpr
394 
395  END SUBROUTINE mercator_grid_area
396 
397 end module ip_mercator_grid_mod
398 
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.
Determine earth radius and shape.
Module containing common constants.
real, parameter pi
PI.
real, parameter dpr
Radians to degrees.
Users derived type grid descriptor objects to abstract away the raw GRIB1 and GRIB2 grid definitions.
Abstract ip_grid type.
Definition: ip_grid_mod.F90:10
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.
Definition: ip_grid_mod.F90:52