NCEPLIBS-ip  4.4.0
ipolatev_mod::ipolatev Interface Reference

Private Member Functions

subroutine ipolatev_grib1 (ip, ipopt, kgdsi, kgdso, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
 This subprogram interpolates vector field from any grid to any grid given a grib1 Grid Descriptor Section. More...
 
subroutine ipolatev_grib1_single_field (ip, ipopt, kgdsi, kgdso, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
 Special case of ipolatev_grib1 when interpolating a single field. More...
 
subroutine ipolatev_grib2 (ip, ipopt, igdtnumi, igdtmpli, igdtleni, igdtnumo, igdtmplo, igdtleno, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
 This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor. More...
 
subroutine ipolatev_grib2_single_field (ip, ipopt, igdtnumi, igdtmpli, igdtleni, igdtnumo, igdtmplo, igdtleno, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
 This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor. More...
 

Detailed Description

Definition at line 21 of file ipolatev.F90.

Member Function/Subroutine Documentation

◆ ipolatev_grib1()

subroutine ipolatev_mod::ipolatev::ipolatev_grib1 ( integer(c_int), intent(in)  ip,
integer(c_int), dimension(20), intent(in)  ipopt,
integer(c_int), dimension(200), intent(inout)  kgdsi,
integer(c_int), dimension(200), intent(inout)  kgdso,
integer(c_int), intent(in)  mi,
integer(c_int), intent(in)  mo,
integer(c_int), intent(in)  km,
integer(c_int), dimension(km), intent(in)  ibi,
logical(c_bool), dimension(mi,km), intent(in)  li,
real(c_double), dimension(mi,km), intent(in)  ui,
real(c_double), dimension(mi,km), intent(in)  vi,
integer(c_int), intent(out)  no,
real(c_double), dimension(mo), intent(inout)  rlat,
real(c_double), dimension(mo), intent(inout)  rlon,
real(c_double), dimension(mo), intent(inout)  crot,
real(c_double), dimension(mo), intent(inout)  srot,
integer(c_int), dimension(km), intent(out)  ibo,
logical(c_bool), dimension(mo,km), intent(out)  lo,
real(c_double), dimension(mo,km), intent(out)  uo,
real(c_double), dimension(mo,km), intent(out)  vo,
integer(c_int), intent(out)  iret 
)
private

This subprogram interpolates vector field from any grid to any grid given a grib1 Grid Descriptor Section.

Only horizontal interpolation is performed. The following interpolation methods are possible:

  • (ip=0) bilinear
  • (ip=1) bicubic
  • (ip=2) neighbor
  • (ip=3) budget
  • (ip=4) spectral
  • (ip=6) neighbor-budget

Some of these methods have interpolation options and/or restrictions on the input or output grids, both of which are documented more fully in their respective subprograms.

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
  • (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. for the budget approach, a subsection of the grid may be output by subtracting kgdso(1) from 255 and passing in the latitudes and longitudes of the points. Input bitmaps will be interpolated to output bitmaps.

Output bitmaps will also be created when the output grid extends outside of the domain of the input grid. the output field is set to 0 where the output bitmap is off.

Parameters
ipInterpolation method
  • ip = BILINEAR_INTERP_ID = 0 for bilinear
  • ip = BICUBIC_INTERP_ID = 1 for bicubic
  • ip = NEIGHBOR_INTERP_ID = 2 for neighbor;
  • ip = BUDGET_INTERP_ID = 3 for budget;
  • ip = SPECTRAL_INTERP_ID = 4 for spectral;
  • ip = NEIGHBOR_BUDGET_INTERP_ID = 6 for neighbor-budget
ipoptInterpolation options
  • ip=0 (bilinear): (No options)
  • ip=1 Cbicubic): constraint option
  • ip=2 (neighbor): (No options)
  • ip=3 (budget): Number in radius, radius weights, search radius
  • ip=4 (spectral): Spectral shape, spectral truncation
  • ip=6 (neighbor-budget): Number in radius, radius weights ...)
[in]kgdsiInput gds parameters as decoded by w3fi63.
[in]kgdsoOutput gds parameters.
[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]ibiInput bitmap flags.
[in]liInput bitmaps (if respective ibi(k)=1).
[in]uiInput u-component fields to interpolate.
[in]viInput v-component fields to interpolate.
[out]noNumber of output points (only if kgdso(1)<0).
[out]rlatOutput latitudes in degrees (if kgdso(1)<0).
[out]rlonOutput longitudes in degrees (if kgdso(1)<0).
[in,out]crotVector rotation cosines (if igdtnumo>=0).
[in,out]srotVector rotation sines (if igdtnumo>=0).
[out]iboOutput bitmap flags.
[out]loOutput bitmaps (always output).
[out]uoOutput u-component fields interpolated.
[out]voOutput v-component fields interpolated.
[out]iretReturn code.
  • 0 Successful interpolation.
  • 1 Unrecognized interpolation method.
  • 2 Unrecognized input grid or no grid overlap.
  • 3 Unrecognized output grid.
  • 1x Invalid bicubic method parameters.
  • 3x Invalid budget method parameters.
  • 4x Invalid spectral method parameters.
Note
Examples demonstrating relative cpu costs. This example is interpolating 12 levels of winds from the 360 x 181 global grid (ncep grid 3) to the 93 x 68 hawaiian mercator grid (ncep grid 204).

The example times are for the c90. As a reference, the cp time for unpacking the global 12 temperature fields is 0.07 seconds.

METHOD IP IPOPT CP SECONDS
BILINEAR 0 0.05
BICUBIC 1 0 0.16
BICUBIC 1 1 0.17
NEIGHBOR 2 0.02
BUDGET 3 -1,-1 0.94
SPECTRAL 4 0,40 0.31
SPECTRAL 4 1,40 0.33
SPECTRAL 4 0,-1 0.59
N-BUDGET 6 -1,-1 0.31

The spectral interpolation is fast for the mercator grid. However, for some grids the spectral interpolation is slow.

The following example is interpolating 12 levels of temperatures from the 360 x 181 global grid (ncep grid 3) to the 93 x 65 conus lambert conformal grid (ncep grid 211).

METHOD IP IPOPT CP SECONDS
BILINEAR 0 0.05
BICUBIC 1 0 0.15
BICUBIC 1 1 0.16
NEIGHBOR 2 0.02
BUDGET 3 -1,-1 0.92
SPECTRAL 4 0,40 4.51
SPECTRAL 4 1,40 5.77
SPECTRAL 4 0,-1 12.60
N-BUDGET 6 -1,-1 0.33
Date
July 2021
Author
Kyle Gerheiser

Definition at line 565 of file ipolatev.F90.

References ipolatev_mod::ipolatev_grid().

◆ ipolatev_grib1_single_field()

subroutine ipolatev_mod::ipolatev::ipolatev_grib1_single_field ( integer(c_int), intent(in)  ip,
integer(c_int), dimension(20), intent(in)  ipopt,
integer(c_int), dimension(200), intent(inout)  kgdsi,
integer(c_int), dimension(200), intent(inout)  kgdso,
integer(c_int), intent(in)  mi,
integer(c_int), intent(in)  mo,
integer(c_int), intent(in)  km,
integer(c_int), intent(in)  ibi,
logical(c_bool), dimension(mi), intent(in)  li,
real(c_double), dimension(mi), intent(in)  ui,
real(c_double), dimension(mi), intent(in)  vi,
integer(c_int), intent(out)  no,
real(c_double), dimension(mo), intent(inout)  rlat,
real(c_double), dimension(mo), intent(inout)  rlon,
real(c_double), dimension(mo), intent(inout)  crot,
real(c_double), dimension(mo), intent(inout)  srot,
integer(c_int), intent(out)  ibo,
logical(c_bool), dimension(mo), intent(out)  lo,
real(c_double), dimension(mo), intent(out)  uo,
real(c_double), dimension(mo), intent(out)  vo,
integer(c_int), intent(out)  iret 
)
private

Special case of ipolatev_grib1 when interpolating a single field.

Removes the km dimension of input arrays so vectors can be passed to ibi/ibo.

Parameters
ipInterpolation method
  • ip = BILINEAR_INTERP_ID = 0 for bilinear
  • ip = BICUBIC_INTERP_ID = 1 for bicubic
  • ip = NEIGHBOR_INTERP_ID = 2 for neighbor;
  • ip = BUDGET_INTERP_ID = 3 for budget;
  • ip = SPECTRAL_INTERP_ID = 4 for spectral;
  • ip = NEIGHBOR_BUDGET_INTERP_ID = 6 for neighbor-budget
ipoptInterpolation options
  • ip=0 (bilinear): (No options)
  • ip=1 Cbicubic): constraint option
  • ip=2 (neighbor): (No options)
  • ip=3 (budget): Number in radius, radius weights, search radius
  • ip=4 (spectral): Spectral shape, spectral truncation
  • ip=6 (neighbor-budget): Number in radius, radius weights ...)
[in]kgdsiInput gds parameters as decoded by w3fi63.
[in]kgdsoOutput gds parameters.
[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]ibiInput bitmap flags.
[in]liInput bitmaps (if respective ibi(k)=1).
[in]uiInput u-component fields to interpolate.
[in]viInput v-component fields to interpolate.
[out]noNumber of output points (only if kgdso(1)<0).
[out]rlatOutput latitudes in degrees (if kgdso(1)<0).
[out]rlonOutput longitudes in degrees (if kgdso(1)<0).
[in,out]crotVector rotation cosines (if igdtnumo>=0).
[in,out]srotVector rotation sines (if igdtnumo>=0).
[out]iboOutput bitmap flags.
[out]loOutput bitmaps (always output).
[out]uoOutput u-component fields interpolated.
[out]voOutput v-component fields interpolated.
[out]iretReturn code.
  • 0 Successful interpolation.
  • 1 Unrecognized interpolation method.
  • 2 Unrecognized input grid or no grid overlap.
  • 3 Unrecognized output grid.
  • 1x Invalid bicubic method parameters.
  • 3x Invalid budget method parameters.
  • 4x Invalid spectral method parameters.
Date
Jan 2022
Author
Kyle Gerheiser

Definition at line 680 of file ipolatev.F90.

References ipolatev_mod::ipolatev_grid().

◆ ipolatev_grib2()

subroutine ipolatev_mod::ipolatev::ipolatev_grib2 ( integer(c_int), intent(in)  ip,
integer(c_int), dimension(20), intent(in)  ipopt,
integer(c_int), intent(in)  igdtnumi,
integer(c_int), dimension(igdtleni), intent(in)  igdtmpli,
integer(c_int), intent(in)  igdtleni,
integer(c_int), intent(in)  igdtnumo,
integer(c_int), dimension(igdtleno), intent(in)  igdtmplo,
integer(c_int), intent(in)  igdtleno,
integer(c_int), intent(in)  mi,
integer(c_int), intent(in)  mo,
integer(c_int), intent(in)  km,
integer(c_int), dimension(km), intent(in)  ibi,
logical(c_bool), dimension(mi,km), intent(in)  li,
real(c_double), dimension(mi,km), intent(in)  ui,
real(c_double), dimension(mi,km), intent(in)  vi,
integer(c_int), intent(out)  no,
real(c_double), dimension(mo), intent(inout)  rlat,
real(c_double), dimension(mo), intent(inout)  rlon,
real(c_double), dimension(mo), intent(inout)  crot,
real(c_double), dimension(mo), intent(inout)  srot,
integer(c_int), dimension(km), intent(out)  ibo,
logical(c_bool), dimension(mo,km), intent(out)  lo,
real(c_double), dimension(mo,km), intent(out)  uo,
real(c_double), dimension(mo,km), intent(out)  vo,
integer(c_int), intent(out)  iret 
)
private

This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor.

This is a wrapper for ipolates_grid which converts a grib1 descriptor into an ip_grid_descriptor, which is used to create an ip_grid. Only horizontal interpolation is performed.

The following interpolation methods are possible:

  • (ip=0) bilinear
  • (ip=1) bicubic
  • (ip=2) neighbor
  • (ip=3) budget
  • (ip=4) spectral
  • (ip=6) neighbor-budget

Some of these methods have interpolation options and/or restrictions on the input or output grids, both of which are documented more fully in their respective subprograms.

Input and output grids are defined by their grib 2 grid definition template as decoded by the ncep g2 library. The current code recognizes the following projections, where "igdtnumi/o" is the grib 2 grid defintion template number for the input and output grids, respectively:

  • (igdtnumi/o=00) equidistant cylindrical
  • (igdtnumi/o=01) rotated equidistant cylindrical. "e" and non-"e" staggered
  • (igdtnumi/o=10) mercator cylindrical
  • (igdtnumi/o=20) polar stereographic azimuthal
  • (igdtnumi/o=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, data may be interpolated to a set of station points if "igdtnumo"<0 (or subtracted from 255 for the budget option), in which case the number of points and their latitudes and longitudes must be input.

Input bitmaps will be interpolated to output bitmaps. Output bitmaps will also be created when the output grid extends outside of the domain of the input grid.

The output field is set to 0 where the output bitmap is off.

Parameters
[in]ipInterpolation method
  • ip=0 for bilinear
  • ip=1 for bicubic
  • ip=2 for neighbor;
  • ip=3 for budget;
  • ip=4 for spectral;
  • ip=6 for neighbor-budget
[in]ipoptInterpolation options
  • ip=0: (No options)
  • ip=1: Constraint option
  • ip=2: (No options)
  • ip=3: Number in radius, radius weights, search radius
  • ip=4: Spectral shape, spectral truncation
  • ip=6: Number in radius, radius weights ...)
[in]igdtnumiGrid definition template number for the input grid. Corresponds to the gfldigdtnum component of the ncep 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]igdtmpliGrid definition template array input grid. Corresponds to the gfldigdtmpl component of the NCEPLIBS-g2 gridmod data structure

