My Project
programmer's documentation
Enumerations | Functions | Variables
cs_math.h File Reference
#include "cs_defs.h"
#include <assert.h>
#include <math.h>
Include dependency graph for cs_math.h:

Go to the source code of this file.

Enumerations

enum  cs_math_sym_tensor_component_t {
  XX, YY, ZZ, XY,
  YZ, XZ
}
 

Functions

static int cs_math_binom (int n, int k)
 Computes the binomial coefficient of n and k. More...
 
static cs_real_t cs_math_fabs (cs_real_t x)
 Compute the absolute value of a real value. More...
 
static cs_real_t cs_math_fmin (cs_real_t x, cs_real_t y)
 Compute the min value of two real values. More...
 
static cs_real_t cs_math_fmax (cs_real_t x, cs_real_t y)
 Compute the max value of two real values. More...
 
static cs_real_t cs_math_sq (cs_real_t x)
 Compute the square of a real value. More...
 
static cs_real_t cs_math_pow2 (cs_real_t x)
 Compute the square of a real value. More...
 
static cs_real_t cs_math_pow3 (cs_real_t x)
 Compute the cube of a real value. More...
 
static cs_real_t cs_math_pow4 (cs_real_t x)
 Compute the 4-th power of a real value. More...
 
static cs_real_t cs_math_3_distance (const cs_real_t xa[3], const cs_real_t xb[3])
 Compute the (euclidean) distance between two points xa and xb in a cartesian coordinate system of dimension 3. More...
 
static cs_real_t cs_math_3_distance_dot_product (const cs_real_t xa[3], const cs_real_t xb[3], const cs_real_t xc[3])
 Compute $ (\vect{x}_b - \vect{x}_a) \cdot \vect{x}_c $. More...
 
static cs_real_t cs_math_3_square_distance (const cs_real_t xa[3], const cs_real_t xb[3])
 Compute the squared distance between two points xa and xb in a cartesian coordinate system of dimension 3. More...
 
static cs_real_t cs_math_3_dot_product (const cs_real_t u[3], const cs_real_t v[3])
 Compute the dot product of two vectors of 3 real values. More...
 
static cs_real_t cs_math_3_33_3_dot_product (const cs_real_t n1[3], const cs_real_t t[3][3], const cs_real_t n2[3])
 Compute the dot product of a tensor t with two vectors n1, and n2 n1 t n2. More...
 
static cs_real_t cs_math_3_norm (const cs_real_t v[3])
 Compute the euclidean norm of a vector of dimension 3. More...
 
static cs_real_t cs_math_3_square_norm (const cs_real_t v[3])
 Compute the square norm of a vector of 3 real values. More...
 
static void cs_math_3_normalise (const cs_real_t vin[3], cs_real_t vout[restrict 3])
 Normalise a vector of 3 real values. More...
 
static void cs_math_3_orthogonal_projection (const cs_real_t n[3], const cs_real_t v[3], cs_real_t vout[restrict 3])
 Orthogonal projection of a vector with respect to a normalised vector. More...
 
static void cs_math_3_normal_scaling (const cs_real_t n[3], cs_real_t factor, cs_real_3_t v)
 Add the dot product with a normal vector to the normal direction to a vector. More...
 
static void cs_math_33_normal_scaling_add (const cs_real_t n[3], cs_real_t factor, cs_real_33_t t)
 Add the dot product with a normal vector to the normal,normal component of a tensor: t += factor * n.t.n n(x)n. More...
 
static void cs_math_33_3_product (const cs_real_t m[3][3], const cs_real_t v[3], cs_real_3_t mv)
 Compute the product of a matrix of 3x3 real values by a vector of 3 real values. More...
 
static void cs_math_33_3_product_add (const cs_real_t m[3][3], const cs_real_t v[3], cs_real_3_t mv)
 Compute the product of a matrix of 3x3 real values by a vector of 3 real values add. More...
 
static void cs_math_33t_3_product (const cs_real_t m[3][3], const cs_real_t v[3], cs_real_3_t mv)
 Compute the product of the transpose of a matrix of 3x3 real values by a vector of 3 real values. More...
 
