19 #include "../general/table.hpp"
26 #if defined(_MSC_VER) && (_MSC_VER < 1800)
28 #define copysign _copysign
47 data =
new double[hw];
49 for (
int i = 0; i < hw; i++)
63 MFEM_ASSERT(s >= 0,
"invalid DenseMatrix size: " << s);
67 data =
new double[capacity]();
77 MFEM_ASSERT(m >= 0 && n >= 0,
78 "invalid DenseMatrix size: " << m <<
" x " << n);
82 data =
new double[capacity]();
91 :
Matrix(mat.width, mat.height)
96 data =
new double[capacity];
98 for (
int i = 0; i <
height; i++)
99 for (
int j = 0; j <
width; j++)
101 (*this)(i,j) = mat(j,i);
112 MFEM_ASSERT(h >= 0 && w >= 0,
113 "invalid DenseMatrix size: " << h <<
" x " << w);
121 if (hw > std::abs(capacity))
128 data =
new double[hw]();
146 for (
int row = 0; row <
height; row++)
152 double *d_col = data;
154 for (
int row = 0; row <
height; row++)
156 y[row] = x_col*d_col[row];
159 for (
int col = 1; col <
width; col++)
162 for (
int row = 0; row <
height; row++)
164 y[row] += x_col*d_col[row];
173 "incompatible dimensions");
175 Mult((
const double *)x, (
double *)y);
181 "incompatible dimensions");
185 for (
int i = 0; i < hw; i++)
187 a += data[i] * m.data[i];
195 double *d_col = data;
196 for (
int col = 0; col <
width; col++)
199 for (
int row = 0; row <
height; row++)
201 y_col += x[row]*d_col[row];
211 "incompatible dimensions");
219 "incompatible dimensions");
221 const double *xp = x;
222 double *d_col = data, *yp = y;
223 for (
int col = 0; col <
width; col++)
225 double x_col = xp[col];
226 for (
int row = 0; row <
height; row++)
228 yp[row] += x_col*d_col[row];
237 "incompatible dimensions");
239 const double *xp = x;
240 double *d_col = data, *yp = y;
241 for (
int col = 0; col <
width; col++)
243 double x_col = a*xp[col];
244 for (
int row = 0; row <
height; row++)
246 yp[row] += x_col*d_col[row];
256 "incompatible dimensions");
258 double *d_col = data;
259 for (
int col = 0; col <
width; col++)
262 for (
int row = 0; row <
height; row++)
264 y_col += x[row]*d_col[row];
275 for (
int i = 0; i <
height; i++)
278 for (
int j = 0; j <
width; j++)
280 Axi += (*this)(i,j) * x[j];
291 double * it_data = data;
292 for (
int j = 0; j <
width; ++j)
293 for (
int i = 0; i <
height; ++i)
295 *(it_data++) *= s(i);
302 double * it_data = data;
303 for (
int j = 0; j <
width; ++j)
304 for (
int i = 0; i <
height; ++i)
306 *(it_data++) /= s(i);
314 double * it_data = data;
315 for (
int j = 0; j <
width; ++j)
318 for (
int i = 0; i <
height; ++i)
329 double * it_data = data;
330 for (
int j = 0; j <
width; ++j)
333 for (
int i = 0; i <
height; ++i)
348 double * ss =
new double[
width];
351 for (
double * end_s = it_s +
width; it_s != end_s; ++it_s)
353 *(it_ss++) = sqrt(*it_s);
356 double * it_data = data;
357 for (
int j = 0; j <
width; ++j)
358 for (
int i = 0; i <
height; ++i)
360 *(it_data++) *= ss[i]*ss[j];
374 double * ss =
new double[
width];
377 for (
double * end_s = it_s +
width; it_s != end_s; ++it_s)
379 *(it_ss++) = 1./sqrt(*it_s);
382 double * it_data = data;
383 for (
int j = 0; j <
width; ++j)
384 for (
int i = 0; i <
height; ++i)
386 *(it_data++) *= ss[i]*ss[j];
397 mfem_error(
"DenseMatrix::Trace() : not a square matrix!");
403 for (
int i = 0; i <
width; i++)
431 return data[0] * data[3] - data[1] * data[2];
435 const double *d = data;
437 d[0] * (d[4] * d[8] - d[5] * d[7]) +
438 d[3] * (d[2] * d[7] - d[1] * d[8]) +
439 d[6] * (d[1] * d[5] - d[2] * d[4]);
454 return sqrt(data[0] * data[0] + data[1] * data[1]);
458 return sqrt(data[0] * data[0] + data[1] * data[1] + data[2] * data[2]);
462 const double *d = data;
463 double E = d[0] * d[0] + d[1] * d[1] + d[2] * d[2];
464 double G = d[3] * d[3] + d[4] * d[4] + d[5] * d[5];
465 double F = d[0] * d[3] + d[1] * d[4] + d[2] * d[5];
466 return sqrt(E * G - F * F);
474 for (
int j = 0; j <
Width(); j++)
475 for (
int i = 0; i <
Height(); i++)
477 (*this)(i,j) += c * A(i,j);
484 for (
int i = 0; i < s; i++)
494 for (
int i = 0; i < s; i++)
506 for (
int i = 0; i < hw; i++)
517 "incompatible matrix sizes.");
519 for (
int j = 0; j <
width; j++)
520 for (
int i = 0; i <
height; i++)
522 (*this)(i, j) += m(i, j);
530 for (
int j = 0; j <
width; j++)
531 for (
int i = 0; i <
height; i++)
533 (*this)(i, j) -= m(i, j);
542 for (
int i = 0; i < s; i++)
552 for (
int i = 0; i < hw; i++)
558 #ifdef MFEM_USE_LAPACK
560 dgetrf_(
int *,
int *,
double *,
int *,
int *,
int *);
562 dgetrs_(
char *,
int *,
int *,
double *,
int *,
int *,
double *,
int *,
int *);
564 dgetri_(
int *N,
double *A,
int *LDA,
int *IPIV,
double *WORK,
565 int *LWORK,
int *INFO);
577 #ifdef MFEM_USE_LAPACK
578 int *ipiv =
new int[
width];
587 mfem_error(
"DenseMatrix::Invert() : Error in DGETRF");
593 work =
new double[lwork];
599 mfem_error(
"DenseMatrix::Invert() : Error in DGETRI");
605 int c, i, j, n =
Width();
609 for (c = 0; c < n; c++)
611 a = fabs((*
this)(c, c));
613 for (j = c + 1; j < n; j++)
615 b = fabs((*
this)(j, c));
624 mfem_error(
"DenseMatrix::Invert() : singular matrix");
627 for (j = 0; j < n; j++)
629 Swap<double>((*this)(c, j), (*
this)(i, j));
632 a = (*this)(c, c) = 1.0 / (*
this)(c, c);
633 for (j = 0; j < c; j++)
637 for (j++; j < n; j++)
641 for (i = 0; i < c; i++)
643 (*this)(i, c) = a * (b = -(*
this)(i, c));
644 for (j = 0; j < c; j++)
646 (*this)(i, j) += b * (*
this)(c, j);
648 for (j++; j < n; j++)
650 (*this)(i, j) += b * (*
this)(c, j);
653 for (i++; i < n; i++)
655 (*this)(i, c) = a * (b = -(*
this)(i, c));
656 for (j = 0; j < c; j++)
658 (*this)(i, j) += b * (*
this)(c, j);
660 for (j++; j < n; j++)
662 (*this)(i, j) += b * (*
this)(c, j);
667 for (c = n - 1; c >= 0; c--)
670 for (i = 0; i < n; i++)
672 Swap<double>((*this)(i, c), (*
this)(i, j));
680 for (
int j = 0; j <
Width(); j++)
683 for (
int i = 0; i <
Height(); i++)
685 v[j] += (*this)(i,j)*(*
this)(i,j);
694 const double *d = data;
695 double norm = 0.0, abs_entry;
697 for (
int i = 0; i < hw; i++)
699 abs_entry = fabs(d[i]);
700 if (norm < abs_entry)
712 double max_norm = 0.0, entry, fnorm2;
714 for (i = 0; i < hw; i++)
716 entry = fabs(data[i]);
717 if (entry > max_norm)
729 for (i = 0; i < hw; i++)
731 entry = data[i] / max_norm;
732 fnorm2 += entry * entry;
735 return max_norm * sqrt(fnorm2);
738 #ifdef MFEM_USE_LAPACK
740 dsyevr_(
char *JOBZ,
char *RANGE,
char *UPLO,
int *N,
double *A,
int *LDA,
741 double *VL,
double *VU,
int *IL,
int *IU,
double *ABSTOL,
int *M,
742 double *W,
double *Z,
int *LDZ,
int *ISUPPZ,
double *WORK,
int *LWORK,
743 int *IWORK,
int *LIWORK,
int *INFO);
745 dsyev_(
char *JOBZ,
char *UPLO,
int *N,
double *A,
int *LDA,
double *W,
746 double *WORK,
int *LWORK,
int *INFO);
748 dgesvd_(
char *JOBU,
char *JOBVT,
int *M,
int *N,
double *A,
int *LDA,
749 double *S,
double *U,
int *LDU,
double *VT,
int *LDVT,
double *WORK,
750 int *LWORK,
int *INFO);
756 #ifdef MFEM_USE_LAPACK
764 double *A =
new double[N*N];
775 int *ISUPPZ =
new int[2*N];
794 double *data = a.
Data();
796 for (
int i = 0; i < hw; i++)
801 dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
802 &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, &QWORK, &LWORK,
803 &QIWORK, &LIWORK, &INFO );
808 WORK =
new double[LWORK];
809 IWORK =
new int[LIWORK];
811 dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
812 &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, WORK, &LWORK,
813 IWORK, &LIWORK, &INFO );
817 cerr <<
"dsyevr_Eigensystem(...): DSYEVR error code: "
825 cerr <<
"dsyevr_Eigensystem(...):\n"
826 <<
" DSYEVR did not find all eigenvalues "
827 << M <<
"/" << N << endl;
830 for (IL = 0; IL < N; IL++)
833 mfem_error(
"dsyevr_Eigensystem(...): !finite value in W");
835 for (IL = 0; IL < N*N; IL++)
838 mfem_error(
"dsyevr_Eigensystem(...): !finite value in Z");
841 for (IL = 0; IL < N; IL++)
842 for (IU = 0; IU <= IL; IU++)
845 for (M = 0; M < N; M++)
847 VL += Z[M+IL*N] * Z[M+IU*N];
863 cerr <<
"dsyevr_Eigensystem(...):"
864 <<
" Z^t Z - I deviation = " << VU
865 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
866 << W[0] <<
", N = " << N << endl;
872 cerr <<
"dsyevr_Eigensystem(...):"
873 <<
" Z^t Z - I deviation = " << VU
874 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
875 << W[0] <<
", N = " << N << endl;
879 mfem_error(
"dsyevr_Eigensystem(...): ERROR: ...");
882 for (IL = 0; IL < N; IL++)
883 for (IU = 0; IU < N; IU++)
886 for (M = 0; M < N; M++)
888 VL += Z[IL+M*N] * W[M] * Z[IU+M*N];
890 VL = fabs(VL-data[IL+N*IU]);
898 cerr <<
"dsyevr_Eigensystem(...):"
899 <<
" max matrix deviation = " << VU
900 <<
"\n W[max] = " << W[N-1] <<
", W[min] = "
901 << W[0] <<
", N = " << N << endl;
905 mfem_error(
"dsyevr_Eigensystem(...): ERROR: ...");
920 #ifdef MFEM_USE_LAPACK
948 double *data = a.
Data();
949 for (
int i = 0; i < hw; i++)
954 dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, &QWORK, &LWORK, &INFO);
957 WORK =
new double[LWORK];
959 dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, WORK, &LWORK, &INFO);
963 cerr <<
"dsyev_Eigensystem: DSYEV error code: " << INFO << endl;
968 if (evect == NULL) {
delete [] A; }
973 void DenseMatrix::Eigensystem(Vector &ev, DenseMatrix *evect)
975 #ifdef MFEM_USE_LAPACK
990 #ifdef MFEM_USE_LAPACK
996 double *a = copy_of_this.data;
1001 double *work = NULL;
1006 dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1007 s, u, &m, vt, &n, &qwork, &lwork, &info);
1009 lwork = (int) qwork;
1010 work =
new double[lwork];
1012 dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1013 s, u, &m, vt, &n, work, &lwork, &info);
1018 cerr <<
"DenseMatrix::SingularValues : info = " << info << endl;
1033 for (
int i=0; i < sv.
Size(); ++i)
1042 static const double sqrt_1_eps = sqrt(1./numeric_limits<double>::epsilon());
1049 double t, zeta = (d2 - d1)/(2*d12);
1050 if (fabs(zeta) < sqrt_1_eps)
1052 t = d12*copysign(1./(fabs(zeta) + sqrt(1. + zeta*zeta)), zeta);
1056 t = d12*copysign(0.5/fabs(zeta), zeta);
1064 double &c,
double &s)
1074 double t, zeta = (d2 - d1)/(2*d12);
1075 if (fabs(zeta) < sqrt_1_eps)
1077 t = copysign(1./(fabs(zeta) + sqrt(1. + zeta*zeta)), zeta);
1081 t = copysign(0.5/fabs(zeta), zeta);
1084 c = sqrt(1./(1. + t*t));
1093 const double &x1,
const double &x2,
const double &x3,
1094 double &n1,
double &n2,
double &n3)
1102 t = sqrt(1./(t + r*r));
1103 n1 = copysign(t, x1);
1110 double &n1,
double &n2,
double &n3)
1114 if (fabs(x1) >= fabs(x2))
1116 if (fabs(x1) >= fabs(x3))
1129 else if (fabs(x2) >= fabs(x3))
1139 double &d1,
double &d12,
double &d21,
double &d2)
1150 double n1 = fabs(d1) + fabs(d21);
1151 double n2 = fabs(d2) + fabs(d12);
1153 bool swap_columns = (n2 > n1);
1165 if (fabs(d1) > fabs(d21))
1173 if (fabs(d1) < fabs(d21))
1185 if (fabs(d12) > fabs(d2))
1198 if (fabs(d12) < fabs(d2))
1211 n1 = hypot(d1, d21);
1217 mu = copysign(n1, d1);
1218 n1 = -d21*(d21/(d1 + mu));
1222 if (fabs(n1) <= fabs(d21))
1226 mu = (2./(1. + n1*n1))*(n1*d12 + d2);
1234 mu = (2./(1. + n2*n2))*(d12 + n2*d2);
1256 n2 = 1./(1. + fabs(mu));
1258 if (fabs(d1) <= n2*fabs(d2))
1280 double &d1,
double &d2,
double &d3,
double &c12,
double &c13,
double &c23,
1281 double &c21,
double &c31,
double &c32)
1284 double mu, n1, n2, n3, s1, s2, s3;
1286 s1 = hypot(c21, c31);
1293 mu = copysign(n1, d1);
1294 n1 = -s1*(s1/(d1 + mu));
1299 if (fabs(n1) >= fabs(c21))
1301 if (fabs(n1) >= fabs(c31))
1306 mu = 2./(1. + s2*s2 + s3*s3);
1307 n2 = mu*(c12 + s2*d2 + s3*c32);
1308 n3 = mu*(c13 + s2*c23 + s3*d3);
1318 else if (fabs(c21) >= fabs(c31))
1323 mu = 2./(1. + s1*s1 + s3*s3);
1324 n2 = mu*(s1*c12 + d2 + s3*c32);
1325 n3 = mu*(s1*c13 + c23 + s3*d3);
1337 mu = 2./(1. + s1*s1 + s2*s2);
1338 n2 = mu*(s1*c12 + s2*d2 + c32);
1339 n3 = mu*(s1*c13 + s2*c23 + d3);
1374 d1 = -(c12*d2 + c13*d3)/d1;
1385 const double &d12,
const double &d13,
const double &d23,
1386 double &d1,
double &d2,
double &d3)
1399 double c12 = d12, c13 = d13, c23 = d23;
1400 double c21, c31, c32;
1404 c32 = fabs(d1) + fabs(c12) + fabs(c13);
1405 c31 = fabs(d2) + fabs(c12) + fabs(c23);
1406 c21 = fabs(d3) + fabs(c13) + fabs(c23);
1411 col = (c32 >= c31) ? 1 : 2;
1415 col = (c31 >= c21) ? 2 : 3;
1447 if (fabs(d1) <= fabs(c13))
1449 row = (fabs(d1) <= fabs(c12)) ? 1 : 2;
1453 row = (fabs(c12) <= fabs(c13)) ? 2 : 3;
1458 if (fabs(d1) >= fabs(c13))
1460 row = (fabs(d1) >= fabs(c12)) ? 1 : 2;
1464 row = (fabs(c12) >= fabs(c13)) ? 2 : 3;
1514 double &d1,
double &d2,
double &d3,
double &d12,
double &d13,
double &d23,
1515 double &z1,
double &z2,
double &z3,
double &v1,
double &v2,
double &v3,
1535 double s, w1, w2, w3;
1541 if (fabs(z1) <= fabs(z3))
1543 k = (fabs(z1) <= fabs(z2)) ? 1 : 2;
1547 k = (fabs(z2) <= fabs(z3)) ? 2 : 3;
1553 if (fabs(z1) >= fabs(z3))
1555 k = (fabs(z1) >= fabs(z2)) ? 1 : 2;
1559 k = (fabs(z2) >= fabs(z3)) ? 2 : 3;
1586 g = copysign(1., z1);
1587 v1 = -s*(s/(z1 + g));
1590 if (fabs(z2) > g) { g = fabs(z2); }
1591 if (fabs(z3) > g) { g = fabs(z3); }
1595 g = 2./(v1*v1 + v2*v2 + v3*v3);
1600 w1 = g*( d1*v1 + d12*v2 + d13*v3);
1601 w2 = g*(d12*v1 + d2*v2 + d23*v3);
1602 w3 = g*(d13*v1 + d23*v2 + d3*v3);
1604 s = (g/2)*(v1*w1 + v2*w2 + v3*w3);
1611 d23 -= v2*w3 + v3*w2;
1616 s = d12 - v1*w2 - v2*w1;
1617 s = d13 - v1*w3 - v3*w1;
1639 mult = frexp(d_max, &d_exp);
1640 if (d_exp == numeric_limits<double>::max_exponent)
1642 mult *= numeric_limits<double>::radix;
1659 mfem_error(
"DenseMatrix::CalcSingularvalue");
1664 const double *d = data;
1672 double d0, d1, d2, d3;
1679 double d_max = fabs(d0);
1680 if (d_max < fabs(d1)) { d_max = fabs(d1); }
1681 if (d_max < fabs(d2)) { d_max = fabs(d2); }
1682 if (d_max < fabs(d3)) { d_max = fabs(d3); }
1696 double t = 0.5*((d0+d2)*(d0-d2)+(d1-d3)*(d1+d3));
1698 double s = d0*d2 + d1*d3;
1699 s = sqrt(0.5*(d0*d0 + d1*d1 + d2*d2 + d3*d3) + sqrt(t*t + s*s));
1704 t = fabs(d0*d3 - d1*d2) / s;
1721 double d0, d1, d2, d3, d4, d5, d6, d7, d8;
1722 d0 = d[0]; d3 = d[3]; d6 = d[6];
1723 d1 = d[1]; d4 = d[4]; d7 = d[7];
1724 d2 = d[2]; d5 = d[5]; d8 = d[8];
1727 double d_max = fabs(d0);
1728 if (d_max < fabs(d1)) { d_max = fabs(d1); }
1729 if (d_max < fabs(d2)) { d_max = fabs(d2); }
1730 if (d_max < fabs(d3)) { d_max = fabs(d3); }
1731 if (d_max < fabs(d4)) { d_max = fabs(d4); }
1732 if (d_max < fabs(d5)) { d_max = fabs(d5); }
1733 if (d_max < fabs(d6)) { d_max = fabs(d6); }
1734 if (d_max < fabs(d7)) { d_max = fabs(d7); }
1735 if (d_max < fabs(d8)) { d_max = fabs(d8); }
1740 d0 /= mult; d1 /= mult; d2 /= mult;
1741 d3 /= mult; d4 /= mult; d5 /= mult;
1742 d6 /= mult; d7 /= mult; d8 /= mult;
1744 double b11 = d0*d0 + d1*d1 + d2*d2;
1745 double b12 = d0*d3 + d1*d4 + d2*d5;
1746 double b13 = d0*d6 + d1*d7 + d2*d8;
1747 double b22 = d3*d3 + d4*d4 + d5*d5;
1748 double b23 = d3*d6 + d4*d7 + d5*d8;
1749 double b33 = d6*d6 + d7*d7 + d8*d8;
1768 double aa = (b11 + b22 + b33)/3;
1774 double b11_b22 = ((d0-d3)*(d0+d3)+(d1-d4)*(d1+d4)+(d2-d5)*(d2+d5));
1775 double b22_b33 = ((d3-d6)*(d3+d6)+(d4-d7)*(d4+d7)+(d5-d8)*(d5+d8));
1776 double b33_b11 = ((d6-d0)*(d6+d0)+(d7-d1)*(d7+d1)+(d8-d2)*(d8+d2));
1777 c1 = (b11_b22 - b33_b11)/3;
1778 c2 = (b22_b33 - b11_b22)/3;
1779 c3 = (b33_b11 - b22_b33)/3;
1782 Q = (2*(b12*b12 + b13*b13 + b23*b23) + c1*c1 + c2*c2 + c3*c3)/6;
1783 R = (c1*(b23*b23 - c2*c3)+ b12*(b12*c3 - 2*b13*b23) +b13*b13*c2)/2;
1816 double sqrtQ = sqrt(Q);
1817 double sqrtQ3 = Q*sqrtQ;
1822 if (fabs(R) >= sqrtQ3)
1844 aa -= 2*sqrtQ*cos(acos(R)/3);
1848 aa -= 2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3);
1852 aa -= 2*sqrtQ*cos((acos(R) - 2.0*M_PI)/3);
1859 r = -2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3);
1868 r = -2*sqrtQ*cos(acos(R)/3);
1915 double v1, v2, v3, g;
1916 Reduce3S(mode, b11, b22, b33, b12, b13, b23,
1917 c1, c2, c3, v1, v2, v3, g);
1928 aa = std::min(std::min(b11, b22), b33);
1934 aa = (b22 <= b33) ? b22 : std::max(b11, b33);
1938 aa = (b11 <= b33) ? b11 : std::max(b33, b22);
1943 aa = std::max(std::max(b11, b22), b33);
1949 return sqrt(fabs(aa))*mult;
1963 const double *d = data;
2003 double d_max = fabs(d11);
2004 if (d_max < fabs(d22)) { d_max = fabs(d22); }
2005 if (d_max < fabs(d33)) { d_max = fabs(d33); }
2006 if (d_max < fabs(d12)) { d_max = fabs(d12); }
2007 if (d_max < fabs(d13)) { d_max = fabs(d13); }
2008 if (d_max < fabs(d23)) { d_max = fabs(d23); }
2013 d11 /= mult; d22 /= mult; d33 /= mult;
2014 d12 /= mult; d13 /= mult; d23 /= mult;
2016 double aa = (d11 + d22 + d33)/3;
2017 double c1 = d11 - aa;
2018 double c2 = d22 - aa;
2019 double c3 = d33 - aa;
2023 Q = (2*(d12*d12 + d13*d13 + d23*d23) + c1*c1 + c2*c2 + c3*c3)/6;
2024 R = (c1*(d23*d23 - c2*c3)+ d12*(d12*c3 - 2*d13*d23) + d13*d13*c2)/2;
2028 lambda[0] = lambda[1] = lambda[2] = aa;
2029 vec[0] = 1.; vec[3] = 0.; vec[6] = 0.;
2030 vec[1] = 0.; vec[4] = 1.; vec[7] = 0.;
2031 vec[2] = 0.; vec[5] = 0.; vec[8] = 1.;
2035 double sqrtQ = sqrt(Q);
2036 double sqrtQ3 = Q*sqrtQ;
2040 if (fabs(R) >= sqrtQ3)
2059 r = -2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3);
2063 r = -2*sqrtQ*cos(acos(R)/3);
2091 lambda[0] = lambda[1] = lambda[2] = aa;
2092 vec[0] = 1.; vec[3] = 0.; vec[6] = 0.;
2093 vec[1] = 0.; vec[4] = 1.; vec[7] = 0.;
2094 vec[2] = 0.; vec[5] = 0.; vec[8] = 1.;
2108 double v1, v2, v3, g;
2109 int k =
Reduce3S(mode, d11, d22, d33, d12, d13, d23,
2110 c1, c2, c3, v1, v2, v3, g);
2121 double *vec_1, *vec_2, *vec_3;
2126 lambda[0] = d11; vec_1 = vec;
2127 lambda[1] = d22; vec_2 = vec + 3;
2128 lambda[2] = d33; vec_3 = vec + 6;
2130 else if (d11 <= d33)
2132 lambda[0] = d11; vec_1 = vec;
2133 lambda[1] = d33; vec_3 = vec + 3;
2134 lambda[2] = d22; vec_2 = vec + 6;
2138 lambda[0] = d33; vec_3 = vec;
2139 lambda[1] = d11; vec_1 = vec + 3;
2140 lambda[2] = d22; vec_2 = vec + 6;
2147 lambda[0] = d22; vec_2 = vec;
2148 lambda[1] = d11; vec_1 = vec + 3;
2149 lambda[2] = d33; vec_3 = vec + 6;
2151 else if (d22 <= d33)
2153 lambda[0] = d22; vec_2 = vec;
2154 lambda[1] = d33; vec_3 = vec + 3;
2155 lambda[2] = d11; vec_1 = vec + 6;
2159 lambda[0] = d33; vec_3 = vec;
2160 lambda[1] = d22; vec_2 = vec + 3;
2161 lambda[2] = d11; vec_1 = vec + 6;
2168 d22 = g*(v2*c - v3*s);
2169 d33 = g*(v2*s + v3*c);
2170 vec_2[0] = - v1*d22; vec_3[0] = - v1*d33;
2171 vec_2[1] = c - v2*d22; vec_3[1] = s - v2*d33;
2172 vec_2[2] = -s - v3*d22; vec_3[2] = c - v3*d33;
2176 Swap(vec_2[0], vec_2[1]);
2177 Swap(vec_3[0], vec_3[1]);
2181 Swap(vec_2[0], vec_2[2]);
2182 Swap(vec_3[0], vec_3[2]);
2203 for (
int i = 0; i < n; i++)
2217 for (
int i = 0; i <
height; ++i)
2219 d(i) = (*this)(i,i);
2233 for (
int j = 0; j <
width; ++j)
2234 for (
int i = 0; i <
height; ++i)
2236 l(i) += fabs((*
this)(i,j));
2243 for (
int i = 0; i <
height; i++)
2246 for (
int j = 0; j <
width; j++)
2259 for (i = 0; i < N; i++)
2263 for (i = 0; i < n; i++)
2274 for (i = 0; i < N; i++)
2278 for (i = 0; i < n; i++)
2280 data[i*(n+1)] = diag[i];
2291 for (i = 0; i <
Height(); i++)
2292 for (j = i+1; j <
Width(); j++)
2295 (*this)(i,j) = (*
this)(j,i);
2310 for (
int i = 0; i <
Height(); i++)
2311 for (
int j = 0; j <
Width(); j++)
2313 (*this)(i,j) = A(j,i);
2322 mfem_error(
"DenseMatrix::Symmetrize() : not a square matrix!");
2326 for (
int i = 0; i <
Height(); i++)
2327 for (
int j = 0; j < i; j++)
2329 double a = 0.5 * ((*this)(i,j) + (*
this)(j,i));
2330 (*this)(j,i) = (*
this)(i,j) = a;
2336 for (
int i = 0; i <
Height(); i++)
2339 for (
int j = 0; j <
Width(); j++)
2342 (*this)(i, j) = 0.0;
2362 for (
int i = 0; i < n; i++)
2365 double x = (*this)(i,0);
2366 double y = (*this)(i,1);
2379 for (
int i = 0; i < n; i++)
2382 double x = (*this)(i,0);
2383 double y = (*this)(i,1);
2384 double z = (*this)(i,2);
2420 double *ddata = div.
GetData();
2422 for (
int i = 0; i < n; i++)
2432 for (
int i = row1; i <= row2; i++)
2433 for (
int j = 0; j <
Width(); j++)
2435 (*this)(i-row1,j) = A(i,j);
2443 for (
int i = 0; i <
Height(); i++)
2444 for (
int j = col1; j <= col2; j++)
2446 (*this)(i,j-col1) = A(i,j);
2456 for (j = 0; j < n; j++)
2457 for (i = 0; i < m; i++)
2459 (*this)(i,j) = A(Aro+i,Aco+j);
2468 for (j = 0; j < A.
Width(); j++)
2469 for (i = 0; i < A.
Height(); i++)
2471 (*this)(row_offset+i,col_offset+j) = *(v++);
2480 for (i = 0; i < A.
Width(); i++)
2481 for (j = 0; j < A.
Height(); j++)
2483 (*this)(row_offset+i,col_offset+j) = *(v++);
2488 int row_offset,
int col_offset)
2492 MFEM_VERIFY(row_offset+m <= this->
Height() && col_offset+n <= this->
Width(),
2493 "this DenseMatrix is too small to accomodate the submatrix.");
2494 MFEM_VERIFY(Aro+m <= A.
Height() && Aco+n <= A.
Width(),
2495 "The A DenseMatrix is too small to accomodate the submatrix.");
2497 for (j = 0; j < n; j++)
2498 for (i = 0; i < m; i++)
2500 (*this)(row_offset+i,col_offset+j) = A(Aro+i,Aco+j);
2508 for (i = 0; i < n; i++)
2509 for (j = i+1; j < n; j++)
2510 (*
this)(row_offset+i,col_offset+j) =
2511 (*
this)(row_offset+j,col_offset+i) = 0.0;
2513 for (i = 0; i < n; i++)
2515 (*this)(row_offset+i,col_offset+i) = c;
2524 for (i = 0; i < n; i++)
2525 for (j = i+1; j < n; j++)
2526 (*
this)(row_offset+i,col_offset+j) =
2527 (*
this)(row_offset+j,col_offset+i) = 0.0;
2529 for (i = 0; i < n; i++)
2531 (*this)(row_offset+i,col_offset+i) = diag[i];
2545 if (co+aw >
Width() || ro+ah > h)
2551 p = data + ro + co * h;
2554 for (
int c = 0; c < aw; c++)
2556 for (
int r = 0; r < ah; r++)
2575 if (co+aw >
Width() || ro+ah > h)
2581 p = data + ro + co * h;
2584 for (
int c = 0; c < aw; c++)
2586 for (
int r = 0; r < ah; r++)
2598 double *vdata = v.
GetData() + offset;
2600 for (i = 0; i < n; i++)
2602 vdata[i] += data[i];
2609 const double *vdata = v.
GetData() + offset;
2611 for (i = 0; i < n; i++)
2624 mfem_error(
"DenseMatrix::AdjustDofDirection(...)");
2629 for (
int i = 0; i < n-1; i++)
2631 int s = (dof[i] < 0) ? (-1) : (1);
2632 for (
int j = i+1; j < n; j++)
2634 int t = (dof[j] < 0) ? (-s) : (s);
2637 (*this)(i,j) = -(*
this)(i,j);
2638 (*this)(j,i) = -(*
this)(j,i);
2646 for (
int j = 0; j <
Width(); j++)
2648 (*this)(row, j) = value;
2654 for (
int i = 0; i <
Height(); i++)
2656 (*this)(i, col) = value;
2662 for (
int col = 0; col <
Width(); col++)
2664 for (
int row = 0; row <
Height(); row++)
2666 if (std::abs(
operator()(row,col)) <= eps)
2677 ios::fmtflags old_flags = out.flags();
2679 out << setiosflags(ios::scientific | ios::showpos);
2680 for (
int i = 0; i <
height; i++)
2682 out <<
"[row " << i <<
"]\n";
2683 for (
int j = 0; j <
width; j++)
2685 out << (*this)(i,j);
2686 if (j+1 == width || (j+1) % width_ == 0)
2697 out.flags(old_flags);
2703 ios::fmtflags old_flags = out.flags();
2705 out << setiosflags(ios::scientific | ios::showpos);
2706 for (
int i = 0; i <
height; i++)
2708 for (
int j = 0; j <
width; j++)
2710 out << (*this)(i,j);
2716 out.flags(old_flags);
2722 ios::fmtflags old_flags = out.flags();
2724 out << setiosflags(ios::scientific | ios::showpos);
2725 for (
int j = 0; j <
width; j++)
2727 out <<
"[col " << j <<
"]\n";
2728 for (
int i = 0; i <
height; i++)
2730 out << (*this)(i,j);
2731 if (i+1 == height || (i+1) % width_ == 0)
2742 out.flags(old_flags);
2751 for (
int i = 0; i <
width; i++)
2755 cout <<
"size = " << width <<
", i_max = " << C.
MaxMaxNorm()
2756 <<
", cond_F = " <<
FNorm()*copy.FNorm() << endl;
2772 for (
int j = 0; j < C.
Width(); j++)
2773 for (
int i = 0; i < C.
Height(); i++)
2775 C(i,j) = A(i,j) + alpha * B(i,j);
2782 for (
int j = 0; j < C.
Width(); j++)
2783 for (
int i = 0; i < C.
Height(); i++)
2785 C(i,j) = alpha * A(i,j) + beta * B(i,j);
2790 #ifdef MFEM_USE_LAPACK
2792 dgemm_(
char *,
char *,
int *,
int *,
int *,
double *,
double *,
2793 int *,
double *,
int *,
double *,
double *,
int *);
2799 b.
Width() == c.
Height(),
"incompatible dimensions");
2801 #ifdef MFEM_USE_LAPACK
2802 static char transa =
'N', transb =
'N';
2803 static double alpha = 1.0, beta = 0.0;
2806 dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.
Data(), &m,
2807 c.
Data(), &k, &beta, a.
Data(), &m);
2809 const int ah = a.
Height();
2810 const int aw = a.
Width();
2811 const int bw = b.
Width();
2812 double *ad = a.
Data();
2813 const double *bd = b.
Data();
2814 const double *cd = c.
Data();
2815 for (
int i = 0; i < ah*aw; i++)
2819 for (
int j = 0; j < aw; j++)
2821 for (
int k = 0; k < bw; k++)
2823 for (
int i = 0; i < ah; i++)
2825 ad[i+j*ah] += bd[i+k*ah] * cd[k+j*bw];
2835 b.
Width() == c.
Height(),
"incompatible dimensions");
2837 #ifdef MFEM_USE_LAPACK
2838 static char transa =
'N', transb =
'N';
2839 static double alpha = 1.0, beta = 1.0;
2842 dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.
Data(), &m,
2843 c.
Data(), &k, &beta, a.
Data(), &m);
2845 const int ah = a.
Height();
2846 const int aw = a.
Width();
2847 const int bw = b.
Width();
2848 double *ad = a.
Data();
2849 const double *bd = b.
Data();
2850 const double *cd = c.
Data();
2851 for (
int j = 0; j < aw; j++)
2853 for (
int k = 0; k < bw; k++)
2855 for (
int i = 0; i < ah; i++)
2857 ad[i+j*ah] += bd[i+k*ah] * cd[k+j*bw];
2879 const double *d = a.
Data();
2880 double *ad = adja.
Data();
2895 e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
2896 g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
2897 f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
2899 ad[0] = d[0]*g - d[3]*f;
2900 ad[1] = d[3]*e - d[0]*f;
2901 ad[2] = d[1]*g - d[4]*f;
2902 ad[3] = d[4]*e - d[1]*f;
2903 ad[4] = d[2]*g - d[5]*f;
2904 ad[5] = d[5]*e - d[2]*f;
2913 else if (a.
Width() == 2)
2916 adja(0,1) = -a(0,1);
2917 adja(1,0) = -a(1,0);
2922 adja(0,0) = a(1,1)*a(2,2)-a(1,2)*a(2,1);
2923 adja(0,1) = a(0,2)*a(2,1)-a(0,1)*a(2,2);
2924 adja(0,2) = a(0,1)*a(1,2)-a(0,2)*a(1,1);
2926 adja(1,0) = a(1,2)*a(2,0)-a(1,0)*a(2,2);
2927 adja(1,1) = a(0,0)*a(2,2)-a(0,2)*a(2,0);
2928 adja(1,2) = a(0,2)*a(1,0)-a(0,0)*a(1,2);
2930 adja(2,0) = a(1,0)*a(2,1)-a(1,1)*a(2,0);
2931 adja(2,1) = a(0,1)*a(2,0)-a(0,0)*a(2,1);
2932 adja(2,2) = a(0,0)*a(1,1)-a(0,1)*a(1,0);
2949 else if (a.
Width() == 2)
2951 adjat(0,0) = a(1,1);
2952 adjat(1,0) = -a(0,1);
2953 adjat(0,1) = -a(1,0);
2954 adjat(1,1) = a(0,0);
2958 adjat(0,0) = a(1,1)*a(2,2)-a(1,2)*a(2,1);
2959 adjat(1,0) = a(0,2)*a(2,1)-a(0,1)*a(2,2);
2960 adjat(2,0) = a(0,1)*a(1,2)-a(0,2)*a(1,1);
2962 adjat(0,1) = a(1,2)*a(2,0)-a(1,0)*a(2,2);
2963 adjat(1,1) = a(0,0)*a(2,2)-a(0,2)*a(2,0);
2964 adjat(2,1) = a(0,2)*a(1,0)-a(0,0)*a(1,2);
2966 adjat(0,2) = a(1,0)*a(2,1)-a(1,1)*a(2,0);
2967 adjat(1,2) = a(0,1)*a(2,0)-a(0,0)*a(2,1);
2968 adjat(2,2) = a(0,0)*a(1,1)-a(0,1)*a(1,0);
2980 MFEM_ASSERT(inva.
Height() == a.
Width(),
"incorrect dimensions");
2981 MFEM_ASSERT(inva.
Width() == a.
Height(),
"incorrect dimensions");
2987 const double *d = a.
Data();
2988 double *
id = inva.
Data();
2991 t = 1.0 / (d[0]*d[0] + d[1]*d[1]);
2999 t = 1.0 / (d[0]*d[0] + d[1]*d[1] + d[2]*d[2]);
3007 e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
3008 g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
3009 f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
3010 t = 1.0 / (e*g - f*f);
3011 e *= t; g *= t; f *= t;
3013 id[0] = d[0]*g - d[3]*f;
3014 id[1] = d[3]*e - d[0]*f;
3015 id[2] = d[1]*g - d[4]*f;
3016 id[3] = d[4]*e - d[1]*f;
3017 id[4] = d[2]*g - d[5]*f;
3018 id[5] = d[5]*e - d[2]*f;
3027 cerr <<
"CalcInverse(...) : singular matrix!"
3040 inva(0,0) = a(1,1) * t ;
3041 inva(0,1) = -a(0,1) * t ;
3042 inva(1,0) = -a(1,0) * t ;
3043 inva(1,1) = a(0,0) * t ;
3046 inva(0,0) = (a(1,1)*a(2,2)-a(1,2)*a(2,1))*t;
3047 inva(0,1) = (a(0,2)*a(2,1)-a(0,1)*a(2,2))*t;
3048 inva(0,2) = (a(0,1)*a(1,2)-a(0,2)*a(1,1))*t;
3050 inva(1,0) = (a(1,2)*a(2,0)-a(1,0)*a(2,2))*t;
3051 inva(1,1) = (a(0,0)*a(2,2)-a(0,2)*a(2,0))*t;
3052 inva(1,2) = (a(0,2)*a(1,0)-a(0,0)*a(1,2))*t;
3054 inva(2,0) = (a(1,0)*a(2,1)-a(1,1)*a(2,0))*t;
3055 inva(2,1) = (a(0,1)*a(2,0)-a(0,0)*a(2,1))*t;
3056 inva(2,2) = (a(0,0)*a(1,1)-a(0,1)*a(1,0))*t;
3071 double t = 1. / a.
Det() ;
3076 inva(0,0) = 1.0 / a(0,0);
3079 inva(0,0) = a(1,1) * t ;
3080 inva(1,0) = -a(0,1) * t ;
3081 inva(0,1) = -a(1,0) * t ;
3082 inva(1,1) = a(0,0) * t ;
3085 inva(0,0) = (a(1,1)*a(2,2)-a(1,2)*a(2,1))*t;
3086 inva(1,0) = (a(0,2)*a(2,1)-a(0,1)*a(2,2))*t;
3087 inva(2,0) = (a(0,1)*a(1,2)-a(0,2)*a(1,1))*t;
3089 inva(0,1) = (a(1,2)*a(2,0)-a(1,0)*a(2,2))*t;
3090 inva(1,1) = (a(0,0)*a(2,2)-a(0,2)*a(2,0))*t;
3091 inva(2,1) = (a(0,2)*a(1,0)-a(0,0)*a(1,2))*t;
3093 inva(0,2) = (a(1,0)*a(2,1)-a(1,1)*a(2,0))*t;
3094 inva(1,2) = (a(0,1)*a(2,0)-a(0,0)*a(2,1))*t;
3095 inva(2,2) = (a(0,0)*a(1,1)-a(0,1)*a(1,0))*t;
3110 const double *d = J.
Data();
3118 n(0) = d[1]*d[5] - d[2]*d[4];
3119 n(1) = d[2]*d[3] - d[0]*d[5];
3120 n(2) = d[0]*d[4] - d[1]*d[3];
3126 for (
int i = 0; i < a.
Height(); i++)
3127 for (
int j = 0; j <= i; j++)
3130 for (
int k = 0; k < a.
Width(); k++)
3132 temp += a(i,k) * a(j,k);
3134 aat(j,i) = aat(i,j) = temp;
3140 for (
int i = 0; i < A.
Height(); i++)
3142 for (
int j = 0; j < i; j++)
3145 for (
int k = 0; k < A.
Width(); k++)
3147 t += D(k) * A(i, k) * A(j, k);
3155 for (
int i = 0; i < A.
Height(); i++)
3158 for (
int k = 0; k < A.
Width(); k++)
3160 t += D(k) * A(i, k) * A(i, k);
3168 for (
int i = 0; i < A.
Height(); i++)
3170 for (
int j = 0; j <= i; j++)
3173 for (
int k = 0; k < A.
Width(); k++)
3175 t += D(k) * A(i, k) * A(j, k);
3177 ADAt(j, i) = ADAt(i, j) = t;
3192 #ifdef MFEM_USE_LAPACK
3193 static char transa =
'N', transb =
'T';
3194 static double alpha = 1.0, beta = 0.0;
3197 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &m,
3198 B.
Data(), &n, &beta, ABt.
Data(), &m);
3200 const int ah = A.
Height();
3201 const int bh = B.
Height();
3202 const int aw = A.
Width();
3203 const double *ad = A.
Data();
3204 const double *bd = B.
Data();
3205 double *cd = ABt.
Data();
3207 for (
int i = 0, s = ah*bh; i < s; i++)
3211 for (
int k = 0; k < aw; k++)
3214 for (
int j = 0; j < bh; j++)
3216 const double bjk = bd[j];
3217 for (
int i = 0; i < ah; i++)
3219 cp[i] += ad[i] * bjk;
3227 const int ah = A.
Height();
3228 const int bh = B.
Height();
3229 const int aw = A.
Width();
3230 const double *ad = A.
Data();
3231 const double *bd = B.
Data();
3232 double *cd = ABt.
Data();
3234 for (
int j = 0; j < bh; j++)
3235 for (
int i = 0; i < ah; i++)
3238 const double *ap = ad + i;
3239 const double *bp = bd + j;
3240 for (
int k = 0; k < aw; k++)
3252 for (i = 0; i < A.
Height(); i++)
3253 for (j = 0; j < B.
Height(); j++)
3256 for (k = 0; k < A.
Width(); k++)
3258 d += A(i, k) * B(j, k);
3276 const int ah = A.
Height();
3277 const int bh = B.
Height();
3278 const int aw = A.
Width();
3279 const double *ad = A.
Data();
3280 const double *bd = B.
Data();
3281 const double *dd = D.
GetData();
3282 double *cd = ADBt.
Data();
3284 for (
int i = 0, s = ah*bh; i < s; i++)
3288 for (
int k = 0; k < aw; k++)
3291 for (
int j = 0; j < bh; j++)
3293 const double dk_bjk = dd[k] * bd[j];
3294 for (
int i = 0; i < ah; i++)
3296 cp[i] += ad[i] * dk_bjk;
3315 #ifdef MFEM_USE_LAPACK
3316 static char transa =
'N', transb =
'T';
3317 static double alpha = 1.0, beta = 1.0;
3320 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &m,
3321 B.
Data(), &n, &beta, ABt.
Data(), &m);
3323 const int ah = A.
Height();
3324 const int bh = B.
Height();
3325 const int aw = A.
Width();
3326 const double *ad = A.
Data();
3327 const double *bd = B.
Data();
3328 double *cd = ABt.
Data();
3330 for (
int k = 0; k < aw; k++)
3333 for (
int j = 0; j < bh; j++)
3335 const double bjk = bd[j];
3336 for (
int i = 0; i < ah; i++)
3338 cp[i] += ad[i] * bjk;
3349 for (i = 0; i < A.
Height(); i++)
3350 for (j = 0; j < B.
Height(); j++)
3353 for (k = 0; k < A.
Width(); k++)
3355 d += A(i, k) * B(j, k);
3372 #ifdef MFEM_USE_LAPACK
3373 static char transa =
'T', transb =
'N';
3374 static double alpha = 1.0, beta = 0.0;
3377 dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.
Data(), &k,
3378 B.
Data(), &k, &beta, AtB.
Data(), &m);
3380 const int ah = A.
Height();
3381 const int aw = A.
Width();
3382 const int bw = B.
Width();
3383 const double *ad = A.
Data();
3384 const double *bd = B.
Data();
3385 double *cd = AtB.
Data();
3387 for (
int j = 0; j < bw; j++)
3389 const double *ap = ad;
3390 for (
int i = 0; i < aw; i++)
3393 for (
int k = 0; k < ah; k++)
3406 for (i = 0; i < A.
Width(); i++)
3407 for (j = 0; j < B.
Width(); j++)
3410 for (k = 0; k < A.
Height(); k++)
3412 d += A(k, i) * B(k, j);
3423 for (
int i = 0; i < A.
Height(); i++)
3425 for (
int j = 0; j < i; j++)
3428 for (
int k = 0; k < A.
Width(); k++)
3430 d += A(i,k) * A(j,k);
3432 AAt(i, j) += (d *= a);
3436 for (
int k = 0; k < A.
Width(); k++)
3438 d += A(i,k) * A(i,k);
3446 for (
int i = 0; i < A.
Height(); i++)
3447 for (
int j = 0; j <= i; j++)
3450 for (
int k = 0; k < A.
Width(); k++)
3452 d += A(i,k) * A(j,k);
3454 AAt(i, j) = AAt(j, i) = a * d;
3460 for (
int i = 0; i < v.
Size(); i++)
3461 for (
int j = 0; j <= i; j++)
3463 vvt(i,j) = vvt(j,i) = v(i) * v(j);
3479 for (i = 0; i < v.
Size(); i++)
3482 for (j = 0; j < w.
Size(); j++)
3484 VWt(i, j) = vi * w(j);
3500 for (
int i = 0; i < m; i++)
3503 for (
int j = 0; j < n; j++)
3505 VWt(i, j) += vi * w(j);
3522 for (
int i = 0; i < m; i++)
3524 double avi = a * v(i);
3525 for (
int j = 0; j < n; j++)
3527 VWt(i, j) += avi * w(j);
3543 for (
int i = 0; i < n; i++)
3545 double avi = a * v(i);
3546 for (
int j = 0; j < i; j++)
3548 double avivj = avi * v(j);
3552 VVt(i, i) += avi * v(i);
3559 #ifdef MFEM_USE_LAPACK
3562 MFEM_VERIFY(!info,
"LAPACK: error in DGETRF")
3566 for (
int i = 0; i < m; i++)
3571 double a = std::abs(data[piv+i*m]);
3572 for (
int j = i+1; j < m; j++)
3574 const double b = std::abs(data[j+i*m]);
3585 for (
int j = 0; j < m; j++)
3587 Swap<double>(data[i+j*m], data[piv+j*m]);
3591 MFEM_ASSERT(data[i+i*m] != 0.0,
"division by zero");
3592 const double a_ii_inv = 1.0/data[i+i*m];
3593 for (
int j = i+1; j < m; j++)
3595 data[j+i*m] *= a_ii_inv;
3597 for (
int k = i+1; k < m; k++)
3599 const double a_ik = data[i+k*m];
3600 for (
int j = i+1; j < m; j++)
3602 data[j+k*m] -= a_ik * data[j+i*m];
3614 for (
int k = 0; k < n; k++)
3617 for (
int i = 0; i < m; i++)
3619 double x_i = x[i] * data[i+i*m];
3620 for (
int j = i+1; j < m; j++)
3622 x_i += x[j] * data[i+j*m];
3627 for (
int i = m-1; i >= 0; i--)
3630 for (
int j = 0; j < i; j++)
3632 x_i += x[j] * data[i+j*m];
3637 for (
int i = m-1; i >= 0; i--)
3639 Swap<double>(x[i], x[ipiv[i]-
ipiv_base]);
3650 for (
int k = 0; k < n; k++)
3653 for (
int i = 0; i < m; i++)
3655 Swap<double>(x[i], x[ipiv[i]-
ipiv_base]);
3658 for (
int j = 0; j < m; j++)
3660 const double x_j = x[j];
3661 for (
int i = j+1; i < m; i++)
3663 x[i] -= data[i+j*m] * x_j;
3675 for (
int k = 0; k < n; k++)
3677 for (
int j = m-1; j >= 0; j--)
3679 const double x_j = ( x[j] /= data[j+j*m] );
3680 for (
int i = 0; i < j; i++)
3682 x[i] -= data[i+j*m] * x_j;
3691 #ifdef MFEM_USE_LAPACK
3695 MFEM_VERIFY(!info,
"LAPACK: error in DGETRS")
3710 for (
int k = 0; k < m; k++)
3712 const double minus_x_k = -( x[k] = 1.0/data[k+k*m] );
3713 for (
int i = 0; i < k; i++)
3715 x[i] = data[i+k*m] * minus_x_k;
3717 for (
int j = k-1; j >= 0; j--)
3719 const double x_j = ( x[j] /= data[j+j*m] );
3720 for (
int i = 0; i < j; i++)
3722 x[i] -= data[i+j*m] * x_j;
3730 for (
int j = 0; j < k; j++)
3732 const double minus_L_kj = -data[k+j*m];
3733 for (
int i = 0; i <= j; i++)
3735 X[i+j*m] += X[i+k*m] * minus_L_kj;
3737 for (
int i = j+1; i < m; i++)
3739 X[i+j*m] = X[i+k*m] * minus_L_kj;
3743 for (
int k = m-2; k >= 0; k--)
3745 for (
int j = 0; j < k; j++)
3747 const double L_kj = data[k+j*m];
3748 for (
int i = 0; i < m; i++)
3750 X[i+j*m] -= X[i+k*m] * L_kj;
3755 for (
int k = m-1; k >= 0; k--)
3760 for (
int i = 0; i < m; i++)
3762 Swap<double>(X[i+k*m], X[i+piv_k*m]);
3769 const double *X1,
double *X2)
3772 for (
int k = 0; k < r; k++)
3774 for (
int j = 0; j < m; j++)
3776 const double x1_jk = X1[j+k*m];
3777 for (
int i = 0; i < n; i++)
3779 X2[i+k*n] -= A21[i+j*n] * x1_jk;
3786 int m,
int n,
double *A12,
double *A21,
double *A22)
const
3792 for (
int j = 0; j < m; j++)
3794 const double u_jj_inv = 1.0/data[j+j*m];
3795 for (
int i = 0; i < n; i++)
3797 A21[i+j*n] *= u_jj_inv;
3799 for (
int k = j+1; k < m; k++)
3801 const double u_jk = data[j+k*m];
3802 for (
int i = 0; i < n; i++)
3804 A21[i+k*n] -= A21[i+j*n] * u_jk;
3809 SubMult(m, n, n, A21, A12, A22);
3813 double *B1,
double *B2)
const
3818 SubMult(m, n, r, L21, B1, B2);
3822 const double *X2,
double *Y1)
const
3825 SubMult(n, m, r, U12, X2, Y1);
3834 MFEM_ASSERT(
height ==
width,
"not a square matrix");
3844 MFEM_ASSERT(
height ==
width,
"not a square matrix");
3852 MFEM_ASSERT(a,
"DenseMatrix is not given");
3853 const double *adata = a->data;
3856 lu.
data[i] = adata[i];
3863 MFEM_VERIFY(mat.
height == mat.
width,
"DenseMatrix is not square!");
3879 MFEM_VERIFY(p != NULL,
"Operator is not a DenseMatrix!");
3899 for (
int i = 0; i <
width; i++)
3903 cout <<
"size = " << width <<
", i_max = " << C.
MaxMaxNorm() << endl;
3921 #ifdef MFEM_USE_LAPACK
3927 &qwork, &lwork, &info);
3929 lwork = (int) qwork;
3930 work =
new double[lwork];
3937 if (mat.
Width() != n)
3939 mfem_error(
"DenseMatrixEigensystem::Eval()");
3943 #ifdef MFEM_USE_LAPACK
3946 work, &lwork, &info);
3950 cerr <<
"DenseMatrixEigensystem::Eval(): DSYEV error code: "
3955 mfem_error(
"DenseMatrixEigensystem::Eval(): Compiled without LAPACK");
3961 #ifdef MFEM_USE_LAPACK
3981 void DenseMatrixSVD::Init()
3983 #ifdef MFEM_USE_LAPACK
3992 NULL, &n, &qwork, &lwork, &info);
3994 lwork = (int) qwork;
3995 work =
new double[lwork];
3997 mfem_error(
"DenseMatrixSVD::Init(): Compiled without LAPACK");
4010 #ifdef MFEM_USE_LAPACK
4012 NULL, &n, work, &lwork, &info);
4016 cerr <<
"DenseMatrixSVD::Eval() : info = " << info << endl;
4020 mfem_error(
"DenseMatrixSVD::Eval(): Compiled without LAPACK");
4026 #ifdef MFEM_USE_LAPACK
4036 const int *I = elem_dof.
GetI(), *J = elem_dof.
GetJ(), *dofs;
4037 double *d_col = tdata, *yp = y, x_col;
4038 const double *xp = x;
4042 for (
int i = 0; i < ne; i++)
4045 for (
int col = 0; col < n; col++)
4047 x_col = xp[dofs[col]];
4048 for (
int row = 0; row < n; row++)
4050 yp[dofs[row]] += x_col*d_col[row];
4059 for (
int i = 0; i < ne; i++)
4062 x_col = xp[dofs[0]];
4063 for (
int row = 0; row < n; row++)
4065 ye(row) = x_col*d_col[row];
4068 for (
int col = 1; col < n; col++)
4070 x_col = xp[dofs[col]];
4071 for (
int row = 0; row < n; row++)
4073 ye(row) += x_col*d_col[row];
4077 for (
int row = 0; row < n; row++)
4079 yp[dofs[row]] += ye(row);
virtual void PrintT(std::ostream &out=std::cout, int width_=4) const
Prints the transpose matrix to stream out.
void GetColumn(int c, Vector &col)
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
Logical size of the array.
void SymmetricScaling(const Vector &s)
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
DenseMatrix & operator*=(double c)
bool KernelVector2G(const int &mode, double &d1, double &d12, double &d21, double &d2)
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
DenseMatrix & operator+=(DenseMatrix &m)
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
void Eigenvalues2S(const double &d12, double &d1, double &d2)
void GetScalingFactor(const double &d_max, double &mult)
void SingularValues(Vector &sv) const
void dsyev_Eigensystem(DenseMatrix &a, Vector &ev, DenseMatrix *evect)
void SetSize(int s)
Resizes the vector if the new size is different.
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
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.
void BlockFactor(int m, int n, double *A12, double *A21, double *A22) const
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 CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
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.
void Eval(DenseMatrix &M)
Abstract data type for matrix inverse.
void dgetrs_(char *, int *, int *, double *, int *, int *, double *, int *, int *)
void SetCol(int col, double value)
Set all entries of a column to the specified value.
void Factor()
Factor the current DenseMatrix, *a.
void GetInverseMatrix(int m, double *X) const
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
void CalcOrtho(const DenseMatrix &J, Vector &n)
DenseMatrix & operator=(double c)
Sets the matrix elements equal to constant c.
void vec_normalize3_aux(const double &x1, const double &x2, const double &x3, double &n1, double &n2, double &n3)
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 Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
static void SubMult(int m, int n, int r, const double *A21, const double *X1, double *X2)
void dgetri_(int *N, double *A, int *LDA, int *IPIV, double *WORK, int *LWORK, int *INFO)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication with the inverse of dense matrix.
int KernelVector3G_aux(const int &mode, double &d1, double &d2, double &d3, double &c12, double &c13, double &c23, double &c21, double &c31, double &c32)
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.
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 dgemm_(char *, char *, int *, int *, int *, double *, double *, int *, double *, int *, double *, double *, int *)
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 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 BlockForwSolve(int m, int n, int r, const double *L21, double *B1, double *B2) const
DenseMatrixSVD(DenseMatrix &M)
Abstract data type matrix.
void Norm2(double *v) const
Take the 2-norm of the columns of A and store in v.
void Invert()
Replaces the current matrix with its inverse.
virtual ~DenseMatrixInverse()
Destroys dense inverse matrix.
void LSolve(int m, int n, double *X) const
void LeftScaling(const Vector &s)
LeftScaling this = diag(s) * this.
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.
void vec_normalize3(const double &x1, const double &x2, const double &x3, double &n1, double &n2, double &n3)
virtual void PrintMatlab(std::ostream &out=std::cout) const
~DenseMatrixEigensystem()
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
void Neg()
(*this) = -(*this)
virtual void Print(std::ostream &out=std::cout, int width_=4) const
Prints matrix to stream out.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void Solve(int m, int n, double *X) const
int KernelVector3S(const int &mode, const double &d12, const double &d13, const double &d23, double &d1, double &d2, double &d3)
void AddToVector(int offset, Vector &v) const
Add the matrix 'data' to the Vector 'v' at the given 'offset'.
void AddMult(const Vector &x, Vector &y) const
y += A.x
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||.
void Eigensystem2S(const double &d12, double &d1, double &d2, double &c, double &s)
double MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
void dsyev_(char *JOBZ, char *UPLO, int *N, double *A, int *LDA, double *W, double *WORK, int *LWORK, int *INFO)
double * Data() const
Returns vector of the elements.
void Swap(Array< T > &, Array< T > &)
void Transpose()
(*this) = (*this)^t
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
double Trace() const
Trace of a square matrix.
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
void mfem_error(const char *msg)
int Reduce3S(const int &mode, double &d1, double &d2, double &d3, double &d12, double &d13, double &d23, double &z1, double &z2, double &z3, double &v1, double &v2, double &v3, double &g)
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
void dgetrf_(int *, int *, double *, int *, int *, int *)
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 CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
void SetDataAndSize(double *d, int s)
DenseMatrix & operator-=(DenseMatrix &m)
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
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 ~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 GetDiag(Vector &d)
Returns the diagonal of the matrix.
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
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
static const int ipiv_base
void GradToCurl(DenseMatrix &curl)
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)
DenseMatrixInverse()
Default constructor.
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.
void CalcEigenvalues(double *lambda, double *vec) const
void Getl1Diag(Vector &l)
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
DenseMatrixEigensystem(DenseMatrix &m)
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 MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
void Mult(const double *x, double *y) const
Matrix vector multiplication.
void AddMultTranspose_a(double a, const Vector &x, Vector &y) const
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 InvLeftScaling(const Vector &s)
InvLeftScaling this = diag(1./s) * this.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void SetRow(int row, double value)
Set all entries of a row to the specified value.
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
void AdjustDofDirection(Array< int > &dofs)
void GradToDiv(Vector &div)
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.