Section 3 Info:

All map projections:

  1. Shape of earth, octet 15.
  2. Scale factor of spherical earth radius, octet 16.
  3. Scaled value of radius of spherical earth, octets 17-20.
  4. Scale factor of major axis of elliptical earth, octet 21.
  5. Scaled value of major axis of elliptical earth, octets 22-25.
  6. Scale factor of minor axis of elliptical earth, octet 26. 7: Scaled value of minor axis of elliptical earth, octets 27-30.

Equidistant Cyclindrical:

  1. Number of points along a parallel, octs 31-34.
  2. Number of points along a meridian, octs 35-38.
  3. Basic angle of initial production domain, octets 39-42.
  4. Subdivisions of basic angle, octets 43-46.
  5. Latitude of first grid point, octets 47-50.
  6. Longitude of first grid point, octets 51-54.
  7. Resolution and component flags, octet 55.
  8. Latitude of last grid point, octets 56-59.
  9. Longitude of last grid point, octets 60-63.
  10. i-direction increment, octets 64-67.
  11. j-direction increment, octets 68-71.
  12. Scanning mode, octet 72.

Mercator Cyclindrical:

  1. Number of points along a parallel, octs 31-34.
  2. Number of points along a meridian, octs 35-38.
  3. Latitude of first point, octets 39-42.
  4. Longitude of first point, octets 43-46.
  5. Resolution and component flags, octet 47.
  6. Tangent latitude, octets 48-51.
  7. Latitude of last point, octets 52-55.
  8. Longitude of last point, octets 56-59.
  9. Scanning mode flags, octet 60.
  10. Orientation of grid, octets 61-64.
  11. Longitudinal grid length, octets 65-68.
  12. Latitudinal grid length, octets 69-72.