static void cs_math_sym_33_3_product (const cs_real_t m[6], const cs_real_t v[3], cs_real_t mv[restrict 3])
 Compute the product of a symmetric matrix of 3x3 real values by a vector of 3 real values. NB: Symmetric matrix are stored as follows (s11, s22, s33, s12, s23, s13) More...
 
static void cs_math_sym_33_3_product_add (const cs_real_t m[6], const cs_real_t v[3], cs_real_t mv[restrict 3])
 Compute the product of a symmetric matrix of 3x3 real values by a vector of 3 real values and add it to the vector. NB: Symmetric matrix are stored as follows (s11, s22, s33, s12, s23, s13) More...
 
static void cs_math_66_6_product (const cs_real_t m[6][6], const cs_real_t v[6], cs_real_t mv[restrict 6])
 Compute the product of a matrix of 6x6 real values by a vector of 6 real values. More...
 
static void cs_math_66_6_product_add (const cs_real_t m[6][6], const cs_real_t v[6], cs_real_t mv[restrict 6])
 Compute the product of a matrix of 6x6 real values by a vector of 6 real values and add it to the vector. More...
 
static cs_real_t cs_math_33_determinant (const cs_real_t m[3][3])
 Compute the determinant of a 3x3 matrix. More...
 
static cs_real_t cs_math_sym_33_determinant (const cs_real_6_t m)
 Compute the determinant of a 3x3 symmetric matrix. More...
 
static void cs_math_3_cross_product (const cs_real_t u[3], const cs_real_t v[3], cs_real_t uv[restrict 3])
 Compute the cross product of two vectors of 3 real values. More...
 
static cs_real_t cs_math_3_triple_product (const cs_real_t u[3], const cs_real_t v[3], const cs_real_t w[3])
 Compute the triple product. More...
 
static void cs_math_33_inv_cramer (const cs_real_t in[3][3], cs_real_t out[3][3])
 Inverse a 3x3 matrix. More...
 
static void cs_math_33_inv_cramer_in_place (cs_real_t a[3][3])
 Inverse a 3x3 matrix in place, using Cramer's rule. More...
 
static void cs_math_33_inv_cramer_sym_in_place (cs_real_t a[3][3])
 Inverse a 3x3 symmetric matrix (with non-symmetric storage) in place, using Cramer's rule. More...
 
static void cs_math_sym_33_inv_cramer (const cs_real_t s[6], cs_real_t sout[restrict 6])
 Compute the inverse of a symmetric matrix using Cramer's rule. More...
 
static void cs_math_33_product (const cs_real_t m1[3][3], const cs_real_t m2[3][3], cs_real_33_t mout)
 Compute the product of a matrix of 3x3 real values by a matrix of 3x3 real values. More...
 
static void cs_math_33_product_add (const cs_real_t m1[3][3], const cs_real_t m2[3][3], cs_real_33_t mout)
 Compute the product of a matrix of 3x3 real values by a matrix of 3x3 real values and add. More...
 
static void cs_math_sym_33_product (const cs_real_t s1[6], const cs_real_t s2[6], cs_real_t sout[restrict 6])
 Compute the product of two symmetric matrices. Warning: this is valid if and only if s1 and s2 commute (otherwise sout is not symmetric). More...
 
static void cs_math_reduce_sym_prod_33_to_66 (const cs_real_t s[3][3], cs_real_t sout[restrict 6][6])
 Compute a 6x6 matrix A, equivalent to a 3x3 matrix s, such as: A*R_6 = R*s^t + s*R. More...
 
static void cs_math_sym_33_double_product (const cs_real_t s1[6], const cs_real_t s2[6], const cs_real_t s3[6], cs_real_t sout[restrict 3][3])
 Compute the product of three symmetric matrices. More...
 
static void cs_nvec3 (const cs_real_3_t v, cs_nvec3_t *qv)
 Define a cs_nvec3_t structure from a cs_real_3_t. More...
 
void cs_math_set_machine_epsilon (void)
 Compute the value related to the machine precision. More...
 
double cs_math_get_machine_epsilon (void)
 Get the value related to the machine precision. More...
 
void cs_math_3_length_unitv (const cs_real_t xa[3], const cs_real_t xb[3], cs_real_t *len, cs_real_3_t unitv)
 Compute the length (euclidien norm) between two points xa and xb in a cartesian coordinate system of dimension 3. More...
 
