-
-
Notifications
You must be signed in to change notification settings - Fork 288
/
solvers_direct_cholesky_band.c
208 lines (171 loc) · 5.53 KB
/
solvers_direct_cholesky_band.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <grass/gis.h>
#include <grass/gmath.h>
#include <grass/glocale.h>
/**
* \brief Cholesky decomposition of a symmetric band matrix
*
* \param A (double**) the input symmetric band matrix
* \param T (double**) the resulting lower tringular symmetric band matrix
* \param rows (int) number of rows
* \param bandwidth (int) the bandwidth of the symmetric band matrix
*
* */
void G_math_cholesky_sband_decomposition(double **A, double **T, int rows, int bandwidth)
{
int i, j, k, end;
double sum;
G_debug(2, "G_math_cholesky_sband_decomposition(): n=%d bandwidth=%d", rows, bandwidth);
for (i = 0; i < rows; i++) {
G_percent(i, rows, 9);
/* For j = 0 */
sum = A[i][0];
end = ((bandwidth - 0) < (i + 1) ? (bandwidth - 0) : (i + 1));
for (k = 1; k < end; k++)
sum -= T[i - k][k] * T[i - k][0 + k];
if (sum <= 0.0)
G_fatal_error(_("Decomposition failed at row %i and col %i"), i, 0);
T[i][0] = sqrt(sum);
#pragma omp parallel for schedule (static) private(j, k, end, sum) shared(A, T, i, bandwidth)
for (j = 1; j < bandwidth; j++) {
sum = A[i][j];
end = ((bandwidth - j) < (i + 1) ? (bandwidth - j) : (i + 1));
for (k = 1; k < end; k++)
sum -= T[i - k][k] * T[i - k][j + k];
T[i][j] = sum / T[i][0];
}
}
G_percent(i, rows, 2);
return;
}
/**
* \brief Cholesky symmetric band matrix solver for linear equation systems of type Ax = b
*
* \param A (double**) the input symmetric band matrix
* \param x (double*) the resulting vector, result is written in here
* \param b (double*) the right hand side of Ax = b
* \param rows (int) number of rows
* \param bandwidth (int) the bandwidth of the symmetric band matrix
*
* */
void G_math_solver_cholesky_sband(double **A, double *x, double *b, int rows, int bandwidth)
{
double **T;
T = G_alloc_matrix(rows, bandwidth);
G_math_cholesky_sband_decomposition(A, T, rows, bandwidth); /* T computation */
G_math_cholesky_sband_substitution(T, x, b, rows, bandwidth);
G_free_matrix(T);
return;
}
/**
* \brief Forward and backward substitution of a lower tringular symmetric band matrix of A from system Ax = b
*
* \param T (double**) the lower triangle symmetric band matrix
* \param x (double*) the resulting vector
* \param b (double*) the right hand side of Ax = b
* \param rows (int) number of rows
* \param bandwidth (int) the bandwidth of the symmetric band matrix
*
* */
void G_math_cholesky_sband_substitution(double **T, double *x, double *b, int rows, int bandwidth)
{
int i, j, start, end;
/* Forward substitution */
x[0] = b[0] / T[0][0];
for (i = 1; i < rows; i++) {
x[i] = b[i];
/* start = 0 or i - bandwidth + 1 */
start = ((i - bandwidth + 1) < 0 ? 0 : (i - bandwidth + 1));
/* end = i */
for (j = start; j < i; j++)
x[i] -= T[j][i - j] * x[j];
x[i] = x[i] / T[i][0];
}
/* Backward substitution */
x[rows - 1] = x[rows - 1] / T[rows - 1][0];
for (i = rows - 2; i >= 0; i--) {
/* start = i + 1 */
/* end = rows or bandwidth + i */
end = (rows < (bandwidth + i) ? rows : (bandwidth + i));
for (j = i + 1; j < end; j++)
x[i] -= T[i][j - i] * x[j];
x[i] = x[i] / T[i][0];
}
return;
}
/*--------------------------------------------------------------------------------------*/
/* Tcholetsky matrix invertion */
void G_math_cholesky_sband_invert(double **A, double *invAdiag, int rows, int bandwidth)
{
double **T = NULL;
double *vect = NULL;
int i, j, k, start;
double sum;
T = G_alloc_matrix(rows, bandwidth);
vect = G_alloc_vector(rows);
/* T computation */
G_math_cholesky_sband_decomposition(A, T, rows, bandwidth);
/* T Diagonal invertion */
for (i = 0; i < rows; i++) {
T[i][0] = 1.0 / T[i][0];
}
/* A Diagonal invertion */
for (i = 0; i < rows; i++) {
vect[0] = T[i][0];
invAdiag[i] = vect[0] * vect[0];
for (j = i + 1; j < rows; j++) {
sum = 0.0;
/* start = i or j - bandwidth + 1 */
start = ((j - bandwidth + 1) < i ? i : (j - bandwidth + 1));
/* end = j */
for (k = start; k < j; k++) {
sum -= vect[k - i] * T[k][j - k];
}
vect[j - i] = sum * T[j][0];
invAdiag[i] += vect[j - i] * vect[j - i];
}
}
G_free_matrix(T);
G_free_vector(vect);
return;
}
/*--------------------------------------------------------------------------------------*/
/* Tcholetsky matrix solution and invertion */
void G_math_solver_cholesky_sband_invert(double **A, double *x, double *b, double *invAdiag,
int rows, int bandwidth)
{
double **T = NULL;
double *vect = NULL;
int i, j, k, start;
double sum;
T = G_alloc_matrix(rows, bandwidth);
vect = G_alloc_vector(rows);
/* T computation */
G_math_cholesky_sband_decomposition(A, T, rows, bandwidth);
G_math_cholesky_sband_substitution(T, x, b, rows, bandwidth);
/* T Diagonal invertion */
for (i = 0; i < rows; i++) {
T[i][0] = 1.0 / T[i][0];
}
/* A Diagonal invertion */
for (i = 0; i < rows; i++) {
vect[0] = T[i][0];
invAdiag[i] = vect[0] * vect[0];
for (j = i + 1; j < rows; j++) {
sum = 0.0;
/* start = i or j - bandwidth + 1 */
start = ((j - bandwidth + 1) < i ? i : (j - bandwidth + 1));
/* end = j */
for (k = start; k < j; k++) {
sum -= vect[k - i] * T[k][j - k];
}
vect[j - i] = sum * T[j][0];
invAdiag[i] += vect[j - i] * vect[j - i];
}
}
G_free_matrix(T);
G_free_vector(vect);
return;
}