Lambert Conformal Conical:

  1. Number of points along x-axis, octs 31-34.
  2. Number of points along y-axis, octs 35-38.
  3. Latitude of first point, octets 39-42.
  4. Longitude of first point, octets 43-46.
  5. Resolution of component flag, octet 47.
  6. Latitude where grid lengths specified,octets 48-51.
  7. Longitude of meridian that is parallel to y-axis, octets 52-55.
  8. x-direction grid length, octets 56-59.
  9. y-direction grid length, octets 60-63.
  10. Projection center flag, octet 64.
  11. Scanning mode, octet 65.
  12. First tangent latitude from pole, octets 66-69.
  13. Second tangent latitude from pole, octets 70-73.
  14. Latitude of south pole of projection, octets 74-77.
  15. Longitude of south pole of projection, octets 78-81.

Gaussian Cylindrical:

  1. Number of points along a parallel, octs 31-34.
  2. Number of points along a meridian, octs 35-38.
  3. Basic angle of initial production domain, octets 39-42.
  4. Subdivisions of basic angle, octets 43-46.
  5. Latitude of first grid point, octets 47-50.
  6. Longitude of first grid point, octets 51-54.
  7. Resolution and component flags, octet 55.
  8. Latitude of last grid point, octets 56-59.
  9. Longitude of last grid point, octets 60-63.
  10. i-direction increment, octets 64-67.
  11. Number of parallels between pole and equator, octets 68-71.
  12. Scanning mode, octet 72.

