GRASS GIS 7 Programmer's Manual  7.9.dev(2021)-e5379bbd7
defs/gmath.h
Go to the documentation of this file.
1 #ifndef GRASS_GMATHDEFS_H
2 #define GRASS_GMATHDEFS_H
3 
4 /* dalloc.c */
5 double *G_alloc_vector(size_t);
6 double **G_alloc_matrix(int, int);
7 float *G_alloc_fvector(size_t);
8 float **G_alloc_fmatrix(int, int);
9 void G_free_vector(double *);
10 void G_free_matrix(double **);
11 void G_free_fvector(float *);
12 void G_free_fmatrix(float **);
13 
14 /* ialloc.c */
15 int *G_alloc_ivector(size_t);
16 int **G_alloc_imatrix(int, int);
17 void G_free_ivector(int *);
18 void G_free_imatrix(int **);
19 
20 /* fft.c */
21 extern int fft(int, double *[2], int, int, int);
22 extern int fft2(int, double (*)[2], int, int, int);
23 
24 /* gauss.c */
25 extern double G_math_rand_gauss(double);
26 
27 /* max_pow2.c */
28 extern long G_math_max_pow2 (long n);
29 extern long G_math_min_pow2 (long n);
30 
31 /* rand1.c */
32 extern void G_math_srand(int);
33 extern int G_math_srand_auto(void);
34 extern float G_math_rand(void);
35 
36 /* del2g.c */
37 extern int del2g(double *[2], int, double);
38 
39 /* getg.c */
40 extern int getg(double, double *[2], int);
41 
42 /* eigen_tools.c */
43 extern int G_math_egvorder(double *, double **, long);
44 
45 /* mult.c */
46 extern int G_math_complex_mult (double *v1[2], int size1, double *v2[2], int size2, double *v3[2], int size3);
47 
48 /* lu.c*/
49 extern int G_ludcmp(double **, int, int *, double *);
50 extern void G_lubksb(double **a, int n, int *indx, double b[]);
51 
52 /* findzc.c */
53 extern int G_math_findzc(double conv[], int size, double zc[], double thresh, int num_orients);
54 
55 
56 /* *************************************************************** */
57 /* ***** WRAPPER FOR CCMATH FUNCTIONS USED IN GRASS ************** */
58 /* *************************************************************** */
59 extern int G_math_solv(double **,double *,int);
60 extern int G_math_solvps(double **,double *,int);
61 extern void G_math_solvtd(double *,double *,double *,double *,int);
62 extern int G_math_solvru(double **,double *,int);
63 extern int G_math_minv(double **,int);
64 extern int G_math_psinv(double **,int);
65 extern int G_math_ruinv(double **,int);
66 extern void G_math_eigval(double **,double *,int);
67 extern void G_math_eigen(double **,double *,int);
68 extern double G_math_evmax(double **,double *,int);
69 extern int G_math_svdval(double *,double **,int,int);
70 extern int G_math_sv2val(double *,double **,int,int);
71 extern int G_math_svduv(double *,double **,double **, int,double **,int);
72 extern int G_math_sv2uv(double *,double **,double **,int,double **,int);
73 extern int G_math_svdu1v(double *,double **,int,double **,int);
74 
75 /* *************************************************************** */
76 /* *************** LINEARE EQUATION SYSTEM PART ****************** */
77 /* *************************************************************** */
78 
79 /* Sparse matrix and sparse vector functions
80  * */
83 extern void G_math_free_spmatrix(G_math_spvector ** , int );
85 extern int G_math_add_spvector(G_math_spvector **, G_math_spvector * , int );
86 extern G_math_spvector **G_math_A_to_Asp(double **, int, double);
87 extern double **G_math_Asp_to_A(G_math_spvector **, int);
88 extern double **G_math_Asp_to_sband_matrix(G_math_spvector **, int, int);
89 extern G_math_spvector **G_math_sband_matrix_to_Asp(double **, int, int, double);
90 extern void G_math_print_spmatrix(G_math_spvector **, int);
91 extern void G_math_Ax_sparse(G_math_spvector **, double *, double *, int );
92 
93 /*Symmetric band matrix handling */
94 extern double **G_math_matrix_to_sband_matrix(double **, int, int);
95 extern double **G_math_sband_matrix_to_matrix(double **, int, int);
96 extern void G_math_Ax_sband(double ** A, double *x, double *y, int rows, int bandwidth);
97 
98 /*linear equation solver, most of them are multithreaded with OpenMP*/
99 extern int G_math_solver_gauss(double **, double *, double *, int );
100 extern int G_math_solver_lu(double **, double *, double *, int );
101 extern int G_math_solver_cholesky(double **, double *, double *, int , int );
102 extern void G_math_solver_cholesky_sband(double **, double *, double *, int, int);
103 extern int G_math_solver_jacobi(double **, double *, double *, int , int , double , double );
104 extern int G_math_solver_gs(double **, double *, double *, int , int , double , double );
105 
106 extern int G_math_solver_pcg(double **, double *, double *, int , int , double , int );
107 extern int G_math_solver_cg(double **, double *, double *, int , int , double );
108 extern int G_math_solver_cg_sband(double **, double *, double *, int, int, int, double);
109 extern int G_math_solver_bicgstab(double **, double *, double *, int , int , double );
110 extern int G_math_solver_sparse_jacobi(G_math_spvector **, double *, double *, int , int , double , double );
111 extern int G_math_solver_sparse_gs(G_math_spvector **, double *, double *, int , int , double , double );
112 extern int G_math_solver_sparse_pcg(G_math_spvector **, double *, double *, int , int , double , int );
113 extern int G_math_solver_sparse_cg(G_math_spvector **, double *, double *, int , int , double );
114 extern int G_math_solver_sparse_bicgstab(G_math_spvector **, double *, double *, int , int , double );
115 
116 /* solver algoithms and helper functions*/
117 extern void G_math_gauss_elimination(double **, double *, int );
118 extern void G_math_lu_decomposition(double **, double *, int );
119 extern int G_math_cholesky_decomposition(double **, int , int );
120 extern void G_math_cholesky_sband_decomposition(double **, double **, int, int);
121 extern void G_math_backward_substitution(double **, double *, double *, int );
122 extern void G_math_forward_substitution(double **, double *, double *, int );
123 extern void G_math_cholesky_sband_substitution(double **, double *, double *, int, int);
124 
125 /*BLAS like level 1,2 and 3 functions*/
126 
127 /*level 1 vector - vector grass implementation with OpenMP thread support*/
128 extern void G_math_d_x_dot_y(double *, double *, double *, int );
129 extern void G_math_d_asum_norm(double *, double *, int );
130 extern void G_math_d_euclid_norm(double *, double *, int );
131 extern void G_math_d_max_norm(double *, double *, int );
132 extern void G_math_d_ax_by(double *, double *, double *, double , double , int );
133 extern void G_math_d_copy(double *, double *, int );
134 
135 extern void G_math_f_x_dot_y(float *, float *, float *, int );
136 extern void G_math_f_asum_norm(float *, float *, int );
137 extern void G_math_f_euclid_norm(float *, float *, int );
138 extern void G_math_f_max_norm(float *, float *, int );
139 extern void G_math_f_ax_by(float *, float *, float *, float , float , int );
140 extern void G_math_f_copy(float *, float *, int );
141 
142 extern void G_math_i_x_dot_y(int *, int *, double *, int );
143 extern void G_math_i_asum_norm(int *, double *, int );
144 extern void G_math_i_euclid_norm(int *, double *,int );
145 extern void G_math_i_max_norm(int *, int *, int );
146 extern void G_math_i_ax_by(int *, int *, int *, int , int , int );
147 extern void G_math_i_copy(int *, int *, int );
148 
149 /*ATLAS blas level 1 wrapper*/
150 extern double G_math_ddot(double *, double *, int );
151 extern float G_math_sdot(float *, float *, int );
152 extern float G_math_sdsdot(float *, float *, float , int );
153 extern double G_math_dnrm2(double *, int );
154 extern double G_math_dasum(double *, int );
155 extern double G_math_idamax(double *, int );
156 extern float G_math_snrm2(float *, int );
157 extern float G_math_sasum(float *, int );
158 extern float G_math_isamax(float *, int );
159 extern void G_math_dscal(double *, double , int );
160 extern void G_math_sscal(float *, float , int );
161 extern void G_math_dcopy(double *, double *, int );
162 extern void G_math_scopy(float *, float *, int );
163 extern void G_math_daxpy(double *, double *, double , int );
164 extern void G_math_saxpy(float *, float *, float , int );
165 
166 /*level 2 matrix - vector grass implementation with OpenMP thread support*/
167 extern void G_math_d_Ax(double **, double *, double *, int , int );
168 extern void G_math_f_Ax(float **, float *, float *, int , int );
169 extern void G_math_d_x_dyad_y(double *, double *, double **, int, int );
170 extern void G_math_f_x_dyad_y(float *, float *, float **, int, int );
171 extern void G_math_d_aAx_by(double **, double *, double *, double , double , double *, int , int );
172 extern void G_math_f_aAx_by(float **, float *, float *, float , float , float *, int , int );
173 extern int G_math_d_A_T(double **A, int rows);
174 extern int G_math_f_A_T(float **A, int rows);
175 
176 /*level 3 matrix - matrix grass implementation with OpenMP thread support*/
177 extern void G_math_d_aA_B(double **, double **, double , double **, int , int );
178 extern void G_math_f_aA_B(float **, float **, float , float **, int , int );
179 extern void G_math_d_AB(double **, double **, double **, int , int , int );
180 extern void G_math_f_AB(float **, float **, float **, int , int , int );
181 
182 #endif /* GRASS_GMATHDEFS_H */
void G_math_solver_cholesky_sband(double **, double *, double *, int, int)
Cholesky symmetric band matrix solver for linear equation systems of type Ax = b. ...
void G_math_eigen(double **, double *, int)
Compute the eigenvalues and eigenvectors of a real symmetric matrix A.
double G_math_dnrm2(double *, int)
Compute the euclidean norm of vector x using the ATLAS routine cblas_dnrm2.
void G_math_i_x_dot_y(int *, int *, double *, int)
Compute the dot product of vector x and y.
Definition: blas_level_1.c:482
double G_math_ddot(double *, double *, int)
Compute the dot product of vector x and y using the ATLAS routine cblas_ddot.
int * G_alloc_ivector(size_t)
Vector matrix memory allocation.
Definition: ialloc.c:41
int ** G_alloc_imatrix(int, int)
Matrix memory allocation.
Definition: ialloc.c:58
double G_math_evmax(double **, double *, int)
void G_math_d_AB(double **, double **, double **, int, int, int)
Matrix multiplication.
Definition: blas_level_3.c:174
int G_math_solver_cholesky(double **, double *, double *, int, int)
The choleksy decomposition solver for quardatic, symmetric positiv definite matrices.
void G_math_backward_substitution(double **, double *, double *, int)
backward substitution
int fft2(int, double(*)[2], int, int, int)
Fast Fourier Transform for two-dimensional array.
Definition: fft.c:70
void G_math_forward_substitution(double **, double *, double *, int)
forward substitution
G_math_spvector * G_math_alloc_spvector(int)
Allocate memory for a sparse vector.
Definition: sparse_matrix.c:76
void G_free_vector(double *)
Vector memory deallocation.
Definition: dalloc.c:129
The row vector of the sparse matrix.
Definition: gmath.h:54
int G_math_solver_bicgstab(double **, double *, double *, int, int, double)
The iterative biconjugate gradients solver with stabilization for unsymmetric non-definite matrices...
int G_math_solver_sparse_gs(G_math_spvector **, double *, double *, int, int, double, double)
The iterative gauss seidel solver for sparse matrices.
int del2g(double *[2], int, double)
Definition: del2g.c:46
double ** G_math_Asp_to_A(G_math_spvector **, int)
Convert a sparse matrix into a quadratic matrix.
int G_math_psinv(double **, int)
Invert (in place) a symmetric real matrix, V -> Inv(V).
void G_math_d_aA_B(double **, double **, double, double **, int, int)
Add two matrices and scale matrix A with the scalar a.
Definition: blas_level_3.c:50
int G_math_solvru(double **, double *, int)
float * G_alloc_fvector(size_t)
Floating point vector memory allocation.
Definition: dalloc.c:85
int G_math_f_A_T(float **A, int rows)
Compute the transposition of matrix A. Matrix A will be overwritten.
Definition: blas_level_2.c:372
int G_math_solver_sparse_pcg(G_math_spvector **, double *, double *, int, int, double, int)
The iterative preconditioned conjugate gradients solver for sparse symmetric positive definite matric...
double G_math_idamax(double *, int)
Compute the maximum norm of vector x using the ATLAS routine cblas_idamax.
void G_math_free_spmatrix(G_math_spvector **, int)
Release the memory of the sparse matrix.
int G_math_solver_pcg(double **, double *, double *, int, int, double, int)
The iterative preconditioned conjugate gradients solver for symmetric positive definite matrices...
void G_math_f_x_dot_y(float *, float *, float *, int)
Compute the dot product of vector x and y.
Definition: blas_level_1.c:264
void G_math_daxpy(double *, double *, double, int)
Scale vector x with scalar a and add it to y.
void G_math_cholesky_sband_substitution(double **, double *, double *, int, int)
Forward and backward substitution of a lower tringular symmetric band matrix of A from system Ax = b...
void G_math_d_max_norm(double *, double *, int)
Compute the maximum norm of vector x.
Definition: blas_level_1.c:142
void G_math_free_spvector(G_math_spvector *)
Release the memory of the sparse vector.
Definition: sparse_matrix.c:98
double ** G_math_Asp_to_sband_matrix(G_math_spvector **, int, int)
Convert a symmetric sparse matrix into a symmetric band matrix.
float G_math_rand(void)
Definition: rand1.c:17
int getg(double, double *[2], int)
Definition: getg.c:16
void G_math_d_copy(double *, double *, int)
Copy the vector x to y.
Definition: blas_level_1.c:237
G_math_spvector ** G_math_A_to_Asp(double **, int, double)
Convert a quadratic matrix into a sparse matrix.
int G_math_solver_lu(double **, double *, double *, int)
The LU solver for quardatic matrices.
double ** G_alloc_matrix(int, int)
Matrix memory allocation.
Definition: dalloc.c:60
void G_math_d_ax_by(double *, double *, double *, double, double, int)
Scales vectors x and y with the scalars a and b and adds them.
Definition: blas_level_1.c:173
void G_math_d_x_dot_y(double *, double *, double *, int)
Compute the dot product of vector x and y.
Definition: blas_level_1.c:48
G_math_spvector ** G_math_sband_matrix_to_Asp(double **, int, int, double)
Convert a symmetric band matrix into a sparse matrix.
void G_math_eigval(double **, double *, int)
Compute the eigenvalues of a real symmetric matrix A.
void G_math_f_ax_by(float *, float *, float *, float, float, int)
Scales vectors x and y with the scalars a and b and adds them.
Definition: blas_level_1.c:392
int G_math_svdu1v(double *, double **, int, double **, int)
Compute the singular value transformation with A overloaded by the partial U-matrix.
float G_math_snrm2(float *, int)
Compute the euclidean norm of vector x using the ATLAS routine cblas_dnrm2.
#define x
float G_math_sasum(float *, int)
Compute the absolute sum norm of vector x using the ATLAS routine cblas_dasum.
void G_math_i_max_norm(int *, int *, int)
Compute the maximum norm of vector x.
Definition: blas_level_1.c:576
void G_math_f_x_dyad_y(float *, float *, float **, int, int)
Compute the dyadic product of two vectors. The result is stored in the matrix A.
Definition: blas_level_2.c:145
void G_math_srand(int)
Seed the pseudo-random number generator.
Definition: rand1.c:28
void G_math_Ax_sparse(G_math_spvector **, double *, double *, int)
Compute the matrix - vector product of sparse matrix **Asp and vector x.
void G_math_i_copy(int *, int *, int)
Copy the vector x to y.
Definition: blas_level_1.c:670
double ** G_math_matrix_to_sband_matrix(double **, int, int)
Convert a symmetrix matrix into a symmetric band matrix.
int G_math_ruinv(double **, int)
Invert an upper right triangular matrix T -> Inv(T).
int G_math_svduv(double *, double **, double **, int, double **, int)
void G_math_f_Ax(float **, float *, float *, int, int)
Compute the matrix - vector product of matrix A and vector x.
Definition: blas_level_2.c:81
int G_math_srand_auto(void)
Seed the pseudo-random number generator from the time and PID.
Definition: rand1.c:39
void G_math_f_copy(float *, float *, int)
Copy the vector x to y.
Definition: blas_level_1.c:455
void G_math_f_aA_B(float **, float **, float, float **, int, int)
Add two matrices and scale matrix A with the scalar a.
Definition: blas_level_3.c:113
void G_math_sscal(float *, float, int)
Scale vector x with scalar a using the ATLAS routine cblas_dscal.
void G_free_imatrix(int **)
Matrix memory deallocation.
Definition: ialloc.c:99
float G_math_isamax(float *, int)
Compute the maximum norm of vector x using the ATLAS routine cblas_idamax.
int G_math_sv2val(double *, double **, int, int)
Compute singular values when m >> n.
int G_math_cholesky_decomposition(double **, int, int)
cholesky decomposition for symmetric, positiv definite matrices with bandwidth optimization ...
double b
Definition: r_raster.c:39
void G_math_f_asum_norm(float *, float *, int)
Compute the asum norm of vector x.
Definition: blas_level_1.c:328
double * G_alloc_vector(size_t)
Vector matrix memory allocation.
Definition: dalloc.c:41
void G_math_gauss_elimination(double **, double *, int)
Gauss elimination.
int G_math_solver_jacobi(double **, double *, double *, int, int, double, double)
The iterative jacobi solver for quadratic matrices.
void G_free_fvector(float *)
Vector memory deallocation.
Definition: dalloc.c:149
double G_math_rand_gauss(double)
Definition: gauss.c:17
void G_math_print_spmatrix(G_math_spvector **, int)
print the sparse matrix Asp to stdout
long G_math_max_pow2(long n)
Finds least power of 2 >= n
Definition: max_pow2.c:16
int G_math_egvorder(double *, double **, long)
Definition: eigen_tools.c:9
long G_math_min_pow2(long n)
Finds largest power of 2 <= n
Definition: max_pow2.c:44
void G_math_dcopy(double *, double *, int)
Copy vector x to vector y.
void G_math_f_euclid_norm(float *, float *, int)
Compute the euclid norm of vector x.
Definition: blas_level_1.c:296
void G_math_saxpy(float *, float *, float, int)
Scale vector x with scalar a and add it to y.
void G_math_f_aAx_by(float **, float *, float *, float, float, float *, int, int)
Compute the scaled matrix - vector product of matrix A and vectors x and y.
Definition: blas_level_2.c:261
void G_math_f_max_norm(float *, float *, int)
Compute the maximum norm of vector x.
Definition: blas_level_1.c:361
int G_ludcmp(double **, int, int *, double *)
LU decomposition.
Definition: lu.c:18
int G_math_solver_cg(double **, double *, double *, int, int, double)
The iterative conjugate gradients solver for symmetric positive definite matrices.
int G_math_solver_cg_sband(double **, double *, double *, int, int, int, double)
The iterative conjugate gradients solver for symmetric positive definite band matrices.
double G_math_dasum(double *, int)
Compute the absolute sum norm of vector x using the ATLAS routine cblas_dasum.
int G_math_sv2uv(double *, double **, double **, int, double **, int)
Compute the singular value transformation when m >> n.
int G_math_solv(double **, double *, int)
Solve a general linear system A*x = b.
G_math_spvector ** G_math_alloc_spmatrix(int)
Allocate memory for a sparse matrix.
Definition: sparse_matrix.c:58
int G_math_findzc(double conv[], int size, double zc[], double thresh, int num_orients)
Finds locations and orientations of zero crossings.
Definition: findzc.c:56
void G_math_f_AB(float **, float **, float **, int, int, int)
Matrix multiplication.
Definition: blas_level_3.c:215
void G_math_d_x_dyad_y(double *, double *, double **, int, int)
Compute the dyadic product of two vectors. The result is stored in the matrix A.
Definition: blas_level_2.c:115
int G_math_svdval(double *, double **, int, int)
Compute the singular values of a real m by n matrix A.
void G_math_lu_decomposition(double **, double *, int)
lu decomposition
void G_math_dscal(double *, double, int)
Scale vector x with scalar a using the ATLAS routine cblas_dscal.
void G_math_i_ax_by(int *, int *, int *, int, int, int)
Scales vectors x and y with the scalars a and b and adds them.
Definition: blas_level_1.c:607
int G_math_d_A_T(double **A, int rows)
Compute the transposition of matrix A. Matrix A will be overwritten.
Definition: blas_level_2.c:339
int G_math_solver_sparse_jacobi(G_math_spvector **, double *, double *, int, int, double, double)
The iterative jacobi solver for sparse matrices.
int G_math_minv(double **, int)
Invert (in place) a general real matrix A -> Inv(A).
int G_math_solver_gauss(double **, double *, double *, int)
The gauss elimination solver for quardatic matrices.
float G_math_sdot(float *, float *, int)
Compute the dot product of vector x and y using the ATLAS routine cblas_sdot.
void G_math_d_euclid_norm(double *, double *, int)
Compute the euclid norm of vector x.
Definition: blas_level_1.c:80
void G_math_scopy(float *, float *, int)
Copy vector x to vector y.
double ** G_math_sband_matrix_to_matrix(double **, int, int)
Convert a symmetric band matrix into a symmetric matrix.
void G_free_ivector(int *)
Vector memory deallocation.
Definition: ialloc.c:81
void G_math_cholesky_sband_decomposition(double **, double **, int, int)
Cholesky decomposition of a symmetric band matrix.
void G_lubksb(double **a, int n, int *indx, double b[])
LU backward substitution.
Definition: lu.c:104
void G_math_d_asum_norm(double *, double *, int)
Compute the asum norm of vector x.
Definition: blas_level_1.c:112
float G_math_sdsdot(float *, float *, float, int)
Compute the dot product of vector x and y using the ATLAS routine cblas_sdsdot.
void G_free_matrix(double **)
Matrix memory deallocation.
Definition: dalloc.c:169
void G_math_i_euclid_norm(int *, double *, int)
Compute the euclid norm of vector x.
Definition: blas_level_1.c:514
int G_math_solvps(double **, double *, int)
Solve a symmetric positive definite linear system S*x = b.
void G_math_d_Ax(double **, double *, double *, int, int)
Compute the matrix - vector product of matrix A and vector x.
Definition: blas_level_2.c:47
void G_math_solvtd(double *, double *, double *, double *, int)
Solve a tridiagonal linear system M*x = y.
void G_free_fmatrix(float **)
Floating point matrix memory deallocation.
Definition: dalloc.c:190
void G_math_i_asum_norm(int *, double *, int)
Compute the asum norm of vector x.
Definition: blas_level_1.c:546
void G_math_d_aAx_by(double **, double *, double *, double, double, double *, int, int)
Compute the scaled matrix - vector product of matrix double **A and vector x and y.
Definition: blas_level_2.c:179
int fft(int, double *[2], int, int, int)
Fast Fourier Transform for two-dimensional array.
Definition: fft.c:127
float ** G_alloc_fmatrix(int, int)
Floating point matrix memory allocation.
Definition: dalloc.c:104
int G_math_solver_gs(double **, double *, double *, int, int, double, double)
The iterative gauss seidel solver for quadratic matrices.
int G_math_solver_sparse_cg(G_math_spvector **, double *, double *, int, int, double)
The iterative conjugate gradients solver for sparse symmetric positive definite matrices.
int G_math_add_spvector(G_math_spvector **, G_math_spvector *, int)
Adds a sparse vector to a sparse matrix at position row.
Definition: sparse_matrix.c:35
int G_math_complex_mult(double *v1[2], int size1, double *v2[2], int size2, double *v3[2], int size3)
Multiply two complex vectors, point by point.
Definition: mult.c:23
int G_math_solver_sparse_bicgstab(G_math_spvector **, double *, double *, int, int, double)
The iterative biconjugate gradients solver with stabilization for unsymmetric non-definite matrices...
void G_math_Ax_sband(double **A, double *x, double *y, int rows, int bandwidth)
Compute the matrix - vector product of symmetric band matrix A and vector x.