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