NCEPLIBS-ip  4.1.0
ip_rot_equid_cylind_grid_mod.F90
Go to the documentation of this file.
1 
7 
20  use iso_fortran_env, only: real64
22  use ip_grid_mod
23  use constants_mod, only: dpr, pi
25  implicit none
26 
27  private
28  public :: ip_rot_equid_cylind_grid
29 
30  integer, parameter :: kd = real64
31 
33  real(kd) :: clat0
34  real(kd) :: dlats
35  real(kd) :: dlons
36  real(kd) :: rlon0
37  real(kd) :: slat0
38  real(kd) :: wbd
39  real(kd) :: sbd
43  integer :: irot
44  contains
47  procedure :: init_grib1
50  procedure :: init_grib2
55 
56  INTEGER :: irot
57  REAL(kind=kd) :: rerth
58  REAL(kind=kd) :: clat0
59  REAL(kind=kd) :: dlats
60  REAL(kind=kd) :: dlons
61  REAL(kind=kd) :: rlon0
62  REAL(kind=kd) :: slat0
63 
64 CONTAINS
65 
73  subroutine init_grib1(self, g1_desc)
74  class(ip_rot_equid_cylind_grid), intent(inout) :: self
75  type(grib1_descriptor), intent(in) :: g1_desc
76 
77  real(kd) :: rlat1, rlon1, rlat0, rlat2, rlon2, nbd, ebd
78  real(kd) :: hs, hs2, slat1, slat2, slatr, clon1, clon2, clat1, clat2, clatr, clonr, rlonr, rlatr
79 
80  associate(kgds => g1_desc%gds)
81  self%rerth = 6.3712e6_kd
82  self%eccen_squared = 0d0
83 
84  rlat1=kgds(4)*1.e-3_kd
85  rlon1=kgds(5)*1.e-3_kd
86  rlat0=kgds(7)*1.e-3_kd
87  self%RLON0=kgds(8)*1.e-3_kd
88  rlat2=kgds(12)*1.e-3_kd
89  rlon2=kgds(13)*1.e-3_kd
90 
91  self%IROT=mod(kgds(6)/8,2)
92  self%IM=kgds(2)
93  self%JM=kgds(3)
94 
95  slat1=sin(rlat1/dpr)
96  clat1=cos(rlat1/dpr)
97  self%SLAT0=sin(rlat0/dpr)
98  self%CLAT0=cos(rlat0/dpr)
99 
100  hs=sign(1._kd,mod(rlon1-self%RLON0+180+3600,360._kd)-180)
101  clon1=cos((rlon1-self%RLON0)/dpr)
102  slatr=self%CLAT0*slat1-self%SLAT0*clat1*clon1
103  clatr=sqrt(1-slatr**2)
104  clonr=(self%CLAT0*clat1*clon1+self%SLAT0*slat1)/clatr
105  rlatr=dpr*asin(slatr)
106  rlonr=hs*dpr*acos(clonr)
107 
108  self%WBD=rlonr
109  self%SBD=rlatr
110  slat2=sin(rlat2/dpr)
111  clat2=cos(rlat2/dpr)
112  hs2=sign(1._kd,mod(rlon2-self%RLON0+180+3600,360._kd)-180)
113  clon2=cos((rlon2-self%RLON0)/dpr)
114  slatr=self%CLAT0*slat2-self%SLAT0*clat2*clon2
115  clatr=sqrt(1-slatr**2)
116  clonr=(self%CLAT0*clat2*clon2+self%SLAT0*slat2)/clatr
117  nbd=dpr*asin(slatr)
118  ebd=hs2*dpr*acos(clonr)
119  self%DLATS=(nbd-self%SBD)/float(self%JM-1)
120  self%DLONS=(ebd-self%WBD)/float(self%IM-1)
121 
122  self%iwrap = 0
123  self%jwrap1 = 0
124  self%jwrap2 = 0
125  self%nscan = mod(kgds(11) / 32, 2)
126  self%nscan_field_pos = self%nscan
127  self%kscan = 0
128  end associate
129 
130  end subroutine init_grib1
131 
139  subroutine init_grib2(self, g2_desc)
140  class(ip_rot_equid_cylind_grid), intent(inout) :: self
141  type(grib2_descriptor), intent(in) :: g2_desc
142 
143  real(kd) :: rlat1, rlon1, rlat0, rlat2, rlon2, nbd, ebd
144  integer :: iscale
145  integer :: i_offset_odd, i_offset_even, j_offset
146 
147  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
148 
149  CALL earth_radius(igdtmpl,igdtlen,self%rerth,self%eccen_squared)
150 
151  i_offset_odd=mod(igdtmpl(19)/8,2)
152  i_offset_even=mod(igdtmpl(19)/4,2)
153  j_offset=mod(igdtmpl(19)/2,2)
154 
155  iscale=igdtmpl(10)*igdtmpl(11)
156  IF(iscale==0) iscale=10**6
157 
158  rlat1=float(igdtmpl(12))/float(iscale)
159  rlon1=float(igdtmpl(13))/float(iscale)
160  rlat0=float(igdtmpl(20))/float(iscale)
161  rlat0=rlat0+90.0_kd
162 
163  self%RLON0=float(igdtmpl(21))/float(iscale)
164 
165  rlat2=float(igdtmpl(15))/float(iscale)
166  rlon2=float(igdtmpl(16))/float(iscale)
167 
168  self%IROT=mod(igdtmpl(14)/8,2)
169  self%IM=igdtmpl(8)
170  self%JM=igdtmpl(9)
171 
172  self%SLAT0=sin(rlat0/dpr)
173  self%CLAT0=cos(rlat0/dpr)
174 
175  self%WBD=rlon1
176  IF (self%WBD > 180.0) self%WBD = self%WBD - 360.0
177  self%SBD=rlat1
178 
179  nbd=rlat2
180  ebd=rlon2
181 
182  self%DLATS=(nbd-self%SBD)/float(self%JM-1)
183  self%DLONS=(ebd-self%WBD)/float(self%IM-1)
184 
185  IF(i_offset_odd==1) self%WBD=self%WBD+(0.5_kd*self%DLONS)
186  IF(j_offset==1) self%SBD=self%SBD+(0.5_kd*self%DLATS)
187 
188  self%iwrap = 0
189  self%jwrap1 = 0
190  self%jwrap2 = 0
191  self%kscan = 0
192  self%nscan = mod(igdtmpl(19) / 32, 2)
193  self%nscan_field_pos = self%nscan
194  end associate
195  end subroutine init_grib2
196 
256  SUBROUTINE gdswzd_rot_equid_cylind(self,IOPT,NPTS, &
257  FILL,XPTS,YPTS,RLON,RLAT,NRET, &
258  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
259  IMPLICIT NONE
260 
261  class(ip_rot_equid_cylind_grid), intent(in) :: self
262  INTEGER, INTENT(IN ) :: IOPT, NPTS
263  INTEGER, INTENT( OUT) :: NRET
264  !
265  REAL, INTENT(IN ) :: FILL
266  REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
267  REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
268  REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
269  REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
270  REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
271  !
272  INTEGER :: IM,JM,N
273  !
274  LOGICAL :: LROT, LMAP, LAREA
275  !
276  REAL(KIND=kd) :: hs
277  REAL(KIND=kd) :: clonr,clatr,slatr
278  REAL(KIND=kd) :: clat,slat,clon
279  REAL(KIND=kd) :: rlatr,rlonr
280  REAL(KIND=kd) :: wbd,sbd
281  REAL :: XMIN,XMAX,YMIN,YMAX
282  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283  IF(PRESENT(crot)) crot=fill
284  IF(PRESENT(srot)) srot=fill
285  IF(PRESENT(xlon)) xlon=fill
286  IF(PRESENT(xlat)) xlat=fill
287  IF(PRESENT(ylon)) ylon=fill
288  IF(PRESENT(ylat)) ylat=fill
289  IF(PRESENT(area)) area=fill
290  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291  ! IS THE EARTH RADIUS DEFINED?
292  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293  ! IS THIS AN "E"-STAGGER GRID? ROUTINE CAN'T PROCESS THOSE.
294  ! I_OFFSET_ODD=MOD(IGDTMPL(19)/8,2)
295  ! I_OFFSET_EVEN=MOD(IGDTMPL(19)/4,2)
296  ! J_OFFSET=MOD(IGDTMPL(19)/2,2)
297  ! IF(I_OFFSET_ODD/=I_OFFSET_EVEN) THEN
298  ! CALL ROT_EQUID_CYLIND_ERROR(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
299  ! RETURN
300  ! ENDIF
301  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
302 
303 
304  rlon0=self%rlon0
305  irot=self%irot
306 
307  im=self%im
308  jm=self%jm
309 
310  slat0=self%slat0
311  clat0=self%clat0
312 
313  wbd=self%wbd
314  sbd=self%sbd
315 
316  dlats=self%dlats
317  dlons=self%dlons
318 
319  xmin=0
320  xmax=im+1
321  ymin=0
322  ymax=jm+1
323  nret=0
324 
325  rerth = self%rerth
326  IF(rerth<0.)THEN
327  CALL rot_equid_cylind_error(iopt,fill,rlat,rlon,xpts,ypts,npts)
328  RETURN
329  ENDIF
330 
331  IF(PRESENT(crot).AND.PRESENT(srot))THEN
332  lrot=.true.
333  ELSE
334  lrot=.false.
335  ENDIF
336  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
337  lmap=.true.
338  ELSE
339  lmap=.false.
340  ENDIF
341  IF(PRESENT(area))THEN
342  larea=.true.
343  ELSE
344  larea=.false.
345  ENDIF
346  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
347  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
348  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
349  !$OMP PARALLEL DO PRIVATE(N,RLONR,RLATR,HS,CLONR,SLATR,CLATR,SLAT,CLAT,CLON) &
350  !$OMP& REDUCTION(+:NRET) SCHEDULE(STATIC)
351  DO n=1,npts
352  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
353  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
354  rlonr=wbd+(xpts(n)-1._kd)*dlons
355  rlatr=sbd+(ypts(n)-1._kd)*dlats
356  IF(rlonr <= 0._kd) THEN
357  hs=-1.0_kd
358  ELSE
359  hs=1.0_kd
360  ENDIF
361  clonr=cos(rlonr/dpr)
362  slatr=sin(rlatr/dpr)
363  clatr=cos(rlatr/dpr)
364  slat=clat0*slatr+slat0*clatr*clonr
365  IF(slat.LE.-1) THEN
366  clat=0.
367  clon=cos(rlon0/dpr)
368  rlon(n)=0.
369  rlat(n)=-90.
370  ELSEIF(slat.GE.1) THEN
371  clat=0.
372  clon=cos(rlon0/dpr)
373  rlon(n)=0.
374  rlat(n)=90.
375  ELSE
376  clat=sqrt(1-slat**2)
377  clon=(clat0*clatr*clonr-slat0*slatr)/clat
378  clon=min(max(clon,-1._kd),1._kd)
379  rlon(n)=real(mod(rlon0+hs*dpr*acos(clon)+3600,360._kd))
380  rlat(n)=real(dpr*asin(slat))
381  ENDIF
382  nret=nret+1
383  IF(lrot) CALL rot_equid_cylind_vect_rot(rlon(n), clatr, slatr, &
384  clat, slat, clon, crot(n), srot(n))
385  IF(lmap) CALL rot_equid_cylind_map_jacob(fill, rlon(n), clatr, &
386  clat, slat, clon, xlon(n), xlat(n), ylon(n), ylat(n))
387  IF(larea) CALL rot_equid_cylind_grid_area(clatr, fill, area(n))
388  ELSE
389  rlon(n)=fill
390  rlat(n)=fill
391  ENDIF
392  ENDDO
393  !$OMP END PARALLEL DO
394  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
395  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
396  ELSEIF(iopt.EQ.-1) THEN
397  !$OMP PARALLEL DO PRIVATE(N,HS,CLON,SLAT,CLAT,SLATR,CLATR,CLONR,RLONR,RLATR) &
398  !$OMP& REDUCTION(+:NRET) SCHEDULE(STATIC)
399  DO n=1,npts
400  IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90) THEN
401  hs=sign(1._kd,mod(rlon(n)-rlon0+180+3600,360._kd)-180)
402  clon=cos((rlon(n)-rlon0)/dpr)
403  slat=sin(rlat(n)/dpr)
404  clat=cos(rlat(n)/dpr)
405  slatr=clat0*slat-slat0*clat*clon
406  IF(slatr.LE.-1) THEN
407  clatr=0._kd
408  rlonr=0.
409  rlatr=-90.
410  ELSEIF(slatr.GE.1) THEN
411  clatr=0._kd
412  rlonr=0.
413  rlatr=90.
414  ELSE
415  clatr=sqrt(1-slatr**2)
416  clonr=(clat0*clat*clon+slat0*slat)/clatr
417  clonr=min(max(clonr,-1._kd),1._kd)
418  rlonr=hs*dpr*acos(clonr)
419  rlatr=dpr*asin(slatr)
420  ENDIF
421  xpts(n)=real((rlonr-wbd)/dlons+1._kd)
422  ypts(n)=real((rlatr-sbd)/dlats+1._kd)
423  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
424  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
425  nret=nret+1
426  IF(lrot) CALL rot_equid_cylind_vect_rot(rlon(n), clatr, slatr, &
427  clat, slat, clon, crot(n), srot(n))
428  IF(lmap) CALL rot_equid_cylind_map_jacob(fill, rlon(n), clatr, &
429  clat, slat, clon, xlon(n), xlat(n), ylon(n), ylat(n))
430  IF(larea) CALL rot_equid_cylind_grid_area(clatr, fill, area(n))
431  ELSE
432  xpts(n)=fill
433  ypts(n)=fill
434  ENDIF
435  ELSE
436  xpts(n)=fill
437  ypts(n)=fill
438  ENDIF
439  ENDDO
440  !$OMP END PARALLEL DO
441  ENDIF
442  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
443  END SUBROUTINE gdswzd_rot_equid_cylind
444 
462  SUBROUTINE rot_equid_cylind_error(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
463  IMPLICIT NONE
464  !
465  INTEGER, INTENT(IN ) :: IOPT, NPTS
466  !
467  REAL, INTENT(IN ) :: FILL
468  REAL, INTENT( OUT) :: RLAT(NPTS),RLON(NPTS)
469  REAL, INTENT( OUT) :: XPTS(NPTS),YPTS(NPTS)
470  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
471  IF(iopt>=0) THEN
472  rlon=fill
473  rlat=fill
474  ENDIF
475  IF(iopt<=0) THEN
476  xpts=fill
477  ypts=fill
478  ENDIF
479  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
480  END SUBROUTINE rot_equid_cylind_error
481 
507  SUBROUTINE rot_equid_cylind_vect_rot(RLON, CLATR, SLATR, CLAT, SLAT, &
508  CLON, CROT, SROT)
509  IMPLICIT NONE
510 
511  REAL(KIND=kd), INTENT(IN ) :: clat, clatr, clon, slat, slatr
512  REAL , INTENT(IN ) :: RLON
513  REAL , INTENT( OUT) :: CROT, SROT
514 
515  REAL(KIND=kd) :: slon
516 
517  IF(irot.EQ.1) THEN
518  IF(clatr.LE.0._kd) THEN
519  crot=real(-sign(1._kd,slatr*slat0))
520  srot=0.
521  ELSE
522  slon=sin((rlon-rlon0)/dpr)
523  crot=real((clat0*clat+slat0*slat*clon)/clatr)
524  srot=real(slat0*slon/clatr)
525  ENDIF
526  ELSE
527  crot=1.
528  srot=0.
529  ENDIF
530 
531  END SUBROUTINE rot_equid_cylind_vect_rot
532 
558  SUBROUTINE rot_equid_cylind_map_jacob(FILL, RLON, CLATR, CLAT, &
559  SLAT, CLON, XLON, XLAT, YLON, YLAT)
560  IMPLICIT NONE
561 
562  REAL(KIND=kd), INTENT(IN ) :: clatr, clat, slat, clon
563  REAL , INTENT(IN ) :: FILL, RLON
564  REAL , INTENT( OUT) :: XLON, XLAT, YLON, YLAT
565 
566  REAL(KIND=kd) :: slon, term1, term2
567 
568  IF(clatr.LE.0._kd) THEN
569  xlon=fill
570  xlat=fill
571  ylon=fill
572  ylat=fill
573  ELSE
574  slon=sin((rlon-rlon0)/dpr)
575  term1=(clat0*clat+slat0*slat*clon)/clatr
576  term2=slat0*slon/clatr
577  xlon=real(term1*clat/(dlons*clatr))
578  xlat=real(-term2/(dlons*clatr))
579  ylon=real(term2*clat/dlats)
580  ylat=real(term1/dlats)
581  ENDIF
582 
583  END SUBROUTINE rot_equid_cylind_map_jacob
584 
603  SUBROUTINE rot_equid_cylind_grid_area(CLATR, FILL, AREA)
604  IMPLICIT NONE
605 
606  REAL(KIND=kd), INTENT(IN ) :: clatr
607  REAL, INTENT(IN ) :: FILL
608  REAL, INTENT( OUT) :: AREA
609 
610  IF(clatr.LE.0._kd) THEN
611  area=fill
612  ELSE
613  area=real(2._kd*(rerth**2)*clatr*(dlons/dpr)*sin(0.5_kd*dlats/dpr))
614  ENDIF
615 
616  END SUBROUTINE rot_equid_cylind_grid_area
617 
619 
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 pi
PI.
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
Rotated equidistant cylindrical GRIB decoder and grid coordinate transformations for Arakawa grids A ...
subroutine rot_equid_cylind_error(IOPT, FILL, RLAT, RLON, XPTS, YPTS, NPTS)
Error handler.
subroutine rot_equid_cylind_vect_rot(RLON, CLATR, SLATR, CLAT, SLAT, CLON, CROT, SROT)
Vector rotation fields for rotated equidistant cylindrical grids - non "e" stagger.
subroutine gdswzd_rot_equid_cylind(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for rotated equidistant cylindrical.
real(kind=kd) rlon0
Local copy of rlon0.
real(kind=kd) slat0
Local copy of slat0.
subroutine rot_equid_cylind_map_jacob(FILL, RLON, CLATR, CLAT, SLAT, CLON, XLON, XLAT, YLON, YLAT)
Map jacobians for rotated equidistant cylindrical grids - non "e" stagger.
real(kind=kd) dlons
Local copy of dlons.
real(kind=kd) dlats
Local copy of dlats.
subroutine init_grib1(self, g1_desc)
Initializes a Rotated equidistant cylindrical grid given a grib1_descriptor object.
subroutine rot_equid_cylind_grid_area(CLATR, FILL, AREA)
Grid box area for rotated equidistant cylindrical grids - non "e" stagger.
subroutine init_grib2(self, g2_desc)
Initializes a Rotated equidistant cylindrical grid given a grib2_descriptor object.
real(kind=kd) rerth
Radius of the Earth.
integer, parameter kd
Fortran kind for reals.
real(kind=kd) clat0
Local copy of clat0.
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