Polar Stereographic Azimuthal:

  1. Number of points along x-axis, octets 31-34.
  2. Number of points along y-axis, octets 35-38.
  3. Latitude of first grid point, octets 39-42.
  4. Longitude of first grid point, octets 43-46.
  5. Resolution and component flags, octet 47.
  6. True latitude, octets 48-51.
  7. Orientation longitude, octets 52-55.
  8. x-direction grid length, octets 56-59.
  9. y-direction grid length, octets 60-63.
  10. Projection center flag, octet 64.
  11. Scanning mode flags, octet 65.

Rotated Equidistant Cyclindrical:

  1. Number of points along a parallel, octs 31-34.
  2. Number of points along a meridian, octs 35-38.
  3. Basic angle of initial production domain, octets 39-42.
  4. Subdivisions of basic angle, octets 43-46.
  5. Latitude of first grid point, octets 47-50.
  6. Longitude of first grid point, octets 51-54.
  7. Resolution and component flags, octet 55.
  8. Latitude of last grid point, octets 56-59.
  9. Longitude of last grid point, octets 60-63.
  10. i-direction increment, octets 64-67.
  11. j-direction increment, octets 68-71.
  12. Scanning mode, octet 72.
  13. Latitude of southern pole of projection, octets 73-76.
  14. Longitude of southern pole of projection, octets 77-80.
  15. Angle of rotation of projection, octs 81-84.
