Detailed description of the procedures for cssgrid


CSSGRID - interpolation on a sphere

CSSGRID is called to interpolate data randomly spaced on a sphere.

CSSGRID uses tension splines to construct an interpolatory surface. As a first step, CSSGRID creates a Delaunay triangulation of the input data points; the companion entry CSSTRI is provided only for those times when you want to calculate the triangulation for independent inspection.

Various aspects of the interpolation algorithm can be changed by using the parmater access routines CSSETI and CSSETR to change values for the cssgrid control parematers.

------------------------------------------------------------------
             Argument | Type                |  Mode  | Dimension
------------------------------------------------------------------
CALL CSSGRID N,       | Integer             | Input  |
             RLAT,    | Real                | Input  | N
             RLON,    | Real                | Input  | N
             F,       | Real                | Input  | N
             NI,      | Integer             | Input  | 
             NJ,      | Integer             | Input  | 
             PLAT,    | Real                | Input  | NI
             PLON,    | Real                | Input  | NJ
             FF,      | Real                | Output | NI x NJ
             IWK,     | Integer             | Input  | 27*N
             RWK,     | Double precision    | Input  | 13*N
             IER)     | Integer             | Output | 
------------------------------------------------------------------
N
The number of input data points (N > 2).
RLAT,RLON
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
F
Array containing data values. F(I) is a functional value at (RLAT(I),RLON(I)) for I = 1 to N.
NI
The number of rows in the uniform output grid. NI can be 1.
NJ
The number of columns in the uniform output grid. NJ can be 1.
PLAT,PLON
Arrays of length NI and NJ, respectively, containing the latitudes and longitudes of the grid lines. The values for PLAT and PLON should be in degrees.
FF
An NI by NJ array containing the desired interpolated values. FF(I,J) is the interpolated value at the coordinate specified by PLAT(I) and PLON(J) for I = 1 to NI and J = 1 to NJ.
IWK
An integer workspace of length 27*N.
RWK
A DOUBLE PRECISION workspace of length at least 13*N. A common source of error in using CSSGRID is failing to type RWK as double precision.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

CSSTRI - calculates a Delaunay triangulation

CSSTRI is called to find a Delaunay triangulation of data randomly positioned on the surface of a sphere.
-------------------------------------------------------------------
             Argument | Type                |  Mode  | Dimension
-------------------------------------------------------------------
CALL CSSTRI (N,       | Integer             | Input  |
             RLAT,    | Real                | Input  | N
             RLON,    | Real                | Input  | N
             NT,      | Integer             | Output | 
             NTRI,    | Integer             | Output | 3 x NT where NT=2*N 
             IWK,     | Integer             | Input  | 27*N
             RWK,     | Double precision    | Input  | 13*N
             IER)     | Integer             | Output | 