void cs_math_sym_33_eigen (const cs_real_t m[6], cs_real_t eig_vals[3])
 Compute all eigenvalues of a 3x3 symmetric matrix with symmetric storage. More...
 
void cs_math_33_eigen (const cs_real_t m[3][3], cs_real_t *eig_ratio, cs_real_t *eig_max)
 Compute max/min eigenvalues ratio and max. eigenvalue of a 3x3 symmetric matrix with non-symmetric storage. More...
 
double cs_math_surftri (const cs_real_t xv[3], const cs_real_t xe[3], const cs_real_t xf[3])
 Compute the area of the convex_hull generated by 3 points. This corresponds to the computation of the surface of a triangle. More...
 
double cs_math_voltet (const cs_real_t xv[3], const cs_real_t xe[3], const cs_real_t xf[3], const cs_real_t xc[3])
 Compute the volume of the convex_hull generated by 4 points. This is equivalent to the computation of the volume of a tetrahedron. More...
 
void cs_math_fact_lu (cs_lnum_t n_blocks, const int b_size, const cs_real_t *a, cs_real_t *a_lu)
 Compute LU factorization of an array of dense matrices of identical size. More...
 
void cs_math_fw_and_bw_lu (const cs_real_t a_lu[], const int n, cs_real_t x[], const cs_real_t b[])
 Block Jacobi utilities. Compute forward and backward to solve an LU P*P system. More...
 

Variables

const cs_real_t cs_math_zero_threshold
 
const cs_real_t cs_math_1ov3
 
const cs_real_t cs_math_2ov3
 
const cs_real_t cs_math_4ov3
 
const cs_real_t cs_math_5ov3
 
const cs_real_t cs_math_1ov6
 
const cs_real_t cs_math_1ov12
 
const cs_real_t cs_math_1ov24
 
const cs_real_t cs_math_epzero
 
const cs_real_t cs_math_infinite_r
 
const cs_real_t cs_math_big_r
 
const cs_real_t cs_math_pi
 

Enumeration Type Documentation

◆ cs_math_sym_tensor_component_t

Enumerator
XX 
YY 
ZZ 
XY 
YZ 
XZ 

Function Documentation

◆ cs_math_33_3_product()

static void cs_math_33_3_product ( const cs_real_t  m[3][3],
const cs_real_t  v[3],
cs_real_3_t  mv 
)
inlinestatic

Compute the product of a matrix of 3x3 real values by a vector of 3 real values.

Parameters
[in]mmatrix of 3x3 real values
[in]vvector of 3 real values
[out]mvvector of 3 real values

◆ cs_math_33_3_product_add()

static void cs_math_33_3_product_add ( const cs_real_t  m[3][3],
const cs_real_t  v[3],
cs_real_3_t  mv 
)
inlinestatic

Compute the product of a matrix of 3x3 real values by a vector of 3 real values add.

Parameters
[in]mmatrix of 3x3 real values
[in]vvector of 3 real values
[in,out]mvvector of 3 real values

◆ cs_math_33_determinant()

static cs_real_t cs_math_33_determinant ( const cs_real_t  m[3][3])
inlinestatic

Compute the determinant of a 3x3 matrix.

Parameters
[in]m3x3 matrix
Returns
the determinant

◆ cs_math_33_eigen()

void cs_math_33_eigen ( const cs_real_t  m[3][3],
cs_real_t eig_ratio,
cs_real_t eig_max 
)

Compute max/min eigenvalues ratio and max. eigenvalue of a 3x3 symmetric matrix with non-symmetric storage.

Based on: Oliver K. Smith "eigenvalues of a symmetric 3x3 matrix", Communication of the ACM (April 1961) (Wikipedia article entitled "Eigenvalue algorithm")

Parameters
[in]m3x3 matrix
[out]eig_ratiomax/min
[out]eig_maxmax. eigenvalue

◆ cs_math_33_inv_cramer()

static void cs_math_33_inv_cramer ( const cs_real_t  in[3][3],
cs_real_t  out[3][3] 
)
inlinestatic

Inverse a 3x3 matrix.

Parameters
[in]inmatrix to inverse
[out]outinversed matrix

◆ cs_math_33_inv_cramer_in_place()

static void cs_math_33_inv_cramer_in_place ( cs_real_t  a[3][3])
inlinestatic

Inverse a 3x3 matrix in place, using Cramer's rule.

