Gwyddion – Free SPM (AFM, SNOM/NSOM, STM, MFM, …) data analysis software

Math

Math — Mathematical utility functions

Functions

#define ROUND()
#define GWY_ROUND()
GwyXY * gwy_xy_new ()
GwyXY * gwy_xy_copy ()
void gwy_xy_free ()
GwyXYZ * gwy_xyz_new ()
GwyXYZ * gwy_xyz_copy ()
void gwy_xyz_free ()
gdouble (*GwyRealFunc) ()
gdouble gwy_math_humanize_numbers ()
gboolean gwy_math_is_in_polygon ()
gint gwy_math_find_nearest_line ()
gint gwy_math_find_nearest_point ()
gdouble * gwy_math_lin_solve ()
gdouble * gwy_math_lin_solve_rewrite ()
gboolean gwy_math_tridiag_solve_rewrite ()
gdouble * gwy_math_fit_polynom ()
gboolean gwy_math_choleski_decompose ()
void gwy_math_choleski_solve ()
gboolean gwy_math_choleski_invert ()
guint gwy_math_curvature ()
guint gwy_math_curvature_at_apex ()
guint gwy_math_curvature_at_origin ()
gboolean gwy_math_refine_maximum_1d ()
gboolean gwy_math_refine_maximum_2d ()
gboolean gwy_math_refine_maximum ()
gdouble gwy_math_find_minimum_1d ()
guint * gwy_check_regular_2d_grid ()
gdouble gwy_math_median ()
gdouble gwy_math_kth_rank ()
void gwy_math_kth_ranks ()
void gwy_math_percentiles ()
void gwy_math_sort ()
void gwy_guint_sort ()
void gwy_math_sort_with_index ()
gdouble gwy_math_trimmed_mean ()
gdouble gwy_math_median_uncertainty ()
gint gwy_compare_double ()
gdouble gwy_xlnx_int ()
gdouble gwy_sinc ()
gdouble gwy_canonicalize_angle ()
guint gwy_math_histogram ()

Types and Values

#define GWY_SQRT3
#define GWY_SQRT_PI
  GwyXY
  GwyXYZ

Includes

#include <libgwyddion/gwyddion.h>

Description

Function gwy_math_humanize_numbers() deals with number representation.

Nearest object finding functions gwy_math_find_nearest_line() and gwy_math_find_nearest_point() can be useful in widget and vector layer implementation.

And gwy_math_lin_solve(), gwy_math_lin_solve_rewrite(), and gwy_math_fit_polynom() are general purpose numeric methods.

Functions

ROUND()

#define ROUND(x) ((gint)floor((x) + 0.5))

Warning

ROUND is deprecated and should not be used in newly-written code.

Rounds a number to nearest integer. Use GWY_ROUND instead.

Parameters

x

A double value.

 

GWY_ROUND()

#define GWY_ROUND(x) ((gint)floor((x) + 0.5))

Rounds a number to nearest integer.

Parameters

x

A double value.

 

Since: 2.5

gwy_xy_new ()

GwyXY *
gwy_xy_new (gdouble x,
            gdouble y);

Creates Cartesian coordinates in plane.

This is mostly useful for language bindings.

Parameters

x

X-coordinate.

 

y

Y-coordinate.

 

Returns

New XY structure. The result must be freed using gwy_xy_free(), not g_free().

Since: 2.47

gwy_xy_copy ()

GwyXY *
gwy_xy_copy (const GwyXY *xy);

Copies Cartesian coordinates in plane.

Parameters

xy

Cartesian coordinates in plane.

 

Returns

A copy of xy . The result must be freed using gwy_xy_free(), not g_free().

Since: 2.45

gwy_xy_free ()

void
gwy_xy_free (GwyXY *xy);

Frees Cartesian coordinates in plane created with gwy_xy_copy().

Parameters

xy

Cartesian coordinates in plane.

 

Since: 2.45

gwy_xyz_new ()

GwyXYZ *
gwy_xyz_new (gdouble x,
             gdouble y,
             gdouble z);

Creates Cartesian coordinates in space.

This is mostly useful for language bindings.

Parameters

x

X-coordinate.

 

y

Y-coordinate.

 

z

Z-coordinate.

 

Returns