-------------------------------------------------------------------
N
The number of input data points (N > 2).
RLAT,RLON
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
NT
The number of triangles in the triangulation, unless IER .NE. 0, in which case NT = 0. Let NB be the number of boundary points on the convex hull of the data. If NB .GE. 3, then NT = 2N-NB-2, otherwise NT=2N-4.
NTRI
A two-dimensional integer array dimensioned for 3 x NT where NT is the number of triangles in the triangulation (NT is at most 2*N). NTRI contains the triangulation data. The vertices of the Kth triangle are:
         (PLAT(NTRI((1,K)),PLON(NTRI(1,K))
         (PLAT(NTRI((2,K)),PLON(NTRI(2,K))
         (PLAT(NTRI((3,K)),PLON(NTRI(3,K))
         
IWK
An integer workspace of length 27*N.
RWK
A DOUBLE PRECISION workspace of length at least 13*N. A common source of error in using CSSTRI is failing to type RWK as double precision.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

CSS2C - convert from lat/lon coordinates to Cartesian coordinates.

CSS2C is called to find the equivalent Cartesian coordinates on a unit sphere to specified latitude and longitude coordinates. The coordinate of 0. latitude and 0. longitude is converted to Cartesian coordinate (1.,0.,0.). Latitudes and longitudes are assumed to be in degrees.
------------------------------------------------------------------
            Argument | Type    |  Mode  | Dimension
------------------------------------------------------------------
CALL CSS2C (N,       | Integer | Input  |
            RLAT,    | Real    | Input  | N
            RLON,    | Real    | Input  | N
            X,       | Real    | Output | N
            Y,       | Real    | Output | N
            Z)       | Real    | Output | N
------------------------------------------------------------------
N
The number of input lat/lon coordinates.
RLAT
An array containing the latitudes of the input coordinates, expressed in degrees.
RLON
An array containing the longitudes of the input coordinates, expressed in degrees.
X,Y,Z
Arrays containing the Cartesian coordinates of the output points. (X(I),Y(I),Z(I)) is the Cartesian coordinate corresponding to the lat/lon coordinate (RLAT(I),RLON(I)) for I=1 to N. X(I)**2 + Y(I)**2 + Z(I)**2 = 1 for I = 1 to N.

CSC2S - convert from Cartesian coordinates on a unit sphere to lat/lon coordinates

CSC2S is called to find an equivalent latitude and longitude coordinates on a sphere to a specified Cartesian coordinate on the unit sphere. The coordinate (1.,0.,0.) is mapped to the latitude/longitude coordinate (0.,0.). The latitude/longitude coordinates are returned in degrees, latitudes between -90. and 90. (inclusive) and longitudes between -180. and 180. (inclusive).
------------------------------------------------------------------
            Argument | Type    |  Mode  | Dimension
------------------------------------------------------------------
CALL CSC2S (N,       | Integer | Input  | 
            X,       | Real    | Input  |  N
            Y,       | Real    | Input  |  N
            Z,       | Real    | Input  |  N
            RLAT,    | Real    | Output |  N
            RLON)    | Real    | Output |  N
------------------------------------------------------------------
N
The number of input (X,Y,Z) coordinates.
X,Y,Z
The Cartesian coordinates of the input points, X(I)**2 + Y(I)**2 + Z(I)**2 = 1 for I = 1 to N.
RLAT
The latitudes of the output coordinates, in degrees.
RLON
The longitudes of the output coordinates, in degrees.

CSVORO - calculate Voronoi polygons.

CSVORO is called if you want to determine the Voronoi polygons for data randomly positioned on a sphere. Each call to CSVORO calculates the vertices for the Voronoi polygon surrounding a specified input point.
------------------------------------------------------------------
             Argument | Type              |  Mode  | Dimension
------------------------------------------------------------------
CALL CSVORO (NPTS,    | Integer           | Input  |
             RLATI,   | Real              | Input  | N
             RLONI,   | Real              | Input  | N
             NI,      | Integer           | Input  | 
             NF,      | Integer           | Input  |
             IWK,     | Integer           | Input  | 27*NPTS
             RWK,     | Double precision  | Input  | 9*NPTS
             NC,      | Integer           | Input  |
             RLATO,   | Real              | Input  | NC
             RLONO,   | Real              | Input  | NC
             RC,      | Real              | Output | NC
             NCA,     | Integer           | Output |
             NUMV,    | Integer           | Output |
             NV,      | Integer           | Output | NPTS
             IER)     | Integer           | Output | 
------------------------------------------------------------------
NPTS
The number of input data points (NPTS > 3).
RLATI
An array, dimensioned for NPTS, containing the latitudes of the input coordinates, in degrees.
RLONI
An array, dimensioned for NPTS, containing the longitudes of the input coordinates, in degrees.
NI
The index of the input coordinate for which you want to determine the Voronoi polygon (1 .LE. NI .LE. NPTS).
NF
Flag indicating if this is the first call to CSVORO to retrieve Voronoi polygons for this dataset (1=yes, 0=no). Calls subsequent to the first call for a given dataset are much faster than the first call.
IWK
Integer work space dimensioned for 27*NPTS.
RWK
DOUBLE PRECISION work space dimensioned for 9*NPTS.
NC
The maximum size of the output arrays RLATO, RLONO, and RC. NC should be 2*NPTS.
RLATO
An array of latitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. RLATO is returned in degrees.
RLONO
An array of longitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. RLONO is returned in degrees.
RC
Array containing circumradii (arc lengths in degrees of the angle between a circumcenter and its associated triangle vertices).
NCA
The actual number of circumcenters returned in RLATO and RLONO. This number may be larger than NPTS if the input dataset has boundary points since certain "pseudo points" may have been added to the original dataset in order to complete the Voronoi polygon set.
NUMV
The number of vertices in the Voronoi polygon enclosing the coordinate (RLATI(NI),RLONI(NI)).
NV
An array (dimensioned for NPTS) containing NUMV indices for the Voronoi polygon enclosing the coordinate (RLATI(NI),RLONI(NI)). The indices returned in this array refer to the coordinates returned in RLATO, RLONO, and RC. For example, if the integer "J" is an element of the NV array, then (RLATO(J),RLONO(J)) is a vertex of the Voronoi polygon enclosing (RLATI(NI),RLONI(NI)). The indices in NV list out the vertices of the Voronoi polygon in counter-clockwise order.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

CSSETI - set INTEGER parameter values

CSSETI is used to set values for any of the cssgrid control parameters that take integer values. The values set by CSSETI remain in effect until changed by subsequent calls to CSSETI.
-----------------------------------------------------
            Argument | Type      |  Mode  | Dimension
-----------------------------------------------------
CALL CSSETI (PNAM,   | Character | Input  |
             IVAL)   | Integer   | Input  | 
-----------------------------------------------------
PNAM
The name of the control parameter to be set.
IVAL
The value to be assigned to the parameter.

CSGETI - retrieve values for INTEGER parameters

CSGETI is a called to obtain current values for any of the INTEGER valued cssgrid control parameters.
-----------------------------------------------------
            Argument | Type      |  Mode  | Dimension
-----------------------------------------------------
CALL CSGETI (PNAM,   | Character | Input  |
             IVAL)   | Integer   | Output | 
-----------------------------------------------------
PNAM
The name of the control parameter whose value is to be retrieved.
IVAL
The current value assigned to the control parameter.

CSSETR - set REAL parameter values

CSSETR is used to set values for any of the cssgrid control parameters that take REAL values. The values set by CSSETR remain in effect until changed by subsequent calls to CSSETR.
-----------------------------------------------------
            Argument | Type      |  Mode  | Dimension
-----------------------------------------------------
CALL CSSETR (PNAM,   | Character | Input  |
             RVAL)   | Real      | Input  | 
-----------------------------------------------------
PNAM
The name of the control parameter to be set.
RVAL
The value to be assigned to the parameter.

CSGETR - retrieve values for REAL parameters

CSGETR is a called to obtain current values for any of the REAL valued cssgrid control parameters.
-----------------------------------------------------
            Argument | Type      |  Mode  | Dimension
-----------------------------------------------------
CALL CSGETR (PNAM,   | Character | Input  |
             RVAL)   | Real      | Output | 
-----------------------------------------------------
PNAM
The name of the control parameter whose value is to be retrieved.
RVAL
The current value assigned to the control parameter.

CSSETD - set DOUBLE PRECISION parameter values

CSSETD is used to set values for any of the cssgrid control parameters that take DOUBLE PRECISION values. The values set by CSSETD remain in effect until changed by subsequent calls to CSSETD.
---------------------------------------------------------------
            Argument | Type               |  Mode  | Dimension
---------------------------------------------------------------
CALL CSSETD (PNAM,   | Character          | Input  |
             DVAL)   | Double precision   | Input  | 
---------------------------------------------------------------
PNAM
The name of the control parameter to be set.
DVAL
The value to be assigned to the parameter.

CSGETD - retrieve values for DOUBLE PRECISION parameters

CSGETD is a called to obtain current values for any of the REAL valued cssgird control parameters.
------------------------------------------------------------------
            Argument | Type                  |  Mode  | Dimension
------------------------------------------------------------------
CALL CSGETD (PNAM,   | Character             | Input  |
             DVAL)   | Double precision      | Output | 
------------------------------------------------------------------
PNAM
The name of the control parameter whose value is to be retrieved.
DVAL
The current value assigned to the control parameter.

CSSGRIDD - interpolation on a sphere

CSSGRIDD is called to interpolate data randomly spaced on a sphere. CSSGRIDD is a double precision version of CSSGRID.

CSSGRIDD uses tension splines to construct an interpolatory surface. As a first step, CSSGRIDD creates a Delaunay triangulation of the input data points; the companion entry CSSTRI is provided only for those times when you want to calculate the triangulation for independent inspection.

Various aspects of the interpolation algorithm can be changed by using the parmater access routines CSSETI and CSSETD to change values for the cssgridd control parematers.

------------------------------------------------------------------
              Argument | Type                |  Mode  | Dimension
------------------------------------------------------------------
CALL CSSGRIDD N,       | Integer             | Input  |
              RLAT,    | Double precision    | Input  | N
              RLON,    | Double precision    | Input  | N
              F,       | Double precision    | Input  | N
              NI,      | Integer             | Input  | 
              NJ,      | Integer             | Input  | 
              PLAT,    | Double precision    | Input  | NI
              PLON,    | Double precision    | Input  | NJ
              FF,      | Double precision    | Output | NI x NJ
              IWK,     | Integer             | Input  | 27*N
              RWK,     | Double precision    | Input  | 13*N
              IER)     | Integer             | Output | 
------------------------------------------------------------------
N
The number of input data points (N > 2).
RLAT,RLON
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
F
Array containing data values. F(I) is a functional value at (X(I),Y(I),Z(I)) for I = 1 to N.
NI
The number of rows in the uniform output grid. NI can be 1.
NJ
The number of columns in the uniform output grid. NJ can be 1.
PLAT,PLON
Arrays of length NI and NJ, respectively, containing the latitudes and longitudes of the grid lines. The values for PLAT and PLON should be in radians.
FF
An NI by NJ array containing the desired interpolated values. FF(I,J) is the interpolated value at the coordinate specified by PLAT(I) and PLON(J) for I = 1 to NI and J = 1 to NJ.
IWK
An integer workspace of length 27*N.
RWK
A DOUBLE PRECISION workspace of length at least 13*N.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

CSSTRID - calculates a Delaunay triangulation

CSSTRID is called to find a Delaunay triangulation of data randomly positioned on the surface of a sphere. CSSTRID is a double precision version of CSSTRI.
-------------------------------------------------------------------
              Argument | Type                |  Mode  | Dimension
-------------------------------------------------------------------
CALL CSSTRID (N,       | Integer             | Input  |
              RLAT,    | Double precision    | Input  | N
              RLON,    | Double precision    | Input  | N
              NT,      | Integer             | Output | 
              NTRI,    | Integer             | Output | 3 x NT where NT=2*N 
              IWK,     | Integer             | Input  | 27*N
              RWK,     | Double precision    | Input  | 13*N
              IER)     | Integer             | Output | 
-------------------------------------------------------------------
N
The number of input data points (N > 2).
RLAT,RLON
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
NT
The number of triangles in the triangulation, unless IER .NE. 0, in which case NT = 0. Let NB be the number of boundary points on the convex hull of the data. If NB .GE. 3, then NT = 2N-NB-2, otherwise NT=2N-4. The input data are considered to be bounded if they all lie in one hemisphere.
NTRI
A two-dimensional integer array dimensioned for 3 x NT where NT is the number of triangles in the triangulation (NT is at most 2*N). NTRI contains the triangulation data. The vertices of the Kth triangle are:
         (PLAT(NTRI((1,K)),PLON(NTRI(1,K))
         (PLAT(NTRI((2,K)),PLON(NTRI(2,K))
         (PLAT(NTRI((3,K)),PLON(NTRI(3,K))
         
IWK
An integer workspace of length 27*N.
RWK
A DOUBLE PRECISION workspace of length at least 13*N.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

CSS2CD - convert from lat/lon coordinates to Cartesian coordinates.

CSS2CD is called to find the equivalent Cartesian coordinates on a unit sphere to specified latitude and longitude coordinates. The coordinate of 0. latitude and 0. longitude is converted to Cartesian coordinate (1.,0.,0.). Latitudes and longitudes are assumed to be in degrees. CSS2CD is a double precision version of CSS2C.
------------------------------------------------------------------
             Argument | Type             |  Mode  | Dimension
------------------------------------------------------------------
CALL CSS2CD (N,       | Integer          | Input  |
             RLAT,    | Double precision | Input  | N
             RLON,    | Double precision | Input  | N
             X,       | Double precision | Output | N
             Y,       | Double precision | Output | N
             Z)       | Double precision | Output | N
------------------------------------------------------------------
N
The number of input lat/lon coordinates.
RLAT
An array containing the latitudes of the input coordinates, expressed in degrees.
RLON
An array containing the longitudes of the input coordinates, expressed in degrees.
X,Y,Z
Arrays containing the Cartesian coordinates of the output points. (X(I),Y(I),Z(I)) is the Cartesian coordinate corresponding to the lat/lon coordinate (RLAT(I),RLON(I)) for I=1 to N. X(I)**2 + Y(I)**2 + Z(I)**2 = 1 for I = 1 to N.

CSC2SD - convert from Cartesian coordinates on a unit sphere to lat/lon coordinates

CSC2SD is called to find an equivalent latitude and longitude coordinates on a sphere to a specified Cartesian coordinate on the unit sphere. The coordinate (1.,0.,0.) is mapped to the latitude/longitude coordinate (0.,0.). The latitude/longitude coordinates are returned in degrees, latitudes between -90. and 90. (inclusive) and longitudes between -180. and 180. (inclusive). CSC2SD is a double precision version of CSC2S.
------------------------------------------------------------------
            Argument | Type             |  Mode  | Dimension
------------------------------------------------------------------
CALL CSC2SD (N,       | Integer          | Input  | 
             X,       | Double precision | Input  |  N
             Y,       | Double precision | Input  |  N
             Z,       | Double precision | Input  |  N
             RLAT,    | Double precision | Output |  N
             RLON)    | Double precision | Output |  N
------------------------------------------------------------------
N
The number of input (X,Y,Z) coordinates.
X,Y,Z
The Cartesian coordinates of the input point, X(I)**2 + Y(I)**2 + Z(I)**2 = 1 for I = 1 to N.
RLAT
The latitude of the output coordinates, in degrees.
RLON
The longitude of the output coordinates, in degrees.

CSVOROD - calculate Voronoi polygons.

CSVOROD is called if you want to determine the Voronoi polygons for data randomly positioned on a sphere. Each call to CSVOROD calculates the vertices for the Voronoi polygon surrounding a specified input point. CSVOROD is a double precision version of CSVORO.
------------------------------------------------------------------
              Argument | Type             |  Mode  | Dimension
------------------------------------------------------------------
CALL CSVOROD (NPTS,    | Integer          | Input  |
              RLATI,   | Double precision | Input  | N
              RLONI,   | Double precision | Input  | N
              NI,      | Integer          | Input  | 
              NF,      | Integer          | Input  |
              IWK,     | Integer          | Input  | 27*NPTS
              RWK,     | Double precision | Input  | 9*NPTS
              NC,      | Integer          | Input  |
              RLATO,   | Double precision | Input  | NC
              RLONO,   | Double precision | Input  | NC
              RC,      | Double precision | Output | NC
              NCA,     | Integer          | Output |
              NUMV,    | Integer          | Output |
              NV,      | Integer          | Output | NPTS
              IER)     | Integer          | Output | 
------------------------------------------------------------------
NPTS
The number of input data points (NPTS > 3).
RLATI
An array, dimensioned for NPTS, containing the latitudes of the input coordinates, in degrees.
RLONI
An array, dimensioned for NPTS, containing the longitudes of the input coordinates, in degrees.
NI
The index of the input coordinate for which you want to determine the Voronoi polygon (1 .LE. NI .LE. NPTS).
NF
Flag indicating if this is the first call to CSVORO to retrieve Voronoi polygons for this dataset (1=yes, 0=no). Calls subsequent to the first call for a given dataset are much faster than the first call.
IWK
Integer work space dimensioned for 27*NPTS.
RWK
DOUBLE PRECISION work space dimensioned for 9*NPTS.
NC
The maximum size of the output arrays RLATO, RLONO, and RC. NC should be 2*NPTS.
RLATO
An array of latitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. RLATO is returned in degrees.
RLONO
An array of longitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. RLONO is returned in degrees.
RC
Array containing circumradii (arc lengths in degrees of the angle between a circumcenter and its associated triangle vertices).
NCA
The actual number of circumcenters returned in RLATO and RLONO. This number may be larger than NPTS if the input dataset has boundary points since certain "pseudo points" may have been added to the original dataset in order to complete the Voronoi polygon set.
NUMV
The number of vertices in the Voronoi polygon enclosing the coordinate (RLATI(NI),RLONI(NI)).
NV
An array (dimensioned for NPTS) containing NUMV indices for the Voronoi polygon enclosing the coordinate (RLATI(NI),RLONI(NI)). The indices returned in this array refer to the coordinates returned in RLATO, RLONO, and RC. For example, if the integer "J" is an element of the NV array, then (RLATO(J),RLONO(J)) is a vertex of the Voronoi polygon enclosing (RLATI(NI),RLONI(NI)). The indices in NV list out the vertices of the Voronoi polygon in counter-clockwise order.
IER
An error return value. If IER is returned as 0, then no errors were detected. If IER is non-zero, then refer to the error list in the error table for details.

c_cssgrid - interpolation on a sphere


c_cssgrid is called to interpolate function values starting with data randomly spaced on a sphere.

c_cssgrid uses tension splines to construct an interpolatory surface. As a first step, c_cssgrid creates a Delaunay triangulation of the input data points; the companion entry c_csstri is provided only for those times when you want to calculate the triangulation for independent inspection.

Various aspects of the interpolation algorithm can be changed by using the parmater access routines c_csseti and c_cssetr to change values for the cssgrid control parematers.

Function prototype:

  float *c_cssgrid(int, float [], float [], float [],
                   int, int, float [], float [], int *);

Usage:

-------------------------------------------------
                  Argument | Type     |  Size
-------------------------------------------------
float *c_cssgrid (n,       | int      |
                  rlat,    | float [] | n
                  rlon,    | float [] | n
                  f,       | float [] | n
                  nlat,    | int      | 
                  nlon,    | int      |  
                  plat,    | float [] | ni
                  plon,    | float [] | nj
                  ier      | int *    |
                 );
-------------------------------------------------
n
The number of input data points, n > 2.
rlat,rlon
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
f
Array containing data values. f[i] is the functional value at (rlat[i],rlon[i]) for i = 0 to n-1.
ni
The number of latitudes in the interpolated grid.
nj
The number of longitudes in the interpolated grid. ni and nj can both be 1, allowing for interpolation at a single point.
plat
An array containing the latitudes of the points where interpolated values are to be computed. The values in plat should be in degrees.
plon
An array containing the longitudes of the points where interpolated values are to be computed. The values in plon should be in degrees.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.
Return value:

c_cssgrid returns a pointer to a linear array of data that contains interpolated values at user-specified lat/lon pairs. The returned array stores its values as if they were a 2-dimensional C array with latitude being the first dimension and longitude the second dimension. That is, if out is declared as

  float *out;
and we set:
  out = c_cssgrid(n, rlat, rlon, f, nlat, nlon, plat, plon, &ier);
then out[i*nlon+j] is the interpolated function value at coordinate point (plat[i], plon[j]) for 0 <= i < nlat and 0 <= j < nlon. The space for out is allocated internal to c_cssgrid and is nlat * nlon floats in size.

c_csstri - calculates a Delaunay triangulation


c_csstri is called to find a Delaunay triangulation of data randomly positioned on the surface of a sphere.

Function prototype:

  int   *c_csstri(int, float [], float [], int *, int *);

Usage:
 
-------------------------------------------------
               Argument | Type     |  Size
-------------------------------------------------
int *c_csstri (n,       | int      |
               rlat,    | float [] | n
               rlon,    | float [] | n
               nt,      | int *    |
               ier      | int *    |
              );
-------------------------------------------------
n
The number of input data points, n > 2.
rlat,rlon
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
nt
*nt is the number of triangles in the triangulation, unless *ier is non-zero, in which case *nt = 0. Where nb is the number of boundary points on the convex hull of the data, if nb is greater than 3, then *nt = 2n-nb-2, otherwise *nt = 2n-4. The input data are considered to be bounded if they all lie in one hemisphere.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.
Return value:

c_csstri returns a pointer to a linear array that contains a sequence of integer triples. The elements of a triple are indices of vertices of a triangle. Each index references an original data point as it occurs in sequence in the input data set (numbering starts at 0). For example, if the triple <5,0,2> were in the list of triples, then (rlat[5],rlon[5]), (rlat[0],rlon[0]), and (rlat[2],rlon[2]) would be vertices of a triangle in the Delaunay triangulation.

c_css2c - converts from lat/lon to Cartesian coordinates


c_css2c is called to find the equivalent Cartesian coordinates on a unit sphere to specified latitude and longitude coordinates on a sphere. The coordinate of 0. latitude and 0. longitude is converted to Cartesian coordinate (1.,0.,0.). Latitudes and longitudes are assumed to be in degrees.

Function prototype:

  void   c_css2c(int, float *, float *, float *, float *, float *);

Usage:
 
-------------------------------------------------
              Argument | Type     |  Size
-------------------------------------------------
void c_css2c (n,       | int      |
              plat,    | float *  | n
              plon,    | float *  | n
              x,       | float *  | n
              y,       | float *  | n
              z        | float *  | n
             );
-------------------------------------------------
n
The number of coordinates to be converted.
plat
Contains the latitudes of the input coordinates.
plon
Contains the longitudes of the input coordinates.
x,y,z
Contains the Cartesian coordinates of the output points. These lie on the unit sphere.

c_csc2s - convert from Cartesian coordinates to lat/lon coordinates

c_csc2s converts Cartesian coordinates on a unit sphere to lat/lon coordinates. The coordinate (1.,0.,0.) is mapped to the latitude/longitude coordinate (0.,0.). The latitude/longitude coordinates are returned in degrees.

Function prototype:

     void   c_csc2s(int, float *, float *, float *, float *, float *);
Usage:
--------------------------------------------
              Argument | Type     | Size
--------------------------------------------
void c_csc2s (n,       | int      |
              x,       | float *  |  
              y,       | float *  |  
              z,       | float *  |  
              plat,    | float *  |  
              plon,    | float *  |  
             )   
--------------------------------------------
n
The number of coordinates to be converted.
x,y,z
Cartesian coordinates on the unit sphere.
plat, plon
(plat[i], plon[i]) is the latitude/longitude coordinate equivalent to (x[i],y[i],z[i]) for i equal 0 to n-1. plat and plon are in degrees.

c_csvoro - calculate Voronoi polygons


c_csvoro is called if you want to determine the Voronoi polygons for data randomly positioned on a sphere. Each call to c_csvoro calculates the vertices for the Voronoi polygon surrounding a specified input point.

Function prototype:

void   c_csvoro(int, float [], float [], int, int,
                float [], float [], float [], int *,
                int *, int [], int *);
Usage:
-------------------------------------------------
                Argument | Type     |  Size
-------------------------------------------------
void *c_csvoro (n,       | int      |
                rlat,    | float [] | n
                rlon,    | float [] | n
                ni,      | int      |
                nf,      | int      | 
                plat,    | float [] | 2*n 
                plon,    | float [] | 2*n
                rc,      | float [] | 2*n
                nca,     | int *    |
                numv,    | int *    |
                nv,      | int []   | n
                ier      | int *    |
               );
-------------------------------------------------
n
The number of input data points (n > 3).
rlat
An array containing the latitudes of the input coordinates, in degrees.
rlon
An array containing the longitudes of the input coordinates, in degrees.
ni
The index of the input coordinate for which you want to determine the Voronoi polygon (0 <= ni < n-1).
nf
Flag indicating if this is the first call to c_csvoro to retrieve Voronoi polygons for this dataset (1=yes, 0=no). Calls subsequent to the first call for a given dataset are much faster than the first call.
plat
An array of latitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. plat is returned in degrees.
plon
An array of longitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. plon is returned in degrees.
rc
Array containing circumradii (arc lengths in degrees of the angle between a circumcenter and its associated triangle vertices).
nca
*nca is the actual number of circumcenters returned in plat and plon. This number may be larger than n if the input dataset has boundary points, since certain "pseudo points" may have been added to the original dataset in order to complete the Voronoi polygon set.
numv
*numv is the number of vertices in the Voronoi polygon enclosing the coordinate (rlat[ni],rlon[ni]).
nv
An array containing numv indices for the Voronoi polygon enclosing the coordinate (rlat[ni],rlon[ni]). The indices returned in this array refer to the coordinates returned in plat and plon. For example, if the integer "j" is an element of the nv array, then (plat[j],plon[j]) is a vertex of the Voronoi polygon enclosing (rlat[ni],rlon[ni]). The indices in nv list out the vertices of the Voronoi polygon in counter-clockwise order.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.

c_csseti - Set int valued parameters

c_csseti is used to set values for any of the cssgrid control parameters that take int values. The values set by c_csseti remain in effect until changed by subsequent calls to c_csseti.

Function prototype:

  void c_csseti(char *, int);
Argument description:
-------------------------------------------
              Argument | Type    |  Size  
-------------------------------------------
void c_csseti (pnam,   | char *  |        
               ival);  | int     |
-------------------------------------------
pnam
The name of the control parameter to be assigned an int value.
ival
The value to be assigned to the control parameter whose name is pointed to by pnam.

c_csgeti - Retrieve an int valued parameter

c_csgeti is a called to obtain current values for any of the int valued cssgrid control parameters.

Function prototype:

  void c_csgeti(char *, int *);
Argument description:
-------------------------------------------
              Argument | Type    |  Size  
-------------------------------------------
void c_csgeti (pnam,   | char *  |        
               ival);  | int  *  |
-------------------------------------------
pnam
The name of the control parameter whose value is to be retrieved.
ival
*ival will be the value currently assigned to the control parameter whose name is pointed to by pnam.

c_cssetr - Set float valued parameters

c_cssetr is used to set values for any of the cssgrid control parameters that take float values. The values set by c_cssetr remain in effect until changed by subsequent calls to c_cssetr.

Function prototype:

  void c_cssetr(char *, float);
Argument description:
-------------------------------------------
              Argument | Type    |  Size  
-------------------------------------------
void c_cssetr (pnam,   | char *  |        
               fval);  | float   |
-------------------------------------------
pnam
The name of the control parameter to be assigned a float value.
fval
The value to be assigned to the control parameter whose name is pointed to by pnam.

c_csgetr - Retrieve a float valued parameter

c_csgetr is a called to obtain current values for any of the float valued cssgrid control parameters.

Function prototype:

  void c_csgetr(char *, float *);
Argument description:
-------------------------------------------
              Argument | Type    |  Size  
-------------------------------------------
void c_csgetr (pnam,   | char *  |        
               fval);  | float * |
-------------------------------------------
pnam
The name of the control parameter whose value is to be retrieved.
fval
*fval will be the value currently assigned to the control parameter whose name is pointed to by pnam.

c_cssetd - Set double precision valued parameters

c_cssetd is used to set values for any of the cssgrid control parameters that take double precision values. The values set by c_cssetd remain in effect until changed by subsequent calls to c_cssetd.

Function prototype:

  void c_cssetd(char *, double);
Argument description:
-------------------------------------------
              Argument | Type     |  Size  
-------------------------------------------
void c_cssetd (pnam,   | char *   |        
               dval);  | double   |
-------------------------------------------
pnam
The name of the control parameter to be assigned a float value.
dval
The value to be assigned to the control parameter whose name is pointed to by pnam.

c_csgetd - Retrieve a float valued parameter

c_csgetd is a called to obtain current values for any of the double precision valued cssgrid control parameters.

Function prototype:

  void c_csgetd(char *, double *);
Argument description:
-------------------------------------------
              Argument | Type     |  Size  
-------------------------------------------
void c_csgetd (pnam,   | char *   |        
               dval);  | double * |
-------------------------------------------
pnam
The name of the control parameter whose value is to be retrieved.
dval
*dval will be the value currently assigned to the control parameter whose name is pointed to by pnam.

c_cssgridd - interpolation on a sphere


c_cssgridd is called to interpolate function values starting with data randomly spaced on a sphere. c_cssgridd is a double precision version of c_cssgrid.

c_cssgridd uses tension splines to construct an interpolatory surface. As a first step, c_cssgridd creates a Delaunay triangulation of the input data points; the companion entry c_csstrid is provided only for those times when you want to calculate the triangulation for independent inspection.

Various aspects of the interpolation algorithm can be changed by using the parmater access routines c_csseti and c_cssetd to change values for the cssgrid control parematers.

Function prototype:

  float *c_cssgridd(int, double [], double [], double [],
                   int, int, double [], double [], int *);

Usage:

-------------------------------------------------
                   Argument | Type      |  Size
-------------------------------------------------
float *c_cssgridd (n,       | int       |
                   rlat,    | double [] | n
                   rlon,    | double [] | n
                   f,       | double [] | n
                   nlat,    | int       | 
                   nlon,    | int       |  
                   plat,    | double [] | ni
                   plon,    | double [] | nj
                   ier      | int *     |
                  );
-------------------------------------------------
n
The number of input data points, n > 2.
rlat,rlon
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
f
Array containing data values. f[i] is the functional value at (rlat[i],rlon[i]) for i = 0 to n-1.
ni
The number of latitudes in the interpolated grid.
nj
The number of longitudes in the interpolated grid. ni and nj can both be 1, allowing for interpolation at a single point.
plat
An array containing the latitudes of the points where interpolated values are to be computed. The values in plat should be in degrees.
plon
An array containing the longitudes of the points where interpolated values are to be computed. The values in plon should be in degrees.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.
Return value:

c_cssgridd returns a pointer to a linear array of data that contains interpolated values at user-specified lat/lon pairs. The returned array stores its values as if they were a 2-dimensional C array with latitude being the first dimension and longitude the second dimension. That is, if out is declared as

  double *out;
and we set:
  out = c_cssgridd(n, rlat, rlon, f, nlat, nlon, plat, plon, &ier);
then out[i*nlon+j] is the interpolated function value at coordinate point (plat[i], plon[j]) for 0 <= i < nlat and 0 <= j < nlon. The space for out is allocated internal to c_cssgridd and is nlat * nlon floats in size.

c_csstrid - calculates a Delaunay triangulation


c_csstrid is called to find a Delaunay triangulation of data randomly positioned on the surface of a sphere. c_csstrid is a double precision version of c_csstri.

Function prototype:

  int   *c_csstrid(int, double [], double [], int *, int *)

Usage:
 
-------------------------------------------------
                Argument | Type      |  Size
-------------------------------------------------
int *c_csstrid (n,       | int       |
                rlat,    | double [] | n
                rlon,    | double [] | n
                nt,      | int *     |
                ier      | int *     |
               );
-------------------------------------------------
n
The number of input data points, n > 2.
rlat,rlon
Arrays containing latitude and longitude coordinates, expressed in degrees, of the input data. The first three points must not be collinear (lie on a common great circle).
nt
*nt is the number of triangles in the triangulation, unless *ier is non-zero, in which case *nt = 0. Where nb is the number of boundary points on the convex hull of the data, if nb is greater than 3, then *nt = 2n-nb-2, otherwise *nt = 2n-4. The input data are considered to be bounded if they all lie in one hemisphere.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.
Return value:

c_csstrid returns a pointer to a linear array that contains a sequence of integer triples. The elements of a triple are indices of vertices of a triangle. Each index references an original data point as it occurs in sequence in the input data set (numbering starts at 0). For example, if the triple <5,0,2> were in the list of triples, then (rlat[5],rlon[5]), (rlat[0],rlon[0]), and (rlat[2],rlon[2]) would be vertices of a triangle in the Delaunay triangulation.

c_css2cd - converts from lat/lon to Cartesian coordinates


c_css2cd is called to find the equivalent Cartesian coordinates on a unit sphere to specified latitude and longitude coordinates on a sphere. The coordinate of 0. latitude and 0. longitude is converted to Cartesian coordinate (1.,0.,0.). Latitudes and longitudes are assumed to be in degrees. c_css2cd is a double precision version of c_css2c.

Function prototype:

  void   c_css2cd(int, double *, double *, double *, double *, double *);

Usage:
 
-------------------------------------------------
               Argument | Type      |  Size
-------------------------------------------------
void c_css2cd (n,       | int       |
               plat,    | double *  | n
               plon,    | double *  | n
               x,       | double *  | n
               y,       | double *  | n
               z        | double *  | n
              );
-------------------------------------------------
n
The number of coordinates to be converted.
plat
Contains the latitudes of the input coordinates.
plon
Contains the longitudes of the input coordinates.
x,y,z
Contains the Cartesian coordinates of the output points. These lie on the unit sphere.

c_csc2sd - convert from a Cartesian to a lat/lon coordinate

c_csc2sd converts Cartesian coordinates on a unit sphere to lat/lon coordinates. The coordinate (1.,0.,0.) is mapped to the latitude/longitude coordinate (0.,0.). The latitude/longitude coordinate is returned in degrees. c_csc2sd is a double precision version of c_csc2s.

Function prototype:

     void   c_csc2sd(int, double *, double *, double *, double *, double *);
Usage:
--------------------------------------------
               Argument | Type      | Size
--------------------------------------------
void c_csc2sd (n,       | int       |
               x,       | double *  |  
               y,       | double *  |  
               z,       | double *  |  
               plat,    | double *  |  
               plon,    | double *  |  
              )   
--------------------------------------------
n
The number of coordinates to be converted.
x,y,z
Cartesian coordinates on the unit sphere.
plat, plon
(plat[i], plon[i]) is the latitude/longitude coordinate equivalent to (x[i],y[i],z[i]) for i equal 0 to n-1. plat and plon are in degrees.

c_csvorod - calculate Voronoi polygons


c_csvorod is called if you want to determine the Voronoi polygons for data randomly positioned on a sphere. Each call to c_csvorod calculates the vertices for the Voronoi polygon surrounding a specified input point. c_csvorod is a double precision version of c_csvoro.

Function prototype:

void   c_csvorod(int, double [], double [], int, int,
                 double [], double [], double [], int *,
                 int *, int [], int *);
Usage:
-------------------------------------------------
                 Argument | Type      |  Size
-------------------------------------------------
void *c_csvorod (n,       | int       |
                 rlat,    | double [] | n
                 rlon,    | double [] | n
                 ni,      | int       |
                 nf,      | int       | 
                 plat,    | double [] | 2*n 
                 plon,    | double [] | 2*n
                 rc,      | double [] | 2*n
                 nca,     | int *     |
                 numv,    | int *     |
                 nv,      | int []    | n
                 ier      | int *     |
                );
-------------------------------------------------
n
The number of input data points (n > 3).
rlat
An array containing the latitudes of the input coordinates, in degrees.
rlon
An array containing the longitudes of the input coordinates, in degrees.
ni
The index of the input coordinate for which you want to determine the Voronoi polygon (0 <= ni < n).
nf
Flag indicating if this is the first call to c_csvorod to retrieve Voronoi polygons for this dataset (1=yes, 0=no). Calls subsequent to the first call for a given dataset are much faster than the first call.
plat
An array of latitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. plat is returned in degrees.
plon
An array of longitude values for the Voronoi indices. These are circumcenters of circles passing through the Delaunay triangles. If a coordinate is a boundary point, then the circle may pass through certain "pseudo points" that have been added to the original dataset in order to complete the Voronoi polygon. RLONO is returned in degrees.
rc
Array containing circumradii (arc lengths in degrees of the angle between a circumcenter and its associated triangle vertices).
nca
*nca is the actual number of circumcenters returned in plat and plon. This number may be larger than n if the input dataset has boundary points, since certain "pseudo points" may have been added to the original dataset in order to complete the Voronoi polygon set.
numv
*numv is the number of vertices in the Voronoi polygon enclosing the coordinate (rlat[ni],rlon[ni]).
nv
An array containing numv indices for the Voronoi polygon enclosing the coordinate (rlat[ni],rlon[ni]). The indices returned in this array refer to the coordinates returned in plat and plon. For example, if the integer "j" is an element of the nv array, then (plat[j],plon[j]) is a vertex of the Voronoi polygon enclosing (rlat[ni],rlon[ni]). The indices in nv list out the vertices of the Voronoi polygon in counter-clockwise order.
ier
An error return value. If *ier is returned as 0, then no errors were detected. If *ier is non-zero, then refer to the error list in the error table for details.

home | contents | defs | procedures | examples | errors