NCEPLIBS-ip 5.2.0
Loading...
Searching...
No Matches
spectral_interp_mod::polates4 Interface Reference

Public Member Functions

subroutine polates4_grib1 (ipopt, kgdsi, kgdso, mi, mo, km, ibi, gi, no, rlat, rlon, ibo, lo, go, iret)
 Interpolate scalar fields (spectral).
 
subroutine polates4_grib2 (ipopt, igdtnumi, igdtmpli, igdtleni, igdtnumo, igdtmplo, igdtleno, mi, mo, km, ibi, gi, no, rlat, rlon, ibo, lo, go, iret)
 Interpolate scalar fields (spectral).
 

Detailed Description

Definition at line 25 of file spectral_interp_mod.F90.

Member Function/Subroutine Documentation

◆ polates4_grib1()

subroutine spectral_interp_mod::polates4::polates4_grib1 ( integer, dimension(20), intent(in)  ipopt,
integer, dimension(200), intent(in)  kgdsi,
integer, dimension(200), intent(in)  kgdso,
integer, intent(in)  mi,
integer, intent(in)  mo,
integer, intent(in)  km,
integer, dimension(km), intent(in)  ibi,
real, dimension(mi,km), intent(in)  gi,
integer  no,
real, dimension(mo), intent(inout)  rlat,
real, dimension(mo), intent(inout)  rlon,
integer, dimension(km), intent(out)  ibo,
logical*1, dimension(mo,km), intent(out)  lo,
real, dimension(mo,km), intent(out)  go,
integer, intent(out)  iret 
)

Interpolate scalar fields (spectral).

This subprogram performs spectral interpolation from any grid to any grid for scalar fields. It requires that the input fields be uniformly global.

Options allow choices between triangular shape (ipopt(1)=0) and rhomboidal shape (ipopt(1)=1) which has no default; a second option is the truncation (ipopt(2)) which defaults to a sensible truncation for the input grid (if opt(2)=-1).

Note
If the output grid is not found in a special list, then the transform back to grid is not very fast. This special list contains global cylindrical grids, polar stereographic grids centered at the pole and mercator grids.

Only horizontal interpolation is performed. The grids are defined by their grid description sections (passed in integer form as decoded by subprogram w3fi63()).

The current code recognizes the following projections:

  • kgds(1) = 000 equidistant cylindrical
  • kgds(1) = 001 mercator cylindrical
  • kgds(1) = 003 lambert conformal conical
  • kgds(1) = 004 gaussian cylindrical (spectral native)
  • kgds(1) = 005 polar stereographic azimuthal
  • kgds(1) = 203 rotated equidistant cylindrical (e-stagger)
  • kgds(1) = 205 rotated equidistant cylindrical (b-stagger)

Where kgds could be either input kgdsi or output kgdso. As an added bonus the number of output grid points and their latitudes and longitudes are also returned. On the other hand, the output can be a set of station points if kgdso(1)<0, in which case the number of points and their latitudes and longitudes must be input. Output bitmaps will not be created.

Program History Log

Date Programmer Comments
96-04-10 Iredell Initial
2001-06-18 Iredell improve detection of special fast transform
2015-01-27 Gayno replace calls to gdswiz() with new merged version of gdswzd().
Parameters
[in]ipopt(20) interpolation options ipopt(1)=0 for triangular, ipopt(1)=1 for rhomboidal; ipopt(2) is truncation number (defaults to sensible if ipopt(2)=-1).
[in]kgdsi(200) input gds parameters as decoded by w3fi63
[in]kgdso(200) output gds parameters (kgdso(1)<0 implies random station points)
[in]miskip number between input grid fields if km>1 or dimension of input grid fields if km=1
[in]moskip number between output grid fields if km>1 or dimension of output grid fields if km=1
[in]kmnumber of fields to interpolate
[in]ibi(km) input bitmap flags (must be all 0)
[in]gi(mi,km) input fields to interpolate
[out]nonumber of output points (only if kgdso(1)<0)
[out]rlat(no) output latitudes in degrees (if kgdso(1)<0)
[out]rlon(no) output longitudes in degrees (if kgdso(1)<0)
[out]ibo(km) output bitmap flags
[out]lo(mo,km) output bitmaps (always output)
[out]go(mo,km) output fields interpolated
[out]iretreturn code
  • 0 successful interpolation
  • 2 unrecognized input grid or no grid overlap
  • 3 unrecognized output grid
  • 41 invalid nonglobal input grid
  • 42 invalid spectral method parameters
Author
Iredell
Date
96-04-10

Definition at line 560 of file spectral_interp_mod.F90.

References sptrun(), sptrung(), sptrunm(), and sptruns().

◆ polates4_grib2()

