21 MFEM_ASSERT(
Op_Real_,
"ComplexDenseMatrix has no real part!");
27 MFEM_ASSERT(
Op_Imag_,
"ComplexDenseMatrix has no imaginary part!");
33 MFEM_ASSERT(
Op_Real_,
"ComplexDenseMatrix has no real part!");
34 return dynamic_cast<const DenseMatrix &
>(*Op_Real_);
39 MFEM_ASSERT(
Op_Imag_,
"ComplexDenseMatrix has no imaginary part!");
40 return dynamic_cast<const DenseMatrix &
>(*Op_Imag_);
57 for (
int j = 0; j<w; j++)
59 for (
int i = 0; i<h; i++)
61 data[i+j*
height] = data_r[i+j*h];
62 data[i+h+(j+h)*
height] = factor*data_r[i+j*h];
69 for (
int j = 0; j<w; j++)
71 for (
int i = 0; i<h; i++)
73 data[i+h+j*
height] = factor*data_i[i+j*h];
74 data[i+(j+h)*
height] = -data_i[i+j*h];
83 MFEM_VERIFY(
height ==
width,
"Matrix has to be square");
88 std::complex<real_t> * data =
new std::complex<real_t>[h*w];
95 for (
int i = 0; i < h*w; i++)
97 data[i] = std::complex<real_t> (data_r[i], data_i[i]);
103 for (
int i = 0; i < h*w; i++)
105 data[i] = std::complex<real_t> (data_r[i], 0.);
111 for (
int i = 0; i < h*w; i++)
113 data[i] = std::complex<real_t> (0., data_i[i]);
118 MFEM_ABORT(
"ComplexDenseMatrix has neither only a real nor an imag part");
121#ifdef MFEM_USE_LAPACK
122 int *ipiv =
new int[w];
124 std::complex<real_t> qwork, *work;
127 MFEM_LAPACK_COMPLEX(getrf_)(&w, &w, data, &w, ipiv, &info);
130 mfem_error(
"DenseMatrix::Invert() : Error in ZGETRF");
133 MFEM_LAPACK_COMPLEX(getri_)(&w, data, &w, ipiv, &qwork, &lwork, &info);
134 lwork = (int) qwork.real();
135 work =
new std::complex<real_t>[lwork];
137 MFEM_LAPACK_COMPLEX(getri_)(&w, data, &w, ipiv, work, &lwork, &info);
140 mfem_error(
"DenseMatrix::Invert() : Error in ZGETRI");
151 std::complex<real_t> ac,bc;
153 for (c = 0; c < n; c++)
155 a = std::abs(data[c+c*h]);
157 for (j = c + 1; j < n; j++)
159 b = std::abs(data[j+c*h]);
168 mfem_error(
"DenseMatrix::Invert() : singular matrix");
171 for (j = 0; j < n; j++)
177 ac = data[c+c*h] = fone / data[c+c*h];
178 for (j = 0; j < c; j++)
182 for (j++; j < n; j++)
186 for (i = 0; i < c; i++)
188 data[i+c*h] = ac * (bc = -data[i+c*h]);
189 for (j = 0; j < c; j++)
191 data[i+j*h] += bc * data[c+j*h];
193 for (j++; j < n; j++)
195 data[i+j*h] += bc * data[c+j*h];
198 for (i++; i < n; i++)
200 data[i+c*h] = ac * (bc = -data[i+c*h]);
201 for (j = 0; j < c; j++)
203 data[i+j*h] += bc * data[c+j*h];
205 for (j++; j < n; j++)
207 data[i+j*h] += bc * data[c+j*h];
212 for (c = n - 1; c >= 0; c--)
215 for (i = 0; i < n; i++)
229 for (
int k = 0; k < h*w; k++)
231 datac_r[k] = data[k].real();
232 datac_i[k] = data[k].imag();
248 MFEM_VERIFY(A.
Width() == B.
Height(),
"Incompatible matrix dimensions");
264 MFEM_VERIFY(C_r || C_i,
"Both real and imag parts are null");
311 MFEM_VERIFY(A.
Height() == B.
Height(),
"Incompatible matrix dimensions");
327 MFEM_VERIFY(C_r || C_i,
"Both real and imag parts are null");
372 std::complex<real_t> * x =
new std::complex<real_t>[m];
375 for (
int i = 0; i<m; i++)
377 x[i] = std::complex<real_t>(x_r[i], x_i[i]);
382 for (
int i = 0; i<m; i++)
384 x[i] = std::complex<real_t>(x_r[i], 0.);
389 for (
int i = 0; i<m; i++)
391 x[i] = std::complex<real_t>(0., x_i[i]);
396 MFEM_ABORT(
"ComplexFactors::RealToComplex:both real and imag part are null");
405 for (
int i = 0; i<m; i++)
407 x_r[i] = x[i].real();
408 x_i[i] = x[i].imag();
415 if (
data) {
return; }
416 MFEM_VERIFY(
data_r ||
data_i,
"ComplexFactors data not set");
424#ifdef MFEM_USE_LAPACK
426 MFEM_VERIFY(
data,
"Matrix data not set");
427 if (m) { MFEM_LAPACK_COMPLEX(getrf_)(&m, &m,
data, &m,
ipiv, &info); }
431 std::complex<real_t> *data_ptr = this->
data;
432 for (
int i = 0; i < m; i++)
437 real_t a = std::abs(data_ptr[piv+i*m]);
438 for (
int j = i+1; j < m; j++)
440 const real_t b = std::abs(data_ptr[j+i*m]);
451 for (
int j = 0; j < m; j++)
458 if (abs(data_ptr[i + i*m]) <= TOL)
464 const std::complex<real_t> a_ii_inv = fone / data_ptr[i+i*m];
465 for (
int j = i+1; j < m; j++)
467 data_ptr[j+i*m] *= a_ii_inv;
469 for (
int k = i+1; k < m; k++)
471 const std::complex<real_t> a_ik = data_ptr[i+k*m];
472 for (
int j = i+1; j < m; j++)
474 data_ptr[j+k*m] -= a_ik * data_ptr[j+i*m];
485 std::complex<real_t> det(1.0,0.);
486 for (
int i=0; i<m; i++)
490 det *= -
data[m * i + i];
494 det *=
data[m * i + i];
503 for (
int k = 0; k < n; k++)
506 for (
int i = 0; i < m; i++)
508 std::complex<real_t> x_i = x[i] *
data[i+i*m];
509 for (
int j = i+1; j < m; j++)
511 x_i += x[j] *
data[i+j*m];
516 for (
int i = m-1; i >= 0; i--)
518 std::complex<real_t> x_i = x[i];
519 for (
int j = 0; j < i; j++)
521 x_i += x[j] *
data[i+j*m];
526 for (
int i = m-1; i >= 0; i--)
539 std::complex<real_t> *x = X;
540 for (
int k = 0; k < n; k++)
543 for (
int i = 0; i < m; i++)
548 for (
int j = 0; j < m; j++)
550 const std::complex<real_t> x_j = x[j];
551 for (
int i = j+1; i < m; i++)
553 x[i] -=
data[i+j*m] * x_j;
565 std::complex<real_t> *x = X;
567 for (
int k = 0; k < n; k++)
569 for (
int j = m-1; j >= 0; j--)
571 const std::complex<real_t> x_j = ( x[j] /=
data[j+j*m] );
572 for (
int i = 0; i < j; i++)
574 x[i] -=
data[i+j*m] * x_j;
585#ifdef MFEM_USE_LAPACK
591 MFEM_LAPACK_COMPLEX(getrs_)(&
trans, &m, &n,
data, &m,
ipiv, x, &m, &info);
593 MFEM_VERIFY(!info,
"LAPACK: error in ZGETRS");
606 std::complex<real_t> * x = X;
607#ifdef MFEM_USE_LAPACK
608 char n_ch =
'N', side =
'R', u_ch =
'U', l_ch =
'L';
611 std::complex<real_t>
alpha(1.0,0.0);
612 MFEM_LAPACK_COMPLEX(trsm_)(&side,&u_ch,&n_ch,&n_ch,&n,&m,&
alpha,
data,&m,X,&n);
613 MFEM_LAPACK_COMPLEX(trsm_)(&side,&l_ch,&n_ch,&u_ch,&n,&m,&
alpha,
data,&m,X,&n);
618 for (
int k = 0; k < n; k++)
620 for (
int j = 0; j < m; j++)
622 const std::complex<real_t> x_j = ( x[j*n] /=
data[j+j*m]);
623 for (
int i = j+1; i < m; i++)
625 x[i*n] -=
data[j + i*m] * x_j;
633 for (
int k = 0; k < n; k++)
635 for (
int j = m-1; j >= 0; j--)
637 const std::complex<real_t> x_j = x[j*n];
638 for (
int i = 0; i < j; i++)
640 x[i*n] -=
data[j + i*m] * x_j;
648 for (
int k = 0; k < n; k++)
650 for (
int i = m-1; i >= 0; --i)
665 std::complex<real_t> * x = X;
667 for (
int k = 0; k < m; k++)
669 const std::complex<real_t> minus_x_k = -( x[k] = fone/
data[k+k*m] );
670 for (
int i = 0; i < k; i++)
672 x[i] =
data[i+k*m] * minus_x_k;
674 for (
int j = k-1; j >= 0; j--)
676 const std::complex<real_t> x_j = ( x[j] /=
data[j+j*m] );
677 for (
int i = 0; i < j; i++)
679 x[i] -=
data[i+j*m] * x_j;
687 for (
int j = 0; j < k; j++)
689 const std::complex<real_t> minus_L_kj = -
data[k+j*m];
690 for (
int i = 0; i <= j; i++)
692 X[i+j*m] += X[i+k*m] * minus_L_kj;
694 for (
int i = j+1; i < m; i++)
696 X[i+j*m] = X[i+k*m] * minus_L_kj;
700 for (
int k = m-2; k >= 0; k--)
702 for (
int j = 0; j < k; j++)
704 const std::complex<real_t> L_kj =
data[k+j*m];
705 for (
int i = 0; i < m; i++)
707 X[i+j*m] -= X[i+k*m] * L_kj;
712 for (
int k = m-1; k >= 0; k--)
717 for (
int i = 0; i < m; i++)
731#ifdef MFEM_USE_LAPACK
734 MFEM_VERIFY(
data,
"Matrix data not set");
735 if (m) { MFEM_LAPACK_COMPLEX(potrf_)(&uplo, &m,
data, &m, &info); }
740 for (
int j = 0; j<m; j++)
742 std::complex<real_t>
a(0.,0.);
743 for (
int k = 0; k<j; k++)
748 MFEM_VERIFY((
data[j+j*m] -
a).real() > 0.,
749 "CholeskyFactors::Factor: The matrix is not SPD");
751 data[j+j*m] = std::sqrt((
data[j+j*m] -
a).real());
753 if (
data[j + j*m].real() <= TOL) {
return false; }
755 for (
int i = j+1; i<m; i++)
757 a = std::complex<real_t>(0.,0.);
758 for (
int k = 0; k<j; k++)
771 std::complex<real_t> det(1.0,0.0);
772 for (
int i=0; i<m; i++)
774 det *=
data[i + i*m];
784 std::complex<real_t> *x = X;
785 for (
int k = 0; k < n; k++)
787 for (
int j = m-1; j >= 0; j--)
789 std::complex<real_t> x_j = x[j] *
data[j+j*m];
790 for (
int i = 0; i < j; i++)
792 x_j += x[i] *
data[j+i*m];
806 std::complex<real_t> *x = X;
807 for (
int k = 0; k < n; k++)
809 for (
int i = 0; i < m; i++)
811 std::complex<real_t> x_i = x[i] *
data[i+i*m];
812 for (
int j = i+1; j < m; j++)
814 x_i += x[j] * std::conj(
data[j+i*m]);
828 std::complex<real_t> *x = X;
829#ifdef MFEM_USE_LAPACK
835 MFEM_LAPACK_COMPLEX(trtrs_)(&uplo, &
trans, &diag, &m, &n,
data, &m, x, &m,
837 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:LSolve:: info");
839 for (
int k = 0; k < n; k++)
842 for (
int j = 0; j < m; j++)
844 const std::complex<real_t> x_j = (x[j] /=
data[j+j*m]);
845 for (
int i = j+1; i < m; i++)
847 x[i] -=
data[i+j*m] * x_j;
861 std::complex<real_t> *x = X;
862#ifdef MFEM_USE_LAPACK
869 MFEM_LAPACK_COMPLEX(trtrs_)(&uplo, &
trans, &diag, &m, &n,
data, &m, x, &m,
871 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:USolve:: info");
874 for (
int k = 0; k < n; k++)
876 for (
int j = m-1; j >= 0; j--)
878 const std::complex<real_t> x_j = ( x[j] /=
data[j+j*m] );
879 for (
int i = 0; i < j; i++)
881 x[i] -= std::conj(
data[j+i*m]) * x_j;
894#ifdef MFEM_USE_LAPACK
898 MFEM_LAPACK_COMPLEX(potrs_)(&uplo, &m, &n,
data, &m, x, &m, &info);
899 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:Solve:: info");
913 std::complex<real_t> *x = X;
914#ifdef MFEM_USE_LAPACK
921 std::complex<real_t>
alpha(1.0,0.0);
924 MFEM_LAPACK_COMPLEX(trsm_)(&side,&uplo,&transt,&diag,&n,&m,&
alpha,
data,&m,x,&n);
925 MFEM_LAPACK_COMPLEX(trsm_)(&side,&uplo,&
trans,&diag,&n,&m,&
alpha,
data,&m,x,&n);
929 for (
int k = 0; k < n; k++)
931 for (
int j = 0; j < m; j++)
933 const std::complex<real_t> x_j = ( x[j*n] /=
data[j+j*m]);
934 for (
int i = j+1; i < m; i++)
936 x[i*n] -= std::conj(
data[i + j*m]) * x_j;
943 for (
int k = 0; k < n; k++)
945 for (
int j = m-1; j >= 0; j--)
947 const std::complex<real_t> x_j = (x[j*n] /=
data[j+j*m]);
948 for (
int i = 0; i < j; i++)
950 x[i*n] -=
data[j + i*m] * x_j;
964 std::complex<real_t> * X =
new std::complex<real_t>[m*m];
965#ifdef MFEM_USE_LAPACK
967 for (
int i = 0; i<m; i++)
969 for (
int j = i; j<m; j++)
971 X[j+i*m] =
data[j+i*m];
976 MFEM_LAPACK_COMPLEX(potri_)(&uplo, &m, X, &m, &info);
977 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:GetInverseMatrix:: info");
979 for (
int i = 0; i<m; i++)
981 for (
int j = i+1; j<m; j++)
983 X[i+j*m] = std::conj(X[j+i*m]);
989 for (
int k = 0; k<m; k++)
991 X[k+k*m] = fone/
data[k+k*m];
992 for (
int i = k+1; i < m; i++)
994 std::complex<real_t> s(0.,0.);
995 for (
int j=k; j<i; j++)
997 s -=
data[i+j*m] * X[j+k*m]/
data[i+i*m];
1002 for (
int i = 0; i < m; i++)
1004 for (
int j = i; j < m; j++)
1006 std::complex<real_t> s(0.,0.);
1007 for (
int k=j; k<m; k++)
1009 s += X[k+i*m] * std::conj(X[k+j*m]);
1012 X[i+j*m] = std::conj(s);
void UMult(int m, int n, real_t *X_r, real_t *X_i) const
void LSolve(int m, int n, real_t *X_r, real_t *X_i) const
void USolve(int m, int n, real_t *X_r, real_t *X_i) const
std::complex< real_t > Det(int m) const override
void RightSolve(int m, int n, real_t *X_r, real_t *X_i) const
bool Factor(int m, real_t TOL=0.0) override
Compute the Cholesky factorization of the current matrix.
void GetInverseMatrix(int m, real_t *X_r, real_t *X_i) const override
Assuming L.L^H = A factored data of size (m x m), compute X <- A^{-1}.
void Solve(int m, int n, real_t *X_r, real_t *X_i) const override
void LMult(int m, int n, real_t *X_r, real_t *X_i) const
Specialization of the ComplexOperator built from a pair of Dense Matrices. The purpose of this specia...
DenseMatrix & imag() override
ComplexDenseMatrix(DenseMatrix *A_Real, DenseMatrix *A_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
DenseMatrix * GetSystemMatrix() const
DenseMatrix & real() override
Real or imaginary part accessor methods.
ComplexDenseMatrix * ComputeInverse()
std::complex< real_t > * data
void ComplexToReal(int m, const std::complex< real_t > *x, real_t *x_r, real_t *x_i) const
void SetComplexData(int m)
std::complex< real_t > * RealToComplex(int m, const real_t *x_r, const real_t *x_i) const
void LSolve(int m, int n, real_t *X_r, real_t *X_i) const
void USolve(int m, int n, real_t *X_r, real_t *X_i) const
std::complex< real_t > Det(int m) const override
void RightSolve(int m, int n, real_t *X_r, real_t *X_i) const
void Solve(int m, int n, real_t *X_r, real_t *X_i) const override
static const int ipiv_base
bool Factor(int m, real_t TOL=0.0) override
Compute the LU factorization of the current matrix.
void GetInverseMatrix(int m, real_t *X_r, real_t *X_i) const override
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
void Mult(int m, int n, real_t *X_r, real_t *X_i) const
bool hasRealPart() const
Check for existence of real or imaginary part of the operator.
@ HERMITIAN
Native convention for Hermitian operators.
Data type dense matrix using column-major storage.
real_t * Data() const
Returns the matrix data array.
int width
Dimension of the input / number of columns in the matrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int height
Dimension of the output / number of rows in the matrix.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void trans(const Vector &u, Vector &x)
void Swap(Array< T > &, Array< T > &)
void mfem_error(const char *msg)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void AddMult_a(real_t alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.