New XYZ structure. The result must be freed using gwy_xyz_free(), not g_free().

Since: 2.47

gwy_xyz_copy ()

GwyXYZ *
gwy_xyz_copy (const GwyXYZ *xyz);

Copies Cartesian coordinates in space.

Parameters

xyz

Cartesian coordinates in space.

 

Returns

A copy of xyz . The result must be freed using gwy_xyz_free(), not g_free().

Since: 2.45

gwy_xyz_free ()

void
gwy_xyz_free (GwyXYZ *xyz);

Frees Cartesian coordinates in space created with gwy_xyz_copy().

Parameters

xyz

Cartesian coordinates in space.

 

Since: 2.45

GwyRealFunc ()

gdouble
(*GwyRealFunc) (gdouble x,
                gpointer user_data);

gwy_math_humanize_numbers ()

gdouble
gwy_math_humanize_numbers (gdouble unit,
                           gdouble maximum,
                           gint *precision);

Finds a human-friendly representation for a range of numbers.

Parameters

unit

The smallest possible step.

 

maximum

The maximum possible value.

 

precision

A location to store printf() precession, if not NULL.

 

Returns

The magnitude i.e., a power of 1000.

gwy_math_is_in_polygon ()

gboolean
gwy_math_is_in_polygon (gdouble x,
                        gdouble y,
                        const gdouble *poly,
                        guint n);

Establishes wether the test point x , y is inside the polygon poly . The polygon can be defined either clockwise or anti-clockwise and can be a concave, convex or self-intersecting polygon.

Warning

Result can be either TRUE or FALSE if the test point is *exactly* on an edge.

Parameters

x

The x coordinate of the test point.

 

y

The y coordinate of the test point.

 

poly

An array of coordinate pairs (points) that define a polygon.

 

n

The number of corners of the polygon.

 

Returns

TRUE if the test point is inside poly and FALSE otherwise.

Since: 2.7

gwy_math_find_nearest_line ()

gint
gwy_math_find_nearest_line (gdouble x,
                            gdouble y,
                            gdouble *d2min,
                            gint n,
                            const gdouble *coords,
                            const gdouble *metric);

Finds the line from coords nearest to the point (x , y ).

Parameters

x

X-coordinate of the point to search.

 

y

Y-coordinate of the point to search.

 

d2min

Where to store the squared minimal distance, or NULL.

 

n

The number of lines (i.e. coords has 4n items).

 

coords

Line coordinates stored as x00, y00, x01, y01, x10, y10, etc.

 

metric

Metric matrix (2x2, but stored sequentially by rows: m11, m12, m21, m22), it must be positive definite. Vector norm is then calculated as m11*x*x + (m12 + m21)*x*y + m22*y*y. It can be NULL, standard Euclidean metric is then used.

 

Returns

The line number. It may return -1 if (x , y ) doesn't lie in the orthogonal stripe of any of the lines.

gwy_math_find_nearest_point ()

gint
gwy_math_find_nearest_point (gdouble x,
                             gdouble y,
                             gdouble *d2min,
                             gint n,
                             const gdouble *coords,
                             const gdouble *metric);

Finds the point from coords nearest to the point (x , y ).

Parameters

x

X-coordinate of the point to search.

 

y

Y-coordinate of the point to search.

 

d2min

Location to store the squared minimal distance to, or NULL.

 

n

The number of points (i.e. coords has 2n items).

 

coords

Point coordinates stored as x0, y0, x1, y1, x2, y2, etc.

 

metric

Metric matrix (2x2, but stored sequentially by rows: m11, m12, m21, m22). Vector norm is then calculated as m11*x*x + (m12 + m21)*x*y + m22*y*y. It can be NULL, standard Euclidean metric is then used.

 

Returns

The point number.

gwy_math_lin_solve ()

gdouble *
gwy_math_lin_solve (gint n,
                    const gdouble *matrix,
                    const gdouble *rhs,
                    gdouble *result);

Solve a regular system of linear equations.

Parameters

n

The size of the system.

 

matrix

The matrix of the system (n times n ), ordered by row, then column.

 

rhs

The right hand side of the sytem.

 

result

Where the result should be stored. May be NULL to allocate a fresh array for the result.

 

Returns

