56 for(
size_t i=0;i<N;i++) {
57 if (O2SCL_IX2(A,i,i)==0.0)
return 1;
94 for (j = 0; j < N - 1; j++) {
97 double ajj, max = fabs(O2SCL_IX2(A,j,j));
100 for (i = j + 1; i < N; i++) {
101 double aij = fabs (O2SCL_IX2(A,i,j));
114 temp=O2SCL_IX2(A,j,k);
115 O2SCL_IX2(A,j,k)=O2SCL_IX2(A,i_pivot,k);
116 O2SCL_IX2(A,i_pivot,k)=temp;
122 ajj = O2SCL_IX2(A,j,j);
125 for (i = j + 1; i < N; i++) {
126 double aij = O2SCL_IX2(A,i,j) / ajj;
127 O2SCL_IX2(A,i,j)=aij;
128 for (k = j + 1; k < N; k++) {
129 double aik = O2SCL_IX2(A,i,k);
130 double ajk = O2SCL_IX2(A,j,k);
131 O2SCL_IX2(A,i,k)=aik - aij * ajk;
147 template<
class mat_t,
class vec_t>
148 int LU_svx(
const size_t N,
const mat_t &LU,
152 O2SCL_ERR(
"LU matrix is singular in LU_svx().",
161 o2scl_cblas::o2cblas_Lower,
162 o2scl_cblas::o2cblas_NoTrans,
163 o2scl_cblas::o2cblas_Unit,
168 o2scl_cblas::o2cblas_Upper,
169 o2scl_cblas::o2cblas_NoTrans,
170 o2scl_cblas::o2cblas_NonUnit,
184 template<
class mat_t,
class mat_row_t>
193 for (j = 0; j < N - 1; j++) {
196 double ajj, max = fabs(O2SCL_IX2(A,j,j));
199 for (i = j + 1; i < N; i++) {
200 double aij = fabs (O2SCL_IX2(A,i,j));
213 mat_row_t r2(A,i_pivot);
216 O2SCL_IX(r1,k)=O2SCL_IX(r2,k);
223 ajj = O2SCL_IX2(A,j,j);
226 for (i = j + 1; i < N; i++) {
227 double aij = O2SCL_IX2(A,i,j) / ajj;
228 O2SCL_IX2(A,i,j)=aij;
229 for (k = j + 1; k < N; k++) {
230 double aik = O2SCL_IX2(A,i,k);
231 double ajk = O2SCL_IX2(A,j,k);
232 O2SCL_IX2(A,i,k)=aik - aij * ajk;
247 template<
class mat_t,
class vec_t,
class vec2_t>
249 const vec_t &b, vec2_t &x) {
252 O2SCL_ERR(
"LU matrix is singular in LU_solve().",
270 template<
class mat_t,
class mat2_t,
class vec_t,
class vec2_t,
class vec3_t>
271 int LU_refine(
const size_t N,
const mat_t &A,
const mat2_t &LU,
276 O2SCL_ERR(
"LU matrix is singular in LU_refine().",
283 o2scl_cblas::o2cblas_NoTrans,
284 N,N,1.0,A,x,-1.0,residual);
288 int status=
LU_svx(N,LU,p,residual);
308 template<
class mat_t,
class mat2_t,
class mat_col_t>
313 O2SCL_ERR(
"LU matrix is singular in LU_invert().",
323 for(
size_t j=0;j<N;j++) {
324 if (i==j) O2SCL_IX2(inverse,i,j)=1.0;
325 else O2SCL_IX2(inverse,i,j)=0.0;
329 for (i = 0; i < N; i++) {
330 mat_col_t c(inverse,i);
331 int status_i=
LU_svx(N,LU,p,c);
347 template<
class mat_t>
348 double LU_det(
const size_t N,
const mat_t &LU,
int signum) {
352 double det=(double)signum;
355 det*=O2SCL_IX2(LU,i,i);
369 template<
class mat_t>
375 for (i = 0; i < N; i++) {
376 lndet+=log(fabs(O2SCL_IX2(LU,i,i)));
388 template<
class mat_t>
389 int LU_sgndet(
const size_t N,
const mat_t &LU,
int signum) {
395 double u=O2SCL_IX2(LU,i,i);
void dtrsv(const enum o2cblas_order order, const enum o2cblas_uplo Uplo, const enum o2cblas_transpose TransA, const enum o2cblas_diag Diag, const size_t M, const size_t N, const mat_t &A, vec_t &X)
Compute .
double LU_lndet(const size_t N, const mat_t &LU)
Compute the logarithm of the absolute value of the determinant of a matrix from its LU decomposition...
int LU_svx(const size_t N, const mat_t &LU, const o2scl::permutation &p, vec_t &x)
Solve a linear system after LU decomposition in place.
A class for representing permutations.
int LU_decomp(const size_t N, mat_t &A, o2scl::permutation &p, int &signum)
Compute the LU decomposition of the matrix A.
double LU_det(const size_t N, const mat_t &LU, int signum)
Compute the determinant of a matrix from its LU decomposition.
void dgemv(const enum o2cblas_order order, const enum o2cblas_transpose TransA, const size_t M, const size_t N, const double alpha, const mat_t &A, const vec_t &X, const double beta, vec2_t &Y)
Compute .
void vector_copy(const vec_t &src, vec2_t &dest)
Simple vector copy.
int LU_sgndet(const size_t N, const mat_t &LU, int signum)
Compute the sign of the determinant of a matrix from its LU decomposition.
The namespace for linear algebra classes and functions.
input domain error, e.g sqrt(-1)
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
int LU_decomp_alt(const size_t N, mat_t &A, o2scl::permutation &p, int &signum)
An alternate form of LU decomposition with matrix row objects.
int init()
Initialize permutation to the identity.
void daxpy(const double alpha, const size_t N, const vec_t &X, vec2_t &Y)
Compute .
int LU_refine(const size_t N, const mat_t &A, const mat2_t &LU, const o2scl::permutation &p, const vec_t &b, vec2_t &x, vec3_t &residual)
Refine the solution of a linear system.
int LU_invert(const size_t N, const mat_t &LU, const o2scl::permutation &p, mat2_t &inverse)
Compute the inverse of a matrix from its LU decomposition.
int LU_solve(const size_t N, const mat_t &LU, const o2scl::permutation &p, const vec_t &b, vec2_t &x)
Solve a linear system after LU decomposition.
int diagonal_has_zero(const size_t N, mat_t &A)
Return 1 if at least one of the elements in the diagonal is zero.
int swap(const size_t i, const size_t j)
Swap two elements of a permutation.
int apply(vec_t &v) const
Apply the permutation to a vector.