NCEPLIBS-ip  4.3.0
ip_equid_cylind_grid_mod.F90
Go to the documentation of this file.
1 
7 
19  use ip_grid_mod
21  implicit none
22 
23  private
24  public :: ip_equid_cylind_grid
25 
26  type, extends(ip_grid) :: ip_equid_cylind_grid
27  real :: hi
28  real :: rlat1
29  real :: rlon1
30  real :: rlat2
31  real :: rlon2
32  real :: dlat
33  real :: dlon
34  contains
35  procedure :: init_grib1
36  procedure :: init_grib2
37  procedure :: gdswzd => gdswzd_equid_cylind
38  end type ip_equid_cylind_grid
39 
40  REAL :: dlat
41  REAL :: dlon
42  REAL :: rerth
43 
44 contains
45 
53  subroutine init_grib1(self, g1_desc)
54  class(ip_equid_cylind_grid), intent(inout) :: self
55  type(grib1_descriptor), intent(in) :: g1_desc
56 
57  integer :: iscan
58 
59  associate(kgds => g1_desc%gds)
60  self%IM=kgds(2)
61  self%JM=kgds(3)
62  self%RLAT1=kgds(4)*1.e-3
63  self%RLON1=kgds(5)*1.e-3
64  self%RLAT2=kgds(7)*1.e-3
65  self%RLON2=kgds(8)*1.e-3
66  iscan=mod(kgds(11)/128,2)
67  self%HI=(-1.)**iscan
68  self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
69  self%DLAT=(self%RLAT2-self%RLAT1)/(self%JM-1)
70 
71  ! defaults
72  self%iwrap = 0
73  self%jwrap1 = 0
74  self%jwrap2 = 0
75  self%nscan = mod(kgds(11) / 32, 2)
76  self%nscan_field_pos = self%nscan
77  self%kscan = 0
78 
79  self%iwrap = nint(360/abs(self%dlon))
80 
81  if(self%im < self%iwrap) self%iwrap=0
82  self%jwrap1 = 0
83  self%jwrap2 = 0
84  if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0) then
85  if(abs(self%rlat1) > 90-0.25*self%dlat) then
86  self%jwrap1 = 2
87  elseif(abs(self%rlat1) > 90-0.75*self%dlat) then
88  self%jwrap1 = 1
89  endif
90  if(abs(self%rlat2) > 90-0.25*self%dlat) then
91  self%jwrap2 = 2 * self%jm
92  elseif(abs(self%rlat2) > 90-0.75*self%dlat) then
93  self%jwrap2 = 2 * self%jm+1
94  endif
95  endif
96 
97  self%rerth = 6.3712e6
98  self%eccen_squared = 0.0
99  end associate
100 
101  end subroutine init_grib1
102 
109  subroutine init_grib2(self, g2_desc)
110  class(ip_equid_cylind_grid), intent(inout) :: self
111  type(grib2_descriptor), intent(in) :: g2_desc
112 
113  integer :: iscale, iscan
114 
115  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
116  self%IM=igdtmpl(8)
117  self%JM=igdtmpl(9)
118  iscale=igdtmpl(10)*igdtmpl(11)
119  IF(iscale==0) iscale=10**6
120  self%RLAT1=float(igdtmpl(12))/float(iscale)
121  self%RLON1=float(igdtmpl(13))/float(iscale)
122  self%RLAT2=float(igdtmpl(15))/float(iscale)
123  self%RLON2=float(igdtmpl(16))/float(iscale)
124  iscan=mod(igdtmpl(19)/128,2)
125  self%HI=(-1.)**iscan
126  self%DLON=self%HI*(mod(self%HI*(self%RLON2-self%RLON1)-1+3600,360.)+1)/(self%IM-1)
127  self%DLAT=(self%RLAT2-self%RLAT1)/(self%JM-1)
128 
129  self%nscan = mod(igdtmpl(19)/32,2)
130  self%nscan_field_pos = self%nscan
131  self%kscan = 0
132  self%iwrap = nint(360/abs(self%DLON))
133 
134  if(self%im.lt.self%iwrap) self%iwrap=0
135  self%jwrap1=0
136  self%jwrap2=0
137 
138  if(self%im < self%iwrap) self%iwrap=0
139  self%jwrap1 = 0
140  self%jwrap2 = 0
141  if(self%iwrap > 0 .and. mod(self%iwrap,2) == 0) then
142  if(abs(self%rlat1) > 90-0.25*self%dlat) then
143  self%jwrap1 = 2
144  elseif(abs(self%rlat1) > 90-0.75*self%dlat) then
145  self%jwrap1 = 1
146  endif
147  if(abs(self%rlat2) > 90-0.25*self%dlat) then
148  self%jwrap2 = 2 * self%jm
149  elseif(abs(self%rlat2) > 90-0.75*self%dlat) then
150  self%jwrap2 = 2 * self%jm+1
151  endif
152  endif
153 
154  call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
155 
156  end associate
157  end subroutine init_grib2
158 
203  SUBROUTINE gdswzd_equid_cylind(self,IOPT,NPTS,FILL, &
204  XPTS,YPTS,RLON,RLAT,NRET, &
205  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
206  IMPLICIT NONE
207  !
208  class(ip_equid_cylind_grid), intent(in) :: self
209  INTEGER, INTENT(IN ) :: IOPT, NPTS
210  INTEGER, INTENT( OUT) :: NRET
211  !
212  REAL, INTENT(IN ) :: FILL
213  REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
214  REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
215  REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
216  REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
217  REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
218  !
219  INTEGER :: IM, JM, N
220  !
221  LOGICAL :: LROT, LMAP, LAREA
222  !
223  REAL :: HI, RLAT1, RLON1, RLAT2, RLON2
224  REAL :: XMAX, XMIN, YMAX, YMIN
225  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226  IF(PRESENT(crot)) crot=fill
227  IF(PRESENT(srot)) srot=fill
228  IF(PRESENT(xlon)) xlon=fill
229  IF(PRESENT(xlat)) xlat=fill
230  IF(PRESENT(ylon)) ylon=fill
231  IF(PRESENT(ylat)) ylat=fill
232  IF(PRESENT(area)) area=fill
233 
234  im=self%im
235  jm=self%jm
236 
237  rlat1=self%rlat1
238  rlon1=self%rlon1
239  rlat2=self%rlat2
240  rlon2=self%rlon2
241 
242  hi=self%hi
243 
244  rerth = self%rerth
245  dlat = self%dlat
246  dlon = self%dlon
247 
248  xmin=0
249  xmax=im+1
250  IF(im.EQ.nint(360/abs(dlon))) xmax=im+2
251  ymin=0
252  ymax=jm+1
253  nret=0
254  IF(PRESENT(crot).AND.PRESENT(srot))THEN
255  lrot=.true.
256  ELSE
257  lrot=.false.
258  ENDIF
259  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
260  lmap=.true.
261  ELSE
262  lmap=.false.
263  ENDIF
264  IF(PRESENT(area))THEN
265  larea=.true.
266  ELSE
267  larea=.false.
268  ENDIF
269  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
271  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
272  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
273  DO n=1,npts
274  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
275  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
276  rlon(n)=mod(rlon1+dlon*(xpts(n)-1)+3600,360.)
277  rlat(n)=min(max(rlat1+dlat*(ypts(n)-1),-90.),90.)
278  nret=nret+1
279  IF(lrot) CALL equid_cylind_vect_rot(crot(n),srot(n))
280  IF(lmap) CALL equid_cylind_map_jacob(xlon(n),xlat(n),ylon(n),ylat(n))
281  IF(larea) CALL equid_cylind_grid_area(rlat(n),area(n))
282  ELSE
283  rlon(n)=fill
284  rlat(n)=fill
285  ENDIF
286  ENDDO
287  !$OMP END PARALLEL DO
288  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
290  ELSEIF(iopt.EQ.-1) THEN
291  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:NRET) SCHEDULE(STATIC)
292  DO n=1,npts
293  IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90) THEN
294  xpts(n)=1+hi*mod(hi*(rlon(n)-rlon1)+3600,360.)/dlon
295  ypts(n)=1+(rlat(n)-rlat1)/dlat
296  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
297  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
298  nret=nret+1
299  IF(lrot) CALL equid_cylind_vect_rot(crot(n),srot(n))
300  IF(lmap) CALL equid_cylind_map_jacob(xlon(n),xlat(n),ylon(n),ylat(n))
301  IF(larea) CALL equid_cylind_grid_area(rlat(n),area(n))
302  ELSE
303  xpts(n)=fill
304  ypts(n)=fill
305  ENDIF
306  ELSE
307  xpts(n)=fill
308  ypts(n)=fill
309  ENDIF
310  ENDDO
311  !$OMP END PARALLEL DO
312  ENDIF
313  END SUBROUTINE gdswzd_equid_cylind
314 
327  SUBROUTINE equid_cylind_vect_rot(CROT,SROT)
328  IMPLICIT NONE
329 
330  REAL, INTENT( OUT) :: CROT, SROT
331 
332  crot=1.0
333  srot=0.0
334 
335  END SUBROUTINE equid_cylind_vect_rot
336 
346  SUBROUTINE equid_cylind_map_jacob(XLON,XLAT,YLON,YLAT)
347  REAL, INTENT( OUT) :: XLON,XLAT,YLON,YLAT
348 
349  xlon=1.0/dlon
350  xlat=0.
351  ylon=0.
352  ylat=1.0/dlat
353 
354  END SUBROUTINE equid_cylind_map_jacob
355 
363  SUBROUTINE equid_cylind_grid_area(RLAT,AREA)
364  IMPLICIT NONE
365 
366  REAL, INTENT(IN ) :: RLAT
367  REAL, INTENT( OUT) :: AREA
368 
369  REAL, PARAMETER :: PI=3.14159265358979
370  REAL, PARAMETER :: DPR=180./pi
371 
372  REAL :: DSLAT, RLATU, RLATD
373 
374  rlatu=min(max(rlat+dlat/2,-90.),90.)
375  rlatd=min(max(rlat-dlat/2,-90.),90.)
376  dslat=sin(rlatu/dpr)-sin(rlatd/dpr)
377  area=rerth**2*abs(dslat*dlon)/dpr
378 
379  END SUBROUTINE equid_cylind_grid_area
380 
381 end module ip_equid_cylind_grid_mod
382 
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.
Equidistant cylindrical grib decoder and grid coordinate transformations.
subroutine equid_cylind_map_jacob(XLON, XLAT, YLON, YLAT)
Computes the map jacobians for a equidistant cylindrical grid.
real dlat
Grid resolution in degrees n/s direction.
subroutine init_grib2(self, g2_desc)
Initializes an equidistant cylindrical grid given a grib2_descriptor object.
subroutine init_grib1(self, g1_desc)
Initializes an equidistant cylindrical grid given a grib1_descriptor object.
real dlon
Grid resolution in degrees e/w direction.
subroutine equid_cylind_grid_area(RLAT, AREA)
Computes the grid box area for a equidistant cylindrical grid.
real rerth
Radius of the Earth.
subroutine equid_cylind_vect_rot(CROT, SROT)
Computes the vector rotation sines and cosines for a equidistant cylindrical grid.
subroutine gdswzd_equid_cylind(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
Calculates Earth coordinates (iopt = 1) or grid coorindates (iopt = -1) for equidistant cylindrical g...
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
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