NCEPLIBS-ip  4.1.0
ip_lambert_conf_grid_mod.F90
Go to the documentation of this file.
1 
5 
16  use ip_grid_mod
18  use constants_mod
19  implicit none
20 
21  private
22  public :: ip_lambert_conf_grid
23 
24  type, extends(ip_grid) :: ip_lambert_conf_grid
25  real :: rlat1
26  real :: rlon1
27  real :: rlati1
28  real :: rlati2
29  real :: orient
30  real :: dxs
31  real :: dys
32  real :: h
33  integer :: irot
34  contains
36  procedure :: init_grib1
38  procedure :: init_grib2
41  procedure :: gdswzd => gdswzd_lambert_conf
42  end type ip_lambert_conf_grid
43 
44  INTEGER :: irot
45  REAL :: an
46  REAL :: dxs
47  REAL :: dys
48  REAL :: h
49  REAL :: rerth
50 
51 contains
52 
59  subroutine init_grib1(self, g1_desc)
60  class(ip_lambert_conf_grid), intent(inout) :: self
61  type(grib1_descriptor), intent(in) :: g1_desc
62 
63  real :: dx, dy, hi, hj
64  integer :: iproj, iscan, jscan
65 
66  associate(kgds => g1_desc%gds)
67  self%rerth = 6.3712e6
68  self%eccen_squared = 0.0
69 
70  self%IM=kgds(2)
71  self%JM=kgds(3)
72 
73  self%RLAT1=kgds(4)*1.e-3
74  self%RLON1=kgds(5)*1.e-3
75 
76  self%IROT=mod(kgds(6)/8,2)
77  self%ORIENT=kgds(7)*1.e-3
78 
79  dx=kgds(8)
80  dy=kgds(9)
81 
82  iproj=mod(kgds(10)/128,2)
83  iscan=mod(kgds(11)/128,2)
84  jscan=mod(kgds(11)/64,2)
85 
86  self%RLATI1=kgds(12)*1.e-3
87  self%RLATI2=kgds(13)*1.e-3
88  self%H=(-1.)**iproj
89 
90  hi=(-1.)**iscan
91  hj=(-1.)**(1-jscan)
92  self%DXS=dx*hi
93  self%DYS=dy*hj
94 
95  self%iwrap = 0
96  self%jwrap1 = 0
97  self%jwrap2 = 0
98  self%nscan = mod(kgds(11) / 32, 2)
99  self%nscan_field_pos = self%nscan
100  self%kscan = 0
101  end associate
102 
103  end subroutine init_grib1
104 
111  subroutine init_grib2(self, g2_desc)
112  class(ip_lambert_conf_grid), intent(inout) :: self
113  type(grib2_descriptor), intent(in) :: g2_desc
114 
115  real :: dx, dy, hi, hj
116  integer :: iproj, iscan, jscan
117 
118 
119  associate(igdtmpl => g2_desc%gdt_tmpl, igdtlen => g2_desc%gdt_len)
120  call earth_radius(igdtmpl, igdtlen, self%rerth, self%eccen_squared)
121 
122  self%IM=igdtmpl(8)
123  self%JM=igdtmpl(9)
124 
125  self%RLAT1=float(igdtmpl(10))*1.0e-6
126  self%RLON1=float(igdtmpl(11))*1.0e-6
127 
128  self%IROT=mod(igdtmpl(12)/8,2)
129  self%ORIENT=float(igdtmpl(14))*1.0e-6
130 
131  dx=float(igdtmpl(15))*1.0e-3
132  dy=float(igdtmpl(16))*1.0e-3
133 
134  iproj=mod(igdtmpl(17)/128,2)
135  iscan=mod(igdtmpl(18)/128,2)
136  jscan=mod(igdtmpl(18)/64,2)
137 
138  self%RLATI1=float(igdtmpl(19))*1.0e-6
139  self%RLATI2=float(igdtmpl(20))*1.0e-6
140 
141  self%H=(-1.)**iproj
142  hi=(-1.)**iscan
143  hj=(-1.)**(1-jscan)
144  self%DXS=dx*hi
145  self%DYS=dy*hj
146 
147  self%nscan = mod(igdtmpl(18) / 32, 2)
148  self%nscan_field_pos = self%nscan
149  self%iwrap = 0
150  self%jwrap1 = 0
151  self%jwrap2 = 0
152  self%kscan = 0
153  end associate
154  end subroutine init_grib2
155 
218  SUBROUTINE gdswzd_lambert_conf(self,IOPT,NPTS,FILL, &
219  XPTS,YPTS,RLON,RLAT,NRET, &
220  CROT,SROT,XLON,XLAT,YLON,YLAT,AREA)
221  IMPLICIT NONE
222  !
223  class(ip_lambert_conf_grid), intent(in) :: self
224  INTEGER, INTENT(IN ) :: IOPT, NPTS
225  INTEGER, INTENT( OUT) :: NRET
226  !
227  REAL, INTENT(IN ) :: FILL
228  REAL, INTENT(INOUT) :: RLON(NPTS),RLAT(NPTS)
229  REAL, INTENT(INOUT) :: XPTS(NPTS),YPTS(NPTS)
230  REAL, OPTIONAL, INTENT( OUT) :: CROT(NPTS),SROT(NPTS)
231  REAL, OPTIONAL, INTENT( OUT) :: XLON(NPTS),XLAT(NPTS)
232  REAL, OPTIONAL, INTENT( OUT) :: YLON(NPTS),YLAT(NPTS),AREA(NPTS)
233  !
234  INTEGER :: IM, JM, N
235  !
236  LOGICAL :: LROT, LMAP, LAREA
237  !
238  REAL :: ANTR, DI, DJ
239  REAL :: DLON1
240  REAL :: DE, DE2, DR2
241  REAL :: ORIENT, RLAT1, RLON1
242  REAL :: RLATI1, RLATI2
243  REAL :: XMAX, XMIN, YMAX, YMIN, XP, YP
244  REAL :: DLON, DR
245  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
246  IF(PRESENT(crot)) crot=fill
247  IF(PRESENT(srot)) srot=fill
248  IF(PRESENT(xlon)) xlon=fill
249  IF(PRESENT(xlat)) xlat=fill
250  IF(PRESENT(ylon)) ylon=fill
251  IF(PRESENT(ylat)) ylat=fill
252  IF(PRESENT(area)) area=fill
253 
254  im=self%im
255  jm=self%jm
256 
257  rlat1=self%rlat1
258  rlon1=self%rlon1
259 
260  irot=self%irot
261  orient=self%orient
262 
263  rlati1=self%rlati1
264  rlati2=self%rlati2
265 
266  h=self%h
267  dxs=self%dxs
268  dys=self%dys
269 
270  rerth = self%rerth
271 
272  IF(rlati1.EQ.rlati2) THEN
273  an=sin(rlati1/dpr)
274  ELSE
275  an=log(cos(rlati1/dpr)/cos(rlati2/dpr))/ &
276  log(tan((90-rlati1)/2/dpr)/tan((90-rlati2)/2/dpr))
277  ENDIF
278  de=rerth*cos(rlati1/dpr)*tan((rlati1+90)/2/dpr)**an/an
279  IF(h*rlat1.EQ.90) THEN
280  xp=1
281  yp=1
282  ELSE
283  dr=de/tan((rlat1+90)/2/dpr)**an
284  dlon1=mod(rlon1-orient+180+3600,360.)-180
285  xp=1-sin(an*dlon1/dpr)*dr/dxs
286  yp=1+cos(an*dlon1/dpr)*dr/dys
287  ENDIF
288  antr=1/(2*an)
289  de2=de**2
290  xmin=0
291  xmax=im+1
292  ymin=0
293  ymax=jm+1
294  nret=0
295  IF(PRESENT(crot).AND.PRESENT(srot))THEN
296  lrot=.true.
297  ELSE
298  lrot=.false.
299  ENDIF
300  IF(PRESENT(xlon).AND.PRESENT(xlat).AND.PRESENT(ylon).AND.PRESENT(ylat))THEN
301  lmap=.true.
302  ELSE
303  lmap=.false.
304  ENDIF
305  IF(PRESENT(area))THEN
306  larea=.true.
307  ELSE
308  larea=.false.
309  ENDIF
310  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311  ! TRANSLATE GRID COORDINATES TO EARTH COORDINATES
312  IF(iopt.EQ.0.OR.iopt.EQ.1) THEN
313  !$OMP PARALLEL DO PRIVATE(N,DI,DJ,DR2,DR,DLON) REDUCTION(+:NRET) SCHEDULE(STATIC)
314  DO n=1,npts
315  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
316  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
317  di=h*(xpts(n)-xp)*dxs
318  dj=h*(ypts(n)-yp)*dys
319  dr2=di**2+dj**2
320  dr=sqrt(dr2)
321  IF(dr2.LT.de2*1.e-6) THEN
322  rlon(n)=0.
323  rlat(n)=h*90.
324  ELSE
325  rlon(n)=mod(orient+1./an*dpr*atan2(di,-dj)+3600,360.)
326  rlat(n)=(2*dpr*atan((de2/dr2)**antr)-90)
327  ENDIF
328  nret=nret+1
329  dlon=mod(rlon(n)-orient+180+3600,360.)-180
330  IF(lrot) CALL lambert_conf_vect_rot(dlon,crot(n),srot(n))
331  IF(lmap) CALL lambert_conf_map_jacob(rlat(n),fill, dlon, dr, &
332  xlon(n),xlat(n),ylon(n),ylat(n))
333  IF(larea) CALL lambert_conf_grid_area(rlat(n),fill,dr,area(n))
334  ELSE
335  rlon(n)=fill
336  rlat(n)=fill
337  ENDIF
338  ENDDO
339  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
340  ! TRANSLATE EARTH COORDINATES TO GRID COORDINATES
341  ELSEIF(iopt.EQ.-1) THEN
342  !$OMP PARALLEL DO PRIVATE(N,DR,DLON) REDUCTION(+:NRET) SCHEDULE(STATIC)
343  DO n=1,npts
344  IF(abs(rlon(n)).LE.360.AND.abs(rlat(n)).LE.90.AND. &
345  h*rlat(n).NE.-90) THEN
346  dr=h*de*tan((90-rlat(n))/2/dpr)**an
347  dlon=mod(rlon(n)-orient+180+3600,360.)-180
348  xpts(n)=xp+h*sin(an*dlon/dpr)*dr/dxs
349  ypts(n)=yp-h*cos(an*dlon/dpr)*dr/dys
350  IF(xpts(n).GE.xmin.AND.xpts(n).LE.xmax.AND. &
351  ypts(n).GE.ymin.AND.ypts(n).LE.ymax) THEN
352  nret=nret+1
353  IF(lrot) CALL lambert_conf_vect_rot(dlon,crot(n),srot(n))
354  IF(lmap) CALL lambert_conf_map_jacob(rlat(n),fill,dlon,dr, &
355  xlon(n),xlat(n),ylon(n),ylat(n))
356  IF(larea) CALL lambert_conf_grid_area(rlat(n),fill,dr,area(n))
357  ELSE
358  xpts(n)=fill
359  ypts(n)=fill
360  ENDIF
361  ELSE
362  xpts(n)=fill
363  ypts(n)=fill
364  ENDIF
365  ENDDO
366  !$OMP END PARALLEL DO
367  ENDIF
368  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
369  END SUBROUTINE gdswzd_lambert_conf
370 
389  SUBROUTINE lambert_conf_vect_rot(DLON,CROT,SROT)
390  IMPLICIT NONE
391  REAL, INTENT( IN) :: DLON
392  REAL, INTENT( OUT) :: CROT, SROT
393 
394  IF(irot.EQ.1) THEN
395  crot=cos(an*dlon/dpr)
396  srot=sin(an*dlon/dpr)
397  ELSE
398  crot=1.
399  srot=0.
400  ENDIF
401 
402  END SUBROUTINE lambert_conf_vect_rot
403 
426  SUBROUTINE lambert_conf_map_jacob(RLAT,FILL,DLON,DR,XLON,XLAT,YLON,YLAT)
427  IMPLICIT NONE
428 
429  REAL, INTENT(IN ) :: RLAT, FILL, DLON, DR
430  REAL, INTENT( OUT) :: XLON, XLAT, YLON, YLAT
431 
432  REAL :: CLAT
433 
434  clat=cos(rlat/dpr)
435  IF(clat.LE.0.OR.dr.LE.0) THEN
436  xlon=fill
437  xlat=fill
438  ylon=fill
439  ylat=fill
440  ELSE
441  xlon=h*cos(an*dlon/dpr)*an/dpr*dr/dxs
442  xlat=-h*sin(an*dlon/dpr)*an/dpr*dr/dxs/clat
443  ylon=h*sin(an*dlon/dpr)*an/dpr*dr/dys
444  ylat=h*cos(an*dlon/dpr)*an/dpr*dr/dys/clat
445  ENDIF
446 
447  END SUBROUTINE lambert_conf_map_jacob
448 
467  SUBROUTINE lambert_conf_grid_area(RLAT,FILL,DR,AREA)
468  IMPLICIT NONE
469 
470  REAL, INTENT(IN ) :: RLAT
471  REAL, INTENT(IN ) :: FILL
472  REAL, INTENT(IN ) :: DR
473  REAL, INTENT( OUT) :: AREA
474 
475  REAL :: CLAT
476 
477  clat=cos(rlat/dpr)
478  IF(clat.LE.0.OR.dr.LE.0) THEN
479  area=fill
480  ELSE
481  area=rerth**2*clat**2*abs(dxs)*abs(dys)/(an*dr)**2
482  ENDIF
483 
484  END SUBROUTINE lambert_conf_grid_area
485 
486 end module ip_lambert_conf_grid_mod
487 
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.
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
Lambert conformal grib decoder and grid coordinate transformations.
subroutine lambert_conf_grid_area(RLAT, FILL, DR, AREA)
Grid box area for lambert conformal conical.
real dxs
x-direction grid length adjusted for scan mode.
subroutine init_grib1(self, g1_desc)
Initializes a Lambert Conformal grid given a grib1_descriptor object.
subroutine lambert_conf_vect_rot(DLON, CROT, SROT)
Vector rotation fields for lambert conformal conical.
subroutine lambert_conf_map_jacob(RLAT, FILL, DLON, DR, XLON, XLAT, YLON, YLAT)
Map jacobians for lambert conformal conical.
subroutine gdswzd_lambert_conf(self, IOPT, NPTS, FILL, XPTS, YPTS, RLON, RLAT, NRET, CROT, SROT, XLON, XLAT, YLON, YLAT, AREA)
GDS wizard for lambert conformal conical.
integer irot
vector rotation flag.
real rerth
Radius of the earth.
real dys
y-direction grid length adjusted for scan model.
subroutine init_grib2(self, g2_desc)
Initializes a Lambert Conformal grid given a grib2_descriptor object.
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