Parameters
[in,out]amatrix to inverse

◆ cs_math_33_inv_cramer_sym_in_place()

static void cs_math_33_inv_cramer_sym_in_place ( cs_real_t  a[3][3])
inlinestatic

Inverse a 3x3 symmetric matrix (with non-symmetric storage) in place, using Cramer's rule.

Parameters
[in,out]amatrix to inverse

◆ cs_math_33_normal_scaling_add()

static void cs_math_33_normal_scaling_add ( const cs_real_t  n[3],
cs_real_t  factor,
cs_real_33_t  t 
)
inlinestatic

Add the dot product with a normal vector to the normal,normal component of a tensor: t += factor * n.t.n n(x)n.

Parameters
[in]nnormalised face normal vector
[in]factorfactor
[in,out]tvector to be scaled

◆ cs_math_33_product()

static void cs_math_33_product ( const cs_real_t  m1[3][3],
const cs_real_t  m2[3][3],
cs_real_33_t  mout 
)
inlinestatic

Compute the product of a matrix of 3x3 real values by a matrix of 3x3 real values.

Parameters
[in]m1matrix of 3x3 real values
[in]m2matrix of 3x3 real values
[out]moutmatrix of 3x3 real values

◆ cs_math_33_product_add()

static void cs_math_33_product_add ( const cs_real_t  m1[3][3],
const cs_real_t  m2[3][3],
cs_real_33_t  mout 
)
inlinestatic

Compute the product of a matrix of 3x3 real values by a matrix of 3x3 real values and add.

Parameters
[in]m1matrix of 3x3 real values
[in]m2matrix of 3x3 real values
[out]moutmatrix of 3x3 real values

◆ cs_math_33t_3_product()

static void cs_math_33t_3_product ( const cs_real_t  m[3][3],
const cs_real_t  v[3],
cs_real_3_t  mv 
)
inlinestatic

Compute the product of the transpose of a matrix of 3x3 real values by a vector of 3 real values.

Parameters
[in]mmatrix of 3x3 real values
[in]vvector of 3 real values
[out]mvvector of 3 real values

◆ cs_math_3_33_3_dot_product()

static cs_real_t cs_math_3_33_3_dot_product ( const cs_real_t  n1[3],
const cs_real_t  t[3][3],
const cs_real_t  n2[3] 
)
inlinestatic

Compute the dot product of a tensor t with two vectors n1, and n2 n1 t n2.

Parameters
[in]n1vector of 3 real values
[in]tvector of 3 real values
[in]n2vector of 3 real values
Returns
the resulting dot product n1.t.n2.

◆ cs_math_3_cross_product()

static void cs_math_3_cross_product ( const cs_real_t  u[3],
const cs_real_t  v[3],
cs_real_t  uv[restrict 3] 
)
inlinestatic

Compute the cross product of two vectors of 3 real values.

Parameters
[in]uvector of 3 real values
[in]vvector of 3 real values
[out]uvvector of 3 real values

◆ cs_math_3_distance()

static cs_real_t cs_math_3_distance ( const cs_real_t  xa[3],
const cs_real_t  xb[3] 
)
inlinestatic

Compute the (euclidean) distance between two points xa and xb in a cartesian coordinate system of dimension 3.

Parameters
[in]xafirst coordinate
[in]xbsecond coordinate
Returns
the length between two points xa and xb

◆ cs_math_3_distance_dot_product()

static cs_real_t cs_math_3_distance_dot_product ( const cs_real_t  xa[3],
const cs_real_t  xb[3],
const cs_real_t  xc[3] 
)
inlinestatic

Compute $ (\vect{x}_b - \vect{x}_a) \cdot \vect{x}_c $.

Parameters
[in]xafirst coordinate
[in]xbsecond coordinate
[in]xcthird coordinate
Returns
the dot product

◆ cs_math_3_dot_product()

static cs_real_t cs_math_3_dot_product ( const cs_real_t  u[3],
const cs_real_t  v[3] 
)
inlinestatic

Compute the dot product of two vectors of 3 real values.

Parameters
[in]uvector of 3 real values
[in]vvector of 3 real values
Returns
the resulting dot product u.v.

◆ cs_math_3_length_unitv()