Parameters
[in]igdtleniNumber of elements of the grid definition template array for the input grid. Corresponds to the gfldigdtlen component of the ncep g2 library gridmod data structure.
[in]igdtnumoGrid definition template number for the output grid. Corresponds to the gfldigdtnum component of the ncep g2 library gridmod data structure. See "igdtnumi" for specific template definitions. Note: igdtnumo<0 means interpolate to random station points.
[in]igdtmploGrid definition template array for the output grid. Corresponds to the gfldigdtmpl component of the ncep g2 library gridmod data structure. See "igdtmpli" for definition of array elements.
[in]igdtlenoNumber of elements of the grid definition template array for the output grid. Corresponds to the gfldigdtlen component of the ncep 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.
[in]kmNumber of fields to interpolate.
[in]ibiInput bitmap flags.
[in]liInput bitmaps (if respective ibi(k)=1).
[in]uiInput u-component fields to interpolate.
[in]viInput v-component fields to interpolate.
[out]noNumber of output points (only if kgdso(1)<0).
[in,out]rlatOutput latitudes in degrees (if kgdso(1)<0).
[in,out]rlonOutput longitudes in degrees (if kgdso(1)<0).
[in,out]crotVector rotation cosines (if igdtnumo>=0).
[in,out]srotVector rotation sines (if igdtnumo>=0).
[out]iboOutput bitmap flags.
[out]loOutput bitmaps (always output).
[out]uoOutput u-component fields interpolated.
[out]voOutput v-component fields interpolated.
[out]iretReturn code.
  • 0 Successful interpolation.
  • 1 Unrecognized interpolation method.
  • 2 Unrecognized input grid or no grid overlap.
  • 3 Unrecognized output grid.
  • 1x Invalid bicubic method parameters.
  • 3x Invalid budget method parameters.
  • 4x Invalid spectral method parameters.
Note
Examples demonstrating relative cpu costs. This example is interpolating 12 levels of winds from the 360 x 181 global grid (ncep grid 3) to the 93 x 68 hawaiian mercator grid (ncep grid 204).

The example times are for the c90. As a reference, the cp time for unpacking the global 12 wind fields is 0.07 seconds.

METHOD IP IPOPT CP SECONDS
BILINEAR 0 0.05
BICUBIC 1 0 0.16
BICUBIC 1 1 0.17
NEIGHBOR 2 0.02
BUDGET 3 -1,-1 0.94
SPECTRAL 4 0,40 0.31
SPECTRAL 4 1,40 0.33
SPECTRAL 4 0,-1 0.59
N-BUDGET 6 -1,-1 0.31

The spectral interpolation is fast for the mercator grid. However, for some grids the spectral interpolation is slow.

The following example is interpolating 12 levels of winds from the 360 x 181 global grid (ncep grid 3) to the 93 x 65 conus lambert conformal grid (ncep grid 211).

METHOD IP IPOPT CP SECONDS
BILINEAR 0 0.05
BICUBIC 1 0 0.15
BICUBIC 1 1 0.16
NEIGHBOR 2 0.02
BUDGET 3 -1,-1 0.92
SPECTRAL 4 0,40 4.51
SPECTRAL 4 1,40 5.77
SPECTRAL 4 0,-1 12.60
N-BUDGET 6 -1,-1 0.33
Author
Kyle Gerheiser
Date
July 2021