The solution (result if it wasn't NULL), may be NULL if the matrix is singular.

gwy_math_lin_solve_rewrite ()

gdouble *
gwy_math_lin_solve_rewrite (gint n,
                            gdouble *matrix,
                            gdouble *rhs,
                            gdouble *result);

Solves a regular system of linear equations.

This is a memory-conservative version of gwy_math_lin_solve() overwriting matrix and rhs with intermediate results.

Parameters

n

The size of the system.

 

matrix

The matrix of the system (n times n ), ordered by row, then column.

 

rhs

The right hand side of the sytem.

 

result

Where the result should be stored. May be NULL to allocate a fresh array for the result.

 

Returns

The solution (result if it wasn't NULL), may be NULL if the matrix is singular.

gwy_math_tridiag_solve_rewrite ()

gboolean
gwy_math_tridiag_solve_rewrite (gint n,
                                gdouble *d,
                                const gdouble *a,
                                const gdouble *b,
                                gdouble *rhs);

Solves a tridiagonal system of linear equations.

Parameters

n

The dimension of d .

 

d

The diagonal of a tridiagonal matrix, its contents will be overwritten.

 

a

The above-diagonal stripe (it has n -1 elements).

 

b

The below-diagonal stripe (it has n -1 elements).

 

rhs

The right hand side of the system, upon return it will contain the solution.

 

Returns

TRUE if the elimination suceeded, FALSE if the system is (numerically) singular. The contents of d and rhs may be overwritten in the case of failure too, but not to any meaningful values.

gwy_math_fit_polynom ()

gdouble *
gwy_math_fit_polynom (gint ndata,
                      const gdouble *xdata,
                      const gdouble *ydata,
                      gint n,
                      gdouble *coeffs);

Fits a polynom through a general (x, y) data set.

Parameters

ndata

The number of items in xdata , ydata .

 

xdata

Independent variable data (of size ndata ).

 

ydata

Dependent variable data (of size ndata ).

 

n

The degree of polynom to fit.

 

coeffs

An array of size n +1 to store the coefficients to, or NULL (a fresh array is allocated then).

 

Returns

The coefficients of the polynom (coeffs when it was not NULL, otherwise a newly allocated array).

gwy_math_choleski_decompose ()

gboolean
gwy_math_choleski_decompose (gint n,
                             gdouble *matrix);

Decomposes a symmetric positive definite matrix in place.

Parameters

n

The dimension of a .

 

matrix

Lower triangular part of a symmetric matrix, stored by rows, i.e., matrix = [a_00 a_10 a_11 a_20 a_21 a_22 a_30 ...].

 

Returns

Whether the matrix was really positive definite. If FALSE, the decomposition failed and a does not contain any meaningful values.

gwy_math_choleski_solve ()

void
gwy_math_choleski_solve (gint n,
                         const gdouble *decomp,
                         gdouble *rhs);

Solves a system of linear equations with predecomposed symmetric positive definite matrix a and right hand side b .

Parameters

n

The dimension of a .

 

decomp

Lower triangular part of Choleski decomposition as computed by gwy_math_choleski_decompose().

 

rhs

Right hand side vector. Is is modified in place, on return it contains the solution.

 

gwy_math_choleski_invert ()

gboolean
gwy_math_choleski_invert (gint n,
                          gdouble *matrix);

Inverts a symmetric positive definite matrix in place.

Parameters

n

Matrix size.

 

matrix

Lower triangular part of a symmetric matrix, stored by rows, i.e., matrix = [a_00 a_10 a_11 a_20 a_21 a_22 a_30 ...].

 

Returns

Whether the matrix was really positive definite. If FALSE, the inversion failed and matrix does not contain any meaningful values.

Since: 2.46

gwy_math_curvature ()

guint
gwy_math_curvature (const gdouble *coeffs,
                    gdouble *kappa1,
                    gdouble *kappa2,
                    gdouble *phi1,
                    gdouble *phi2,
                    gdouble *xc,
                    gdouble *yc,
                    gdouble *zc);

Calculates curvature parameters at the apex from two-dimensional quadratic polynomial coefficients.

This is an old name for gwy_math_curvature_at_apex(). See the description there.

Parameters

coeffs

Array of the six polynomial coefficients of a quadratic surface in the following order: 1, x, y, x², xy, y².

 

kappa1

Location to store the smaller curvature to.

 

kappa2

Location to store the larger curvature to.

 

phi1

Location to store the direction of the smaller curvature to.

 

phi2

Location to store the direction of the larger curvature to.

 

xc

Location to store x-coordinate of the centre of the quadratic surface.

 

yc

Location to store y-coordinate of the centre of the quadratic surface.

 

zc

Location to store value at the centre of the quadratic surface.

 

Returns

The number of curved dimensions (0 to 2).

Since: 2.22

gwy_math_curvature_at_apex ()

guint
gwy_math_curvature_at_apex (const gdouble *coeffs,
                            gdouble *kappa1,
                            gdouble *kappa2,
                            gdouble *phi1,
                            gdouble *phi2,
                            gdouble *xc,
                            gdouble *yc,
                            gdouble *zc);

Calculates curvature parameters at the apex from two-dimensional quadratic polynomial coefficients.

See also gwy_math_curvature_at_origin() which computes the local surface curvature at x =0 and y =0.

If the quadratic surface was obtained by fitting the dimensions of the fitted area should not differ, in the lateral coordinates used, by many orders from 1. Otherwise the recognition of flat surfaces might not work.

Curvatures have signs, positive mean a concave (cup-like) surface, negative mean a convex (cap-like) surface. They are ordered including the sign.

Directions are angles from the interval (-π/2, π/2].

If the quadratic surface is degenerate, i.e. flat in at least one direction, the centre is undefined. The centre is then chosen as the closest point the origin of coordinates. For flat surfaces this means the origin is simply returned as the centre position. Consequently, you should use Cartesian coordinates with origin in a natural centre, for instance centre of image or grain.

Parameters

coeffs

Array of the six polynomial coefficients of a quadratic surface in the following order: 1, x, y, x², xy, y².

 

kappa1

Location to store the smaller curvature to.

 

kappa2

Location to store the larger curvature to.

 

phi1

Location to store the direction of the smaller curvature to.

 

phi2

Location to store the direction of the larger curvature to.

 

xc

Location to store x-coordinate of the centre of the quadratic surface.

 

yc

Location to store y-coordinate of the centre of the quadratic surface.

 

zc

Location to store value at the centre of the quadratic surface.

 

Returns

The number of curved dimensions (0 to 2).

Since: 2.61

gwy_math_curvature_at_origin ()

guint
gwy_math_curvature_at_origin (const gdouble *coeffs,
                              gdouble *kappa1,
                              gdouble *kappa2,
                              gdouble *phi1,
                              gdouble *phi2);

Calculates curvature parameters at origin from two-dimensional quadratic polynomial coefficients.

See gwy_math_curvature() for discussion of scaling and sign convenrions. This function function differs from it by computing the local surface curvature at x =0 and y =0, whereas gwy_math_curvature() computes the curvature at the apex of the parabolic surface.

The array coeffs is consistent with gwy_math_curvature_at_apex(), even though here the constant term is not used.

Parameters

coeffs

Array of the six polynomial coefficients of a quadratic surface in the following order: 1, x, y, x², xy, y².

 

kappa1

Location to store the smaller curvature to.

 

kappa2

Location to store the larger curvature to.

 

phi1

Location to store the direction of the smaller curvature to.

 

phi2

Location to store the direction of the larger curvature to.

 

Returns

The number of curved dimensions (0 to 2).

Since: 2.61

gwy_math_refine_maximum_1d ()

gboolean
gwy_math_refine_maximum_1d (const gdouble *y,
                            gdouble *x);

Performs subpixel refinement of parabolic a one-dimensional maximum.

The central value corresponds to x-coordinate 0, distances between values are unity. The refinement is based by fitting a parabola through the maximum. If it fails or the calculated maximum lies farther than the surrounding values the function sets the refined maximum to the origin and returns FALSE.

Parameters

y

Array of length 3, containing the neighbourhood values with the maximum in the centre.

 

x

Location to store the refined x -coordinate.

 

Returns

TRUE if the refinement succeeded, FALSE if it failed. The value of x is usable regardless of the return value.

Since: 2.49

gwy_math_refine_maximum_2d ()

gboolean
gwy_math_refine_maximum_2d (const gdouble *z,
                            gdouble *x,
                            gdouble *y);

Performs subpixel refinement of parabolic a two-dimensional maximum.

The central value corresponds to coordinates (0,0), distances between values are unity. The refinement is based by fitting a two-dimensional parabola through the maximum. If it fails or the calculated maximum lies farther than the surrounding values the function sets the refined maximum to the origin and returns FALSE.

Parameters

z

Array of length 9, containing the square 3x3 neighbourhood values in matrix order and with the maximum in the centre.

 

x

Location to store the refined x -coordinate.

 

y

Location to store the refined y -coordinate.

 

Returns

TRUE if the refinement succeeded, FALSE if it failed. The values of x and y are usable regardless of the return value.

Since: 2.49

gwy_math_refine_maximum ()

gboolean
gwy_math_refine_maximum (const gdouble *z,
                         gdouble *x,
                         gdouble *y);

Performs subpixel refinement of parabolic a two-dimensional maximum.

An alias for gwy_math_refine_maximum_2d().

Parameters

z

Array of length 9, containing the square 3x3 neighbourhood values in matrix order and with the maximum in the centre.

 

x

Location to store the refined x -coordinate.

 

y

Location to store the refined y -coordinate.

 

Returns

TRUE if the refinement succeeded, FALSE if it failed. The values of x and y are usable regardless of the return value.

Since: 2.42

gwy_math_find_minimum_1d ()

gdouble
gwy_math_find_minimum_1d (GwyRealFunc function,
                          gdouble a,
                          gdouble b,
                          gpointer user_data);

Finds a minimum of a real function in a finite interval.

The function simply does what it says on the tin. If there are multiple minima in [a,b] any of them can be returned, even though some effort to scan the interval is made. There is no requiement for the minimum to lie inside [a,b]; if it occurrs at one of the endpoints, the endpoint is returned.

Parameters

function

Function to minimize.

 

a

First interval endpoint.

 

b

Second interval endpoint.

 

user_data

User data passed to function .

 

Since: 2.51

gwy_check_regular_2d_grid ()

guint *
gwy_check_regular_2d_grid (const gdouble *coords,
                           guint stride,
                           guint n,
                           gdouble tolerance,
                           guint *xres,
                           guint *yres,
                           GwyXY *xymin,
                           GwyXY *xystep);

Detects if points in plane form a regular rectangular grid oriented along the Cartesian axes.

Points lying in one straight line are not considered to form a rectangle.

When the function fails, i.e. the points do not form a regular grid, the values of output arguments are undefined.

Parameters

coords

Array of n coordinate pairs in plane. You can also typecast GwyXY or GwyXYZ to doubles.

 

stride

Actual number of double values in one block. It must be at least 2 if coords contains just alternating x and y . If you pass an typecast GwyXYZ array give stride as 3, etc.

 

n

Number of items in coords .

 

tolerance

Relative distance from pixel center which is still considered OK. Pass a negative value for some reasonable default. The maximum meaningful value is 0.5, beyond that the point would end up in a different pixel.

 

xres

Location where to store the number of columns.

 

yres

Location where to store the number of rows.

 

xymin

Location where to store the minimum coordinates (top left corner).

 

xystep

Location where to store the pixel size.

 

Returns

On success, a newly allocated array mapping grid indices (i *xres +j ) to indices in coords . NULL is returned on failure.

Since: 2.48

gwy_math_median ()

gdouble
gwy_math_median (gsize n,
                 gdouble *array);

Finds median of an array of values using Quick select algorithm.

See gwy_math_kth_rank() for details of how the values are shuffled.

Parameters

n

Number of items in array .

 

array

Array of doubles. It is shuffled by this function.

 

Returns

The median value of array .

gwy_math_kth_rank ()

gdouble
gwy_math_kth_rank (gsize n,
                   gdouble *array,
                   gsize k);

Finds k-th item of an array of values using Quick select algorithm.

The value positions change as follows. The returned value is guaranteed to be at k -th position in the array (i.e. correctly ranked). All other values are correctly ordered with respect to this value: preceeding values are smaller (or equal) and following values are larger (or equal).

Parameters

n

Number of items in array .

 

array

Array of doubles. It is shuffled by this function.

 

k

Rank of the value to find (from lowest to highest).

 

Returns

The k -th value of array if it was sorted.

Since: 2.50

gwy_math_kth_ranks ()

void
gwy_math_kth_ranks (gsize n,
                    gdouble *array,
                    guint nk,
                    const guint *k,
                    gdouble *values);

Finds simultaneously several k-th items of an array of values.

The values are shuffled similarly to gwy_math_kth_rank(), except that the guarantee holds for all given ranks simultaneously. All values with explicitly requested ranks are at their correct positions and all values lying between them in the array are also between them numerically.

Parameters

n

Number of items in array .

 

array

Array of doubles. It is shuffled by this function.

 

nk

Number of ranked values to find (sizes of arrays k and values ).

 

k

Ranks of the value to find.

 

values

Array where to store values with ranks (from smallest to highest) given in k .

 

Since: 2.50

gwy_math_percentiles ()

void
gwy_math_percentiles (gsize n,
                      gdouble *array,
                      GwyPercentileInterpolationType interp,
                      guint np,
                      const gdouble *p,
                      gdouble *values);

Finds simultaneously several percentiles of an array of values.

The values in array are shuffled similarly to gwy_math_kth_ranks(). However, it is difficult to state how exactly p translates to the values that become correctly ranked (and it depends on interp ). Hence you can only assume the set of values is preserved.

Parameters

n

Number of items in array .

 

array

Array of doubles. It is shuffled by this function.

 

interp

Interpolation method to use for percentiles that do not correspond exactly to an integer rank.

 

np

Number of percentiles to find.

 

p

Array of size np with the percentiles to compute. The values are in percents, i.e. from the range [0,100].

 

values

Array where to store values with percentiles given in p .

 

Since: 2.50

gwy_math_sort ()

void
gwy_math_sort (gsize n,
               gdouble *array);

Sorts an array of doubles using a quicksort algorithm.

This is usually about twice as fast as the generic quicksort function thanks to specialization for doubles.

Parameters

n

Number of items in array .

 

array

Array of doubles to sort in place.

 

gwy_guint_sort ()

void
gwy_guint_sort (gsize n,
                guint *array);

Sorts an array of unsigned integers using a quicksort algorithm.

This is usually about twice as fast as the generic quicksort function thanks to specialization for integers.

Parameters

n

Number of items in array .

 

array

Array of guint values to sort in place.

 

Since: 2.50

gwy_math_sort_with_index ()

void
gwy_math_sort_with_index (gsize n,
                          gdouble *array,
                          guint *index_array);

Sorts an array of doubles using a quicksort algorithm, remembering the permutation.

The simplest and probably most common use of index_array is to fill it with numbers 0 to n -1 before calling gwy_math_sort(). After sorting, index_array [i ] then contains the original position of the i -th item of the sorted array.

Parameters

n

Number of items in array .

 

array

Array of doubles to sort in place.

 

index_array

Array of integer identifiers of the items that are permuted simultaneously with array .

 

Since: 2.50

gwy_math_trimmed_mean ()

gdouble
gwy_math_trimmed_mean (gsize n,
                       gdouble *array,
                       guint nlowest,
                       guint nhighest);

Finds trimmed mean of an array of values.

At least one value must remain after the trimming, i.e. nlowest + nhighest must be smaller than n . Usually one passes the same number as both nlowest and nhighest , but it is not a requirement.

The function can be also used to calculate normal mean values as it implements efficiently the cases when no trimming is done at either end.

Parameters

n

Number of items in array .

 

array

Array of doubles. It is shuffled by this function.

 

nlowest

The number of lowest values to discard.

 

nhighest

The number of highest values to discard.

 

Returns

The trimmed mean.

Since: 2.50

gwy_math_median_uncertainty ()

gdouble
gwy_math_median_uncertainty (gsize n,
                             gdouble *array,
                             gdouble *uarray);

Find the uncertainty value corresponding to data median.

Note that this is not the uncertainty arising from the calculation of the median. It is just the uncertainty of the single value that happens to be the data median. As such, the function is not very useful.

Parameters

n

Number of items in array .

 

array

Array of doubles. It is modified by this function. All values are kept, but their positions in the array change.

 

uarray

Array of value unvertainries. It is modified by this function. All values are kept, but their positions in the array change.

 

Returns

The uncertainty of the median value.

Since: 2.23

gwy_compare_double ()

gint
gwy_compare_double (gconstpointer a,
                    gconstpointer b);

Compares two double values, given as pointers.

This function is suitable as GCompareFunc and can be also used with plain qsort(). The typical usage is sorting of arrays containing structs where the first item is a floating point value/coordinate, followed by additional data. For sorting of plain arrays of doubles use gwy_math_sort().

It should only be used to sort normal numbers. The behaviour for NaNs is undefined.

Parameters

a

Pointer to a double.

 

b

Pointer to a double.

 

Since: 2.62

gwy_xlnx_int ()

gdouble
gwy_xlnx_int (guint x);

Calculates natural logarithm multiplied by the argument for integers.

The value for zero x is taken as the limit, i.e. zero.

This function is useful for entropy calculations where values of n *log(n ) can be evaulated a lot for small n . Therefore, values for small arguments are tabulated. For large arguments the function is evaluated using the standard log() function which is of course slower.

Parameters

x

Value to calculate x *log(x ) of.

 

Returns

Value of x *log(x ).

Since: 2.44

gwy_sinc ()

gdouble
gwy_sinc (gdouble x);

Calculates the sinc function.

The sinc function is equal to sin(x )/x for non-zero x , and defined to the limit 1 for zero x .

Parameters

x

Value to calculate sinc (cardinal sine) of.

 

Returns

Value of sinc(x ).

Since: 2.51

gwy_canonicalize_angle ()

gdouble
gwy_canonicalize_angle (gdouble phi,
                        gboolean positive,
                        gboolean oriented);

Canonicalizes an angle to requested interval.

For positive =FALSE, oriented =FALSE the output interval is [-π/2,π/2].

For positive =FALSE, oriented =TRUE the output interval is [-π,π].

For positive =TRUE, oriented =FALSE the output interval is [0,π).

For positive =TRUE, oriented =TRUE the output interval is [0,2π).

Parameters

phi

Angle to canonicalize, in radians.

 

positive

TRUE if a positive angle is requested, FALSE for outputs symmetrical around zero.

 

oriented

TRUE for direction of a vector, FALSE for the direction of a line (i.e. with no distinction between forward and backward direction).

 

Returns

Canonicalized angle, equivalent (in given sense) to phi .

Since: 2.50

gwy_math_histogram ()

guint
gwy_math_histogram (const gdouble *values,
                    guint n,
                    gdouble min,
                    gdouble max,
                    guint nbins,
                    guint *counts);

Counts the numbers of values falling into equal-sized bins.

The value of min must not be larger than max . The values may lie outside [min ,max ]. They are not counted in the histogram, nor the returned total.

Rounding rules for values exactly at the edge of two bins are arbitrary and must not be relied upon.

Parameters

values

Values to make histogram from.

 

n

Number of values in values .

 

min

Minimum value to consider (left edge of histogram).

 

max

Maximum value to consider (right edge of histogram).

 

nbins

Number of histogram bins (number of counts items), a positive number.

 

counts

Array where to store the counts.

 

Returns

The number of values inside the entire histogram, i.e. at most n but possibly a reduced count.

Since: 2.49

Types and Values

GWY_SQRT3

#define GWY_SQRT3 1.73205080756887729352744634150587236694280525381038

The square root of 3.

GWY_SQRT_PI

#define GWY_SQRT_PI 1.77245385090551602729816748334114518279754945612237

The square root of pi.

GwyXY

typedef struct {
    gdouble x;
    gdouble y;
} GwyXY;

Representation of Cartesian coordinates in plane.

Members

gdouble x;

X-coordinate.

 

gdouble y;

Y-coordinate.

 

Since: 2.45

GwyXYZ

typedef struct {
    gdouble x;
    gdouble y;
    gdouble z;
} GwyXYZ;

Representation of Cartesian coordinates in space.

Members

gdouble x;

X-coordinate.

 

gdouble y;

Y-coordinate.

 

gdouble z;

Z-coordinate.

 

Since: 2.45

See Also

GwyNLFitter, non-linear least square fitter;

Math Fallback,

fallback mathematical functions

© David Nečas and Petr Klapetek

Home Download News Features Screenshots Documentation Communicate Participate Resources Publications Applications Site Map

Valid XHTML 1.0 Valid CSS