GRASS GIS 7 Programmer's Manual
7.8.2(2019)-exported
|
Go to the documentation of this file.
23 #include <grass/gis.h>
24 #include <grass/gmath.h>
25 #include <grass/glocale.h>
28 #define COMP_PIVOT 100
45 G_message(_(
"Starting direct gauss elimination solver"));
73 G_message(_(
"Starting direct lu decomposition solver"));
83 #pragma omp for schedule (static) private(i)
84 for (i = 0; i < rows; i++) {
94 #pragma omp for schedule (static) private(i)
95 for (i = 0; i < rows; i++) {
131 G_message(_(
"Starting cholesky decomposition solver"));
134 G_warning(_(
"Unable to solve the linear equation system"));
162 for (k = 0; k < rows - 1; k++) {
163 #pragma omp parallel for schedule (static) private(i, j, tmpval) shared(k, A, b, rows)
164 for (i = k + 1; i < rows; i++) {
165 tmpval = A[i][k] / A[k][k];
166 b[i] =
b[i] - tmpval *
b[k];
167 for (j = k + 1; j < rows; j++) {
168 A[i][j] = A[i][j] - tmpval * A[k][j];
194 for (k = 0; k < rows - 1; k++) {
195 #pragma omp parallel for schedule (static) private(i, j) shared(k, A, rows)
196 for (i = k + 1; i < rows; i++) {
197 A[i][k] = A[i][k] / A[k][k];
198 for (j = k + 1; j < rows; j++) {
199 A[i][j] = A[i][j] - A[i][k] * A[k][j];
223 int i = 0, j = 0, k = 0;
236 for (k = 0; k < rows; k++) {
237 #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k) reduction(+:sum_1)
238 for (j = 0; j < k; j++) {
239 sum_1 += A[k][j] * A[k][j];
242 if (0 > (A[k][k] - sum_1)) {
243 G_warning(
"Matrix is not positive definite. break.");
246 A[k][k] = sqrt(A[k][k] - sum_1);
249 if ((k + bandwidth) > rows) {
253 colsize = k + bandwidth;
256 #pragma omp parallel for schedule (static) private(i, j, sum_2) shared(A, k, sum_1, colsize)
258 for (i = k + 1; i < colsize; i++) {
260 for (j = 0; j < k; j++) {
261 sum_2 += A[i][j] * A[k][j];
263 A[i][k] = (A[i][k] - sum_2) / A[k][k];
268 #pragma omp parallel for schedule (static) private(i, k) shared(A, rows)
269 for (k = 0; k < rows; k++) {
270 for (i = k + 1; i < rows; i++) {
293 for (i = rows - 1; i >= 0; i--) {
294 for (j = i + 1; j < rows; j++) {
295 b[i] =
b[i] - A[i][j] *
x[j];
297 x[i] = (
b[i]) / A[i][i];
319 for (i = 0; i < rows; i++) {
321 for (j = 0; j < i; j++) {
322 tmpval += A[i][j] *
x[j];
324 x[i] = (
b[i] - tmpval) / A[i][i];
void G_math_lu_decomposition(double **A, double *b, int rows)
lu decomposition
int G_math_cholesky_decomposition(double **A, int rows, int bandwidth)
cholesky decomposition for symmetric, positiv definite matrices with bandwidth optimization
void G_math_backward_substitution(double **A, double *x, double *b, int rows)
backward substitution
void G_math_gauss_elimination(double **A, double *b, int rows)
Gauss elimination.
void G_message(const char *msg,...)
Print a message to stderr.
void G_free(void *buf)
Free allocated memory.
int G_math_solver_gauss(double **A, double *x, double *b, int rows)
The gauss elimination solver for quardatic matrices.
double * G_alloc_vector(size_t n)
Vector matrix memory allocation.
void G_math_forward_substitution(double **A, double *x, double *b, int rows)
forward substitution
void G_warning(const char *msg,...)
Print a warning message to stderr.
int G_math_solver_lu(double **A, double *x, double *b, int rows)
The LU solver for quardatic matrices.
int G_math_solver_cholesky(double **A, double *x, double *b, int bandwidth, int rows)
The choleksy decomposition solver for quardatic, symmetric positiv definite matrices.