subroutine spectral_interp_mod::polates4::polates4_grib2 ( integer, dimension(20), intent(in)  ipopt,
integer, intent(in)  igdtnumi,
integer, dimension(igdtleni), intent(in)  igdtmpli,
integer, intent(in)  igdtleni,
integer, intent(in)  igdtnumo,
integer, dimension(igdtleno), intent(in)  igdtmplo,
integer, intent(in)  igdtleno,
integer, intent(in)  mi,
integer, intent(in)  mo,
integer, intent(in)  km,
integer, dimension(km), intent(in)  ibi,
real, dimension(mi,km), intent(in)  gi,
integer, intent(out)  no,
real, dimension(mo), intent(inout)  rlat,
real, dimension(mo), intent(inout)  rlon,
integer, dimension(km), intent(out)  ibo,
logical*1, dimension(mo,km), intent(out)  lo,
real, dimension(mo,km), intent(out)  go,
integer, intent(out)  iret 
)

Interpolate scalar fields (spectral).

This subprogram performs spectral interpolation from any grid to any grid for scalar fields. It requires that the input fields be uniformly global. Options allow choices between triangular shape (ipopt(1)=0) and rhomboidal shape (ipopt(1)=1) which has no default; a second option is the truncation (ipopt(2)) which defaults to a sensible truncation for the input grid (if opt(2)=-1).

Note
If the output grid is not found in a special list, then the transform back to grid is not very fast. This special list contains global cylindrical grids, polar stereographic grids centered at the pole and mercator grids.

Only horizontal interpolation is performed.

The code recognizes the following projections, where "igdtnumi/o" is the GRIB2 grid defintion template number for the input and onutput grids, respectively:

  • igdtnumi/o = 00 equidistant cylindrical
  • igdtnumo = 01 rotated equidistant cylindrical. "e" and non-"e" staggered
  • igdtnumo = 10 mercator cylindrical
  • igdtnumo = 20 polar stereographic azimuthal
  • igdtnumo = 30 lambert conformal conical
  • igdtnumi/o = 40 gaussian cylindrical

As an added bonus the number of output grid points and their latitudes and longitudes are also returned. On the other hand, the output can be a set of station points if igdtnumo < 0, in which case the number of points and their latitudes and longitudes must be input. Output bitmaps will not be created.

Program History Log

Date Programmer Comments
96-04-10 Iredell initial
2001-06-18 Iredell improve detection of special fast transform
2015-01-27 Gayno replace calls to gdswiz with new merged version of gdswzd.
2015-07-13 Gayno convert to grib 2. replace grib 1 kgds arrays with grib 2 grid definition template arrays.
Parameters
[in]ipopt(20) interpolation options; ipopt(1)=0 for triangular, ipopt(1)=1 for rhomboidal; ipopt(2) is truncation number (defaults to sensible if ipopt(2)=-1).
[in]igdtnumigrid definition template number - input grid. Corresponds to the gfldigdtnum component of the NCEPLIBS-g2 library gridmod data structure.
  • 00 - equidistant cylindrical
  • 01 - rotated equidistant cylindrical. "e" and non-"e" staggered
  • 10 - mercator cyclindrical
  • 20 - polar stereographic azimuthal
  • 30 - lambert conformal conical
  • 40 - gaussian equidistant cyclindrical
[in]igdtmpli(igdtleni) grid definition template array - input grid. corresponds to the gfldigdtmpl component of the ncep g2 library gridmod data structure: (section 3 info). See comments in routine ipolates() for complete definition.
[in]igdtleninumber of elements of the grid definition template array - input grid. corresponds to the gfldigdtlen component of the ncep g2 library gridmod data structure.
[in]igdtnumogrid definition template number - output grid. Corresponds to the gfldigdtnum component of the NCEPLIBS-g2 library gridmod data structure. igdtnumo<0 means interpolate to random station points. Otherwise, same definition as igdtnumi.
[in]igdtmplo(igdtleno) grid definition template array - output grid. Corresponds to the gfldigdtmpl component of the ncep g2 library gridmod data structure (section 3 info). See comments in routine ipolates() for complete definition.
[in]igdtlenonumber of elements of the grid definition template array - output grid. Corresponds to the gfldigdtlen component of the [NCEPLIBS-g2](https://github.com/NOAA-EMC/NCEPLIBS-g2) library gridmod data structure.
[in]miskip number between input grid fields if km>1 or dimension of input grid fields if km=1
[in]moskip number between output grid fields if km>1 or dimension of output grid fields if km=1
[out]kmnumber of fields to interpolate
[out]ibi(km) input bitmap flags (must be all 0)
[out]gi(mi,km) input fields to interpolate
[out]nonumber of output points (only if igdtnumo>=0)
[out]rlat(mo) output latitudes in degrees (if igdtnumo<0)
[out]rlon(mo) output longitudes in degrees (if igdtnumo<0)
[out]ibo(km) output bitmap flags
[out]lo(mo,km) output bitmaps (always output)
[out]go(mo,km) output fields interpolated
[out]iretreturn code
  • 0 successful interpolation
  • 2 unrecognized input grid or no grid overlap
  • 3 unrecognized output grid
  • 41 invalid nonglobal input grid
  • 42 invalid spectral method parameters

!

Author
Mark Iredell
Date
96-04-10

Definition at line 255 of file spectral_interp_mod.F90.

References earth_radius_mod::earth_radius(), sptrun(), sptrung(), sptrunm(), and sptruns().


The documentation for this interface was generated from the following file: