20 bool ownReal,
bool ownImag,
22 :
Operator(2*((Op_Real)?Op_Real->Height():Op_Imag->Height()),
23 2*((Op_Real)?Op_Real->Width():Op_Imag->Width()))
28 , convention_(convention)
47 MFEM_ASSERT(
Op_Real_,
"ComplexOperator has no real part!");
53 MFEM_ASSERT(
Op_Imag_,
"ComplexOperator has no imaginary part!");
59 MFEM_ASSERT(
Op_Real_,
"ComplexOperator has no real part!");
65 MFEM_ASSERT(
Op_Imag_,
"ComplexOperator has no imaginary part!");
170 MFEM_ASSERT(
Op_Real_,
"ComplexSparseMatrix has no real part!");
176 MFEM_ASSERT(
Op_Imag_,
"ComplexSparseMatrix has no imaginary part!");
182 MFEM_ASSERT(
Op_Real_,
"ComplexSparseMatrix has no real part!");
188 MFEM_ASSERT(
Op_Imag_,
"ComplexSparseMatrix has no imaginary part!");
197 const int nrows_r = (A_r)?A_r->
Height():0;
198 const int nrows_i = (A_i)?A_i->
Height():0;
199 const int nrows = std::max(nrows_r, nrows_i);
201 const int *I_r = (A_r)?A_r->
GetI():NULL;
202 const int *I_i = (A_i)?A_i->
GetI():NULL;
204 const int *J_r = (A_r)?A_r->
GetJ():NULL;
205 const int *J_i = (A_i)?A_i->
GetJ():NULL;
207 const double *D_r = (A_r)?A_r->
GetData():NULL;
208 const double *D_i = (A_i)?A_i->
GetData():NULL;
210 const int nnz_r = (I_r)?I_r[nrows]:0;
211 const int nnz_i = (I_i)?I_i[nrows]:0;
212 const int nnz = 2 * (nnz_r + nnz_i);
221 I[nrows] = nnz_r + nnz_i;
222 for (
int i=0; i<nrows; i++)
224 I[i + 1] = ((I_r)?I_r[i+1]:0) + ((I_i)?I_i[i+1]:0);
225 I[i + nrows + 1] = I[i+1] + nnz_r + nnz_i;
229 const int off_i = (I_i)?(I_i[i+1] - I_i[i]):0;
230 for (
int j=0; j<I_r[i+1] - I_r[i]; j++)
232 J[I[i] + j] = J_r[I_r[i] + j];
233 D[I[i] + j] = D_r[I_r[i] + j];
235 J[I[i+nrows] + off_i + j] = J_r[I_r[i] + j] + nrows;
236 D[I[i+nrows] + off_i + j] = factor*D_r[I_r[i] + j];
241 const int off_r = (I_r)?(I_r[i+1] - I_r[i]):0;
242 for (
int j=0; j<I_i[i+1] - I_i[i]; j++)
244 J[I[i] + off_r + j] = J_i[I_i[i] + j] + nrows;
245 D[I[i] + off_r + j] = -D_i[I_i[i] + j];
247 J[I[i+nrows] + j] = J_i[I_i[i] + j];
248 D[I[i+nrows] + j] = factor*D_i[I_i[i] + j];
257 #ifdef MFEM_USE_SUITESPARSE
282 umfpack_zi_free_numeric(&
Numeric);
286 umfpack_zl_free_numeric(&
Numeric);
292 MFEM_VERIFY(
mat,
"not a ComplexSparseMatrix");
295 "Real and imag Sparsity pattern mismatch: Try setting Assemble (skip_zeros = 0)");
306 MFEM_VERIFY(
width ==
height,
"not a square matrix");
316 int status = umfpack_zi_symbolic(
width,
width,Ap,Ai,Ax,Az,&Symbolic,
321 umfpack_zi_report_status(
Control, status);
322 mfem_error(
"ComplexUMFPackSolver::SetOperator :"
323 " umfpack_zi_symbolic() failed!");
326 status = umfpack_zi_numeric(Ap, Ai, Ax, Az, Symbolic, &
Numeric,
331 umfpack_zi_report_status(
Control, status);
332 mfem_error(
"ComplexUMFPackSolver::SetOperator :"
333 " umfpack_zi_numeric() failed!");
335 umfpack_zi_free_symbolic(&Symbolic);
339 SuiteSparse_long status;
343 AI =
new SuiteSparse_long[
width + 1];
344 AJ =
new SuiteSparse_long[Ap[
width]];
345 for (
int i = 0; i <=
width; i++)
347 AI[i] = (SuiteSparse_long)(Ap[i]);
349 for (
int i = 0; i < Ap[
width]; i++)
351 AJ[i] = (SuiteSparse_long)(Ai[i]);
354 status = umfpack_zl_symbolic(
width,
width,
AI,
AJ, Ax, Az, &Symbolic,
359 umfpack_zl_report_status(
Control, status);
360 mfem_error(
"ComplexUMFPackSolver::SetOperator :"
361 " umfpack_zl_symbolic() failed!");
364 status = umfpack_zl_numeric(
AI,
AJ, Ax, Az, Symbolic, &
Numeric,
369 umfpack_zl_report_status(
Control, status);
370 mfem_error(
"ComplexUMFPackSolver::SetOperator :"
371 " umfpack_zl_numeric() failed!");
373 umfpack_zl_free_symbolic(&Symbolic);
380 mfem_error(
"ComplexUMFPackSolver::Mult : matrix is not set!"
381 " Call SetOperator first!");
394 if (conv == ComplexOperator::Convention::BLOCK_SYMMETRIC)
407 umfpack_zi_report_info(Control,
Info);
410 umfpack_zi_report_status(Control, status);
411 mfem_error(
"ComplexUMFPackSolver::Mult : umfpack_zi_solve() failed!");
416 SuiteSparse_long status =
421 umfpack_zl_report_info(Control,
Info);
424 umfpack_zl_report_status(Control, status);
425 mfem_error(
"ComplexUMFPackSolver::Mult : umfpack_zl_solve() failed!");
428 if (conv == ComplexOperator::Convention::BLOCK_SYMMETRIC)
437 mfem_error(
"ComplexUMFPackSolver::Mult : matrix is not set!"
438 " Call SetOperator first!");
463 umfpack_zi_report_info(Control,
Info);
466 umfpack_zi_report_status(Control, status);
467 mfem_error(
"ComplexUMFPackSolver::Mult : umfpack_zi_solve() failed!");
472 SuiteSparse_long status =
477 umfpack_zl_report_info(Control,
Info);
480 umfpack_zl_report_status(Control, status);
481 mfem_error(
"ComplexUMFPackSolver::Mult : umfpack_zl_solve() failed!");
505 umfpack_zi_free_numeric(&
Numeric);
509 umfpack_zl_free_numeric(&
Numeric);
520 bool ownReal,
bool ownImag,
524 comm_ = (A_Real) ? A_Real->
GetComm() :
525 ((A_Imag) ? A_Imag->
GetComm() : MPI_COMM_WORLD);
527 MPI_Comm_rank(comm_, &myid_);
528 MPI_Comm_size(comm_, &nranks_);
533 MFEM_ASSERT(
Op_Real_,
"ComplexHypreParMatrix has no real part!");
539 MFEM_ASSERT(
Op_Imag_,
"ComplexHypreParMatrix has no imaginary part!");
545 MFEM_ASSERT(
Op_Real_,
"ComplexHypreParMatrix has no real part!");
551 MFEM_ASSERT(
Op_Imag_,
"ComplexHypreParMatrix has no imaginary part!");
560 if ( A_r == NULL && A_i == NULL ) {
return NULL; }
564 HYPRE_BigInt global_num_rows = std::max(global_num_rows_r,
569 HYPRE_BigInt global_num_cols = std::max(global_num_cols_r,
572 int row_starts_size = (HYPRE_AssumedPartitionCheck()) ? 2 : nranks_ + 1;
579 ((A_i) ? A_i->
RowPart() : NULL);
581 ((A_i) ? A_i->
ColPart() : NULL);
583 for (
int i = 0; i < row_starts_size; i++)
585 row_starts[i] = 2 * row_starts_z[i];
586 col_starts[i] = 2 * col_starts_z[i];
592 int nrows_r = 0, nrows_i = 0, ncols_r = 0, ncols_i = 0;
593 int ncols_offd_r = 0, ncols_offd_i = 0;
598 nrows_r = diag_r.
Height();
599 ncols_r = diag_r.
Width();
600 ncols_offd_r = offd_r.
Width();
606 nrows_i = diag_i.
Height();
607 ncols_i = diag_i.
Width();
608 ncols_offd_i = offd_i.
Width();
610 int nrows = std::max(nrows_r, nrows_i);
611 int ncols = std::max(ncols_r, ncols_i);
614 std::set<HYPRE_BigInt> cset;
615 for (
int i=0; i<ncols_offd_r; i++)
617 cset.insert(cmap_r[i]);
619 for (
int i=0; i<ncols_offd_i; i++)
621 cset.insert(cmap_i[i]);
623 int num_cols_offd = (int)cset.size();
626 const int * diag_r_I = (A_r) ? diag_r.
GetI() : NULL;
627 const int * diag_i_I = (A_i) ? diag_i.
GetI() : NULL;
629 const int * diag_r_J = (A_r) ? diag_r.
GetJ() : NULL;
630 const int * diag_i_J = (A_i) ? diag_i.
GetJ() : NULL;
632 const double * diag_r_D = (A_r) ? diag_r.
GetData() : NULL;
633 const double * diag_i_D = (A_i) ? diag_i.
GetData() : NULL;
635 int diag_r_nnz = (diag_r_I) ? diag_r_I[nrows] : 0;
636 int diag_i_nnz = (diag_i_I) ? diag_i_I[nrows] : 0;
637 int diag_nnz = 2 * (diag_r_nnz + diag_i_nnz);
640 const int * offd_r_I = (A_r) ? offd_r.
GetI() : NULL;
641 const int * offd_i_I = (A_i) ? offd_i.
GetI() : NULL;
643 const int * offd_r_J = (A_r) ? offd_r.
GetJ() : NULL;
644 const int * offd_i_J = (A_i) ? offd_i.
GetJ() : NULL;
646 const double * offd_r_D = (A_r) ? offd_r.
GetData() : NULL;
647 const double * offd_i_D = (A_i) ? offd_i.
GetData() : NULL;
649 int offd_r_nnz = (offd_r_I) ? offd_r_I[nrows] : 0;
650 int offd_i_nnz = (offd_i_I) ? offd_i_I[nrows] : 0;
651 int offd_nnz = 2 * (offd_r_nnz + offd_i_nnz);
654 HYPRE_Int * diag_I = mfem_hypre_CTAlloc_host(HYPRE_Int, 2 * nrows + 1);
655 HYPRE_Int * diag_J = mfem_hypre_CTAlloc_host(HYPRE_Int, diag_nnz);
656 double * diag_D = mfem_hypre_CTAlloc_host(
double, diag_nnz);
658 HYPRE_Int * offd_I = mfem_hypre_CTAlloc_host(HYPRE_Int, 2 * nrows + 1);
659 HYPRE_Int * offd_J = mfem_hypre_CTAlloc_host(HYPRE_Int, offd_nnz);
660 double * offd_D = mfem_hypre_CTAlloc_host(
double, offd_nnz);
668 diag_I[nrows] = diag_r_nnz + diag_i_nnz;
669 for (
int i=0; i<nrows; i++)
671 diag_I[i + 1] = ((diag_r_I)?diag_r_I[i+1]:0) +
672 ((diag_i_I)?diag_i_I[i+1]:0);
673 diag_I[i + nrows + 1] = diag_I[i+1] + diag_r_nnz + diag_i_nnz;
677 for (
int j=0; j<diag_r_I[i+1] - diag_r_I[i]; j++)
679 diag_J[diag_I[i] + j] = diag_r_J[diag_r_I[i] + j];
680 diag_D[diag_I[i] + j] = diag_r_D[diag_r_I[i] + j];
682 diag_J[diag_I[i+nrows] + j] =
683 diag_r_J[diag_r_I[i] + j] + ncols;
684 diag_D[diag_I[i+nrows] + j] =
685 factor * diag_r_D[diag_r_I[i] + j];
690 const int off_r = (diag_r_I)?(diag_r_I[i+1] - diag_r_I[i]):0;
691 for (
int j=0; j<diag_i_I[i+1] - diag_i_I[i]; j++)
693 diag_J[diag_I[i] + off_r + j] = diag_i_J[diag_i_I[i] + j] + ncols;
694 diag_D[diag_I[i] + off_r + j] = -diag_i_D[diag_i_I[i] + j];
696 diag_J[diag_I[i+nrows] + off_r + j] = diag_i_J[diag_i_I[i] + j];
697 diag_D[diag_I[i+nrows] + off_r + j] =
698 factor * diag_i_D[diag_i_I[i] + j];
704 int num_recv_procs = 0;
706 this->getColStartStop(A_r, A_i, num_recv_procs, offd_col_start_stop);
708 std::set<HYPRE_BigInt>::iterator sit;
709 std::map<HYPRE_BigInt,HYPRE_BigInt> cmapa, cmapb, cinvmap;
710 for (sit=cset.begin(); sit!=cset.end(); sit++)
715 for (
int i=0; i<num_recv_procs; i++)
717 if (offd_col_start_stop[2*i] <= col_orig &&
718 col_orig < offd_col_start_stop[2*i+1])
720 col_2x2 = offd_col_start_stop[2*i] + col_orig;
721 col_size = offd_col_start_stop[2*i+1] - offd_col_start_stop[2*i];
725 cmapa[*sit] = col_2x2;
726 cmapb[*sit] = col_2x2 + col_size;
727 cinvmap[col_2x2] = -1;
728 cinvmap[col_2x2 + col_size] = -1;
730 delete [] offd_col_start_stop;
732 std::map<HYPRE_BigInt, HYPRE_BigInt>::iterator mit;
734 for (mit=cinvmap.begin(); mit!=cinvmap.end(); mit++, i++)
737 cmap[i] = mit->first;
742 offd_I[nrows] = offd_r_nnz + offd_i_nnz;
743 for (
int i=0; i<nrows; i++)
745 offd_I[i + 1] = ((offd_r_I)?offd_r_I[i+1]:0) +
746 ((offd_i_I)?offd_i_I[i+1]:0);
747 offd_I[i + nrows + 1] = offd_I[i+1] + offd_r_nnz + offd_i_nnz;
751 const int off_i = (offd_i_I)?(offd_i_I[i+1] - offd_i_I[i]):0;
752 for (
int j=0; j<offd_r_I[i+1] - offd_r_I[i]; j++)
754 offd_J[offd_I[i] + j] =
755 cinvmap[cmapa[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
756 offd_D[offd_I[i] + j] = offd_r_D[offd_r_I[i] + j];
758 offd_J[offd_I[i+nrows] + off_i + j] =
759 cinvmap[cmapb[cmap_r[offd_r_J[offd_r_I[i] + j]]]];
760 offd_D[offd_I[i+nrows] + off_i + j] =
761 factor * offd_r_D[offd_r_I[i] + j];
766 const int off_r = (offd_r_I)?(offd_r_I[i+1] - offd_r_I[i]):0;
767 for (
int j=0; j<offd_i_I[i+1] - offd_i_I[i]; j++)
769 offd_J[offd_I[i] + off_r + j] =
770 cinvmap[cmapb[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
771 offd_D[offd_I[i] + off_r + j] = -offd_i_D[offd_i_I[i] + j];
773 offd_J[offd_I[i+nrows] + j] =
774 cinvmap[cmapa[cmap_i[offd_i_J[offd_i_I[i] + j]]]];
775 offd_D[offd_I[i+nrows] + j] = factor * offd_i_D[offd_i_I[i] + j];
784 row_starts, col_starts,
785 diag_I, diag_J, diag_D,
786 offd_I, offd_J, offd_D,
787 2 * num_cols_offd, cmap,
791 hypre_ParCSRMatrix *hA = (hypre_ParCSRMatrix*)(*A);
792 hypre_ParCSRMatrixSetRowStartsOwner(hA,1);
793 hypre_ParCSRMatrixSetColStartsOwner(hA,1);
799 ComplexHypreParMatrix::getColStartStop(
const HypreParMatrix * A_r,
801 int & num_recv_procs,
805 hypre_ParCSRCommPkg * comm_pkg_r =
806 (A_r) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_r)) : NULL;
807 hypre_ParCSRCommPkg * comm_pkg_i =
808 (A_i) ? hypre_ParCSRMatrixCommPkg((hypre_ParCSRMatrix*)(*A_i)) : NULL;
810 std::set<HYPRE_Int> send_procs, recv_procs;
813 for (HYPRE_Int i=0; i<comm_pkg_r->num_sends; i++)
815 send_procs.insert(comm_pkg_r->send_procs[i]);
817 for (HYPRE_Int i=0; i<comm_pkg_r->num_recvs; i++)
819 recv_procs.insert(comm_pkg_r->recv_procs[i]);
824 for (HYPRE_Int i=0; i<comm_pkg_i->num_sends; i++)
826 send_procs.insert(comm_pkg_i->send_procs[i]);
828 for (HYPRE_Int i=0; i<comm_pkg_i->num_recvs; i++)
830 recv_procs.insert(comm_pkg_i->recv_procs[i]);
834 num_recv_procs = (int)recv_procs.size();
837 offd_col_start_stop =
new HYPRE_BigInt[2 * num_recv_procs];
840 ((A_i) ? A_i->
RowPart() : NULL);
842 int row_part_ind = (HYPRE_AssumedPartitionCheck()) ? 0 : myid_;
843 loc_start_stop[0] = row_part[row_part_ind];
844 loc_start_stop[1] = row_part[row_part_ind+1];
846 MPI_Request * req =
new MPI_Request[send_procs.size()+recv_procs.size()];
847 MPI_Status * stat =
new MPI_Status[send_procs.size()+recv_procs.size()];
852 std::set<HYPRE_Int>::iterator sit;
853 for (sit=send_procs.begin(); sit!=send_procs.end(); sit++)
855 MPI_Isend(loc_start_stop, 2, HYPRE_MPI_BIG_INT,
856 *sit, tag, comm_, &req[send_count]);
859 for (sit=recv_procs.begin(); sit!=recv_procs.end(); sit++)
861 MPI_Irecv(&offd_col_start_stop[2*recv_count], 2, HYPRE_MPI_BIG_INT,
862 *sit, tag, comm_, &req[send_count+recv_count]);
866 MPI_Waitall(send_count+recv_count, req, stat);
872 #endif // MFEM_USE_MPI
double Control[UMFPACK_CONTROL]
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
MPI_Comm GetComm() const
MPI communicator.
virtual Operator & real()
Real or imaginary part accessor methods.
virtual ~ComplexOperator()
ComplexSparseMatrix * mat
virtual ~ComplexUMFPackSolver()
ComplexHypreParMatrix(HypreParMatrix *A_Real, HypreParMatrix *A_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
void SetSize(int s)
Resize the vector to size s.
Mimic the action of a complex operator using two real operators.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
int * GetJ()
Return the array J.
int Size() const
Returns the size of the vector.
int * GetI()
Return the array I.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Alternate convention for damping operators.
double * GetData() const
Return a pointer to the beginning of the Vector data.
ComplexOperator(Operator *Op_Real, Operator *Op_Imag, bool ownReal, bool ownImag, Convention convention=HERMITIAN)
Constructs complex operator object.
void SortColumnIndices()
Sort the column indices corresponding to each row.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
double * GetData()
Return the element data, i.e. the array A.
void GetOffd(SparseMatrix &offd, HYPRE_BigInt *&cmap) const
Get the local off-diagonal block. NOTE: 'offd' will not own any data.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
const double * HostReadData() const
const int * HostReadJ() const
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
virtual SparseMatrix & real()
Real or imaginary part accessor methods.
SparseMatrix * GetSystemMatrix() const
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
virtual Operator & imag()
Convention GetConvention() const
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
int height
Dimension of the output / number of rows in the matrix.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a ComplexSparseMatrix.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
HYPRE_BigInt * ColPart()
Returns the column partitioning.
virtual void Mult(const Vector &b, Vector &x) const
This is solving the system A x = b.
HYPRE_BigInt * RowPart()
Returns the row partitioning.
virtual SparseMatrix & imag()
Native convention for Hermitian operators.
virtual void MultTranspose(const Vector &b, Vector &x) const
This is solving the system: A^H x = b (when transa = false) This is equivalent to solving the transpo...
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const int * HostReadI() const
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
HypreParMatrix * GetSystemMatrix() const
virtual HypreParMatrix & imag()
Wrapper for hypre's ParCSR matrix class.
double Info[UMFPACK_INFO]
virtual HypreParMatrix & real()
Real or imaginary part accessor methods.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
int width
Dimension of the input / number of columns in the matrix.