Definition at line 382 of file ipolatev.F90.

References ipolatev_mod::ipolatev_grid().

◆ ipolatev_grib2_single_field()

subroutine ipolatev_mod::ipolatev::ipolatev_grib2_single_field ( integer(c_int), intent(in)  ip,
integer(c_int), dimension(20), intent(in)  ipopt,
integer(c_int), intent(in)  igdtnumi,
integer(c_int), dimension(igdtleni), intent(in)  igdtmpli,
integer(c_int), intent(in)  igdtleni,
integer(c_int), intent(in)  igdtnumo,
integer(c_int), dimension(igdtleno), intent(in)  igdtmplo,
integer(c_int), intent(in)  igdtleno,
integer(c_int), intent(in)  mi,
integer(c_int), intent(in)  mo,
integer(c_int), intent(in)  km,
integer(c_int), intent(in)  ibi,
logical(c_bool), dimension(mi), intent(in)  li,
real(c_double), dimension(mi), intent(in)  ui,
real(c_double), dimension(mi), intent(in)  vi,
integer(c_int), intent(out)  no,
real(c_double), dimension(mo), intent(inout)  rlat,
real(c_double), dimension(mo), intent(inout)  rlon,
real(c_double), dimension(mo), intent(inout)  crot,
real(c_double), dimension(mo), intent(inout)  srot,
integer(c_int), intent(out)  ibo,
logical(c_bool), dimension(mo), intent(out)  lo,
real(c_double), dimension(mo), intent(out)  uo,
real(c_double), dimension(mo), intent(out)  vo,
integer(c_int), intent(out)  iret 
)
private

This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor.

Parameters
[in]ipInterpolation method
  • ip=0 for bilinear
  • ip=1 for bicubic
  • ip=2 for neighbor;
  • ip=3 for budget;
  • ip=4 for spectral;
  • ip=6 for neighbor-budget
[in]ipoptInterpolation options
  • ip=0: (No options)
  • ip=1: Constraint option
  • ip=2: (No options)
  • ip=3: Number in radius, radius weights, search radius
  • ip=4: Spectral shape, spectral truncation
  • ip=6: Number in radius, radius weights ...)
[in]igdtnumiGrid definition template number for the input grid. Corresponds to the gfldigdtnum component of the ncep 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]igdtmpliGrid definition template array input grid. Corresponds to the gfldigdtmpl component of the NCEPLIBS-g2 gridmod data structure
[in]igdtleniNumber of elements of the grid definition template array for the input grid. Corresponds to the gfldigdtlen component of the ncep g2 library gridmod data structure.
[in]igdtnumoGrid definition template number for the output grid. Corresponds to the gfldigdtnum component of the ncep g2 library gridmod data structure. See "igdtnumi" for specific template definitions. Note: igdtnumo<0 means interpolate to random station points.
[in]igdtmploGrid definition template array for the output grid. Corresponds to the gfldigdtmpl component of the ncep g2 library gridmod data structure. See "igdtmpli" for definition of array elements.
[in]igdtlenoNumber of elements of the grid definition template array for the output grid. Corresponds to the gfldigdtlen component of the ncep 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.
[in]kmNumber of fields to interpolate.
[in]ibiInput bitmap flags.
[in]liInput bitmaps (if respective ibi(k)=1).
[in]uiInput u-component fields to interpolate.
[in]viInput v-component fields to interpolate.
[out]noNumber of output points (only if kgdso(1)<0).
[in,out]rlatOutput latitudes in degrees (if kgdso(1)<0).
[in,out]rlonOutput longitudes in degrees (if kgdso(1)<0).
[in,out]crotVector rotation cosines (if igdtnumo>=0).
[in,out]srotVector rotation sines (if igdtnumo>=0).
[out]iboOutput bitmap flags.
[out]loOutput bitmaps (always output).
[out]uoOutput u-component fields interpolated.
[out]voOutput v-component fields interpolated.
[out]iretReturn code.
  • 0 Successful interpolation.
  • 1 Unrecognized interpolation method.
  • 2 Unrecognized input grid or no grid overlap.
  • 3 Unrecognized output grid.
  • 1x Invalid bicubic method parameters.
  • 3x Invalid budget method parameters.
  • 4x Invalid spectral method parameters.
Author
Eric Engle
Date
November 2022

Definition at line 832 of file ipolatev.F90.

References ipolatev_mod::ipolatev_grid().


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