NCEPLIBS-ip  4.4.0
ip_polar_stereo_grid_mod.F90
Go to the documentation of this file.
1 
5 
15  use ip_grid_mod
18  implicit none
19 
20  private
21  public :: ip_polar_stereo_grid
22 
23  type, extends(ip_grid) :: ip_polar_stereo_grid
24  logical :: elliptical
25  real :: rlat1
26  real :: rlon1
27  real :: orient
28  real :: h
29  real :: dxs
30  real :: dys
31  real :: slatr
35  integer :: irot
36  contains
37  procedure :: init_grib1
38  procedure :: init_grib2
41  procedure :: gdswzd => gdswzd_polar_stereo
42  end type ip_polar_stereo_grid
43 
44  INTEGER :: irot
45  REAL :: de2
46  REAL :: dxs
47  REAL :: dys
48  REAL :: e2
49  REAL :: rerth
50  REAL :: h
51  REAL :: orient
52  REAL :: tinyreal=tiny(1.0)
53 
54 CONTAINS
55 
63  subroutine init_grib1(self, g1_desc)
64  class(ip_polar_stereo_grid), intent(inout) :: self
65  type(grib1_descriptor), intent(in) :: g1_desc
66 
67  REAL, PARAMETER :: SLAT=60.0 ! standard latitude according grib1 standard
68 
69  real :: dx, dy, hi, hj
70  integer :: iproj, iscan, jscan
71 
72  associate(kgds => g1_desc%gds)
73  self%ELLIPTICAL=mod(kgds(6)/64,2).EQ.1
74 
75  if (.not. self%elliptical) then
76  self%rerth = 6.3712e6
77  self%eccen_squared = 0d0
78  else
79  self%rerth = rerth_wgs84
80  self%eccen_squared = e2_wgs84 !wgs84 datum
81  end if
82 
83  self%IM=kgds(2)
84  self%JM=kgds(3)
85 
86  self%RLAT1=kgds(4)*1.e-3
87  self%RLON1=kgds(5)*1.e-3
88 
89  self%IROT=mod(kgds(6)/8,2)
90 
91  self%SLATR=slat/dpr
92 
93  self%ORIENT=kgds(7)*1.e-3
94 
95  dx=kgds(8)
96  dy=kgds(9)
97 
98  iproj=mod(kgds(10)/128,2)
99  iscan=mod(kgds(11)/128,2)
100  jscan=mod(kgds(11)/64,2)
101 
102  self%H=(-1.)**iproj
103  hi=(-1.)**iscan
104  hj=(-1.)**(1-jscan)
105 
106  IF(abs(self%H+1.).LT.tinyreal) self%ORIENT=self%ORIENT+180.
107 
108  self%DXS=dx*hi
109  self%DYS=dy*hj
110 
111  self%iwrap= 0
112  self%jwrap1 = 0
113  self%jwrap2 = 0
114  self%nscan = mod(kgds(11) / 32, 2)
115  self%nscan_field_pos = self%nscan
116  self%kscan = 0
117  end associate
118 
119  end subroutine init_grib1
120 
128  subroutine init_grib2(self, g2_desc)
129  class(ip_polar_stereo_grid), intent(inout) :: self
130  type(grib2_descriptor), intent(in) :: g2_desc
131 
132  real :: slat, dx, dy, hi, hj
133  integer :: iproj, iscan, jscan
134 
135  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
136  call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
137 
138  self%ELLIPTICAL = self%eccen_squared > 0.0
139 
140  self%IM=igdtmpl(8)
141  self%JM=igdtmpl(9)
142 
143  self%RLAT1=float(igdtmpl(10))*1.e-6
144  self%RLON1=float(igdtmpl(11))*1.e-6
145 
146  self%IROT=mod(igdtmpl(12)/8,2)
147 
148  slat=float(abs(igdtmpl(13)))*1.e-6
149  self%SLATR=slat/dpr
150 
151  self%ORIENT=float(igdtmpl(14))*1.e-6
152 
153  dx=float(igdtmpl(15))*1.e-3
154  dy=float(igdtmpl(16))*1.e-3
155 
156  iproj=mod(igdtmpl(17)/128,2)
157  iscan=mod(igdtmpl(18)/128,2)
158  jscan=mod(igdtmpl(18)/64,2)
159 
160  self%H=(-1.)**iproj
161  hi=(-1.)**iscan
162  hj=(-1.)**(1-jscan)
163 
164  self%DXS=dx*hi
165  self%DYS=dy*hj
166 
167  self%nscan = mod(igdtmpl(18) / 32, 2)
168  self%nscan_field_pos = self%nscan
169  self%iwrap = 0
170  self%jwrap1 = 0
171  self%jwrap2 = 0
172  self%kscan = 0
173  end associate
174  end subroutine init_grib2
175 
239  SUBROUTINE gdswzd_polar_stereo(self,IOPT,NPTS, &
240  FILL,XPTS,YPTS,RLON,RLAT,NRET, &
241  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
242  IMPLICIT NONE
243  !
244 
245  class(ip_polar_stereo_grid), intent(in) :: self
246  INTEGER, INTENT(IN ) :: IOPT, NPTS
247  INTEGER, INTENT( OUT) :: NRET
248  !
249  REAL, INTENT(IN ) :: FILL
250  REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
251  REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
252  REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
253  REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
254  REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
255  !
256  INTEGER :: IM, JM
257  INTEGER :: ITER, N
258  !
259  LOGICAL :: ELLIPTICAL, LROT, LMAP, LAREA
260  !
261  REAL :: ALAT, ALAT1, ALONG, DIFF
262  REAL :: DI, DJ, DE
263  REAL :: DR, E, E_OVER_2
264  REAL :: MC, SLATR
265  REAL :: RLAT1, RLON1, RHO, T, TC
266  REAL :: XMAX, XMIN, YMAX, YMIN
267  REAL :: XP, YP, DR2
268  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269  IF(PRESENT(crot)) crot=fill
270  IF(PRESENT(srot)) srot=fill
271  IF(PRESENT(xlon)) xlon=fill
272  IF(PRESENT(xlat)) xlat=fill
273  IF(PRESENT(ylon)) ylon=fill
274  IF(PRESENT(ylat)) ylat=fill
275  IF(PRESENT(area)) area=fill
276 
277  elliptical = self%elliptical
278  im=self%im
279  jm=self%jm
280 
281  rlat1=self%rlat1
282  rlon1=self%rlon1
283 
284  irot=self%irot
285  slatr=self%slatr
286  orient=self%orient
287 
288  h=self%h
289  dxs=self%dxs
290  dys=self%dys
291 
292  rerth = self%rerth
293  e2 = self%eccen_squared
294  !
295  ! FIND X/Y OF POLE
296  IF (.NOT.elliptical) THEN
297  de=(1.+sin(slatr))*rerth
298  dr=de*cos(rlat1/dpr)/(1+h*sin(rlat1/dpr))
299  xp=1-h*sin((rlon1-orient)/dpr)*dr/dxs
300  yp=1+cos((rlon1-orient)/dpr)*dr/dys
301  de2=de**2
302  ELSE
303  e=sqrt(e2)
304  e_over_2=e*0.5
305  alat=h*rlat1/dpr
306  along = (rlon1-orient)/dpr
307  t=tan(pi4-alat/2.)/((1.-e*sin(alat))/ &
308  (1.+e*sin(alat)))**(e_over_2)
309  tc=tan(pi4-slatr/2.)/((1.-e*sin(slatr))/ &
310  (1.+e*sin(slatr)))**(e_over_2)
311  mc=cos(slatr)/sqrt(1.0-e2*(sin(slatr)**2))
312  rho=rerth*mc*t/tc
313  yp = 1.0 + rho*cos(h*along)/dys
314  xp = 1.0 - rho*sin(h*along)/dxs
315  ENDIF ! ELLIPTICAL
316  xmin=0
317  xmax=im+1
318  ymin=0
319  ymax=jm+1
320  nret=0
321  IF(PRESENT(crot).AND.PRESENT(srot))THEN
322  lrot=.true.
323  ELSE
324  lrot=.false.
325  ENDIF
326  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
327  lmap=.true.
328  ELSE
329  lmap=.false.
330  ENDIF
331  IF(PRESENT(area))THEN
332  larea=.true.
333  ELSE
334  larea=.false.
335  ENDIF
336  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
337  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
338  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
339  IF(.NOT.elliptical)THEN
340  !$OMP PARALLEL DO PRIVATE(N,DI,DJ,DR2) REDUCTION(+:NRET) SCHEDULE(STATIC)
341  DO n=1,npts
342  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
343  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
344  di=(xpts(n)-xp)*dxs
345  dj=(ypts(n)-yp)*dys
346  dr2=di**2+dj**2
347  IF(dr2.LT.de2*1.e-6) THEN
348  rlon(n)=0.
349  rlat(n)=h*90.
350  ELSE
351  rlon(n)=mod(orient+h*dpr*atan2(di,-dj)+3600,360.)
352  rlat(n)=h*dpr*asin((de2-dr2)/(de2+dr2))
353  ENDIF
354  nret=nret+1
355  IF(lrot) CALL polar_stereo_vect_rot(rlon(n),crot(n),srot(n))
356  IF(lmap) CALL polar_stereo_map_jacob(rlon(n),rlat(n),dr2, &
357  xlon(n),xlat(n),ylon(n),ylat(n))
358  IF(larea) CALL polar_stereo_grid_area(rlat(n),dr2,area(n))
359  ELSE
360  rlon(n)=fill
361  rlat(n)=fill
362  ENDIF
363  ENDDO
364  !$OMP END PARALLEL DO
365  ELSE ! ELLIPTICAL
366  !$OMP PARALLEL DO PRIVATE(N,DI,DJ,RHO,T,ALONG,ALAT1,ALAT,DIFF) &
367  !$OMP& REDUCTION(+:NRET) SCHEDULE(STATIC)
368  DO n=1,npts
369  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
370  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
371  di=(xpts(n)-xp)*dxs
372  dj=(ypts(n)-yp)*dys
373  rho=sqrt(di*di+dj*dj)
374  t=(rho*tc)/(rerth*mc)
375  IF(abs(ypts(n)-yp)<0.01)THEN
376  IF(di>0.0) along=orient+h*90.0
377  IF(di<=0.0) along=orient-h*90.0
378  ELSE
379  along=orient+h*atan(di/(-dj))*dpr
380  IF(dj>0) along=along+180.
381  END IF
382  alat1=pi2-2.0*atan(t)
383  DO iter=1,10
384  alat = pi2 - 2.0*atan(t*(((1.0-e*sin(alat1))/ &
385  (1.0+e*sin(alat1)))**(e_over_2)))
386  diff = abs(alat-alat1)*dpr
387  IF (diff < 0.000001) EXIT
388  alat1=alat
389  ENDDO
390  rlat(n)=h*alat*dpr
391  rlon(n)=along
392  IF(rlon(n)<0.0) rlon(n)=rlon(n)+360.
393  IF(rlon(n)>360.0) rlon(n)=rlon(n)-360.0
394  nret=nret+1
395  IF(lrot) CALL polar_stereo_vect_rot(rlon(n),crot(n),srot(n))
396  ELSE
397  rlon(n)=fill
398  rlat(n)=fill
399  ENDIF
400  ENDDO
401  !$OMP END PARALLEL DO
402  ENDIF ! ELLIPTICAL
403  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
404  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
405  ELSEIF(iopt.EQ.-1) THEN
406  IF(.NOT.elliptical)THEN
407  !$OMP PARALLEL DO PRIVATE(N,DR,DR2) REDUCTION(+:NRET) SCHEDULE(STATIC)
408  DO n=1,npts
409  IF(abs(rlon(n)).LT.(360.+tinyreal).AND.abs(rlat(n)).LT.(90.+tinyreal).AND. &
410  abs(h*rlat(n)+90).GT.tinyreal) THEN
411  dr=de*tan((90-h*rlat(n))/2/dpr)
412  dr2=dr**2
413  xpts(n)=xp+h*sin((rlon(n)-orient)/dpr)*dr/dxs
414  ypts(n)=yp-cos((rlon(n)-orient)/dpr)*dr/dys
415  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
416  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
417  nret=nret+1
418  IF(lrot) CALL polar_stereo_vect_rot(rlon(n),crot(n),srot(n))
419  IF(lmap) CALL polar_stereo_map_jacob(rlon(n),rlat(n),dr2, &
420  xlon(n),xlat(n),ylon(n),ylat(n))
421  IF(larea) CALL polar_stereo_grid_area(rlat(n),dr2,area(n))
422  ELSE
423  xpts(n)=fill
424  ypts(n)=fill
425  ENDIF
426  ELSE
427  xpts(n)=fill
428  ypts(n)=fill
429  ENDIF
430  ENDDO
431  !$OMP END PARALLEL DO
432  ELSE ! ELLIPTICAL CASE
433  !$OMP PARALLEL DO PRIVATE(N,ALAT,ALONG,T,RHO) REDUCTION(+:NRET) SCHEDULE(STATIC)
434  DO n=1,npts
435  IF(abs(rlon(n)).LT.(360+tinyreal).AND.abs(rlat(n)).LT.(90+tinyreal).AND. &
436  abs(h*rlat(n)+90).GT.tinyreal) THEN
437  alat = h*rlat(n)/dpr
438  along = (rlon(n)-orient)/dpr
439  t=tan(pi4-alat*0.5)/((1.-e*sin(alat))/ &
440  (1.+e*sin(alat)))**(e_over_2)
441  rho=rerth*mc*t/tc
442  xpts(n)= xp + rho*sin(h*along) / dxs
443  ypts(n)= yp - rho*cos(h*along) / dys
444  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
445  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
446  nret=nret+1
447  IF(lrot) CALL polar_stereo_vect_rot(rlon(n),crot(n),srot(n))
448  ELSE
449  xpts(n)=fill
450  ypts(n)=fill
451  ENDIF
452  ELSE
453  xpts(n)=fill
454  ypts(n)=fill
455  ENDIF
456  ENDDO
457  !$OMP END PARALLEL DO
458  ENDIF
459  ENDIF
460  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
461  END SUBROUTINE gdswzd_polar_stereo
462 
482  SUBROUTINE polar_stereo_vect_rot(RLON, CROT, SROT)
483  IMPLICIT NONE
484 
485  REAL, INTENT(IN ) :: RLON
486  REAL, INTENT( OUT) :: CROT, SROT
487 
488  IF(irot.EQ.1) THEN
489  crot=h*cos((rlon-orient)/dpr)
490  srot=sin((rlon-orient)/dpr)
491  ELSE
492  crot=1.
493  srot=0.
494  ENDIF
495 
496  END SUBROUTINE polar_stereo_vect_rot
497 
520  SUBROUTINE polar_stereo_map_jacob(RLON,RLAT,DR2,XLON,XLAT,YLON,YLAT)
521  IMPLICIT NONE
522 
523  REAL, INTENT(IN ) :: RLON, RLAT, DR2
524  REAL, INTENT( OUT) :: XLON, XLAT, YLON, YLAT
525 
526  REAL :: CLAT, DE, DR
527 
528  IF(dr2.LT.de2*1.e-6) THEN
529  de=sqrt(de2)
530  xlon=0.
531  xlat=-sin((rlon-orient)/dpr)/dpr*de/dxs/2
532  ylon=0.
533  ylat=h*cos((rlon-orient)/dpr)/dpr*de/dys/2
534  ELSE
535  dr=sqrt(dr2)
536  clat=cos(rlat/dpr)
537  xlon=h*cos((rlon-orient)/dpr)/dpr*dr/dxs
538  xlat=-sin((rlon-orient)/dpr)/dpr*dr/dxs/clat
539  ylon=sin((rlon-orient)/dpr)/dpr*dr/dys
540  ylat=h*cos((rlon-orient)/dpr)/dpr*dr/dys/clat
541  ENDIF
542 
543  END SUBROUTINE polar_stereo_map_jacob
544 
563  SUBROUTINE polar_stereo_grid_area(RLAT, DR2, AREA)
564  IMPLICIT NONE
565 
566  REAL, INTENT(IN ) :: RLAT, DR2
567  REAL, INTENT( OUT) :: AREA
568 
569  REAL :: CLAT
570 
571  IF(dr2.LT.de2*1.e-6) THEN
572  area=rerth**2*abs(dxs)*abs(dys)*4/de2
573  ELSE
574  clat=cos(rlat/dpr)
575  area=rerth**2*clat**2*abs(dxs)*abs(dys)/dr2
576  ENDIF
577 
578  END SUBROUTINE polar_stereo_grid_area
579 
580 end module ip_polar_stereo_grid_mod
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 pi2
PI / 2.0.
real, parameter pi
PI.
real, parameter pi4
PI / 4.0.
real, parameter rerth_wgs84
Radius of the Earth defined by WGS-84.
real, parameter dpr
Radians to degrees.
real, parameter e2_wgs84
Eccentricity squared of Earth defined by WGS-84.
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 polar stereographic azimuthal.
integer irot
Local copy of irot.
subroutine polar_stereo_grid_area(RLAT, DR2, AREA)
Grid box area for polar stereographic grids.
real rerth
Radius of the Earth.
subroutine polar_stereo_map_jacob(RLON, RLAT, DR2, XLON, XLAT, YLON, YLAT)
Map jacobians for polar stereographic grids.
subroutine init_grib2(self, g2_desc)
Initializes a polar stereographic grid given a grib2_descriptor object.
subroutine init_grib1(self, g1_desc)
Initializes a polar stereographic grid given a grib1_descriptor object.
real orient
Local copy of orient.
subroutine polar_stereo_vect_rot(RLON, CROT, SROT)
Vector rotation fields for polar stereographic grids.
real e2
Eccentricity squared.
subroutine gdswzd_polar_stereo(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for polar stereographic azimuthal.
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