void cs_math_3_length_unitv ( const cs_real_t  xa[3],
const cs_real_t  xb[3],
cs_real_t len,
cs_real_3_t  unitv 
)
inline

Compute the length (euclidien norm) between two points xa and xb in a cartesian coordinate system of dimension 3.

Parameters
[in]xacoordinate of the first extremity
[in]xbcoordinate of the second extremity
[out]lenpointer to the length of the vector va -> vb
[out]unitvunitary vector along xa -> xb
[in]xacoordinate of the first extremity
[in]xbcoordinate of the second extremity
[out]lenpointer to the length of the vector va -> vb
[out]unitvunitary vector anlong va -> vb

◆ cs_math_3_norm()

static cs_real_t cs_math_3_norm ( const cs_real_t  v[3])
inlinestatic

Compute the euclidean norm of a vector of dimension 3.

Parameters
[in]v
Returns
the value of the norm

◆ cs_math_3_normal_scaling()

static void cs_math_3_normal_scaling ( const cs_real_t  n[3],
cs_real_t  factor,
cs_real_3_t  v 
)
inlinestatic

Add the dot product with a normal vector to the normal direction to a vector.

Parameters
[in]nnormalised face normal vector
[in]factorfactor
[in,out]vvector to be scaled

◆ cs_math_3_normalise()

static void cs_math_3_normalise ( const cs_real_t  vin[3],
cs_real_t  vout[restrict 3] 
)
inlinestatic

Normalise a vector of 3 real values.

Parameters
[in]vinvector
[out]voutnormalized vector

◆ cs_math_3_orthogonal_projection()

static void cs_math_3_orthogonal_projection ( const cs_real_t  n[3],
const cs_real_t  v[3],
cs_real_t  vout[restrict 3] 
)
inlinestatic

Orthogonal projection of a vector with respect to a normalised vector.

Parameters
[in]nnormal vector direction
[in]vvector to be projected
[out]voutprojection

◆ cs_math_3_square_distance()

static cs_real_t cs_math_3_square_distance ( const cs_real_t  xa[3],
const cs_real_t  xb[3] 
)
inlinestatic

Compute the squared distance between two points xa and xb in a cartesian coordinate system of dimension 3.

Parameters
[in]xafirst coordinate
[in]xbsecond coordinate
Returns
the square length between two points xa and xb

◆ cs_math_3_square_norm()

static cs_real_t cs_math_3_square_norm ( const cs_real_t  v[3])
inlinestatic

Compute the square norm of a vector of 3 real values.

Parameters
[in]vvector of 3 real values
Returns
square norm of v.

◆ cs_math_3_triple_product()

static cs_real_t cs_math_3_triple_product ( const cs_real_t  u[3],
const cs_real_t  v[3],
const cs_real_t  w[3] 
)
inlinestatic

Compute the triple product.

Parameters
[in]uvector of 3 real values
[in]vvector of 3 real values
[in]wvector of 3 real values
Returns
the scalar triple product

◆ cs_math_66_6_product()

static void cs_math_66_6_product ( const cs_real_t  m[6][6],
const cs_real_t  v[6],
cs_real_t  mv[restrict 6] 
)
inlinestatic

Compute the product of a matrix of 6x6 real values by a vector of 6 real values.

Parameters
[in]mmatrix of 6x6 real values
[in]vvector of 6 real values
[out]mvvector of 6 real values

◆ cs_math_66_6_product_add()

static void cs_math_66_6_product_add ( const cs_real_t  m[6][6],
const cs_real_t  v[6],
cs_real_t  mv[restrict 6] 
)
inlinestatic

Compute the product of a matrix of 6x6 real values by a vector of 6 real values and add it to the vector.

Parameters
[in]mmatrix of 6x6 real values
[in]vvector of 6 real values
[out]mvvector of 6 real values

◆ cs_math_binom()

static int cs_math_binom ( int  n,
int  k 
)
inlinestatic

Computes the binomial coefficient of n and k.

Parameters
[in]nfirst argument
[in]ksecond argument
Returns
the value of binom(n)(k)

◆ cs_math_fabs()

static cs_real_t cs_math_fabs ( cs_real_t  x)
inlinestatic

Compute the absolute value of a real value.

Parameters
[in]xvalue
Returns
the real of the given value

◆ cs_math_fact_lu()

