15 #ifdef MFEM_USE_LAPACK 17 zgetrf_(
int *,
int *, std::complex<double> *,
int *,
int *,
int *);
19 zgetrs_(
char *,
int *,
int *, std::complex<double> *,
int *,
int *,
20 std::complex<double> *,
int *,
int *);
22 zgetri_(
int *, std::complex<double> *,
int *,
int *,
23 std::complex<double> *,
int *,
int *);
25 ztrsm_(
char *,
char *,
char *,
char *,
int *,
int *, std::complex<double> *,
26 std::complex<double> *,
int *, std::complex<double> *,
int *);
28 zpotrf_(
char *,
int *, std::complex<double> *,
int *,
int *);
31 ztrtrs_(
char *,
char*,
char *,
int *,
int *, std::complex<double> *,
int *,
32 std::complex<double> *,
int *,
int *);
34 zpotri_(
char *,
int *, std::complex<double> *,
int*,
int *);
37 zpotrs_(
char *,
int *,
int *, std::complex<double> *,
int *,
38 std::complex<double> *,
int *,
int *);
46 MFEM_ASSERT(
Op_Real_,
"ComplexDenseMatrix has no real part!");
52 MFEM_ASSERT(
Op_Imag_,
"ComplexDenseMatrix has no imaginary part!");
58 MFEM_ASSERT(
Op_Real_,
"ComplexDenseMatrix has no real part!");
59 return dynamic_cast<const DenseMatrix &
>(*Op_Real_);
64 MFEM_ASSERT(
Op_Imag_,
"ComplexDenseMatrix has no imaginary part!");
65 return dynamic_cast<const DenseMatrix &
>(*Op_Imag_);
73 double * data = A->
Data();
74 double * data_r =
nullptr;
75 double * data_i =
nullptr;
82 for (
int j = 0; j<w; j++)
84 for (
int i = 0; i<h; i++)
86 data[i+j*
height] = data_r[i+j*h];
87 data[i+h+(j+h)*
height] = factor*data_r[i+j*h];
94 for (
int j = 0; j<w; j++)
96 for (
int i = 0; i<h; i++)
98 data[i+h+j*
height] = factor*data_i[i+j*h];
99 data[i+(j+h)*
height] = -data_i[i+j*h];
108 MFEM_VERIFY(
height ==
width,
"Matrix has to be square");
113 std::complex<double> * data =
new std::complex<double>[h*w];
120 for (
int i = 0; i < h*w; i++)
122 data[i] = std::complex<double> (data_r[i], data_i[i]);
128 for (
int i = 0; i < h*w; i++)
130 data[i] = std::complex<double> (data_r[i], 0.);
136 for (
int i = 0; i < h*w; i++)
138 data[i] = std::complex<double> (0., data_i[i]);
143 MFEM_ABORT(
"ComplexDenseMatrix has neither only a real nor an imag part");
146 #ifdef MFEM_USE_LAPACK 147 int *ipiv =
new int[w];
149 std::complex<double> qwork, *work;
152 zgetrf_(&w, &w, data, &w, ipiv, &info);
155 mfem_error(
"DenseMatrix::Invert() : Error in ZGETRF");
158 zgetri_(&w, data, &w, ipiv, &qwork, &lwork, &info);
160 lwork = (int) qwork.real();
161 work =
new std::complex<double>[lwork];
163 zgetri_(&w, data, &w, ipiv, work, &lwork, &info);
167 mfem_error(
"DenseMatrix::Invert() : Error in ZGETRI");
178 std::complex<double> ac,bc;
180 for (c = 0; c < n; c++)
182 a = std::abs(data[c+c*h]);
184 for (j = c + 1; j < n; j++)
186 b = std::abs(data[j+c*h]);
195 mfem_error(
"DenseMatrix::Invert() : singular matrix");
198 for (j = 0; j < n; j++)
200 mfem::Swap<std::complex<double>>(data[c+j*h], data[i+j*h]);
203 ac = data[c+c*h] = 1.0 / data[c+c*h];
204 for (j = 0; j < c; j++)
208 for (j++; j < n; j++)
212 for (i = 0; i < c; i++)
214 data[i+c*h] = ac * (bc = -data[i+c*h]);
215 for (j = 0; j < c; j++)
217 data[i+j*h] += bc * data[c+j*h];
219 for (j++; j < n; j++)
221 data[i+j*h] += bc * data[c+j*h];
224 for (i++; i < n; i++)
226 data[i+c*h] = ac * (bc = -data[i+c*h]);
227 for (j = 0; j < c; j++)
229 data[i+j*h] += bc * data[c+j*h];
231 for (j++; j < n; j++)
233 data[i+j*h] += bc * data[c+j*h];
238 for (c = n - 1; c >= 0; c--)
241 for (i = 0; i < n; i++)
243 mfem::Swap<std::complex<double>>(data[i+c*h], data[i+j*h]);
252 double * datac_r = C_r->
Data();
253 double * datac_i = C_i->
Data();
255 for (
int k = 0; k < h*w; k++)
257 datac_r[k] = data[k].real();
258 datac_i[k] = data[k].imag();
274 MFEM_VERIFY(A.
Width() == B.
Height(),
"Incompatible matrix dimensions");
290 MFEM_VERIFY(C_r || C_i,
"Both real and imag parts are null");
337 MFEM_VERIFY(A.
Height() == B.
Height(),
"Incompatible matrix dimensions");
353 MFEM_VERIFY(C_r || C_i,
"Both real and imag parts are null");
396 (
int m,
const double * x_r,
const double * x_i)
const 398 std::complex<double> * x =
new std::complex<double>[m];
401 for (
int i = 0; i<m; i++)
403 x[i] = std::complex<double>(x_r[i], x_i[i]);
408 for (
int i = 0; i<m; i++)
410 x[i] = std::complex<double>(x_r[i], 0.);
415 for (
int i = 0; i<m; i++)
417 x[i] = std::complex<double>(0., x_i[i]);
422 MFEM_ABORT(
"ComplexFactors::RealToComplex:both real and imag part are null");
429 double * x_r,
double * x_i)
const 431 for (
int i = 0; i<m; i++)
433 x_r[i] = x[i].real();
434 x_i[i] = x[i].imag();
441 if (
data) {
return; }
442 MFEM_VERIFY(
data_r ||
data_i,
"ComplexFactors data not set");
450 #ifdef MFEM_USE_LAPACK 452 MFEM_VERIFY(
data,
"Matrix data not set");
457 std::complex<double> *data_ptr = this->
data;
458 for (
int i = 0; i < m; i++)
463 double a = std::abs(data_ptr[piv+i*m]);
464 for (
int j = i+1; j < m; j++)
466 const double b = std::abs(data_ptr[j+i*m]);
477 for (
int j = 0; j < m; j++)
479 mfem::Swap<std::complex<double>>(data_ptr[i+j*m], data_ptr[piv+j*m]);
484 if (abs(data_ptr[i + i*m]) <= TOL)
489 const std::complex<double> a_ii_inv = 1.0 / data_ptr[i+i*m];
490 for (
int j = i+1; j < m; j++)
492 data_ptr[j+i*m] *= a_ii_inv;
494 for (
int k = i+1; k < m; k++)
496 const std::complex<double> a_ik = data_ptr[i+k*m];
497 for (
int j = i+1; j < m; j++)
499 data_ptr[j+k*m] -= a_ik * data_ptr[j+i*m];
510 std::complex<double> det(1.0,0.);
511 for (
int i=0; i<m; i++)
515 det *= -
data[m * i + i];
519 det *=
data[m * i + i];
528 for (
int k = 0; k < n; k++)
531 for (
int i = 0; i < m; i++)
533 std::complex<double> x_i = x[i] *
data[i+i*m];
534 for (
int j = i+1; j < m; j++)
536 x_i += x[j] *
data[i+j*m];
541 for (
int i = m-1; i >= 0; i--)
543 std::complex<double> x_i = x[i];
544 for (
int j = 0; j < i; j++)
546 x_i += x[j] *
data[i+j*m];
551 for (
int i = m-1; i >= 0; i--)
553 mfem::Swap<std::complex<double>>(x[i], x[
ipiv[i]-
ipiv_base]);
564 std::complex<double> *x = X;
565 for (
int k = 0; k < n; k++)
568 for (
int i = 0; i < m; i++)
570 mfem::Swap<std::complex<double>>(x[i], x[
ipiv[i]-
ipiv_base]);
573 for (
int j = 0; j < m; j++)
575 const std::complex<double> x_j = x[j];
576 for (
int i = j+1; i < m; i++)
578 x[i] -=
data[i+j*m] * x_j;
590 std::complex<double> *x = X;
592 for (
int k = 0; k < n; k++)
594 for (
int j = m-1; j >= 0; j--)
596 const std::complex<double> x_j = ( x[j] /=
data[j+j*m] );
597 for (
int i = 0; i < j; i++)
599 x[i] -=
data[i+j*m] * x_j;
610 #ifdef MFEM_USE_LAPACK 615 MFEM_VERIFY(!info,
"LAPACK: error in ZGETRS");
628 std::complex<double> * x = X;
629 #ifdef MFEM_USE_LAPACK 630 char n_ch =
'N', side =
'R', u_ch =
'U', l_ch =
'L';
633 std::complex<double>
alpha(1.0,0.0);
634 ztrsm_(&side,&u_ch,&n_ch,&n_ch,&n,&m,&
alpha,
data,&m,X,&n);
635 ztrsm_(&side,&l_ch,&n_ch,&u_ch,&n,&m,&
alpha,
data,&m,X,&n);
640 for (
int k = 0; k < n; k++)
642 for (
int j = 0; j < m; j++)
644 const std::complex<double> x_j = ( x[j*n] /=
data[j+j*m]);
645 for (
int i = j+1; i < m; i++)
647 x[i*n] -=
data[j + i*m] * x_j;
655 for (
int k = 0; k < n; k++)
657 for (
int j = m-1; j >= 0; j--)
659 const std::complex<double> x_j = x[j*n];
660 for (
int i = 0; i < j; i++)
662 x[i*n] -=
data[j + i*m] * x_j;
670 for (
int k = 0; k < n; k++)
672 for (
int i = m-1; i >= 0; --i)
674 mfem::Swap<std::complex<double>>(x[i*n], x[(
ipiv[i]-
ipiv_base)*n]);
687 std::complex<double> * x = X;
688 for (
int k = 0; k < m; k++)
690 const std::complex<double> minus_x_k = -( x[k] = 1.0/
data[k+k*m] );
691 for (
int i = 0; i < k; i++)
693 x[i] =
data[i+k*m] * minus_x_k;
695 for (
int j = k-1; j >= 0; j--)
697 const std::complex<double> x_j = ( x[j] /=
data[j+j*m] );
698 for (
int i = 0; i < j; i++)
700 x[i] -=
data[i+j*m] * x_j;
708 for (
int j = 0; j < k; j++)
710 const std::complex<double> minus_L_kj = -
data[k+j*m];
711 for (
int i = 0; i <= j; i++)
713 X[i+j*m] += X[i+k*m] * minus_L_kj;
715 for (
int i = j+1; i < m; i++)
717 X[i+j*m] = X[i+k*m] * minus_L_kj;
721 for (
int k = m-2; k >= 0; k--)
723 for (
int j = 0; j < k; j++)
725 const std::complex<double> L_kj =
data[k+j*m];
726 for (
int i = 0; i < m; i++)
728 X[i+j*m] -= X[i+k*m] * L_kj;
733 for (
int k = m-1; k >= 0; k--)
738 for (
int i = 0; i < m; i++)
740 Swap<std::complex<double>>(X[i+k*m], X[i+piv_k*m]);
752 #ifdef MFEM_USE_LAPACK 755 MFEM_VERIFY(
data,
"Matrix data not set");
760 for (
int j = 0; j<m; j++)
762 std::complex<double>
a(0.,0.);
763 for (
int k = 0; k<j; k++)
768 MFEM_VERIFY((
data[j+j*m] -
a).real() > 0.,
769 "CholeskyFactors::Factor: The matrix is not SPD");
771 data[j+j*m] = std::sqrt((
data[j+j*m] -
a).real());
773 if (
data[j + j*m].real() <= TOL) {
return false; }
775 for (
int i = j+1; i<m; i++)
777 a = std::complex<double>(0.,0.);
778 for (
int k = 0; k<j; k++)
791 std::complex<double> det(1.0,0.0);
792 for (
int i=0; i<m; i++)
794 det *=
data[i + i*m];
804 std::complex<double> *x = X;
805 for (
int k = 0; k < n; k++)
807 for (
int j = m-1; j >= 0; j--)
809 std::complex<double> x_j = x[j] *
data[j+j*m];
810 for (
int i = 0; i < j; i++)
812 x_j += x[i] *
data[j+i*m];
826 std::complex<double> *x = X;
827 for (
int k = 0; k < n; k++)
829 for (
int i = 0; i < m; i++)
831 std::complex<double> x_i = x[i] *
data[i+i*m];
832 for (
int j = i+1; j < m; j++)
834 x_i += x[j] * std::conj(
data[j+i*m]);
848 std::complex<double> *x = X;
849 #ifdef MFEM_USE_LAPACK 856 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:LSolve:: info");
858 for (
int k = 0; k < n; k++)
861 for (
int j = 0; j < m; j++)
863 const std::complex<double> x_j = (x[j] /=
data[j+j*m]);
864 for (
int i = j+1; i < m; i++)
866 x[i] -=
data[i+j*m] * x_j;
880 std::complex<double> *x = X;
881 #ifdef MFEM_USE_LAPACK 889 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:USolve:: info");
892 for (
int k = 0; k < n; k++)
894 for (
int j = m-1; j >= 0; j--)
896 const std::complex<double> x_j = ( x[j] /=
data[j+j*m] );
897 for (
int i = 0; i < j; i++)
899 x[i] -= std::conj(
data[j+i*m]) * x_j;
912 #ifdef MFEM_USE_LAPACK 917 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:Solve:: info");
931 std::complex<double> *x = X;
932 #ifdef MFEM_USE_LAPACK 939 std::complex<double>
alpha(1.0,0.0);
942 ztrsm_(&side,&uplo,&transt,&diag,&n,&m,&
alpha,
data,&m,x,&n);
943 ztrsm_(&side,&uplo,&
trans,&diag,&n,&m,&
alpha,
data,&m,x,&n);
947 for (
int k = 0; k < n; k++)
949 for (
int j = 0; j < m; j++)
951 const std::complex<double> x_j = ( x[j*n] /=
data[j+j*m]);
952 for (
int i = j+1; i < m; i++)
954 x[i*n] -= std::conj(
data[i + j*m]) * x_j;
961 for (
int k = 0; k < n; k++)
963 for (
int j = m-1; j >= 0; j--)
965 const std::complex<double> x_j = (x[j*n] /=
data[j+j*m]);
966 for (
int i = 0; i < j; i++)
968 x[i*n] -=
data[j + i*m] * x_j;
982 std::complex<double> * X =
new std::complex<double>[m*m];
983 #ifdef MFEM_USE_LAPACK 985 for (
int i = 0; i<m; i++)
987 for (
int j = i; j<m; j++)
989 X[j+i*m] =
data[j+i*m];
994 zpotri_(&uplo, &m, X, &m, &info);
995 MFEM_VERIFY(!info,
"ComplexCholeskyFactors:GetInverseMatrix:: info");
997 for (
int i = 0; i<m; i++)
999 for (
int j = i+1; j<m; j++)
1001 X[i+j*m] = std::conj(X[j+i*m]);
1006 for (
int k = 0; k<m; k++)
1008 X[k+k*m] = 1./
data[k+k*m];
1009 for (
int i = k+1; i < m; i++)
1011 std::complex<double>
s(0.,0.);
1012 for (
int j=k; j<i; j++)
1014 s -=
data[i+j*m] * X[j+k*m]/
data[i+i*m];
1019 for (
int i = 0; i < m; i++)
1021 for (
int j = i; j < m; j++)
1023 std::complex<double>
s(0.,0.);
1024 for (
int k=j; k<m; k++)
1026 s += X[k+i*m] * std::conj(X[k+j*m]);
1029 X[i+j*m] = std::conj(
s);
ComplexDenseMatrix * ComputeInverse()
void trans(const Vector &u, Vector &x)
virtual DenseMatrix & real()
Real or imaginary part accessor methods.
void USolve(int m, int n, double *X_r, double *X_i) const
virtual bool Factor(int m, double TOL=0.0)
Compute the Cholesky factorization of the current matrix.
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}.
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 bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Data type dense matrix using column-major storage.
virtual std::complex< double > Det(int m) const
double * Data() const
Returns the matrix data array.
bool hasRealPart() const
Check for existence of real or imaginary part of the operator.
std::complex< double > * data
virtual void Solve(int m, int n, double *X_r, double *X_i) const
void SetComplexData(int m)
void UMult(int m, int n, double *X_r, double *X_i) const
virtual void Solve(int m, int n, double *X_r, double *X_i) const
void ztrtrs_(char *, char *, char *, int *, int *, std::complex< double > *, int *, std::complex< double > *, int *, int *)
DenseMatrix * GetSystemMatrix() const
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
virtual std::complex< double > Det(int m) const
void USolve(int m, int n, double *X_r, double *X_i) const
void ztrsm_(char *, char *, char *, char *, int *, int *, std::complex< double > *, std::complex< double > *, int *, std::complex< double > *, int *)
void RightSolve(int m, int n, double *X_r, double *X_i) const
void LSolve(int m, int n, double *X_r, double *X_i) const
void zgetrs_(char *, int *, int *, std::complex< double > *, int *, int *, std::complex< double > *, int *, int *)
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.
void zgetrf_(int *, int *, std::complex< double > *, int *, int *, int *)
void LMult(int m, int n, double *X_r, double *X_i) const
void ComplexToReal(int m, const std::complex< double > *x, double *x_r, double *x_i) const
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.
void RightSolve(int m, int n, double *X_r, double *X_i) const
std::complex< double > * RealToComplex(int m, const double *x_r, const double *x_i) const
static const int ipiv_base
Native convention for Hermitian operators.
void LSolve(int m, int n, double *X_r, double *X_i) const
void zpotri_(char *, int *, std::complex< double > *, int *, int *)
void zpotrs_(char *, int *, int *, std::complex< double > *, int *, std::complex< double > *, int *, int *)
void Mult(int m, int n, double *X_r, double *X_i) const
void zgetri_(int *, std::complex< double > *, int *, int *, std::complex< double > *, int *, int *)
void zpotrf_(char *, int *, std::complex< double > *, int *, int *)
void AddMult_a(double alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
virtual DenseMatrix & imag()
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}.
int width
Dimension of the input / number of columns in the matrix.
Specialization of the ComplexOperator built from a pair of Dense Matrices. The purpose of this specia...