NCEPLIBS-ip  5.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 ip_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 
140  subroutine init_grib2(self, g2_desc)
141  class(ip_rot_equid_cylind_grid), intent(inout) :: self
142  type(grib2_descriptor), intent(in) :: g2_desc
143 
144  if (ncep_post_arakawa.and.(g2_desc%gdt_num.eq.32769)) then
145  call init_grib2_ncep_post(self, g2_desc)
146  else
147  call init_grib2_default(self, g2_desc)
148  endif
149 
150  end subroutine init_grib2
151 
160  subroutine init_grib2_default(self, g2_desc)
161  class(ip_rot_equid_cylind_grid), intent(inout) :: self
162  type(grib2_descriptor), intent(in) :: g2_desc
163 
164  real(kd) :: rlat1, rlon1, rlat0, rlat2, rlon2, nbd, ebd
165  integer :: iscale
166  integer :: i_offset_odd, i_offset_even, j_offset
167 
168  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
169 
170  CALL earth_radius(igdtmpl,igdtlen,self%rerth,self%eccen_squared)
171 
172  i_offset_odd=mod(igdtmpl(19)/8,2)
173  i_offset_even=mod(igdtmpl(19)/4,2)
174  j_offset=mod(igdtmpl(19)/2,2)
175 
176  iscale=igdtmpl(10)*igdtmpl(11)
177  IF(iscale==0) iscale=10**6
178 
179  rlat1=float(igdtmpl(12))/float(iscale)
180  rlon1=float(igdtmpl(13))/float(iscale)
181  rlat0=float(igdtmpl(20))/float(iscale)
182  rlat0=rlat0+90.0_kd
183 
184  self%RLON0=float(igdtmpl(21))/float(iscale)
185 
186  rlat2=float(igdtmpl(15))/float(iscale)
187  rlon2=float(igdtmpl(16))/float(iscale)
188 
189  self%IROT=mod(igdtmpl(14)/8,2)
190  self%IM=igdtmpl(8)
191  self%JM=igdtmpl(9)
192 
193  self%SLAT0=sin(rlat0/dpr)
194  self%CLAT0=cos(rlat0/dpr)
195 
196  self%WBD=rlon1
197  IF (self%WBD > 180.0) self%WBD = self%WBD - 360.0
198  self%SBD=rlat1
199 
200  nbd=rlat2
201  ebd=rlon2
202 
203  self%DLATS=(nbd-self%SBD)/float(self%JM-1)
204  self%DLONS=(ebd-self%WBD)/float(self%IM-1)
205 
206  IF(i_offset_odd==1) self%WBD=self%WBD+(0.5_kd*self%DLONS)
207  IF(j_offset==1) self%SBD=self%SBD+(0.5_kd*self%DLATS)
208 
209  self%iwrap = 0
210  self%jwrap1 = 0
211  self%jwrap2 = 0
212  self%kscan = 0
213  self%nscan = mod(igdtmpl(19) / 32, 2)
214  self%nscan_field_pos = self%nscan
215  end associate
216  end subroutine init_grib2_default
217 
227  subroutine init_grib2_ncep_post(self, g2_desc)
228  class(ip_rot_equid_cylind_grid), intent(inout) :: self
229  type(grib2_descriptor), intent(in) :: g2_desc
230 
231  real(kd) :: rlat1, rlon1, rlat0, rlat2, rlon2, nbd, ebd
232  integer :: iscale
233  integer :: i_offset_odd, i_offset_even, j_offset
234  real(kd) :: hs, hs2, slat1, slat2, slatr, clon1, clon2, clat1, clat2, clatr, clonr, rlonr, rlatr
235 
236  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
237 
238  CALL earth_radius(igdtmpl,igdtlen,self%rerth,self%eccen_squared)
239 
240  i_offset_odd=mod(igdtmpl(19)/8,2)
241  i_offset_even=mod(igdtmpl(19)/4,2)
242  j_offset=mod(igdtmpl(19)/2,2)
243 
244  iscale=igdtmpl(10)*igdtmpl(11)
245  IF(iscale==0) iscale=10**6
246 
247  rlat1=float(igdtmpl(12))/float(iscale)
248  rlon1=float(igdtmpl(13))/float(iscale)
249  rlat0=float(igdtmpl(15))/float(iscale)
250 
251  self%RLON0=float(igdtmpl(16))/float(iscale)
252 
253  rlat2=float(igdtmpl(20))/float(iscale)
254  rlon2=float(igdtmpl(21))/float(iscale)
255 
256  self%IROT=mod(igdtmpl(14)/8,2)
257  self%IM=igdtmpl(8)
258  self%JM=igdtmpl(9)
259 
260  slat1=sin(rlat1/dpr)
261  clat1=cos(rlat1/dpr)
262 
263  self%SLAT0=sin(rlat0/dpr)
264  self%CLAT0=cos(rlat0/dpr)
265 
266  hs=sign(1._kd,mod(rlon1-self%RLON0+180+3600,360._kd)-180)
267 
268  clon1=cos((rlon1-self%RLON0)/dpr)
269  slatr=self%CLAT0*slat1-self%SLAT0*clat1*clon1
270  clatr=sqrt(1-slatr**2)
271  clonr=(self%CLAT0*clat1*clon1+self%SLAT0*slat1)/clatr
272  rlatr=dpr*asin(slatr)
273  rlonr=hs*dpr*acos(clonr)
274 
275  self%WBD=rlonr
276  IF (self%WBD > 180.0) self%WBD = self%WBD - 360.0
277  self%SBD=rlatr
278 
279  slat2=sin(rlat2/dpr)
280  clat2=cos(rlat2/dpr)
281  hs2=sign(1._kd,mod(rlon2-self%RLON0+180+3600,360._kd)-180)
282  clon2=cos((rlon2-self%RLON0)/dpr)
283  slatr=self%CLAT0*slat2-self%SLAT0*clat2*clon2
284  clatr=sqrt(1-slatr**2)
285  clonr=(self%CLAT0*clat2*clon2+self%SLAT0*slat2)/clatr
286  nbd=dpr*asin(slatr)
287  ebd=hs2*dpr*acos(clonr)
288  self%DLATS=(nbd-self%SBD)/float(self%JM-1)
289  self%DLONS=(ebd-self%WBD)/float(self%IM-1)
290 
291  IF(i_offset_odd==1) self%WBD=self%WBD+(0.5_kd*self%DLONS)
292  IF(j_offset==1) self%SBD=self%SBD+(0.5_kd*self%DLATS)
293 
294  self%iwrap = 0
295  self%jwrap1 = 0
296  self%jwrap2 = 0
297  self%kscan = 0
298  self%nscan = mod(igdtmpl(19) / 32, 2)
299  self%nscan_field_pos = self%nscan
300  end associate
301  end subroutine init_grib2_ncep_post
302 
362  SUBROUTINE gdswzd_rot_equid_cylind(self,IOPT,NPTS, &
363  FILL,XPTS,YPTS,RLON,RLAT,NRET, &
364  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
365  IMPLICIT NONE
366 
367  class(ip_rot_equid_cylind_grid), intent(in) :: self
368  INTEGER, INTENT(IN ) :: IOPT, NPTS
369  INTEGER, INTENT( OUT) :: NRET
370  !
371  REAL, INTENT(IN ) :: FILL
372  REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
373  REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
374  REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
375  REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
376  REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
377  !
378  INTEGER :: IM,JM,N
379  !
380  LOGICAL :: LROT, LMAP, LAREA
381  !
382  REAL(KIND=kd) :: hs
383  REAL(KIND=kd) :: clonr,clatr,slatr
384  REAL(KIND=kd) :: clat,slat,clon
385  REAL(KIND=kd) :: rlatr,rlonr
386  REAL(KIND=kd) :: wbd,sbd
387  REAL :: XMIN,XMAX,YMIN,YMAX
388  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
389  IF(PRESENT(crot)) crot=fill
390  IF(PRESENT(srot)) srot=fill
391  IF(PRESENT(xlon)) xlon=fill
392  IF(PRESENT(xlat)) xlat=fill
393  IF(PRESENT(ylon)) ylon=fill
394  IF(PRESENT(ylat)) ylat=fill
395  IF(PRESENT(area)) area=fill
396  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
397  ! IS THE EARTH RADIUS DEFINED?
398  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
399  ! IS THIS AN "E"-STAGGER GRID? ROUTINE CAN'T PROCESS THOSE.
400  ! I_OFFSET_ODD=MOD(IGDTMPL(19)/8,2)
401  ! I_OFFSET_EVEN=MOD(IGDTMPL(19)/4,2)
402  ! J_OFFSET=MOD(IGDTMPL(19)/2,2)
403  ! IF(I_OFFSET_ODD/=I_OFFSET_EVEN) THEN
404  ! CALL ROT_EQUID_CYLIND_ERROR(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
405  ! RETURN
406  ! ENDIF
407  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
408 
409 
410  rlon0=self%rlon0
411  irot=self%irot
412 
413  im=self%im
414  jm=self%jm
415 
416  slat0=self%slat0
417  clat0=self%clat0
418 
419  wbd=self%wbd
420  sbd=self%sbd
421 
422  dlats=self%dlats
423  dlons=self%dlons
424 
425  xmin=0
426  xmax=im+1
427  ymin=0
428  ymax=jm+1
429  nret=0
430 
431  rerth = self%rerth
432  IF(rerth<0.)THEN
433  CALL rot_equid_cylind_error(iopt,fill,rlat,rlon,xpts,ypts,npts)
434  RETURN
435  ENDIF
436 
437  IF(PRESENT(crot).AND.PRESENT(srot))THEN
438  lrot=.true.
439  ELSE
440  lrot=.false.
441  ENDIF
442  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
443  lmap=.true.
444  ELSE
445  lmap=.false.
446  ENDIF
447  IF(PRESENT(area))THEN
448  larea=.true.
449  ELSE
450  larea=.false.
451  ENDIF
452  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
453  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
454  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
455  !$OMP PARALLEL DO PRIVATE(N,RLONR,RLATR,HS,CLONR,SLATR,CLATR,SLAT,CLAT,CLON) &
456  !$OMP& REDUCTION(+:NRET) SCHEDULE(STATIC)
457  DO n=1,npts
458  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
459  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
460  rlonr=wbd+(xpts(n)-1._kd)*dlons
461  rlatr=sbd+(ypts(n)-1._kd)*dlats
462  IF(rlonr <= 0._kd) THEN
463  hs=-1.0_kd
464  ELSE
465  hs=1.0_kd
466  ENDIF
467  clonr=cos(rlonr/dpr)
468  slatr=sin(rlatr/dpr)
469  clatr=cos(rlatr/dpr)
470  slat=clat0*slatr+slat0*clatr*clonr
471  IF(slat.LE.-1) THEN
472  clat=0.
473  clon=cos(rlon0/dpr)
474  rlon(n)=0.
475  rlat(n)=-90.
476  ELSEIF(slat.GE.1) THEN
477  clat=0.
478  clon=cos(rlon0/dpr)
479  rlon(n)=0.
480  rlat(n)=90.
481  ELSE
482  clat=sqrt(1-slat**2)
483  clon=(clat0*clatr*clonr-slat0*slatr)/clat
484  clon=min(max(clon,-1._kd),1._kd)
485  rlon(n)=real(mod(rlon0+hs*dpr*acos(clon)+3600,360._kd))
486  rlat(n)=real(dpr*asin(slat))
487  ENDIF
488  nret=nret+1
489  IF(lrot) CALL rot_equid_cylind_vect_rot(rlon(n), clatr, slatr, &
490  clat, slat, clon, crot(n), srot(n))
491  IF(lmap) CALL rot_equid_cylind_map_jacob(fill, rlon(n), clatr, &
492  clat, slat, clon, xlon(n), xlat(n), ylon(n), ylat(n))
493  IF(larea) CALL rot_equid_cylind_grid_area(clatr, fill, area(n))
494  ELSE
495  rlon(n)=fill
496  rlat(n)=fill
497  ENDIF
498  ENDDO
499  !$OMP END PARALLEL DO
500  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
501  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
502  ELSEIF(iopt.EQ.-1) THEN
503  !$OMP PARALLEL DO PRIVATE(N,HS,CLON,SLAT,CLAT,SLATR,CLATR,CLONR,RLONR,RLATR) &
504  !$OMP& REDUCTION(+:NRET) SCHEDULE(STATIC)
505  DO n=1,npts
506  IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90) THEN
507  hs=sign(1._kd,mod(rlon(n)-rlon0+180+3600,360._kd)-180)
508  clon=cos((rlon(n)-rlon0)/dpr)
509  slat=sin(rlat(n)/dpr)
510  clat=cos(rlat(n)/dpr)
511  slatr=clat0*slat-slat0*clat*clon
512  IF(slatr.LE.-1) THEN
513  clatr=0._kd
514  rlonr=0.
515  rlatr=-90.
516  ELSEIF(slatr.GE.1) THEN
517  clatr=0._kd
518  rlonr=0.
519  rlatr=90.
520  ELSE
521  clatr=sqrt(1-slatr**2)
522  clonr=(clat0*clat*clon+slat0*slat)/clatr
523  clonr=min(max(clonr,-1._kd),1._kd)
524  rlonr=hs*dpr*acos(clonr)
525  rlatr=dpr*asin(slatr)
526  ENDIF
527  xpts(n)=real((rlonr-wbd)/dlons+1._kd)
528  ypts(n)=real((rlatr-sbd)/dlats+1._kd)
529  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
530  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
531  nret=nret+1
532  IF(lrot) CALL rot_equid_cylind_vect_rot(rlon(n), clatr, slatr, &
533  clat, slat, clon, crot(n), srot(n))
534  IF(lmap) CALL rot_equid_cylind_map_jacob(fill, rlon(n), clatr, &
535  clat, slat, clon, xlon(n), xlat(n), ylon(n), ylat(n))
536  IF(larea) CALL rot_equid_cylind_grid_area(clatr, fill, area(n))
537  ELSE
538  xpts(n)=fill
539  ypts(n)=fill
540  ENDIF
541  ELSE
542  xpts(n)=fill
543  ypts(n)=fill
544  ENDIF
545  ENDDO
546  !$OMP END PARALLEL DO
547  ENDIF
548  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
549  END SUBROUTINE gdswzd_rot_equid_cylind
550 
568  SUBROUTINE rot_equid_cylind_error(IOPT,FILL,RLAT,RLON,XPTS,YPTS,NPTS)
569  IMPLICIT NONE
570  !
571  INTEGER, INTENT(IN ) :: IOPT, NPTS
572  !
573  REAL, INTENT(IN ) :: FILL
574  REAL, INTENT( OUT) :: RLAT(NPTS),RLON(NPTS)
575  REAL, INTENT( OUT) :: XPTS(NPTS),YPTS(NPTS)
576  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
577  IF(iopt>=0) THEN
578  rlon=fill
579  rlat=fill
580  ENDIF
581  IF(iopt<=0) THEN
582  xpts=fill
583  ypts=fill
584  ENDIF
585  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
586  END SUBROUTINE rot_equid_cylind_error
587 
613  SUBROUTINE rot_equid_cylind_vect_rot(RLON, CLATR, SLATR, CLAT, SLAT, &
614  CLON, CROT, SROT)
615  IMPLICIT NONE
616 
617  REAL(KIND=kd), INTENT(IN ) :: clat, clatr, clon, slat, slatr
618  REAL , INTENT(IN ) :: RLON
619  REAL , INTENT( OUT) :: CROT, SROT
620 
621  REAL(KIND=kd) :: slon
622 
623  IF(irot.EQ.1) THEN
624  IF(clatr.LE.0._kd) THEN
625  crot=real(-sign(1._kd,slatr*slat0))
626  srot=0.
627  ELSE
628  slon=sin((rlon-rlon0)/dpr)
629  crot=real((clat0*clat+slat0*slat*clon)/clatr)
630  srot=real(slat0*slon/clatr)
631  ENDIF
632  ELSE
633  crot=1.
634  srot=0.
635  ENDIF
636 
637  END SUBROUTINE rot_equid_cylind_vect_rot
638 
664  SUBROUTINE rot_equid_cylind_map_jacob(FILL, RLON, CLATR, CLAT, &
665  SLAT, CLON, XLON, XLAT, YLON, YLAT)
666  IMPLICIT NONE
667 
668  REAL(KIND=kd), INTENT(IN ) :: clatr, clat, slat, clon
669  REAL , INTENT(IN ) :: FILL, RLON
670  REAL , INTENT( OUT) :: XLON, XLAT, YLON, YLAT
671 
672  REAL(KIND=kd) :: slon, term1, term2
673 
674  IF(clatr.LE.0._kd) THEN
675  xlon=fill
676  xlat=fill
677  ylon=fill
678  ylat=fill
679  ELSE
680  slon=sin((rlon-rlon0)/dpr)
681  term1=(clat0*clat+slat0*slat*clon)/clatr
682  term2=slat0*slon/clatr
683  xlon=real(term1*clat/(dlons*clatr))
684  xlat=real(-term2/(dlons*clatr))
685  ylon=real(term2*clat/dlats)
686  ylat=real(term1/dlats)
687  ENDIF
688 
689  END SUBROUTINE rot_equid_cylind_map_jacob
690 
709  SUBROUTINE rot_equid_cylind_grid_area(CLATR, FILL, AREA)
710  IMPLICIT NONE
711 
712  REAL(KIND=kd), INTENT(IN ) :: clatr
713  REAL, INTENT(IN ) :: FILL
714  REAL, INTENT( OUT) :: AREA
715 
716  IF(clatr.LE.0._kd) THEN
717  area=fill
718  ELSE
719  area=real(2._kd*(rerth**2)*clatr*(dlons/dpr)*sin(0.5_kd*dlats/dpr))
720  ENDIF
721 
722  END SUBROUTINE rot_equid_cylind_grid_area
723 
725 
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
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 init_grib2_ncep_post(self, g2_desc)
Initializes a Rotated equidistant cylindrical grid given a grib2_descriptor object.
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.
subroutine init_grib2_default(self, g2_desc)
Initializes a Rotated equidistant cylindrical grid given a grib2_descriptor object.
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:58