16 #include "kernels.hpp"
20 #include "kernels.hpp"
21 #include "../general/forall.hpp"
22 #include "../general/table.hpp"
23 #include "../general/globals.hpp"
30 #if defined(_MSC_VER) && (_MSC_VER < 1800)
32 #define copysign _copysign
36 #ifdef MFEM_USE_LAPACK
38 dgemm_(
char *,
char *,
int *,
int *,
int *,
double *,
double *,
39 int *,
double *,
int *,
double *,
double *,
int *);
41 dgetrf_(
int *,
int *,
double *,
int *,
int *,
int *);
43 dgetrs_(
char *,
int *,
int *,
double *,
int *,
int *,
double *,
int *,
int *);
45 dgetri_(
int *N,
double *A,
int *LDA,
int *IPIV,
double *WORK,
46 int *LWORK,
int *INFO);
48 dsyevr_(
char *JOBZ,
char *RANGE,
char *UPLO,
int *N,
double *A,
int *LDA,
49 double *VL,
double *VU,
int *IL,
int *IU,
double *ABSTOL,
int *M,
50 double *W,
double *Z,
int *LDZ,
int *ISUPPZ,
double *WORK,
int *LWORK,
51 int *IWORK,
int *LIWORK,
int *INFO);
53 dsyev_(
char *JOBZ,
char *UPLO,
int *N,
double *A,
int *LDA,
double *W,
54 double *WORK,
int *LWORK,
int *INFO);
56 dsygv_ (
int *ITYPE,
char *JOBZ,
char *UPLO,
int * N,
double *A,
int *LDA,
57 double *B,
int *LDB,
double *W,
double *WORK,
int *LWORK,
int *INFO);
59 dgesvd_(
char *JOBU,
char *JOBVT,
int *M,
int *N,
double *A,
int *LDA,
60 double *S,
double *U,
int *LDU,
double *VT,
int *LDVT,
double *WORK,
61 int *LWORK,
int *INFO);
63 dtrsm_(
char *side,
char *uplo,
char *transa,
char *diag,
int *m,
int *n,
64 double *
alpha,
double *
a,
int *lda,
double *
b,
int *ldb);
66 dggev_(
char *jobvl,
char *jobvr,
int *n,
double *
a,
int *lda,
double *B,
67 int *ldb,
double *alphar,
double *alphai,
double *beta,
double *vl,
68 int * ldvl,
double * vr,
int * ldvr,
double * work,
int * lwork,
int* info);
72 dpotrf_(
char *,
int *,
double *,
int *,
int *);
75 dpotrs_(
char *,
int *,
int *,
double *,
int *,
double *,
int *,
int *);
78 dtrtrs_(
char *,
char*,
char *,
int *,
int *,
double *,
int *,
double *,
int *,
81 dpotri_(
char *,
int *,
double *,
int*,
int *);
97 MFEM_ASSERT(m.data,
"invalid source matrix");
99 std::memcpy(data, m.data,
sizeof(
double)*hw);
105 MFEM_ASSERT(s >= 0,
"invalid DenseMatrix size: " << s);
115 MFEM_ASSERT(m >= 0 && n >= 0,
116 "invalid DenseMatrix size: " << m <<
" x " << n);
117 const int capacity = m*n;
126 :
Matrix(mat.width, mat.height)
128 MFEM_CONTRACT_VAR(ch);
134 for (
int i = 0; i <
height; i++)
136 for (
int j = 0; j <
width; j++)
138 (*this)(i,j) = mat(j,i);
146 MFEM_ASSERT(h >= 0 && w >= 0,
147 "invalid DenseMatrix size: " << h <<
" x " << w);
181 "incompatible dimensions");
183 Mult((
const double *)x, (
double *)y);
189 "incompatible dimensions");
193 for (
int i = 0; i < hw; i++)
195 a += data[i] * m.data[i];
203 double *d_col =
Data();
204 for (
int col = 0; col <
width; col++)
207 for (
int row = 0; row <
height; row++)
209 y_col += x[row]*d_col[row];
219 "incompatible dimensions");
227 "incompatible dimensions");
229 const double *xp = x, *d_col = data;
231 for (
int col = 0; col <
width; col++)
233 double x_col = xp[col];
234 for (
int row = 0; row <
height; row++)
236 yp[row] += x_col*d_col[row];
245 "incompatible dimensions");
247 const double *d_col = data;
248 for (
int col = 0; col <
width; col++)
251 for (
int row = 0; row <
height; row++)
253 y_col += x[row]*d_col[row];
263 "incompatible dimensions");
265 const double *xp = x, *d_col = data;
267 for (
int col = 0; col <
width; col++)
269 const double x_col = a*xp[col];
270 for (
int row = 0; row <
height; row++)
272 yp[row] += x_col*d_col[row];
282 "incompatible dimensions");
284 const double *d_col = data;
285 for (
int col = 0; col <
width; col++)
288 for (
int row = 0; row <
height; row++)
290 y_col += x[row]*d_col[row];
301 for (
int i = 0; i <
height; i++)
304 for (
int j = 0; j <
width; j++)
306 Axi += (*this)(i,j) * x[j];
317 double * it_data = data;
318 for (
int j = 0; j <
width; ++j)
320 for (
int i = 0; i <
height; ++i)
322 *(it_data++) *=
s(i);
330 double * it_data = data;
331 for (
int j = 0; j <
width; ++j)
333 for (
int i = 0; i <
height; ++i)
335 *(it_data++) /=
s(i);
344 double * it_data = data;
345 for (
int j = 0; j <
width; ++j)
348 for (
int i = 0; i <
height; ++i)
358 double * it_data = data;
359 for (
int j = 0; j <
width; ++j)
361 const double sj = 1./
s(j);
362 for (
int i = 0; i <
height; ++i)
374 mfem_error(
"DenseMatrix::SymmetricScaling: dimension mismatch");
377 double * ss =
new double[
width];
380 for (
double * end_s = it_s +
width; it_s != end_s; ++it_s)
382 *(it_ss++) = sqrt(*it_s);
385 double * it_data = data;
386 for (
int j = 0; j <
width; ++j)
388 for (
int i = 0; i <
height; ++i)
390 *(it_data++) *= ss[i]*ss[j];
402 mfem_error(
"DenseMatrix::InvSymmetricScaling: dimension mismatch");
405 double * ss =
new double[
width];
408 for (
double * end_s = it_s +
width; it_s != end_s; ++it_s)
410 *(it_ss++) = 1./sqrt(*it_s);
413 double * it_data = data;
414 for (
int j = 0; j <
width; ++j)
416 for (
int i = 0; i <
height; ++i)
418 *(it_data++) *= ss[i]*ss[j];
430 mfem_error(
"DenseMatrix::Trace() : not a square matrix!");
436 for (
int i = 0; i <
width; i++)
452 "The matrix must be square and "
453 <<
"sized larger than zero to compute the determinant."
454 <<
" Height() = " <<
Height()
455 <<
", Width() = " <<
Width());
463 return data[0] * data[3] - data[1] * data[2];
467 const double *d = data;
469 d[0] * (d[4] * d[8] - d[5] * d[7]) +
470 d[3] * (d[2] * d[7] - d[1] * d[8]) +
471 d[6] * (d[1] * d[5] - d[2] * d[4]);
475 const double *d = data;
477 d[ 0] * (d[ 5] * (d[10] * d[15] - d[11] * d[14]) -
478 d[ 9] * (d[ 6] * d[15] - d[ 7] * d[14]) +
479 d[13] * (d[ 6] * d[11] - d[ 7] * d[10])
481 d[ 4] * (d[ 1] * (d[10] * d[15] - d[11] * d[14]) -
482 d[ 9] * (d[ 2] * d[15] - d[ 3] * d[14]) +
483 d[13] * (d[ 2] * d[11] - d[ 3] * d[10])
485 d[ 8] * (d[ 1] * (d[ 6] * d[15] - d[ 7] * d[14]) -
486 d[ 5] * (d[ 2] * d[15] - d[ 3] * d[14]) +
487 d[13] * (d[ 2] * d[ 7] - d[ 3] * d[ 6])
489 d[12] * (d[ 1] * (d[ 6] * d[11] - d[ 7] * d[10]) -
490 d[ 5] * (d[ 2] * d[11] - d[ 3] * d[10]) +
491 d[ 9] * (d[ 2] * d[ 7] - d[ 3] * d[ 6])
500 return lu_factors.
Det();
515 return sqrt(data[0] * data[0] + data[1] * data[1]);
519 return sqrt(data[0] * data[0] + data[1] * data[1] + data[2] * data[2]);
523 const double *d = data;
524 double E = d[0] * d[0] + d[1] * d[1] + d[2] * d[2];
525 double G = d[3] * d[3] + d[4] * d[4] + d[5] * d[5];
526 double F = d[0] * d[3] + d[1] * d[4] + d[2] * d[5];
527 return sqrt(E * G - F * F);
529 mfem_error(
"DenseMatrix::Weight(): mismatched or unsupported dimensions");
536 for (
int i = 0; i <
s; i++)
538 data[i] = alpha*A[i];
544 for (
int j = 0; j <
Width(); j++)
546 for (
int i = 0; i <
Height(); i++)
548 (*this)(i,j) += c * A(i,j);
556 for (
int i = 0; i <
s; i++)
565 for (
int i = 0; i <
s; i++)
575 for (
int i = 0; i <
s; i++)
587 for (
int i = 0; i < hw; i++)
604 "incompatible matrix sizes.");
610 for (
int j = 0; j <
width; j++)
612 for (
int i = 0; i <
height; i++)
614 (*this)(i, j) -= m(i, j);
624 for (
int i = 0; i <
s; i++)
634 for (
int i = 0; i < hw; i++)
645 mfem_error(
"DenseMatrix::Invert(): dimension mismatch");
649 #ifdef MFEM_USE_LAPACK
650 int *ipiv =
new int[
width];
659 mfem_error(
"DenseMatrix::Invert() : Error in DGETRF");
665 work =
new double[lwork];
671 mfem_error(
"DenseMatrix::Invert() : Error in DGETRI");
677 int c, i, j, n =
Width();
681 for (c = 0; c < n; c++)
683 a = fabs((*
this)(c, c));
685 for (j = c + 1; j < n; j++)
687 b = fabs((*
this)(j, c));
696 mfem_error(
"DenseMatrix::Invert() : singular matrix");
699 for (j = 0; j < n; j++)
701 mfem::Swap<double>((*this)(c, j), (*
this)(i, j));
704 a = (*this)(c, c) = 1.0 / (*
this)(c, c);
705 for (j = 0; j < c; j++)
709 for (j++; j < n; j++)
713 for (i = 0; i < c; i++)
715 (*this)(i, c) = a * (b = -(*
this)(i, c));
716 for (j = 0; j < c; j++)
718 (*this)(i, j) += b * (*
this)(c, j);
720 for (j++; j < n; j++)
722 (*this)(i, j) += b * (*
this)(c, j);
725 for (i++; i < n; i++)
727 (*this)(i, c) = a * (b = -(*
this)(i, c));
728 for (j = 0; j < c; j++)
730 (*this)(i, j) += b * (*
this)(c, j);
732 for (j++; j < n; j++)
734 (*this)(i, j) += b * (*
this)(c, j);
739 for (c = n - 1; c >= 0; c--)
742 for (i = 0; i < n; i++)
744 mfem::Swap<double>((*this)(i, c), (*
this)(i, j));
756 mfem_error(
"DenseMatrix::SquareRootInverse() matrix not square.");
766 for (
int v = 0; v <
Height() ; v++) { (*this)(v,v) = 1.0; }
768 for (
int j = 0; j < 10; j++)
770 for (
int i = 0; i < 10; i++)
785 for (
int v = 0; v <
Height() ; v++) { tmp2(v,v) -= 1.0; }
786 if (tmp2.FNorm() < 1e-10) {
break; }
789 if (tmp2.FNorm() > 1e-10)
791 mfem_error(
"DenseMatrix::SquareRootInverse not converged");
797 for (
int j = 0; j <
Width(); j++)
800 for (
int i = 0; i <
Height(); i++)
802 v[j] += (*this)(i,j)*(*
this)(i,j);
811 const double *d = data;
812 double norm = 0.0, abs_entry;
814 for (
int i = 0; i < hw; i++)
816 abs_entry = fabs(d[i]);
817 if (norm < abs_entry)
829 double max_norm = 0.0, entry, fnorm2;
831 for (i = 0; i < hw; i++)
833 entry = fabs(data[i]);
834 if (entry > max_norm)
842 scale_factor = scaled_fnorm2 = 0.0;
847 for (i = 0; i < hw; i++)
849 entry = data[i] / max_norm;
850 fnorm2 += entry * entry;
853 scale_factor = max_norm;
854 scaled_fnorm2 = fnorm2;
859 #ifdef MFEM_USE_LAPACK
866 double *A =
new double[N*N];
877 int *ISUPPZ =
new int[2*N];
896 double *data = a.
Data();
898 for (
int i = 0; i < hw; i++)
903 dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
904 &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, &QWORK, &LWORK,
905 &QIWORK, &LIWORK, &INFO );
910 WORK =
new double[LWORK];
911 IWORK =
new int[LIWORK];
913 dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
914 &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, WORK, &LWORK,
915 IWORK, &LIWORK, &INFO );
919 mfem::err <<
"dsyevr_Eigensystem(...): DSYEVR error code: "
927 mfem::err <<
"dsyevr_Eigensystem(...):\n"
928 <<
" DSYEVR did not find all eigenvalues "
929 << M <<
"/" << N << endl;
934 mfem_error(
"dsyevr_Eigensystem(...): inf/nan values in W");
938 mfem_error(
"dsyevr_Eigensystem(...): inf/nan values in Z");
941 for (IL = 0; IL < N; IL++)
942 for (IU = 0; IU <= IL; IU++)
945 for (M = 0; M < N; M++)
947 VL += Z[M+IL*N] * Z[M+IU*N];
964 <<
" Z^t Z - I deviation = " << VU
965 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
966 << W[0] <<
", N = " << N << endl;
973 <<
" Z^t Z - I deviation = " << VU
974 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
975 << W[0] <<
", N = " << N << endl;
979 mfem_error(
"dsyevr_Eigensystem(...): ERROR: ...");
982 for (IL = 0; IL < N; IL++)
983 for (IU = 0; IU < N; IU++)
986 for (M = 0; M < N; M++)
988 VL += Z[IL+M*N] * W[M] * Z[IU+M*N];
990 VL = fabs(VL-data[IL+N*IU]);
999 <<
" max matrix deviation = " << VU
1000 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
1001 << W[0] <<
", N = " << N << endl;
1005 mfem_error(
"dsyevr_Eigensystem(...): ERROR: ...");
1014 MFEM_CONTRACT_VAR(a);
1015 MFEM_CONTRACT_VAR(ev);
1016 MFEM_CONTRACT_VAR(evect);
1022 #ifdef MFEM_USE_LAPACK
1034 double *WORK = NULL;
1045 A =
new double[N*N];
1049 double *data = a.
Data();
1050 for (
int i = 0; i < hw; i++)
1055 dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, &QWORK, &LWORK, &INFO);
1057 LWORK = (int) QWORK;
1058 WORK =
new double[LWORK];
1060 dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, WORK, &LWORK, &INFO);
1064 mfem::err <<
"dsyev_Eigensystem: DSYEV error code: " << INFO << endl;
1069 if (evect == NULL) {
delete [] A; }
1071 MFEM_CONTRACT_VAR(a);
1072 MFEM_CONTRACT_VAR(ev);
1073 MFEM_CONTRACT_VAR(evect);
1077 void DenseMatrix::Eigensystem(Vector &ev, DenseMatrix *evect)
1079 #ifdef MFEM_USE_LAPACK
1087 MFEM_CONTRACT_VAR(ev);
1088 MFEM_CONTRACT_VAR(evect);
1089 mfem_error(
"DenseMatrix::Eigensystem: Compiled without LAPACK");
1097 #ifdef MFEM_USE_LAPACK
1110 double *B =
new double[N*N];
1112 double *WORK = NULL;
1123 A =
new double[N*N];
1127 double *a_data = a.
Data();
1128 double *b_data = b.
Data();
1129 for (
int i = 0; i < hw; i++)
1135 dsygv_(&ITYPE, &JOBZ, &UPLO, &N, A, &LDA, B, &LDB, W, &QWORK, &LWORK, &INFO);
1137 LWORK = (int) QWORK;
1138 WORK =
new double[LWORK];
1140 dsygv_(&ITYPE, &JOBZ, &UPLO, &N, A, &LDA, B, &LDB, W, WORK, &LWORK, &INFO);
1144 mfem::err <<
"dsygv_Eigensystem: DSYGV error code: " << INFO << endl;
1150 if (evect == NULL) {
delete [] A; }
1152 MFEM_CONTRACT_VAR(a);
1153 MFEM_CONTRACT_VAR(b);
1154 MFEM_CONTRACT_VAR(ev);
1155 MFEM_CONTRACT_VAR(evect);
1159 void DenseMatrix::Eigensystem(DenseMatrix &
b, Vector &ev,
1162 #ifdef MFEM_USE_LAPACK
1167 MFEM_CONTRACT_VAR(b);
1168 MFEM_CONTRACT_VAR(ev);
1169 MFEM_CONTRACT_VAR(evect);
1170 mfem_error(
"DenseMatrix::Eigensystem(generalized): Compiled without LAPACK");
1176 #ifdef MFEM_USE_LAPACK
1182 double *
a = copy_of_this.data;
1187 double *work = NULL;
1192 dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1193 s, u, &m, vt, &n, &qwork, &lwork, &info);
1195 lwork = (int) qwork;
1196 work =
new double[lwork];
1198 dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1199 s, u, &m, vt, &n, work, &lwork, &info);
1204 mfem::err <<
"DenseMatrix::SingularValues : info = " << info << endl;
1208 MFEM_CONTRACT_VAR(sv);
1210 mfem_error(
"DenseMatrix::SingularValues: Compiled without LAPACK");
1220 for (
int i=0; i < sv.
Size(); ++i)
1232 "The matrix must be square and sized 1, 2, or 3 to compute the"
1234 <<
" Height() = " <<
Height()
1235 <<
", Width() = " <<
Width());
1238 const double *d = data;
1264 const double *d = data;
1282 const double* rp = data + r;
1285 for (
int i = 0; i < n; i++)
1297 double *cp =
Data() + c * m;
1300 for (
int i = 0; i < m; i++)
1314 for (
int i = 0; i <
height; ++i)
1316 d(i) = (*this)(i,i);
1330 for (
int j = 0; j <
width; ++j)
1331 for (
int i = 0; i <
height; ++i)
1333 l(i) += fabs((*
this)(i,j));
1340 for (
int i = 0; i <
height; i++)
1343 for (
int j = 0; j <
width; j++)
1356 for (
int i = 0; i < N; i++)
1360 for (
int i = 0; i < n; i++)
1371 for (i = 0; i < N; i++)
1375 for (i = 0; i < n; i++)
1377 data[i*(n+1)] = diag[i];
1388 for (i = 0; i <
Height(); i++)
1389 for (j = i+1; j <
Width(); j++)
1392 (*this)(i,j) = (*
this)(j,i);
1407 for (
int i = 0; i <
Height(); i++)
1408 for (
int j = 0; j <
Width(); j++)
1410 (*this)(i,j) = A(j,i);
1419 mfem_error(
"DenseMatrix::Symmetrize() : not a square matrix!");
1427 for (
int i = 0; i <
Height(); i++)
1430 for (
int j = 0; j <
Width(); j++)
1433 (*this)(i, j) = 0.0;
1447 mfem_error(
"DenseMatrix::GradToCurl(...): dimension mismatch");
1453 for (
int i = 0; i < n; i++)
1456 double x = (*this)(i,0);
1457 double y = (*this)(i,1);
1470 for (
int i = 0; i < n; i++)
1473 double x = (*this)(i,0);
1474 double y = (*this)(i,1);
1475 double z = (*this)(i,2);
1500 MFEM_VERIFY(
Width() == 2,
1501 "DenseMatrix::GradToVectorCurl2D(...): dimension must be 2")
1505 for (
int i = 0; i < n; i++)
1507 curl(i,0) = (*this)(i,1);
1508 curl(i,1) = -(*this)(i,0);
1514 MFEM_ASSERT(
Width()*
Height() == div.
Size(),
"incompatible Vector 'div'!");
1519 double *ddata = div.
GetData();
1521 for (
int i = 0; i < n; i++)
1531 for (
int j = 0; j <
Width(); j++)
1533 for (
int i = row1; i <= row2; i++)
1535 (*this)(i-row1,j) = A(i,j);
1544 for (
int j = col1; j <= col2; j++)
1546 for (
int i = 0; i <
Height(); i++)
1548 (*this)(i,j-col1) = A(i,j);
1557 for (
int j = 0; j < n; j++)
1559 for (
int i = 0; i < m; i++)
1561 (*this)(i,j) = A(Aro+i,Aco+j);
1568 double *v = A.
Data();
1570 for (
int j = 0; j < A.
Width(); j++)
1572 for (
int i = 0; i < A.
Height(); i++)
1574 (*this)(row_offset+i,col_offset+j) = *(v++);
1581 double *v = A.
Data();
1583 for (
int i = 0; i < A.
Width(); i++)
1585 for (
int j = 0; j < A.
Height(); j++)
1587 (*this)(row_offset+i,col_offset+j) = *(v++);
1593 int row_offset,
int col_offset)
1595 MFEM_VERIFY(row_offset+m <= this->
Height() && col_offset+n <= this->
Width(),
1596 "this DenseMatrix is too small to accommodate the submatrix. "
1597 <<
"row_offset = " << row_offset
1599 <<
", this->Height() = " << this->
Height()
1600 <<
", col_offset = " << col_offset
1602 <<
", this->Width() = " << this->
Width()
1604 MFEM_VERIFY(Aro+m <= A.
Height() && Aco+n <= A.
Width(),
1605 "The A DenseMatrix is too small to accommodate the submatrix. "
1608 <<
", A.Height() = " << A.
Height()
1609 <<
", Aco = " << Aco
1611 <<
", A.Width() = " << A.
Width()
1614 for (
int j = 0; j < n; j++)
1616 for (
int i = 0; i < m; i++)
1618 (*this)(row_offset+i,col_offset+j) = A(Aro+i,Aco+j);
1625 for (
int i = 0; i < n; i++)
1627 for (
int j = i+1; j < n; j++)
1629 (*this)(row_offset+i,col_offset+j) =
1630 (*
this)(row_offset+j,col_offset+i) = 0.0;
1634 for (
int i = 0; i < n; i++)
1636 (*this)(row_offset+i,col_offset+i) = c;
1643 for (
int i = 0; i < n; i++)
1645 for (
int j = i+1; j < n; j++)
1647 (*this)(row_offset+i,col_offset+j) =
1648 (*
this)(row_offset+j,col_offset+i) = 0.0;
1652 for (
int i = 0; i < n; i++)
1654 (*this)(row_offset+i,col_offset+i) = diag[i];
1662 int i, j, i_off = 0, j_off = 0;
1664 for (j = 0; j < A.
Width(); j++)
1671 for (i = 0; i < A.
Height(); i++)
1678 (*this)(i-i_off,j-j_off) = A(i,j);
1694 if (co+aw >
Width() || ro+ah > h)
1696 mfem_error(
"DenseMatrix::AddMatrix(...) 1 : dimension mismatch");
1700 p = data + ro + co * h;
1703 for (
int c = 0; c < aw; c++)
1705 for (
int r = 0; r < ah; r++)
1724 if (co+aw >
Width() || ro+ah > h)
1726 mfem_error(
"DenseMatrix::AddMatrix(...) 2 : dimension mismatch");
1730 p = data + ro + co * h;
1733 for (
int c = 0; c < aw; c++)
1735 for (
int r = 0; r < ah; r++)
1747 int idx_max = idx.
Max();
1748 MFEM_VERIFY(idx.
Min() >=0 && idx_max < this->
height && idx_max < this->
width,
1749 "DenseMatrix::GetSubMatrix: Index out of bounds");
1751 double * adata = A.
Data();
1754 for (
int i = 0; i<k; i++)
1757 for (
int j = 0; j<k; j++)
1760 adata[i+j*k] = this->data[ii+jj*
height];
1768 int k = idx_i.
Size();
1769 int l = idx_j.
Size();
1771 MFEM_VERIFY(idx_i.
Min() >=0 && idx_i.
Max() < this->
height,
1772 "DenseMatrix::GetSubMatrix: Row index out of bounds");
1773 MFEM_VERIFY(idx_j.
Min() >=0 && idx_j.
Max() < this->
width,
1774 "DenseMatrix::GetSubMatrix: Col index out of bounds");
1777 double * adata = A.
Data();
1780 for (
int i = 0; i<k; i++)
1783 for (
int j = 0; j<l; j++)
1786 adata[i+j*k] = this->data[ii+jj*
height];
1793 MFEM_VERIFY(iend >= ibeg,
"DenseMatrix::GetSubMatrix: Inconsistent range");
1794 MFEM_VERIFY(ibeg >=0,
1795 "DenseMatrix::GetSubMatrix: Negative index");
1796 MFEM_VERIFY(iend <= this->
height && iend <= this->
width,
1797 "DenseMatrix::GetSubMatrix: Index bigger than upper bound");
1799 int k = iend - ibeg;
1801 double * adata = A.
Data();
1804 for (
int i = 0; i<k; i++)
1807 for (
int j = 0; j<k; j++)
1810 adata[i+j*k] = this->data[ii+jj*
height];
1818 MFEM_VERIFY(iend >= ibeg,
1819 "DenseMatrix::GetSubMatrix: Inconsistent row range");
1820 MFEM_VERIFY(jend >= jbeg,
1821 "DenseMatrix::GetSubMatrix: Inconsistent col range");
1822 MFEM_VERIFY(ibeg >=0,
1823 "DenseMatrix::GetSubMatrix: Negative row index");
1824 MFEM_VERIFY(jbeg >=0,
1825 "DenseMatrix::GetSubMatrix: Negative row index");
1826 MFEM_VERIFY(iend <= this->
height,
1827 "DenseMatrix::GetSubMatrix: Index bigger than row upper bound");
1828 MFEM_VERIFY(jend <= this->
width,
1829 "DenseMatrix::GetSubMatrix: Index bigger than col upper bound");
1831 int k = iend - ibeg;
1832 int l = jend - jbeg;
1834 double * adata = A.
Data();
1837 for (
int i = 0; i<k; i++)
1840 for (
int j = 0; j<l; j++)
1843 adata[i+j*k] = this->data[ii+jj*
height];
1852 "DenseMatrix::SetSubMatrix:Inconsistent matrix dimensions");
1854 int idx_max = idx.
Max();
1856 MFEM_VERIFY(idx.
Min() >=0,
1857 "DenseMatrix::SetSubMatrix: Negative index");
1858 MFEM_VERIFY(idx_max < this->
height,
1859 "DenseMatrix::SetSubMatrix: Index bigger than row upper bound");
1860 MFEM_VERIFY(idx_max < this->
width,
1861 "DenseMatrix::SetSubMatrix: Index bigger than col upper bound");
1863 double * adata = A.
Data();
1866 for (
int i = 0; i<k; i++)
1869 for (
int j = 0; j<k; j++)
1872 this->data[ii+jj*
height] = adata[i+j*k];
1880 int k = idx_i.
Size();
1881 int l = idx_j.
Size();
1883 "DenseMatrix::SetSubMatrix:Inconsistent matrix dimensions");
1884 MFEM_VERIFY(idx_i.
Min() >=0,
1885 "DenseMatrix::SetSubMatrix: Negative row index");
1886 MFEM_VERIFY(idx_j.
Min() >=0,
1887 "DenseMatrix::SetSubMatrix: Negative col index");
1889 "DenseMatrix::SetSubMatrix: Index bigger than row upper bound");
1890 MFEM_VERIFY(idx_j.
Max() < this->
width,
1891 "DenseMatrix::SetSubMatrix: Index bigger than col upper bound");
1893 double * adata = A.
Data();
1896 for (
int i = 0; i<k; i++)
1899 for (
int j = 0; j<l; j++)
1902 this->data[ii+jj*
height] = adata[i+j*k];
1911 MFEM_VERIFY(A.
Width() == k,
"DenseMatrix::SetSubmatrix: A is not square");
1912 MFEM_VERIFY(ibeg >=0,
1913 "DenseMatrix::SetSubmatrix: Negative index");
1914 MFEM_VERIFY(ibeg + k <= this->
height,
1915 "DenseMatrix::SetSubmatrix: index bigger than row upper bound");
1916 MFEM_VERIFY(ibeg + k <= this->
width,
1917 "DenseMatrix::SetSubmatrix: index bigger than col upper bound");
1919 double * adata = A.
Data();
1922 for (
int i = 0; i<k; i++)
1925 for (
int j = 0; j<k; j++)
1928 this->data[ii+jj*
height] = adata[i+j*k];
1938 MFEM_VERIFY(ibeg>=0,
1939 "DenseMatrix::SetSubmatrix: Negative row index");
1940 MFEM_VERIFY(jbeg>=0,
1941 "DenseMatrix::SetSubmatrix: Negative col index");
1942 MFEM_VERIFY(ibeg + k <= this->
height,
1943 "DenseMatrix::SetSubmatrix: Index bigger than row upper bound");
1944 MFEM_VERIFY(jbeg + l <= this->
width,
1945 "DenseMatrix::SetSubmatrix: Index bigger than col upper bound");
1947 double * adata = A.
Data();
1950 for (
int i = 0; i<k; i++)
1953 for (
int j = 0; j<l; j++)
1956 this->data[ii+jj*
height] = adata[i+j*k];
1965 "DenseMatrix::AddSubMatrix:Inconsistent matrix dimensions");
1967 int idx_max = idx.
Max();
1969 MFEM_VERIFY(idx.
Min() >=0,
"DenseMatrix::AddSubMatrix: Negative index");
1970 MFEM_VERIFY(idx_max < this->
height,
1971 "DenseMatrix::AddSubMatrix: Index bigger than row upper bound");
1972 MFEM_VERIFY(idx_max < this->
width,
1973 "DenseMatrix::AddSubMatrix: Index bigger than col upper bound");
1975 double * adata = A.
Data();
1978 for (
int i = 0; i<k; i++)
1981 for (
int j = 0; j<k; j++)
1984 this->data[ii+jj*
height] += adata[i+j*k];
1992 int k = idx_i.
Size();
1993 int l = idx_j.
Size();
1995 "DenseMatrix::AddSubMatrix:Inconsistent matrix dimensions");
1997 MFEM_VERIFY(idx_i.
Min() >=0,
1998 "DenseMatrix::AddSubMatrix: Negative row index");
1999 MFEM_VERIFY(idx_j.
Min() >=0,
2000 "DenseMatrix::AddSubMatrix: Negative col index");
2002 "DenseMatrix::AddSubMatrix: Index bigger than row upper bound");
2003 MFEM_VERIFY(idx_j.
Max() < this->
width,
2004 "DenseMatrix::AddSubMatrix: Index bigger than col upper bound");
2006 double * adata = A.
Data();
2009 for (
int i = 0; i<k; i++)
2012 for (
int j = 0; j<l; j++)
2015 this->data[ii+jj*
height] += adata[i+j*k];
2023 MFEM_VERIFY(A.
Width() == k,
"DenseMatrix::AddSubmatrix: A is not square");
2025 MFEM_VERIFY(ibeg>=0,
2026 "DenseMatrix::AddSubmatrix: Negative index");
2027 MFEM_VERIFY(ibeg + k <= this->
Height(),
2028 "DenseMatrix::AddSubmatrix: Index bigger than row upper bound");
2029 MFEM_VERIFY(ibeg + k <= this->
Width(),
2030 "DenseMatrix::AddSubmatrix: Index bigger than col upper bound");
2032 double * adata = A.
Data();
2035 for (
int i = 0; i<k; i++)
2038 for (
int j = 0; j<k; j++)
2041 this->data[ii+jj*
height] += adata[i+j*k];
2051 MFEM_VERIFY(ibeg>=0,
2052 "DenseMatrix::AddSubmatrix: Negative row index");
2053 MFEM_VERIFY(jbeg>=0,
2054 "DenseMatrix::AddSubmatrix: Negative col index");
2055 MFEM_VERIFY(ibeg + k <= this->
height,
2056 "DenseMatrix::AddSubmatrix: Index bigger than row upper bound");
2057 MFEM_VERIFY(jbeg + l <= this->
width,
2058 "DenseMatrix::AddSubmatrix: Index bigger than col upper bound");
2060 double * adata = A.
Data();
2063 for (
int i = 0; i<k; i++)
2066 for (
int j = 0; j<l; j++)
2069 this->data[ii+jj*
height] += adata[i+j*k];
2077 double *vdata = v.
GetData() + offset;
2079 for (
int i = 0; i < n; i++)
2081 vdata[i] += data[i];
2088 const double *vdata = v.
GetData() + offset;
2090 for (
int i = 0; i < n; i++)
2103 mfem_error(
"DenseMatrix::AdjustDofDirection(...): dimension mismatch");
2108 for (
int i = 0; i < n-1; i++)
2110 const int s = (dof[i] < 0) ? (-1) : (1);
2111 for (
int j = i+1; j < n; j++)
2113 const int t = (dof[j] < 0) ? (-s) : (
s);
2116 (*this)(i,j) = -(*
this)(i,j);
2117 (*this)(j,i) = -(*
this)(j,i);
2125 for (
int j = 0; j <
Width(); j++)
2127 (*this)(row, j) = value;
2133 for (
int i = 0; i <
Height(); i++)
2135 (*this)(i, col) = value;
2141 MFEM_ASSERT(row !=
nullptr,
"supplied row pointer is null");
2142 for (
int j = 0; j <
Width(); j++)
2144 (*this)(r, j) = row[j];
2156 MFEM_ASSERT(col !=
nullptr,
"supplied column pointer is null");
2157 for (
int i = 0; i <
Height(); i++)
2159 (*this)(i, c) = col[i];
2171 for (
int col = 0; col <
Width(); col++)
2173 for (
int row = 0; row <
Height(); row++)
2175 if (std::abs(
operator()(row,col)) <= eps)
2186 ios::fmtflags old_flags = os.flags();
2188 os << setiosflags(ios::scientific | ios::showpos);
2189 for (
int i = 0; i <
height; i++)
2191 os <<
"[row " << i <<
"]\n";
2192 for (
int j = 0; j <
width; j++)
2195 if (j+1 == width || (j+1) % width_ == 0)
2206 os.flags(old_flags);
2212 ios::fmtflags old_flags = os.flags();
2214 os << setiosflags(ios::scientific | ios::showpos);
2215 for (
int i = 0; i <
height; i++)
2217 for (
int j = 0; j <
width; j++)
2225 os.flags(old_flags);
2231 ios::fmtflags old_flags = os.flags();
2233 os << setiosflags(ios::scientific | ios::showpos);
2234 for (
int j = 0; j <
width; j++)
2236 os <<
"[col " << j <<
"]\n";
2237 for (
int i = 0; i <
height; i++)
2240 if (i+1 == height || (i+1) % width_ == 0)
2251 os.flags(old_flags);
2260 for (
int i = 0; i <
width; i++)
2265 <<
", cond_F = " <<
FNorm()*copy.FNorm() << endl;
2306 MFEM_VERIFY(A.
IsSquare(),
"A must be a square matrix!");
2307 MFEM_ASSERT(A.
NumCols() > 0,
"supplied matrix, A, is empty!");
2308 MFEM_ASSERT(X !=
nullptr,
"supplied vector, X, is null!");
2316 double det = A(0,0);
2317 if (std::abs(det) <= TOL) {
return false; }
2324 double det = A.
Det();
2325 if (std::abs(det) <= TOL) {
return false; }
2327 double invdet = 1. / det;
2332 X[0] = ( A(1,1)*b0 - A(0,1)*b1) * invdet;
2333 X[1] = (-A(1,0)*b0 + A(0,0)*b1) * invdet;
2342 if (!lu.Factor(N,TOL)) {
return false; }
2355 b.
Width() == c.
Height(),
"incompatible dimensions");
2357 #ifdef MFEM_USE_LAPACK
2358 static char transa =
'N', transb =
'N';
2359 static double alpha = 1.0, beta = 0.0;
2362 dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.
Data(), &m,
2363 c.
Data(), &k, &beta, a.
Data(), &m);
2365 const int ah = a.
Height();
2366 const int aw = a.
Width();
2367 const int bw = b.
Width();
2368 double *ad = a.
Data();
2369 const double *bd = b.
Data();
2370 const double *cd = c.
Data();
2379 b.
Width() == c.
Height(),
"incompatible dimensions");
2381 #ifdef MFEM_USE_LAPACK
2382 static char transa =
'N', transb =
'N';
2383 static double beta = 1.0;
2386 dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.
Data(), &m,
2387 c.
Data(), &k, &beta, a.
Data(), &m);
2389 const int ah = a.
Height();
2390 const int aw = a.
Width();
2391 const int bw = b.
Width();
2392 double *ad = a.
Data();
2393 const double *bd = b.
Data();
2394 const double *cd = c.
Data();
2395 for (
int j = 0; j < aw; j++)
2397 for (
int k = 0; k < bw; k++)
2399 for (
int i = 0; i < ah; i++)
2401 ad[i+j*ah] += alpha * bd[i+k*ah] * cd[k+j*bw];
2411 b.
Width() == c.
Height(),
"incompatible dimensions");
2413 #ifdef MFEM_USE_LAPACK
2414 static char transa =
'N', transb =
'N';
2415 static double alpha = 1.0, beta = 1.0;
2418 dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.
Data(), &m,
2419 c.
Data(), &k, &beta, a.
Data(), &m);
2421 const int ah = a.
Height();
2422 const int aw = a.
Width();
2423 const int bw = b.
Width();
2424 double *ad = a.
Data();
2425 const double *bd = b.
Data();
2426 const double *cd = c.
Data();
2427 for (
int j = 0; j < aw; j++)
2429 for (
int k = 0; k < bw; k++)
2431 for (
int i = 0; i < ah; i++)
2433 ad[i+j*ah] += bd[i+k*ah] * cd[k+j*bw];
2445 mfem_error(
"CalcAdjugate(...): unsupported dimensions");
2449 mfem_error(
"CalcAdjugate(...): dimension mismatch");
2455 const double *d = a.
Data();
2456 double *ad = adja.
Data();
2471 e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
2472 g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
2473 f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
2475 ad[0] = d[0]*g - d[3]*
f;
2476 ad[1] = d[3]*e - d[0]*
f;
2477 ad[2] = d[1]*g - d[4]*
f;
2478 ad[3] = d[4]*e - d[1]*
f;
2479 ad[4] = d[2]*g - d[5]*
f;
2480 ad[5] = d[5]*e - d[2]*
f;
2489 else if (a.
Width() == 2)
2492 adja(0,1) = -
a(0,1);
2493 adja(1,0) = -
a(1,0);
2498 adja(0,0) =
a(1,1)*
a(2,2)-
a(1,2)*
a(2,1);
2499 adja(0,1) =
a(0,2)*
a(2,1)-
a(0,1)*
a(2,2);
2500 adja(0,2) =
a(0,1)*
a(1,2)-
a(0,2)*
a(1,1);
2502 adja(1,0) =
a(1,2)*
a(2,0)-
a(1,0)*
a(2,2);
2503 adja(1,1) =
a(0,0)*
a(2,2)-
a(0,2)*
a(2,0);
2504 adja(1,2) =
a(0,2)*
a(1,0)-
a(0,0)*
a(1,2);
2506 adja(2,0) =
a(1,0)*
a(2,1)-
a(1,1)*
a(2,0);
2507 adja(2,1) =
a(0,1)*
a(2,0)-
a(0,0)*
a(2,1);
2508 adja(2,2) =
a(0,0)*
a(1,1)-
a(0,1)*
a(1,0);
2518 mfem_error(
"CalcAdjugateTranspose(...): dimension mismatch");
2525 else if (a.
Width() == 2)
2527 adjat(0,0) =
a(1,1);
2528 adjat(1,0) = -
a(0,1);
2529 adjat(0,1) = -
a(1,0);
2530 adjat(1,1) =
a(0,0);
2534 adjat(0,0) =
a(1,1)*
a(2,2)-
a(1,2)*
a(2,1);
2535 adjat(1,0) =
a(0,2)*
a(2,1)-
a(0,1)*
a(2,2);
2536 adjat(2,0) =
a(0,1)*
a(1,2)-
a(0,2)*
a(1,1);
2538 adjat(0,1) =
a(1,2)*
a(2,0)-
a(1,0)*
a(2,2);
2539 adjat(1,1) =
a(0,0)*
a(2,2)-
a(0,2)*
a(2,0);
2540 adjat(2,1) =
a(0,2)*
a(1,0)-
a(0,0)*
a(1,2);
2542 adjat(0,2) =
a(1,0)*
a(2,1)-
a(1,1)*
a(2,0);
2543 adjat(1,2) =
a(0,1)*
a(2,0)-
a(0,0)*
a(2,1);
2544 adjat(2,2) =
a(0,0)*
a(1,1)-
a(0,1)*
a(1,0);
2551 MFEM_ASSERT(inva.
Height() == a.
Width(),
"incorrect dimensions");
2552 MFEM_ASSERT(inva.
Width() == a.
Height(),
"incorrect dimensions");
2558 const double *d = a.
Data();
2559 double *
id = inva.
Data();
2562 t = 1.0 / (d[0]*d[0] + d[1]*d[1]);
2570 t = 1.0 / (d[0]*d[0] + d[1]*d[1] + d[2]*d[2]);
2578 e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
2579 g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
2580 f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
2581 t = 1.0 / (e*g - f*
f);
2582 e *=
t; g *=
t; f *=
t;
2584 id[0] = d[0]*g - d[3]*
f;
2585 id[1] = d[3]*e - d[0]*
f;
2586 id[2] = d[1]*g - d[4]*
f;
2587 id[3] = d[4]*e - d[1]*
f;
2588 id[4] = d[2]*g - d[5]*
f;
2589 id[5] = d[5]*e - d[2]*
f;
2597 MFEM_ASSERT(std::abs(t) > 1.0e-14 * pow(a.FNorm()/a.
Width(), a.
Width()),
2598 "singular matrix!");
2604 inva(0,0) = 1.0 / a.
Det();
2607 kernels::CalcInverse<2>(a.
Data(), inva.
Data());
2610 kernels::CalcInverse<3>(a.
Data(), inva.
Data());
2621 mfem_error(
"CalcInverseTranspose(...): dimension mismatch");
2625 double t = 1. / a.
Det() ;
2630 inva(0,0) = 1.0 /
a(0,0);
2633 inva(0,0) =
a(1,1) *
t ;
2634 inva(1,0) = -
a(0,1) *
t ;
2635 inva(0,1) = -
a(1,0) *
t ;
2636 inva(1,1) =
a(0,0) *
t ;
2639 inva(0,0) = (
a(1,1)*
a(2,2)-
a(1,2)*
a(2,1))*t;
2640 inva(1,0) = (
a(0,2)*
a(2,1)-
a(0,1)*
a(2,2))*t;
2641 inva(2,0) = (
a(0,1)*
a(1,2)-
a(0,2)*
a(1,1))*t;
2643 inva(0,1) = (
a(1,2)*
a(2,0)-
a(1,0)*
a(2,2))*t;
2644 inva(1,1) = (
a(0,0)*
a(2,2)-
a(0,2)*
a(2,0))*t;
2645 inva(2,1) = (
a(0,2)*
a(1,0)-
a(0,0)*
a(1,2))*t;
2647 inva(0,2) = (
a(1,0)*
a(2,1)-
a(1,1)*
a(2,0))*t;
2648 inva(1,2) = (
a(0,1)*
a(2,0)-
a(0,0)*
a(2,1))*t;
2649 inva(2,2) = (
a(0,0)*
a(1,1)-
a(0,1)*
a(1,0))*t;
2659 "Matrix must be 3x2 or 2x1, "
2660 <<
"and the Vector must be sized with the rows. "
2661 <<
" J.Height() = " << J.
Height()
2662 <<
", J.Width() = " << J.
Width()
2663 <<
", n.Size() = " << n.
Size()
2666 const double *d = J.
Data();
2674 n(0) = d[1]*d[5] - d[2]*d[4];
2675 n(1) = d[2]*d[3] - d[0]*d[5];
2676 n(2) = d[0]*d[4] - d[1]*d[3];
2682 const int height = a.
Height();
2683 const int width = a.
Width();
2684 for (
int i = 0; i < height; i++)
2686 for (
int j = 0; j <= i; j++)
2689 for (
int k = 0; k < width; k++)
2691 temp +=
a(i,k) *
a(j,k);
2693 aat(j,i) = aat(i,j) = temp;
2700 for (
int i = 0; i < A.
Height(); i++)
2702 for (
int j = 0; j < i; j++)
2705 for (
int k = 0; k < A.
Width(); k++)
2707 t += D(k) * A(i, k) * A(j, k);
2715 for (
int i = 0; i < A.
Height(); i++)
2718 for (
int k = 0; k < A.
Width(); k++)
2720 t += D(k) * A(i, k) * A(i, k);
2728 for (
int i = 0; i < A.
Height(); i++)
2730 for (
int j = 0; j <= i; j++)
2733 for (
int k = 0; k < A.
Width(); k++)
2735 t += D(k) * A(i, k) * A(j, k);
2737 ADAt(j, i) = ADAt(i, j) =
t;
2748 mfem_error(
"MultABt(...): dimension mismatch");
2752 #ifdef MFEM_USE_LAPACK
2753 static char transa =
'N', transb =
'T';
2754 static double alpha = 1.0, beta = 0.0;
2757 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &m,
2758 B.
Data(), &n, &beta, ABt.
Data(), &m);
2760 const int ah = A.
Height();
2761 const int bh = B.
Height();
2762 const int aw = A.
Width();
2763 const double *ad = A.
Data();
2764 const double *bd = B.
Data();
2765 double *cd = ABt.
Data();
2769 const int ah = A.
Height();
2770 const int bh = B.
Height();
2771 const int aw = A.
Width();
2772 const double *ad = A.
Data();
2773 const double *bd = B.
Data();
2774 double *cd = ABt.
Data();
2776 for (
int j = 0; j < bh; j++)
2777 for (
int i = 0; i < ah; i++)
2780 const double *ap = ad + i;
2781 const double *bp = bd + j;
2782 for (
int k = 0; k < aw; k++)
2794 for (i = 0; i < A.
Height(); i++)
2795 for (j = 0; j < B.
Height(); j++)
2798 for (k = 0; k < A.
Width(); k++)
2800 d += A(i, k) * B(j, k);
2814 mfem_error(
"MultADBt(...): dimension mismatch");
2818 const int ah = A.
Height();
2819 const int bh = B.
Height();
2820 const int aw = A.
Width();
2821 const double *ad = A.
Data();
2822 const double *bd = B.
Data();
2823 const double *dd = D.
GetData();
2824 double *cd = ADBt.
Data();
2826 for (
int i = 0,
s = ah*bh; i <
s; i++)
2830 for (
int k = 0; k < aw; k++)
2833 for (
int j = 0; j < bh; j++)
2835 const double dk_bjk = dd[k] * bd[j];
2836 for (
int i = 0; i < ah; i++)
2838 cp[i] += ad[i] * dk_bjk;
2853 mfem_error(
"AddMultABt(...): dimension mismatch");
2857 #ifdef MFEM_USE_LAPACK
2858 static char transa =
'N', transb =
'T';
2859 static double alpha = 1.0, beta = 1.0;
2862 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &m,
2863 B.
Data(), &n, &beta, ABt.
Data(), &m);
2865 const int ah = A.
Height();
2866 const int bh = B.
Height();
2867 const int aw = A.
Width();
2868 const double *ad = A.
Data();
2869 const double *bd = B.
Data();
2870 double *cd = ABt.
Data();
2872 for (
int k = 0; k < aw; k++)
2875 for (
int j = 0; j < bh; j++)
2877 const double bjk = bd[j];
2878 for (
int i = 0; i < ah; i++)
2880 cp[i] += ad[i] * bjk;
2891 for (i = 0; i < A.
Height(); i++)
2892 for (j = 0; j < B.
Height(); j++)
2895 for (k = 0; k < A.
Width(); k++)
2897 d += A(i, k) * B(j, k);
2911 mfem_error(
"AddMultADBt(...): dimension mismatch");
2915 const int ah = A.
Height();
2916 const int bh = B.
Height();
2917 const int aw = A.
Width();
2918 const double *ad = A.
Data();
2919 const double *bd = B.
Data();
2920 const double *dd = D.
GetData();
2921 double *cd = ADBt.
Data();
2923 for (
int k = 0; k < aw; k++)
2926 for (
int j = 0; j < bh; j++)
2928 const double dk_bjk = dd[k] * bd[j];
2929 for (
int i = 0; i < ah; i++)
2931 cp[i] += ad[i] * dk_bjk;
2947 mfem_error(
"AddMult_a_ABt(...): dimension mismatch");
2951 #ifdef MFEM_USE_LAPACK
2952 static char transa =
'N', transb =
'T';
2954 static double beta = 1.0;
2957 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &m,
2958 B.
Data(), &n, &beta, ABt.
Data(), &m);
2960 const int ah = A.
Height();
2961 const int bh = B.
Height();
2962 const int aw = A.
Width();
2963 const double *ad = A.
Data();
2964 const double *bd = B.
Data();
2965 double *cd = ABt.
Data();
2967 for (
int k = 0; k < aw; k++)
2970 for (
int j = 0; j < bh; j++)
2972 const double bjk = a * bd[j];
2973 for (
int i = 0; i < ah; i++)
2975 cp[i] += ad[i] * bjk;
2986 for (i = 0; i < A.
Height(); i++)
2987 for (j = 0; j < B.
Height(); j++)
2990 for (k = 0; k < A.
Width(); k++)
2992 d += A(i, k) * B(j, k);
3005 mfem_error(
"MultAtB(...): dimension mismatch");
3009 #ifdef MFEM_USE_LAPACK
3010 static char transa =
'T', transb =
'N';
3011 static double alpha = 1.0, beta = 0.0;
3014 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &k,
3015 B.
Data(), &k, &beta, AtB.
Data(), &m);
3017 const int ah = A.
Height();
3018 const int aw = A.
Width();
3019 const int bw = B.
Width();
3020 const double *ad = A.
Data();
3021 const double *bd = B.
Data();
3022 double *cd = AtB.
Data();
3024 for (
int j = 0; j < bw; j++)
3026 const double *ap = ad;
3027 for (
int i = 0; i < aw; i++)
3030 for (
int k = 0; k < ah; k++)
3043 for (i = 0; i < A.
Width(); i++)
3044 for (j = 0; j < B.
Width(); j++)
3047 for (k = 0; k < A.
Height(); k++)
3049 d += A(k, i) * B(k, j);
3060 for (
int i = 0; i < A.
Height(); i++)
3062 for (
int j = 0; j < i; j++)
3065 for (
int k = 0; k < A.
Width(); k++)
3067 d += A(i,k) * A(j,k);
3069 AAt(i, j) += (d *=
a);
3073 for (
int k = 0; k < A.
Width(); k++)
3075 d += A(i,k) * A(i,k);
3083 for (
int i = 0; i < A.
Height(); i++)
3085 for (
int j = 0; j <= i; j++)
3088 for (
int k = 0; k < A.
Width(); k++)
3090 d += A(i,k) * A(j,k);
3092 AAt(i, j) = AAt(j, i) = a * d;
3099 for (
int i = 0; i < v.
Size(); i++)
3101 for (
int j = 0; j <= i; j++)
3103 vvt(i,j) = vvt(j,i) = v(i) * v(j);
3113 mfem_error(
"MultVWt(...): dimension mismatch");
3117 for (
int i = 0; i < v.
Size(); i++)
3119 const double vi = v(i);
3120 for (
int j = 0; j < w.
Size(); j++)
3122 VWt(i, j) = vi * w(j);
3129 const int m = v.
Size(), n = w.
Size();
3134 mfem_error(
"AddMultVWt(...): dimension mismatch");
3138 for (
int i = 0; i < m; i++)
3140 const double vi = v(i);
3141 for (
int j = 0; j < n; j++)
3143 VWt(i, j) += vi * w(j);
3150 const int n = v.
Size();
3155 mfem_error(
"AddMultVVt(...): dimension mismatch");
3159 for (
int i = 0; i < n; i++)
3161 const double vi = v(i);
3162 for (
int j = 0; j < i; j++)
3164 const double vivj = vi * v(j);
3168 VVt(i, i) += vi * vi;
3175 const int m = v.
Size(), n = w.
Size();
3180 mfem_error(
"AddMult_a_VWt(...): dimension mismatch");
3184 for (
int j = 0; j < n; j++)
3186 const double awj = a * w(j);
3187 for (
int i = 0; i < m; i++)
3189 VWt(i, j) += v(i) * awj;
3197 "incompatible dimensions!");
3199 const int n = v.
Size();
3200 for (
int i = 0; i < n; i++)
3202 double avi = a * v(i);
3203 for (
int j = 0; j < i; j++)
3205 const double avivj = avi * v(j);
3209 VVt(i, i) += avi * v(i);
3232 #ifdef MFEM_USE_LAPACK
3238 double *data_ptr = this->
data;
3239 for (
int i = 0; i < m; i++)
3244 double a = std::abs(data_ptr[piv+i*m]);
3245 for (
int j = i+1; j < m; j++)
3247 const double b = std::abs(data_ptr[j+i*m]);
3258 for (
int j = 0; j < m; j++)
3260 mfem::Swap<double>(data_ptr[i+j*m], data_ptr[piv+j*m]);
3265 if (abs(data_ptr[i + i*m]) <= TOL)
3270 const double a_ii_inv = 1.0 / data_ptr[i+i*m];
3271 for (
int j = i+1; j < m; j++)
3273 data_ptr[j+i*m] *= a_ii_inv;
3275 for (
int k = i+1; k < m; k++)
3277 const double a_ik = data_ptr[i+k*m];
3278 for (
int j = i+1; j < m; j++)
3280 data_ptr[j+k*m] -= a_ik * data_ptr[j+i*m];
3292 for (
int i=0; i<m; i++)
3296 det *= -
data[m * i + i];
3300 det *=
data[m * i + i];
3309 for (
int k = 0; k < n; k++)
3312 for (
int i = 0; i < m; i++)
3314 double x_i = x[i] *
data[i+i*m];
3315 for (
int j = i+1; j < m; j++)
3317 x_i += x[j] * data[i+j*m];
3322 for (
int i = m-1; i >= 0; i--)
3325 for (
int j = 0; j < i; j++)
3327 x_i += x[j] *
data[i+j*m];
3332 for (
int i = m-1; i >= 0; i--)
3343 for (
int k = 0; k < n; k++)
3346 for (
int i = 0; i < m; i++)
3351 for (
int j = 0; j < m; j++)
3353 const double x_j = x[j];
3354 for (
int i = j+1; i < m; i++)
3356 x[i] -=
data[i+j*m] * x_j;
3367 for (
int k = 0; k < n; k++)
3369 for (
int j = m-1; j >= 0; j--)
3371 const double x_j = ( x[j] /=
data[j+j*m] );
3372 for (
int i = 0; i < j; i++)
3374 x[i] -=
data[i+j*m] * x_j;
3383 #ifdef MFEM_USE_LAPACK
3386 if (m > 0 && n > 0) {
dgetrs_(&trans, &m, &n,
data, &m,
ipiv, X, &m, &info); }
3387 MFEM_VERIFY(!info,
"LAPACK: error in DGETRS");
3398 #ifdef MFEM_USE_LAPACK
3399 char n_ch =
'N', side =
'R', u_ch =
'U', l_ch =
'L';
3403 dtrsm_(&side,&u_ch,&n_ch,&n_ch,&n,&m,&alpha,
data,&m,X,&n);
3404 dtrsm_(&side,&l_ch,&n_ch,&u_ch,&n,&m,&alpha,
data,&m,X,&n);
3410 for (
int k = 0; k < n; k++)
3412 for (
int j = 0; j < m; j++)
3414 const double x_j = ( x[j*n] /=
data[j+j*m]);
3415 for (
int i = j+1; i < m; i++)
3417 x[i*n] -=
data[j + i*m] * x_j;
3425 for (
int k = 0; k < n; k++)
3427 for (
int j = m-1; j >= 0; j--)
3429 const double x_j = x[j*n];
3430 for (
int i = 0; i < j; i++)
3432 x[i*n] -=
data[j + i*m] * x_j;
3440 for (
int k = 0; k < n; k++)
3442 for (
int i = m-1; i >= 0; --i)
3455 for (
int k = 0; k < m; k++)
3457 const double minus_x_k = -( x[k] = 1.0/
data[k+k*m] );
3458 for (
int i = 0; i < k; i++)
3460 x[i] =
data[i+k*m] * minus_x_k;
3462 for (
int j = k-1; j >= 0; j--)
3464 const double x_j = ( x[j] /=
data[j+j*m] );
3465 for (
int i = 0; i < j; i++)
3467 x[i] -=
data[i+j*m] * x_j;
3475 for (
int j = 0; j < k; j++)
3477 const double minus_L_kj = -
data[k+j*m];
3478 for (
int i = 0; i <= j; i++)
3480 X[i+j*m] += X[i+k*m] * minus_L_kj;
3482 for (
int i = j+1; i < m; i++)
3484 X[i+j*m] = X[i+k*m] * minus_L_kj;
3488 for (
int k = m-2; k >= 0; k--)
3490 for (
int j = 0; j < k; j++)
3492 const double L_kj =
data[k+j*m];
3493 for (
int i = 0; i < m; i++)
3495 X[i+j*m] -= X[i+k*m] * L_kj;
3500 for (
int k = m-1; k >= 0; k--)
3505 for (
int i = 0; i < m; i++)
3507 Swap<double>(X[i+k*m], X[i+piv_k*m]);
3514 const double *X1,
double *X2)
3517 for (
int k = 0; k < r; k++)
3519 for (
int j = 0; j < m; j++)
3521 const double x1_jk = X1[j+k*m];
3522 for (
int i = 0; i < n; i++)
3524 X2[i+k*n] -= A21[i+j*n] * x1_jk;
3531 int m,
int n,
double *A12,
double *A21,
double *A22)
const
3536 for (
int j = 0; j < m; j++)
3538 const double u_jj_inv = 1.0/
data[j+j*m];
3539 for (
int i = 0; i < n; i++)
3541 A21[i+j*n] *= u_jj_inv;
3543 for (
int k = j+1; k < m; k++)
3545 const double u_jk =
data[j+k*m];
3546 for (
int i = 0; i < n; i++)
3548 A21[i+k*n] -= A21[i+j*n] * u_jk;
3553 SubMult(m, n, n, A21, A12, A22);
3557 double *B1,
double *B2)
const
3562 SubMult(m, n, r, L21, B1, B2);
3566 const double *X2,
double *Y1)
const
3569 SubMult(n, m, r, U12, X2, Y1);
3577 #ifdef MFEM_USE_LAPACK
3580 MFEM_VERIFY(
data,
"Matrix data not set");
3585 for (
int j = 0; j<m; j++)
3588 for (
int k = 0; k<j; k++)
3593 MFEM_VERIFY(
data[j+j*m] - a > 0.,
3594 "CholeskyFactors::Factor: The matrix is not SPD");
3596 data[j+j*m] = std::sqrt(
data[j+j*m] - a);
3598 if (
data[j + j*m] <= TOL)
3603 for (
int i = j+1; i<m; i++)
3606 for (
int k = 0; k<j; k++)
3620 for (
int i=0; i<m; i++)
3622 det *=
data[i + i*m];
3631 for (
int k = 0; k < n; k++)
3633 for (
int j = m-1; j >= 0; j--)
3635 double x_j = x[j] *
data[j+j*m];
3636 for (
int i = 0; i < j; i++)
3638 x_j += x[i] * data[j+i*m];
3649 for (
int k = 0; k < n; k++)
3651 for (
int i = 0; i < m; i++)
3653 double x_i = x[i] *
data[i+i*m];
3654 for (
int j = i+1; j < m; j++)
3656 x_i += x[j] * data[j+i*m];
3667 #ifdef MFEM_USE_LAPACK
3673 dtrtrs_(&uplo, &trans, &diag, &m, &n,
data, &m, X, &m, &info);
3674 MFEM_VERIFY(!info,
"CholeskyFactors:LSolve:: info");
3678 for (
int k = 0; k < n; k++)
3681 for (
int j = 0; j < m; j++)
3683 const double x_j = (x[j] /=
data[j+j*m]);
3684 for (
int i = j+1; i < m; i++)
3686 x[i] -=
data[i+j*m] * x_j;
3696 #ifdef MFEM_USE_LAPACK
3703 dtrtrs_(&uplo, &trans, &diag, &m, &n,
data, &m, X, &m, &info);
3704 MFEM_VERIFY(!info,
"CholeskyFactors:USolve:: info");
3709 for (
int k = 0; k < n; k++)
3711 for (
int j = m-1; j >= 0; j--)
3713 const double x_j = ( x[j] /=
data[j+j*m] );
3714 for (
int i = 0; i < j; i++)
3716 x[i] -=
data[j+i*m] * x_j;
3726 #ifdef MFEM_USE_LAPACK
3730 MFEM_VERIFY(!info,
"CholeskyFactors:Solve:: info");
3740 #ifdef MFEM_USE_LAPACK
3750 dtrsm_(&side,&uplo,&transt,&diag,&n,&m,&alpha,
data,&m,X,&n);
3751 dtrsm_(&side,&uplo,&trans,&diag,&n,&m,&alpha,
data,&m,X,&n);
3756 for (
int k = 0; k < n; k++)
3758 for (
int j = 0; j < m; j++)
3760 const double x_j = ( x[j*n] /=
data[j+j*m]);
3761 for (
int i = j+1; i < m; i++)
3763 x[i*n] -=
data[i + j*m] * x_j;
3770 for (
int k = 0; k < n; k++)
3772 for (
int j = m-1; j >= 0; j--)
3774 const double x_j = (x[j*n] /=
data[j+j*m]);
3775 for (
int i = 0; i < j; i++)
3777 x[i*n] -=
data[j + i*m] * x_j;
3788 #ifdef MFEM_USE_LAPACK
3790 for (
int i = 0; i<m; i++)
3792 for (
int j = i; j<m; j++)
3794 X[j+i*m] =
data[j+i*m];
3799 dpotri_(&uplo, &m, X, &m, &info);
3800 MFEM_VERIFY(!info,
"CholeskyFactors:GetInverseMatrix:: info");
3802 for (
int i = 0; i<m; i++)
3804 for (
int j = i+1; j<m; j++)
3806 X[i+j*m] = X[j+i*m];
3811 for (
int k = 0; k<m; k++)
3813 X[k+k*m] = 1./
data[k+k*m];
3814 for (
int i = k+1; i < m; i++)
3817 for (
int j=k; j<i; j++)
3819 s -=
data[i+j*m] * X[j+k*m]/
data[i+i*m];
3824 for (
int i = 0; i < m; i++)
3826 for (
int j = i; j < m; j++)
3829 for (
int k=j; k<m; k++)
3831 s += X[k+i*m] * X[k+j*m];
3833 X[i+j*m] = X[j+i*m] =
s;
3840 void DenseMatrixInverse::Init(
int m)
3848 factors =
new LUFactors();
3852 factors->
data =
new double[m*m];
3855 dynamic_cast<LUFactors *
>(factors)->ipiv =
new int[m];
3864 MFEM_ASSERT(
height ==
width,
"not a square matrix");
3873 MFEM_ASSERT(
height ==
width,
"not a square matrix");
3880 MFEM_ASSERT(a,
"DenseMatrix is not given");
3881 const double *adata = a->data;
3883 for (
int i = 0; i <
s; i++)
3885 factors->
data[i] = adata[i];
3898 MFEM_VERIFY(mat.
height == mat.
width,
"DenseMatrix is not square!");
3902 if (own_data) {
delete [] factors->
data; }
3908 if (own_data) {
delete [] lu->
ipiv; }
3920 MFEM_VERIFY(p != NULL,
"Operator is not a DenseMatrix!");
3926 for (
int row = 0; row <
height; row++)
3949 for (
int i = 0; i <
width; i++)
3960 delete [] factors->
data;
3963 delete []
dynamic_cast<LUFactors *
>(factors)->ipiv;
3969 #ifdef MFEM_USE_LAPACK
3984 &qwork, &lwork, &info);
3986 lwork = (int) qwork;
3987 work =
new double[lwork];
3992 : mat(other.mat), EVal(other.EVal), EVect(other.EVect), ev(NULL, other.n),
3997 lwork = other.lwork;
3999 work =
new double[lwork];
4005 if (mat.
Width() != n)
4007 mfem_error(
"DenseMatrixEigensystem::Eval(): dimension mismatch");
4013 work, &lwork, &info);
4017 mfem::err <<
"DenseMatrixEigensystem::Eval(): DSYEV error code: "
4031 bool left_eigen_vectors,
4032 bool right_eigen_vectors)
4035 MFEM_VERIFY(A.
Height() == A.
Width(),
"A has to be a square matrix");
4036 MFEM_VERIFY(B.
Height() == B.
Width(),
"B has to be a square matrix");
4038 MFEM_VERIFY(B.
Height() == n,
"A and B dimension mismatch");
4044 if (left_eigen_vectors)
4049 if (right_eigen_vectors)
4058 alphar =
new double[n];
4059 alphai =
new double[n];
4060 beta =
new double[n];
4062 int nl = max(1,Vl.
Height());
4063 int nr = max(1,Vr.
Height());
4065 dggev_(&jobvl,&jobvr,&n,A_copy.
Data(),&n,B_copy.
Data(),&n,alphar,
4066 alphai, beta, Vl.
Data(), &nl, Vr.
Data(), &nr,
4067 &qwork, &lwork, &info);
4069 lwork = (int) qwork;
4070 work =
new double[lwork];
4075 int nl = max(1,Vl.
Height());
4076 int nr = max(1,Vr.
Height());
4080 dggev_(&jobvl,&jobvr,&n,A_copy.Data(),&n,B_copy.
Data(),&n,alphar,
4081 alphai, beta, Vl.
Data(), &nl, Vr.
Data(), &nr,
4082 work, &lwork, &info);
4086 mfem::err <<
"DenseMatrixGeneralizedEigensystem::Eval(): DGGEV error code: "
4092 for (
int i = 0; i<n; i++)
4096 evalues_r(i) = alphar[i]/beta[i];
4097 evalues_i(i) = alphai[i]/beta[i];
4116 bool left_singular_vectors,
4117 bool right_singular_vectors)
4121 jobu = (left_singular_vectors)?
'S' :
'N';
4122 jobvt = (right_singular_vectors)?
'S' :
'N';
4127 bool left_singular_vectors,
4128 bool right_singular_vectors)
4132 jobu = (left_singular_vectors)?
'S' :
'N';
4133 jobvt = (right_singular_vectors)?
'S' :
'N';
4137 void DenseMatrixSVD::Init()
4143 NULL, &n, &qwork, &lwork, &info);
4145 lwork = (int) qwork;
4146 work =
new double[lwork];
4157 double * datau =
nullptr;
4158 double * datavt =
nullptr;
4171 datavt, &n, work, &lwork, &info);
4175 mfem::err <<
"DenseMatrixSVD::Eval() : info = " << info << endl;
4185 #endif // if MFEM_USE_LAPACK
4192 const int *I = elem_dof.
GetI(), *J = elem_dof.
GetJ(), *dofs;
4193 const double *d_col = tdata;
4196 const double *xp = x;
4200 for (
int i = 0; i < ne; i++)
4203 for (
int col = 0; col < n; col++)
4205 x_col = xp[dofs[col]];
4206 for (
int row = 0; row < n; row++)
4208 yp[dofs[row]] += x_col*d_col[row];
4217 for (
int i = 0; i < ne; i++)
4220 x_col = xp[dofs[0]];
4221 for (
int row = 0; row < n; row++)
4223 ye(row) = x_col*d_col[row];
4226 for (
int col = 1; col < n; col++)
4228 x_col = xp[dofs[col]];
4229 for (
int row = 0; row < n; row++)
4231 ye(row) += x_col*d_col[row];
4235 for (
int row = 0; row < n; row++)
4237 yp[dofs[row]] += ye(row);
4246 for (
int i=0; i<
s; i++)
4262 const int m = Mlu.
SizeI();
4263 const int NE = Mlu.
SizeK();
4269 pivot_flag[0] =
true;
4270 bool *d_pivot_flag = pivot_flag.ReadWrite();
4274 for (
int i = 0; i < m; i++)
4279 double a = fabs(data_all(piv,i,e));
4280 for (
int j = i+1; j < m; j++)
4282 const double b = fabs(data_all(j,i,e));
4289 ipiv_all(i,e) = piv;
4293 for (
int j = 0; j < m; j++)
4295 mfem::kernels::internal::Swap<double>(data_all(i,j,e), data_all(piv,j,e));
4300 if (abs(data_all(i,i,e)) <= TOL)
4302 d_pivot_flag[0] =
false;
4305 const double a_ii_inv = 1.0 / data_all(i,i,e);
4306 for (
int j = i+1; j < m; j++)
4308 data_all(j,i,e) *= a_ii_inv;
4311 for (
int k = i+1; k < m; k++)
4313 const double a_ik = data_all(i,k,e);
4314 for (
int j = i+1; j < m; j++)
4316 data_all(j,k,e) -= a_ik * data_all(j,i,e);
4324 MFEM_ASSERT(pivot_flag.HostRead()[0],
"Batch LU factorization failed \n");
4330 const int m = Mlu.
SizeI();
4331 const int NE = Mlu.
SizeK();
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
void dsyevr_Eigensystem(DenseMatrix &a, Vector &ev, DenseMatrix *evect)
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
int Size() const
Return the logical size of the array.
DenseMatrix & operator-=(const DenseMatrix &m)
void SymmetricScaling(const Vector &s)
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
void SquareRootInverse()
Replaces the current matrix with its square root inverse.
int CheckFinite(const double *v, const int n)
void trans(const Vector &u, Vector &x)
void BatchLUFactor(DenseTensor &Mlu, Array< int > &P, const double TOL)
Compute the LU factorization of a batch of matrices.
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
DenseMatrix & operator*=(double c)
void GetDiag(Vector &d) const
Returns the diagonal of the matrix.
void SetCol(int c, const double *col)
DenseTensor & operator=(double c)
Sets the tensor elements equal to constant c.
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
void SetRow(int r, const double *row)
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
MFEM_HOST_DEVICE void MultABt(const int Aheight, const int Awidth, const int Bheight, const TA *Adata, const TB *Bdata, TC *ABtdata)
Multiply a matrix of size Aheight x Awidth and data Adata with the transpose of a matrix of size Bhei...
void dtrtrs_(char *, char *, char *, int *, int *, double *, int *, double *, int *, int *)
void SingularValues(Vector &sv) const
void dsyev_Eigensystem(DenseMatrix &a, Vector &ev, DenseMatrix *evect)
void SetSize(int s)
Resize the vector to size s.
void Delete()
Delete the owned pointers and reset the Memory object.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void GetInverseMatrix(int m, double *X) const
Assuming L.L^t = A factored data of size (m x m), compute X <- A^{-1}.
void BlockFactor(int m, int n, double *A12, double *A21, double *A22) const
MFEM_HOST_DEVICE void Add(const int height, const int width, const TALPHA alpha, const TA *Adata, const TB *Bdata, TC *Cdata)
Compute C = A + alpha*B, where the matrices A, B and C are of size height x width with data Adata...
void BlockBackSolve(int m, int n, int r, const double *U12, const double *X2, double *Y1) const
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
void dpotrf_(char *, int *, double *, int *, int *)
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
DenseMatrixSVD(DenseMatrix &M, bool left_singular_vectors=false, bool right_singlular_vectors=false)
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
void dgetri_(int *N, double *A, int *LDA, int *IPIV, double *WORK, int *LWORK, int *INFO)
void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const
void TestInversion()
Invert and print the numerical conditioning of the inversion.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
void CopyRows(const DenseMatrix &A, int row1, int row2)
Copy rows row1 through row2 from A to *this.
virtual double Det(int m) const
void Swap(DenseTensor &t)
virtual void GetInverseMatrix(int m, double *X) const
void Eval(DenseMatrix &M)
Abstract data type for matrix inverse.
void AddMult_a_ABt(double a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
void LMult(int m, int n, double *X) const
void GetInverseMatrix(DenseMatrix &Ainv) const
Compute and return the inverse matrix in Ainv.
void Factor()
Factor the current DenseMatrix, *a.
MFEM_HOST_DEVICE double CalcSingularvalue< 3 >(const double *data, const int i)
Return the i'th singular value of the matrix of size 3 with given data.
void GradToVectorCurl2D(DenseMatrix &curl)
virtual void GetInverseMatrix(int m, double *X) const
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
bool IsSquare() const
Returns whether the matrix is a square matrix.
double * GetData() const
Returns the matrix data array.
double * GetData() const
Return a pointer to the beginning of the Vector data.
T Min() const
Find the minimal element in the array, using the comparison operator < for class T.
void CalcOrtho(const DenseMatrix &J, Vector &n)
DenseMatrix & operator=(double c)
Sets the matrix elements equal to constant c.
void Set(double alpha, const double *A)
Set the matrix to alpha * A, assuming that A has the same dimensions as the matrix and uses column-ma...
int Capacity() const
Return the size of the allocated memory.
void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
void dtrsm_(char *side, char *uplo, char *transa, char *diag, int *m, int *n, double *alpha, double *a, int *lda, double *b, int *ldb)
void dgesvd_(char *JOBU, char *JOBVT, int *M, int *N, double *A, int *LDA, double *S, double *U, int *LDU, double *VT, int *LDVT, double *WORK, int *LWORK, int *INFO)
void dsygv_Eigensystem(DenseMatrix &a, DenseMatrix &b, Vector &ev, DenseMatrix *evect)
static void SubMult(int m, int n, int r, const double *A21, const double *X1, double *X2)
void USolve(int m, int n, double *X) const
virtual void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
void dggev_(char *jobvl, char *jobvr, int *n, double *a, int *lda, double *B, int *ldb, double *alphar, double *alphai, double *beta, double *vl, int *ldvl, double *vr, int *ldvr, double *work, int *lwork, int *info)
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
double & operator()(int i, int j)
Returns reference to a_{ij}.
void USolve(int m, int n, double *X) const
double FNorm() const
Compute the Frobenius norm of the matrix.
MFEM_HOST_DEVICE void CalcEigenvalues< 2 >(const double *data, double *lambda, double *vec)
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
friend class DenseMatrixInverse
void RightSolve(int m, int n, double *X) const
double f(const Vector &xvec)
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.
double operator*(const DenseMatrix &m) const
Matrix inner product: tr(A^t B)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
void InvSymmetricScaling(const Vector &s)
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
void GetRow(int r, Vector &row) const
void BlockForwSolve(int m, int n, int r, const double *L21, double *B1, double *B2) const
Abstract data type matrix.
void Norm2(double *v) const
Take the 2-norm of the columns of A and store in v.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
void Invert()
Replaces the current matrix with its inverse.
bool LinearSolve(DenseMatrix &A, double *X, double TOL)
Solves the dense linear system, A * X = B for X
virtual ~DenseMatrixInverse()
Destroys dense inverse matrix.
void LSolve(int m, int n, double *X) const
void dpotrs_(char *, int *, int *, double *, int *, double *, int *, int *)
DenseMatrixGeneralizedEigensystem(DenseMatrix &a, DenseMatrix &b, bool left_eigen_vectors=false, bool right_eigen_vectors=false)
void LeftScaling(const Vector &s)
LeftScaling this = diag(s) * this.
double Det() const
Compute the determinant of the original DenseMatrix using the LU factors.
void AddMultVVt(const Vector &v, DenseMatrix &VVt)
VVt += v v^t.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
void CopyMNDiag(double c, int n, int row_offset, int col_offset)
Copy c on the diagonal of size n to *this at row_offset, col_offset.
virtual double Det(int m) const
void dgetrf_(int *, int *, double *, int *, int *, int *)
~DenseMatrixEigensystem()
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void Neg()
(*this) = -(*this)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
MFEM_HOST_DEVICE void LUSolve(const double *data, const int m, const int *ipiv, double *x)
Assuming L.U = P.A for a factored matrix (m x m),.
virtual void Solve(int m, int n, double *X) const
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), on_dev).
void Getl1Diag(Vector &l) const
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
void AddToVector(int offset, Vector &v) const
Add the matrix 'data' to the Vector 'v' at the given 'offset'.
void GetColumn(int c, Vector &col) const
void AddMult(const Vector &x, Vector &y) const
y += A.x
void dgemm_(char *, char *, int *, int *, int *, double *, double *, int *, double *, int *, double *, double *, int *)
void Threshold(double eps)
Replace small entries, abs(a_ij) <= eps, with zero.
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
void TestInversion()
Print the numerical conditioning of the inversion: ||A^{-1} A - I||.
MFEM_HOST_DEVICE void Mult(const int height, const int width, const TA *data, const TX *x, TY *y)
Matrix vector multiplication: y = A x, where the matrix A is of size height x width with given data...
double MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
double * Data() const
Returns the matrix data array.
void Swap(Array< T > &, Array< T > &)
double p(const Vector &x, double t)
void LSolve(int m, int n, double *X) const
void Transpose()
(*this) = (*this)^t
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
double Trace() const
Trace of a square matrix.
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
void dgetrs_(char *, int *, int *, double *, int *, int *, double *, int *, int *)
void BatchLUSolve(const DenseTensor &Mlu, const Array< int > &P, Vector &X)
Solve batch linear systems.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
virtual bool Factor(int m, double TOL=0.0)
Compute the Cholesky factorization of the current matrix.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
int NumCols() const
Get the number of columns (size of input) of the Operator. Synonym with Width().
void Swap(DenseMatrix &other)
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width.
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
A class to initialize the size of a Tensor.
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
MFEM_HOST_DEVICE double CalcSingularvalue< 2 >(const double *data, const int i)
Return the i'th singular value of the matrix of size 2 with given data.
void SetSubMatrix(const Array< int > &idx, const DenseMatrix &A)
Set (*this)(idx[i],idx[j]) = A(i,j)
void dsygv_(int *ITYPE, char *JOBZ, char *UPLO, int *N, double *A, int *LDA, double *B, int *LDB, double *W, double *WORK, int *LWORK, int *INFO)
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
int height
Dimension of the output / number of rows in the matrix.
virtual void PrintT(std::ostream &out=mfem::out, int width_=4) const
Prints the transpose matrix to stream out.
void CopyCols(const DenseMatrix &A, int col1, int col2)
Copy columns col1 through col2 from A to *this.
virtual MatrixInverse * Inverse() const
Returns a pointer to the inverse matrix.
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
void AddMultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
virtual ~DenseMatrix()
Destroys dense matrix.
void CopyMNt(const DenseMatrix &A, int row_offset, int col_offset)
Copy matrix A^t to the location in *this at row_offset, col_offset.
void AddMultTranspose(const Vector &x, Vector &y) const
y += A^t x
void CopyExceptMN(const DenseMatrix &A, int m, int n)
Copy All rows and columns except m and n from A.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
void Mult(int m, int n, double *X) const
MFEM_HOST_DEVICE void CalcEigenvalues< 3 >(const double *data, double *lambda, double *vec)
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
DenseMatrixInverse(bool spd_=false)
Default constructor.
static const int ipiv_base
void GradToCurl(DenseMatrix &curl)
virtual void Solve(int m, int n, double *X) const
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
void GetRowSums(Vector &l) const
Compute the row sums of the DenseMatrix.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
void CalcEigenvalues(double *lambda, double *vec) const
virtual bool Factor(int m, double TOL=0.0)
void RightSolve(int m, int n, double *X) const
DenseMatrixEigensystem(DenseMatrix &m)
void dsyevr_(char *JOBZ, char *RANGE, char *UPLO, int *N, double *A, int *LDA, double *VL, double *VU, int *IL, int *IU, double *ABSTOL, int *M, double *W, double *Z, int *LDZ, int *ISUPPZ, double *WORK, int *LWORK, int *IWORK, int *LIWORK, int *INFO)
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
int Rank(double tol) const
void AddMult_a(double a, const Vector &x, Vector &y) const
y += a * A.x
void RightScaling(const Vector &s)
RightScaling: this = this * diag(s);.
void Mult(const double *x, double *y) const
Matrix vector multiplication.
void dpotri_(char *, int *, double *, int *, int *)
void AddMultTranspose_a(double a, const Vector &x, Vector &y) const
y += a * A^t x
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
void GetFromVector(int offset, const Vector &v)
Get the matrix 'data' from the Vector 'v' at the given 'offset'.
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
void GetSubMatrix(const Array< int > &idx, DenseMatrix &A) const
~DenseMatrixGeneralizedEigensystem()
double u(const Vector &xvec)
void InvLeftScaling(const Vector &s)
InvLeftScaling this = diag(1./s) * this.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void AddMult_a(double alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
Rank 3 tensor (array of matrices)
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
void AdjustDofDirection(Array< int > &dofs)
MFEM_HOST_DEVICE void Symmetrize(const int size, T *data)
Symmetrize a square matrix with given size and data: A -> (A+A^T)/2.
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void GradToDiv(Vector &div)
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
int width
Dimension of the input / number of columns in the matrix.
void AddSubMatrix(const Array< int > &idx, const DenseMatrix &A)
(*this)(idx[i],idx[j]) += A(i,j)
void dsyev_(char *JOBZ, char *UPLO, int *N, double *A, int *LDA, double *W, double *WORK, int *LWORK, int *INFO)
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints operator in Matlab format.
virtual void Solve(int m, int n, double *X) const
DenseMatrix & operator+=(const double *m)
void UMult(int m, int n, double *X) const