void cs_math_fact_lu ( cs_lnum_t  n_blocks,
int  b_size,
const cs_real_t a,
cs_real_t a_lu 
)

Compute LU factorization of an array of dense matrices of identical size.

Parameters
[in]n_blocksnumber of blocks
[in]b_sizeblock size
[in]amatrix blocks
[out]a_luLU factorizations of matrix blocks

◆ cs_math_fmax()

static cs_real_t cs_math_fmax ( cs_real_t  x,
cs_real_t  y 
)
inlinestatic

Compute the max value of two real values.

Parameters
[in]x,yvalues
Returns
the max value

◆ cs_math_fmin()

static cs_real_t cs_math_fmin ( cs_real_t  x,
cs_real_t  y 
)
inlinestatic

Compute the min value of two real values.

Parameters
[in]x,yvalues
Returns
the min value

◆ cs_math_fw_and_bw_lu()

void cs_math_fw_and_bw_lu ( const cs_real_t  a_lu[],
int  n,
cs_real_t  x[],
const cs_real_t  b[] 
)

Block Jacobi utilities. Compute forward and backward to solve an LU P*P system.

Parameters
[in]a_lumatrix LU factorization
[in]nmatrix size
[out]xsolution
[out]bright hand side

◆ cs_math_get_machine_epsilon()

double cs_math_get_machine_epsilon ( void  )

Get the value related to the machine precision.

◆ cs_math_pow2()

static cs_real_t cs_math_pow2 ( cs_real_t  x)
inlinestatic

Compute the square of a real value.

Parameters
[in]xvalue
Returns
the square of the given value

◆ cs_math_pow3()

static cs_real_t cs_math_pow3 ( cs_real_t  x)
inlinestatic

Compute the cube of a real value.

Parameters
[in]xvalue
Returns
the cube of the given value

◆ cs_math_pow4()

static cs_real_t cs_math_pow4 ( cs_real_t  x)
inlinestatic

Compute the 4-th power of a real value.

Parameters
[in]xvalue
Returns
the 4th power of the given value

◆ cs_math_reduce_sym_prod_33_to_66()

static void cs_math_reduce_sym_prod_33_to_66 ( const cs_real_t  s[3][3],
cs_real_t  sout[restrict 6][6] 
)
inlinestatic

Compute a 6x6 matrix A, equivalent to a 3x3 matrix s, such as: A*R_6 = R*s^t + s*R.

Parameters
[in]s3x3 matrix
[out]sout6x6 matrix

◆ cs_math_set_machine_epsilon()

void cs_math_set_machine_epsilon ( void  )

Compute the value related to the machine precision.

◆ cs_math_sq()

static cs_real_t cs_math_sq ( cs_real_t  x)
inlinestatic

Compute the square of a real value.

Parameters
[in]xvalue
Returns
the square of the given value

◆ cs_math_surftri()

double cs_math_surftri ( const cs_real_t  xv[3],
const cs_real_t  xe[3],
const cs_real_t  xf[3] 
)
inline

Compute the area of the convex_hull generated by 3 points. This corresponds to the computation of the surface of a triangle.

Parameters
[in]xvcoordinates of the first vertex
[in]xecoordinates of the second vertex
[in]xfcoordinates of the third vertex
Returns
the surface of a triangle

◆ cs_math_sym_33_3_product()

static void cs_math_sym_33_3_product ( const cs_real_t  m[6],
const cs_real_t  v[3],
cs_real_t  mv[restrict 3] 
)
inlinestatic

Compute the product of a symmetric matrix of 3x3 real values by a vector of 3 real values. NB: Symmetric matrix are stored as follows (s11, s22, s33, s12, s23, s13)

Parameters
[in]mmatrix of 3x3 real values
[in]vvector of 3 real values
[out]mvvector of 3 real values

◆ cs_math_sym_33_3_product_add()

static void cs_math_sym_33_3_product_add ( const cs_real_t  m[6],
const cs_real_t  v[3],
cs_real_t  mv[restrict 3] 
)
inlinestatic

Compute the product of a symmetric matrix of 3x3 real values by a vector of 3 real values and add it to the vector. NB: Symmetric matrix are stored as follows (s11, s22, s33, s12, s23, s13)

Parameters
[in]mmatrix of 3x3 real values
[in]vvector of 3 real values
[out]mvvector of 3 real values

