12 #ifndef MFEM_COMPLEX_DENSEMAT
13 #define MFEM_COMPLEX_DENSEMAT
53 ComplexDenseMatrix *
Mult(
const ComplexDenseMatrix &B,
54 const ComplexDenseMatrix &C);
57 ComplexDenseMatrix *
MultAtB(
const ComplexDenseMatrix &A,
58 const ComplexDenseMatrix &B);
68 std::complex<double> *
RealToComplex(
int m,
const double * x_r,
69 const double * x_i)
const;
71 void ComplexToReal(
int m,
const std::complex<double> * x,
double * x_r,
78 std::complex<double> *
data =
nullptr;
93 virtual bool Factor(
int m,
double TOL = 0.0)
95 mfem_error(
"ComplexFactors::ComplexFactors(...)");
99 virtual std::complex<double>
Det(
int m)
const
105 virtual void Solve(
int m,
int n,
double *X_r,
double * X_i)
const
127 #ifdef MFEM_USE_LAPACK
151 virtual bool Factor(
int m,
double TOL = 0.0);
155 virtual std::complex<double>
Det(
int m)
const;
159 void Mult(
int m,
int n,
double *X_r,
double * X_i)
const;
161 void Mult(
int m,
int n, std::complex<double> *X)
const;
165 void LSolve(
int m,
int n,
double *X_r,
double *X_i)
const;
169 void USolve(
int m,
int n,
double *X_r,
double *X_i)
const;
173 virtual void Solve(
int m,
int n,
double *X_r,
double *X_i)
const;
177 void RightSolve(
int m,
int n,
double *X_r,
double *X_i)
const;
209 virtual bool Factor(
int m,
double TOL = 0.0);
213 virtual std::complex<double>
Det(
int m)
const;
217 void LMult(
int m,
int n,
double *X_r,
double * X_i)
const;
221 void UMult(
int m,
int n,
double *X_r,
double *X_i)
const;
225 void LSolve(
int m,
int n,
double *X_r,
double * X_i)
const;
229 void USolve(
int m,
int n,
double *X_r,
double *X_i)
const;
233 virtual void Solve(
int m,
int n,
double *X_r,
double * X_i)
const;
237 void RightSolve(
int m,
int n,
double *X_r,
double *X_i)
const;
246 #endif // MFEM_COMPLEX_DENSEMAT
ComplexDenseMatrix * ComputeInverse()
void LMult(int m, int n, double *X_r, double *X_i) const
virtual DenseMatrix & real()
Real or imaginary part accessor methods.
void LSolve(int m, int n, double *X_r, double *X_i) const
void ResetComplexData(int m)
virtual bool Factor(int m, double TOL=0.0)
Compute the Cholesky factorization of the current matrix.
virtual ~ComplexFactors()
Mimic the action of a complex operator using two real operators.
virtual void Solve(int m, int n, double *X_r, double *X_i) const
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Data type dense matrix using column-major storage.
std::complex< double > * RealToComplex(int m, const double *x_r, const double *x_i) const
void UMult(int m, int n, double *X_r, double *X_i) const
virtual void GetInverseMatrix(int m, double *X_r, double *X_i) const
std::complex< double > * data
void RightSolve(int m, int n, double *X_r, double *X_i) const
ComplexLUFactors(double *data_r_, double *data_i, int *ipiv_)
void RightSolve(int m, int n, double *X_r, double *X_i) const
void SetComplexData(int m)
ID for class ComplexDenseMatrix.
ComplexCholeskyFactors(double *data_r_, double *data_i_)
void USolve(int m, int n, double *X_r, double *X_i) const
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void Mult(int m, int n, double *X_r, double *X_i) const
virtual std::complex< double > Det(int m) const
virtual Type GetType() const
virtual bool Factor(int m, double TOL=0.0)
DenseMatrix * GetSystemMatrix() const
void LSolve(int m, int n, double *X_r, double *X_i) const
Type
Enumeration defining IDs for some classes derived from Operator.
ComplexDenseMatrix(DenseMatrix *A_Real, DenseMatrix *A_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
virtual void GetInverseMatrix(int m, double *X_r, double *X_i) const
Assuming L.L^H = A factored data of size (m x m), compute X <- A^{-1}.
virtual void Solve(int m, int n, double *X_r, double *X_i) const
virtual void GetInverseMatrix(int m, double *X_r, double *X_i) const
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
ComplexFactors(double *data_r_, double *data_i_)
virtual std::complex< double > Det(int m) const
void ComplexToReal(int m, const std::complex< double > *x, double *x_r, double *x_i) const
void USolve(int m, int n, double *X_r, double *X_i) const
static const int ipiv_base
virtual void Solve(int m, int n, double *X_r, double *X_i) const
Native convention for Hermitian operators.
virtual std::complex< double > Det(int m) const
virtual DenseMatrix & imag()
Specialization of the ComplexOperator built from a pair of Dense Matrices. The purpose of this specia...