23 #ifdef MFEM_DEBUG_HYBRIDIZATION_CP 32 : fes(fespace), c_fes(c_fespace), c_bfi(NULL), Ct(NULL), H(NULL),
33 Af_data(NULL), Af_ipiv(NULL)
60 int c_num_face_nbr_dofs = 0;
64 HYPRE_BigInt num_shared_slave_faces = 0, glob_num_shared_slave_faces = 0;
72 MPI_Allreduce(&num_shared_slave_faces, &glob_num_shared_slave_faces, 1,
73 HYPRE_MPI_BIG_INT, MPI_SUM, pmesh->
GetComm());
74 MFEM_ASSERT(glob_num_shared_slave_faces%2 == 0,
"");
75 glob_num_shared_slave_faces /= 2;
76 if (glob_num_shared_slave_faces)
81 #ifdef MFEM_DEBUG_HERE 82 MFEM_WARNING(
'[' << c_pfes->
GetMyRank() <<
83 "] num_shared_slave_faces = " << num_shared_slave_faces
84 <<
", glob_num_shared_slave_faces = " 85 << glob_num_shared_slave_faces
86 <<
"\n num_face_nbr_dofs = " << c_num_face_nbr_dofs
88 #undef MFEM_DEBUG_HERE 95 Ct =
new SparseMatrix(num_hat_dofs, c_vsize + c_num_face_nbr_dofs);
99 const int skip_zeros = 1;
104 for (
int i = 0; i < num_faces; i++)
107 if (!FTr) {
continue; }
114 for (
int j = 0; j < s1; j++)
118 for (
int j = 0; j < s2; j++)
120 vdofs[s1+j] = o2 + j;
136 for (
int i = 0; i < num_shared_faces; i++)
139 const bool ghost_sface = (face_no >= num_faces);
144 MFEM_ASSERT(FTr->
Elem2No < 0,
"");
150 const int fill2 =
false;
155 for (
int j = 0; j < c_vdofs.
Size(); j++)
157 c_vdofs[j] += c_vsize;
163 for (
int j = 0; j < s1; j++)
173 if (glob_num_shared_slave_faces)
186 for (
int i = 0; i < Ct_J.Size(); i++)
188 Ct_J[i] = J[i] < c_vsize ?
189 J[i] + c_ldof_offset :
190 c_face_nbr_glob_ldof[J[i] - c_vsize];
211 MFEM_ABORT(
"TODO: algebraic definition of C");
222 int num_hat_dofs = 0;
225 for (
int i = 0; i < NE; i++)
228 num_hat_dofs += vdofs.
Size();
235 #ifdef MFEM_DEBUG_HYBRIDIZATION_CP 238 std::ofstream C_file(
"C_matrix.txt");
246 std::ofstream P_file(
"P_matrix.txt");
264 free_tdof_marker = 1;
265 for (
int i = 0; i < ess_tdof_list.
Size(); i++)
267 free_tdof_marker[ess_tdof_list[i]] = 0;
276 free_vdofs_marker.
MakeRef(free_tdof_marker);
281 cP->
BooleanMult(free_tdof_marker, free_vdofs_marker);
288 P->
BooleanMult(1, free_tdof_marker, 0, free_vdofs_marker);
294 free_vdofs_marker.
MakeRef(free_tdof_marker);
299 cP->
BooleanMult(free_tdof_marker, free_vdofs_marker);
302 for (
int i = 0; i < NE; i++)
306 for (
int j = 0; j < vdofs.
Size(); j++)
320 for (
int i = 0; i < num_hat_dofs; i++)
335 #ifdef MFEM_DEBUG_HERE 338 for (
int i = 0; i < NE; i++)
344 #ifdef MFEM_DEBUG_HERE 352 #ifdef MFEM_DEBUG_HERE 356 int myid = pmesh ? pmesh->GetMyRank() : 0;
359 int e_size = num_hat_dofs - (i_size + b_size);
361 <<
" [" << myid <<
"] hat dofs - \"internal\": " << i_size
362 <<
", \"boundary\": " << b_size
363 <<
", \"essential\": " << e_size <<
'\n' << std::endl;
364 #undef MFEM_DEBUG_HERE 376 for (
int i = 0; i < NE; i++)
380 for (
int j = 0; j < vdofs.
Size(); j++)
384 vdof_marker[vdofs[j]] = 1;
388 for (
int tdof = 0; tdof < R->
Height(); tdof++)
390 if (free_tdof_marker[tdof]) {
continue; }
392 const int ncols = R->
RowSize(tdof);
395 for (
int j = 0; j < ncols; j++)
397 if (std::abs(vals[j]) != 0.0 && vdof_marker[cols[j]] == 0)
399 MFEM_ABORT(
"Ref != 0");
418 for (
int i = h_start; i < h_end; i++)
421 if (mark == 0) { i_dofs.
Append(i-h_start); }
422 else if (mark == -1) { b_dofs.
Append(i-h_start); }
434 for (
int i = h_start; i < h_end; i++)
437 if (mark == 0) { num_idofs++; }
438 else if (mark == -1) { b_dofs.
Append(i); }
456 for (
int j = 0; j < i_dofs.
Size(); j++)
458 int j_dof = i_dofs[j];
459 for (
int i = 0; i < i_dofs.
Size(); i++)
461 A_ii(i,j) = A(i_dofs[i],j_dof);
463 for (
int i = 0; i < b_dofs.
Size(); i++)
465 A_bi(i,j) = A(b_dofs[i],j_dof);
468 for (
int j = 0; j < b_dofs.
Size(); j++)
470 int j_dof = b_dofs[j];
471 for (
int i = 0; i < i_dofs.
Size(); i++)
473 A_ib(i,j) = A(i_dofs[i],j_dof);
475 for (
int i = 0; i < b_dofs.
Size(); i++)
477 A_bb(i,j) = A(b_dofs[i],j_dof);
504 Ordering::DofsToVDofs<Ordering::byNODES>(e2f.
Size()/vdim, vdim, lvdofs);
505 MFEM_ASSERT(lvdofs.
Size() == A.
Height(),
"internal error");
508 for (
int i = 0; i < lvdofs.
Size(); i++)
511 bd = (bd >= 0) ? vdofs[bd] : -1-vdofs[-1-bd];
512 MFEM_ASSERT(bvdofs[i] == bd,
"internal error");
518 for (
int i = 0; i < lvdofs.
Size(); i++)
534 for (
int j = 0; j < i_dofs.
Size(); j++)
536 int j_f = e2f[i_dofs[j]];
537 if (j_f == -1) {
continue; }
538 for (
int i = 0; i < i_dofs.Size(); i++)
540 int i_f = e2f[i_dofs[i]];
541 if (i_f == -1) {
continue; }
542 A_ii(i,j) += B(i_f,j_f);
544 for (
int i = 0; i < b_dofs.
Size(); i++)
546 int i_f = e2f[b_dofs[i]];
547 if (i_f == -1) {
continue; }
548 A_bi(i,j) += B(i_f,j_f);
551 for (
int j = 0; j < b_dofs.
Size(); j++)
553 int j_f = e2f[b_dofs[j]];
554 if (j_f == -1) {
continue; }
555 for (
int i = 0; i < i_dofs.
Size(); i++)
557 int i_f = e2f[i_dofs[i]];
558 if (i_f == -1) {
continue; }
559 A_ib(i,j) += B(i_f,j_f);
561 for (
int i = 0; i < b_dofs.Size(); i++)
563 int i_f = e2f[b_dofs[i]];
564 if (i_f == -1) {
continue; }
565 A_bb(i,j) += B(i_f,j_f);
572 const int skip_zeros = 1;
586 int c_mark_start = 0;
587 for (
int el = 0; el < NE; el++)
593 double *A_ib_data = LU_ii.
data + i_dofs_size*i_dofs_size;
594 double *A_bi_data = A_ib_data + i_dofs_size*b_dofs.
Size();
596 LU_ii.
ipiv + i_dofs_size);
598 LU_ii.
Factor(i_dofs_size);
600 A_ib_data, A_bi_data, LU_bb.data);
601 LU_bb.Factor(b_dofs.
Size());
605 for (
int i = 0; i < b_dofs.
Size(); i++)
607 const int row = b_dofs[i];
610 for (
int j = 0; j < ncols; j++)
612 const int c_dof = cols[j];
613 if (c_dof_marker[c_dof] < c_mark_start)
615 c_dof_marker[c_dof] = c_mark_start + c_dofs.
Size();
622 for (
int i = 0; i < b_dofs.
Size(); i++)
624 const int row = b_dofs[i];
628 for (
int j = 0; j < ncols; j++)
630 const int loc_j = c_dof_marker[cols[j]] - c_mark_start;
631 Cb_t(i,loc_j) = vals[j];
643 MultAtB(Cb_t, Sb_inv_Cb_t, Hb);
651 V->AddSubMatrix(b_dofs, c_dofs, Sb_inv_Cb_t, skip_zeros);
655 c_mark_start += c_dofs.
Size();
656 MFEM_VERIFY(c_mark_start >= 0,
"overflow");
658 const bool fix_empty_rows =
true;
666 if (!c_pfes) {
return; }
682 MFEM_ASSERT(c_pfes,
"");
687 for (
int i = 0; i < V_J.Size(); i++)
689 V_J[i] = J[i] < c_vsize ?
690 J[i] + c_ldof_offset :
691 c_face_nbr_glob_ldof[J[i] - c_vsize];
698 V->GetI(), V_J.GetData(), V->GetData(),
709 plpH.ConvertFrom(lpH);
744 Vector el_vals, bf_i, i_vals, b_vals;
768 for (
int i = 0; i < NE; i++)
772 for (
int j = 0; j < vdofs.
Size(); j++)
775 if (vdof < 0) { vdof = -1 - vdof; }
776 if (vdof_marker[vdof]) { el_vals(j) = 0.0; }
777 else { vdof_marker[vdof] =
true; }
789 double *U_ib = LU_ii.
data + i_dofs.
Size()*i_dofs.
Size();
790 double *L_bi = U_ib + i_dofs.
Size()*b_dofs.
Size();
854 s.MakeRef(
fes, sol, 0);
863 for (
int i = 0; i < NE; i++)
870 if (vdof >= 0) {
s(vdof) = bf(j); }
871 else {
s(-1-vdof) = -bf(j); }
Abstract class for all finite elements.
void GetBDofs(int el, int &num_idofs, Array< int > &b_dofs) const
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
void SubDofOrder(Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
const FiniteElement * GetFaceNbrFaceFE(int i) const
Class for grid function - Vector with associated FE space.
int GetConformingVSize() const
FiniteElementSpace * c_fes
Array< int > Af_f_offsets
Geometry::Type GetElementBaseGeometry(int i) const
int TrueVSize() const
Obsolete, kept for backward compatibility.
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
void SetSize(int s)
Resize the vector to size s.
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
int RowSize(const int i) const
Returns the number of elements in row i.
Lists all edges/faces in the nonconforming mesh.
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i'th face element (2D and 3D).
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Pointer to an Operator of a specified type.
int * GetJ()
Return the array J.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
double MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
bool Nonconforming() const
int Size() const
Returns the size of the vector.
Data type dense matrix using column-major storage.
double * Data() const
Returns the matrix data array.
int * GetI()
Return the array I.
HYPRE_BigInt GetMyDofOffset() const
Abstract parallel finite element space.
int GetFaceNbrVSize() const
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
void Reset()
Destroy the current hybridization matrix while preserving the computed constraint matrix and the set ...
void BlockForwSolve(int m, int n, int r, const double *L21, double *B1, double *B2) const
void DeleteAll()
Delete the whole array.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
const NCList & GetSharedList(int entity)
Helper to get shared vertices/edges/faces ('entity' == 0/1/2 resp.).
double * GetData()
Return the element data, i.e. the array A.
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
void ReduceRHS(const Vector &b, Vector &b_r) const
void Finalize()
Finalize the construction of the hybridized matrix.
void ExchangeFaceNbrData()
void AssembleBdrMatrix(int bdr_el, const DenseMatrix &A)
Assemble the boundary element matrix A into the hybridized system matrix.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
const FiniteElementCollection * FEColl() const
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31) override
const HYPRE_BigInt * GetFaceNbrGlobalDofMap()
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
ParMesh * GetParMesh() const
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
int GetNE() const
Returns number of elements in the mesh.
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
double * GetData() const
Return a pointer to the beginning of the Vector data.
void Threshold(double eps)
Replace small entries, abs(a_ij) <= eps, with zero.
BilinearFormIntegrator * c_bfi
int GetVDim() const
Returns vector dimension.
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
void BooleanMult(int alpha, const int *x, int beta, int *y)
The "Boolean" analog of y = alpha * A * x + beta * y, where elements in the sparsity pattern of the m...
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Mesh * GetMesh() const
Returns the mesh.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
void GetIBDofs(int el, Array< int > &i_dofs, Array< int > &b_dofs) const
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 AssembleMatrix(int el, const DenseMatrix &A)
Assemble the element matrix A into the hybridized system matrix.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
HYPRE_BigInt GlobalVSize() const
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
MFEM_HOST_DEVICE void MultTranspose(const int height, const int width, const TA *data, const TX *x, TY *y)
Matrix transpose vector multiplication: y = At x, where the matrix A is of size height x width with g...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
~Hybridization()
Destructor.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
HypreParMatrix * Transpose() const
Returns the transpose of *this.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
HYPRE_BigInt * GetDofOffsets() const
int GetNE() const
Returns number of elements.
Operator::Type Type() const
Get the currently set operator type id.
HYPRE_BigInt * ColPart()
Returns the column partitioning.
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i'th boundary element.
HYPRE_BigInt * RowPart()
Returns the row partitioning.
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
void MultAfInv(const Vector &b, const Vector &lambda, Vector &bf, int mode) const
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Operator * Ptr() const
Access the underlying Operator pointer.
Array< int > hat_dofs_marker
int Size() const
Return the logical size of the array.
T & Last()
Return the last element in the array.
void MakeRef(T *, int)
Make this Array a reference to a pointer.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
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
ID for class HypreParMatrix.
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
void Init(const Array< int > &ess_tdof_list)
Prepare the Hybridization object for assembly.
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A^t * x + beta * y.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Wrapper for hypre's ParCSR matrix class.
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
Hybridization(FiniteElementSpace *fespace, FiniteElementSpace *c_fespace)
Constructor.
Class for parallel meshes.
void AdjustDofDirection(Array< int > &dofs)
ID for class PetscParMatrix, MATIS format.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
static void AdjustVDofs(Array< int > &vdofs)
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
void GetBdrElementAdjacentElement(int bdr_el, int &el, int &info) const
For the given boundary element, bdr_el, return its adjacent element and its info, i...