47 const int NE = fes->
GetNE();
49 for (
int i = 0; i < NE; i++)
55 for (
int i = 0; i < NE; i++)
58 const int nsd = vdofs.
Size()/vdim;
60 const int *dofs = vdofs.
GetData();
61 for (
int vd = 0; vd < vdim; vd++)
64 for (
int j = 0; j < nspd; j++)
66 MFEM_ASSERT(dofs[nsd-nspd+j] >= 0,
"");
77 "incompatible volume and trace FE spaces");
81 for (
int i = 0; i < NE; i++)
85 const int nsd = vdofs.
Size()/vdim;
86 const int nsrd = rvdofs.
Size()/vdim;
87 for (
int vd = 0; vd < vdim; vd++)
89 for (
int j = 0; j < nsrd; j++)
91 int rvdof = rvdofs[j+nsrd*vd];
92 int vdof = vdofs[j+nsd*vd];
98 MFEM_ASSERT(vdof >= 0,
"incompatible volume and trace FE spaces");
99 rdof_edof[rvdof] = vdof;
136 const int NE = fes->
GetNE();
140 A_offsets[0] = A_ipiv_offsets[0] = 0;
142 for (
int i = 0; i < NE; i++)
145 const int ned = rvdofs.
Size();
146 const int npd = elem_pdof.
RowSize(i);
147 A_offsets[i+1] = A_offsets[i] + npd*(npd + (symm ? 1 : 2)*ned);
148 A_ipiv_offsets[i+1] = A_ipiv_offsets[i] + npd;
150 A_data =
new double[A_offsets[NE]];
151 A_ipiv =
new int[A_ipiv_offsets[NE]];
152 const int nedofs = tr_fes->
GetVSize();
158 Table elem_rdof, rdof_elem;
160 for (
int i = 0; i < NE; i++)
166 for (
int i = 0; i < NE; i++)
177 nedofs, nedofs,
true,
true,
false);
193 const int vdim = fes->
GetVDim();
194 const int nvpd = elem_pdof.
RowSize(el);
195 const int nved = rvdofs.
Size();
196 DenseMatrix A_pp(A_data + A_offsets[el], nvpd, nvpd);
199 if (symm) { A_ep.
SetSize(nved, nvpd); }
200 else { A_ep.UseExternalData(A_pe.Data() + nvpd*nved, nved, nvpd); }
203 const int npd = nvpd/vdim;
204 const int ned = nved/vdim;
205 const int nd = npd + ned;
207 for (
int i = 0; i < vdim; i++)
209 for (
int j = 0; j < vdim; j++)
211 A_pp.
CopyMN(elmat, npd, npd, i*nd+ned, j*nd+ned, i*npd, j*npd);
212 A_pe.CopyMN(elmat, npd, ned, i*nd+ned, j*nd, i*npd, j*ned);
213 A_ep.CopyMN(elmat, ned, npd, i*nd, j*nd+ned, i*ned, j*npd);
214 A_ee.
CopyMN(elmat, ned, ned, i*nd, j*nd, i*ned, j*ned);
220 lu.BlockFactor(nvpd, nved, A_pe.Data(), A_ep.Data(), A_ee.
Data());
223 const int skip_zeros = 0;
231 const int skip_zeros = 0;
237 const int skip_zeros = 0;
241 if (S_e) { S_e->
Finalize(skip_zeros); }
264 if (S_e) { S_e->
Finalize(skip_zeros); }
289 const Array<int> &ess_rtdof_list,
int keep_diagonal)
291 if (!Parallel() || S)
297 for (
int i = 0; i < ess_rtdof_list.
Size(); i++)
305 MFEM_ASSERT(pS_e.
Ptr() == NULL,
"essential b.c. already eliminated");
315 MFEM_ASSERT(b.
Size() == fes->
GetVSize(),
"'b' has incorrect size");
317 const int NE = fes->
GetNE();
318 const int nedofs = tr_fes->
GetVSize();
330 for (
int i = 0; i < nedofs; i++)
332 b_r(i) = b(rdof_edof[i]);
338 for (
int i = 0; i < NE; i++)
341 const int ned = rvdofs.
Size();
342 const int *rd = rvdofs.
GetData();
343 const int npd = elem_pdof.
RowSize(i);
344 const int *pd = elem_pdof.
GetRow(i);
347 for (
int j = 0; j < npd; j++)
352 LUFactors lu(A_data + A_offsets[i], A_ipiv + A_ipiv_offsets[i]);
364 L_ep.
Mult(b_p, b_ep);
366 for (
int j = 0; j < ned; j++)
368 if (rd[j] >= 0) { b_r(rd[j]) -= b_ep(j); }
369 else { b_r(-1-rd[j]) += b_ep(j); }
392 MFEM_ASSERT(sol.
Size() == fes->
GetVSize(),
"'sol' has incorrect size");
394 const int nedofs = tr_fes->
GetVSize();
406 for (
int i = 0; i < nedofs; i++)
408 sol_r(i) = sol(rdof_edof[i]);
413 tr_R->
Mult(sol_r, sc_sol);
418 Vector &B,
int copy_interior)
const
430 MFEM_ASSERT(pS.
Type() == pS_e.
Type(),
"type id mismatch");
443 const int nedofs = tr_fes->
GetVSize();
448 ess_dof_marker.
MakeRef(ess_tdof_marker);
459 ess_rtdof_marker.
SetSize(nedofs);
460 ess_rdof_marker.
MakeRef(ess_rtdof_marker);
464 ess_rdof_marker.
SetSize(nedofs);
466 for (
int i = 0; i < nedofs; i++)
468 ess_rdof_marker[i] = ess_dof_marker[rdof_edof[i]];
473 tr_R->
BooleanMult(ess_rdof_marker, ess_rtdof_marker);
483 MFEM_ASSERT(b.
Size() == fes->
GetVSize(),
"'b' has incorrect size");
485 const int nedofs = tr_fes->
GetVSize();
497 tr_cP->
Mult(sc_sol, sol_r);
508 for (
int i = 0; i < nedofs; i++)
510 sol(rdof_edof[i]) = sol_r(i);
512 const int NE = fes->
GetNE();
515 for (
int i = 0; i < NE; i++)
518 const int ned = rvdofs.
Size();
519 const int npd = elem_pdof.
RowSize(i);
520 const int *pd = elem_pdof.
GetRow(i);
523 for (
int j = 0; j < npd; j++)
529 LUFactors lu(A_data + A_offsets[i], A_ipiv + A_ipiv_offsets[i]);
533 for (
int j = 0; j < npd; j++)
virtual const Operator * GetProlongationMatrix()
Ordering::Type GetOrdering() const
Return the ordering method.
int Size() const
Logical size of the array.
virtual FiniteElementCollection * GetTraceCollection() const
const SparseMatrix * GetConformingProlongation() const
~StaticCondensation()
Destroy a StaticCondensation object.
void AddColumnsInRow(int r, int ncol)
bool ReducesTrueVSize() const
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
void ReduceRHS(const Vector &b, Vector &sc_b) const
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
HYPRE_Int * GetDofOffsets() const
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
void SetSize(int s)
Resize the vector to size s.
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
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().
Pointer to an Operator of a specified type.
Operator::Type Type() const
Get the currently set operator type id.
void BlockBackSolve(int m, int n, int r, const double *U12, const double *X2, double *Y1) const
void EliminateBC(const OperatorHandle &A_e, const Array< int > &ess_dof_list, const Vector &X, Vector &B) const
Eliminate essential dofs from the solution X into the r.h.s. B.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
T * GetData()
Returns the data.
void EliminateRowCol(int rc, const double sol, Vector &rhs, int d=0)
Data type dense matrix using column-major storage.
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
int Size() const
Returns the size of the vector.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
int GetNumElementInteriorDofs(int i) const
void ReduceSystem(Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0) const
Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r...
void LoseData()
Call this if data has been stolen.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
StaticCondensation(FiniteElementSpace *fespace)
Construct a StaticCondensation object.
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
int Size_of_connections() const
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, but treat all elements as booleans (zero=false, nonzero=true).
void AddConnections(int r, const int *c, int nc)
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 ...
int GetNE() const
Returns number of elements in the mesh.
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Mesh * GetMesh() const
Returns the mesh.
void EliminateRowsCols(OperatorHandle &A, const Array< int > &ess_dof_list)
Reset the OperatorHandle to be the eliminated part of A after elimination of the essential dofs ess_d...
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
HYPRE_Int GlobalTrueVSize() const
void LSolve(int m, int n, double *X) const
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
int GetVDim() const
Returns vector dimension.
virtual const SparseMatrix * GetRestrictionMatrix()
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Operator * Ptr() const
Access the underlying Operator pointer.
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the 'dofs' array to the given 'val'.
HYPRE_Int GlobalVSize() const
double * Data() const
Returns the matrix data array.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void Finalize()
Finalize the construction of the Schur complement matrix.
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
void AssembleBdrMatrix(int el, const DenseMatrix &elmat)
void EliminateReducedTrueDofs(const Array< int > &ess_rtdof_list, int keep_diagonal)
Eliminate the given reduced true dofs from the Schur complement matrix S.
void MakeSquareBlockDiag(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel square block-diagonal matrix using the currently set type...
void AssembleMatrix(int el, const DenseMatrix &elmat)
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
void MakeRef(T *, int)
Make this Array a reference to a pointer.
const FiniteElementCollection * FEColl() const
void Mult(const double *x, double *y) const
Matrix vector multiplication.
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
void ReduceSolution(const Vector &sol, Vector &sc_sol) const
void UseExternalData(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i'th boundary element.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void Init(bool symmetric, bool block_diagonal)
void ConvertMarkerToReducedTrueDofs(const Array< int > &ess_tdof_marker, Array< int > &ess_rtdof_marker) const
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
static void AdjustVDofs(Array< int > &vdofs)
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const