◆ cs_math_sym_33_determinant()

static cs_real_t cs_math_sym_33_determinant ( const cs_real_6_t  m)
inlinestatic

Compute the determinant of a 3x3 symmetric matrix.

Parameters
[in]m3x3 symmetric matrix
Returns
the determinant

◆ cs_math_sym_33_double_product()

static void cs_math_sym_33_double_product ( const cs_real_t  s1[6],
const cs_real_t  s2[6],
const cs_real_t  s3[6],
cs_real_t  sout[restrict 3][3] 
)
inlinestatic

Compute the product of three symmetric matrices.

Remarks
Symmetric matrix coefficients are stored as follows: (s11, s22, s33, s12, s23, s13)
Parameters
[in]s1symmetric matrix
[in]s2symmetric matrix
[in]s3symmetric matrix
[out]soutsout = s1 * s2 * s3

◆ cs_math_sym_33_eigen()

void cs_math_sym_33_eigen ( const cs_real_t  m[6],
cs_real_t  eig_vals[3] 
)

Compute all eigenvalues of a 3x3 symmetric matrix with symmetric storage.

Based on: Oliver K. Smith "eigenvalues of a symmetric 3x3 matrix", Communication of the ACM (April 1961) (Wikipedia article entitled "Eigenvalue algorithm")

Parameters
[in]m3x3 symmetric matrix (m11, m22, m33, m12, m23, m13)
[out]eig_valssize 3 vector

◆ cs_math_sym_33_inv_cramer()

static void cs_math_sym_33_inv_cramer ( const cs_real_t  s[6],
cs_real_t  sout[restrict 6] 
)
inlinestatic

Compute the inverse of a symmetric matrix using Cramer's rule.

Remarks
Symmetric matrix coefficients are stored as follows: (s11, s22, s33, s12, s23, s13)
Parameters
[in]ssymmetric matrix
[out]soutsout = 1/s1

◆ cs_math_sym_33_product()

static void cs_math_sym_33_product ( const cs_real_t  s1[6],
const cs_real_t  s2[6],
cs_real_t  sout[restrict 6] 
)
inlinestatic

Compute the product of two symmetric matrices. Warning: this is valid if and only if s1 and s2 commute (otherwise sout is not symmetric).

Remarks
Symmetric matrix coefficients are stored as follows: (s11, s22, s33, s12, s23, s13)
Parameters
[in]s1symmetric matrix
[in]s2symmetric matrix
[out]soutsout = s1 * s2

◆ cs_math_voltet()

double cs_math_voltet ( const cs_real_t  xv[3],
const cs_real_t  xe[3],
const cs_real_t  xf[3],
const cs_real_t  xc[3] 
)

Compute the volume of the convex_hull generated by 4 points. This is equivalent to the computation of the volume of a tetrahedron.

Parameters
[in]xvcoordinates of the first vertex
[in]xecoordinates of the second vertex
[in]xfcoordinates of the third vertex
[in]xccoordinates of the fourth vertex
Returns
the volume of the tetrahedron.

◆ cs_nvec3()

static void cs_nvec3 ( const cs_real_3_t  v,
cs_nvec3_t qv 
)
inlinestatic

Define a cs_nvec3_t structure from a cs_real_3_t.

Parameters
[in]vvector of size 3
[out]qvpointer to a cs_nvec3_t structure

Variable Documentation

◆ cs_math_1ov12

const cs_real_t cs_math_1ov12

◆ cs_math_1ov24

const cs_real_t cs_math_1ov24

◆ cs_math_1ov3

const cs_real_t cs_math_1ov3

◆ cs_math_1ov6

const cs_real_t cs_math_1ov6

◆ cs_math_2ov3

const cs_real_t cs_math_2ov3

◆ cs_math_4ov3

const cs_real_t cs_math_4ov3

◆ cs_math_5ov3

const cs_real_t cs_math_5ov3

◆ cs_math_big_r

const cs_real_t cs_math_big_r

◆ cs_math_epzero

const cs_real_t cs_math_epzero

◆ cs_math_infinite_r

const cs_real_t cs_math_infinite_r

◆ cs_math_pi

const cs_real_t cs_math_pi

◆ cs_math_zero_threshold

const cs_real_t cs_math_zero_threshold