MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pbilinearform.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "fem.hpp"
17 #include "../general/sort_pairs.hpp"
18 
19 namespace mfem
20 {
21 
23 {
24  int nbr_size = pfes->GetFaceNbrVSize();
25 
26  if (precompute_sparsity == 0 || fes->GetVDim() > 1)
27  {
28  if (keep_nbr_block)
29  {
30  mat = new SparseMatrix(height + nbr_size, width + nbr_size);
31  }
32  else
33  {
34  mat = new SparseMatrix(height, width + nbr_size);
35  }
36  return;
37  }
38 
39  // the sparsity pattern is defined from the map: face->element->dof
40  const Table &lelem_ldof = fes->GetElementToDofTable(); // <-- dofs
41  const Table &nelem_ndof = pfes->face_nbr_element_dof; // <-- vdofs
42  Table elem_dof; // element + nbr-element <---> dof
43  if (nbr_size > 0)
44  {
45  // merge lelem_ldof and nelem_ndof into elem_dof
46  int s1 = lelem_ldof.Size(), s2 = nelem_ndof.Size();
47  const int *I1 = lelem_ldof.GetI(), *J1 = lelem_ldof.GetJ();
48  const int *I2 = nelem_ndof.GetI(), *J2 = nelem_ndof.GetJ();
49  const int nnz1 = I1[s1], nnz2 = I2[s2];
50 
51  elem_dof.SetDims(s1 + s2, nnz1 + nnz2);
52 
53  int *I = elem_dof.GetI(), *J = elem_dof.GetJ();
54  for (int i = 0; i <= s1; i++)
55  {
56  I[i] = I1[i];
57  }
58  for (int j = 0; j < nnz1; j++)
59  {
60  J[j] = J1[j];
61  }
62  for (int i = 0; i <= s2; i++)
63  {
64  I[s1+i] = I2[i] + nnz1;
65  }
66  for (int j = 0; j < nnz2; j++)
67  {
68  J[nnz1+j] = J2[j] + height;
69  }
70  }
71  // dof_elem x elem_face x face_elem x elem_dof (keep_nbr_block = true)
72  // ldof_lelem x lelem_face x face_elem x elem_dof (keep_nbr_block = false)
73  Table dof_dof;
74  {
75  Table face_dof; // face_elem x elem_dof
76  {
77  Table *face_elem = pfes->GetParMesh()->GetFaceToAllElementTable();
78  if (nbr_size > 0)
79  {
80  mfem::Mult(*face_elem, elem_dof, face_dof);
81  }
82  else
83  {
84  mfem::Mult(*face_elem, lelem_ldof, face_dof);
85  }
86  delete face_elem;
87  if (nbr_size > 0)
88  {
89  elem_dof.Clear();
90  }
91  }
92 
93  if (keep_nbr_block)
94  {
95  Table dof_face;
96  Transpose(face_dof, dof_face, height + nbr_size);
97  mfem::Mult(dof_face, face_dof, dof_dof);
98  }
99  else
100  {
101  Table ldof_face;
102  {
103  Table face_ldof;
104  Table *face_lelem = fes->GetMesh()->GetFaceToElementTable();
105  mfem::Mult(*face_lelem, lelem_ldof, face_ldof);
106  delete face_lelem;
107  Transpose(face_ldof, ldof_face, height);
108  }
109  mfem::Mult(ldof_face, face_dof, dof_dof);
110  }
111  }
112 
113  int *I = dof_dof.GetI();
114  int *J = dof_dof.GetJ();
115  int nrows = dof_dof.Size();
116  double *data = new double[I[nrows]];
117 
118  mat = new SparseMatrix(I, J, data, nrows, height + nbr_size);
119  *mat = 0.0;
120 
121  dof_dof.LoseData();
122 }
123 
125 {
126  A.Clear();
127 
128  if (A_local == NULL) { return; }
129  MFEM_VERIFY(A_local->Finalized(), "the local matrix must be finalized");
130 
131  OperatorHandle dA(A.Type()), Ph(A.Type()), hdA;
132 
133  if (fbfi.Size() == 0)
134  {
135  // construct a parallel block-diagonal matrix 'A' based on 'a'
137  pfes->GetDofOffsets(), A_local);
138  }
139  else
140  {
141  // handle the case when 'a' contains offdiagonal
142  int lvsize = pfes->GetVSize();
143  const HYPRE_Int *face_nbr_glob_ldof = pfes->GetFaceNbrGlobalDofMap();
144  HYPRE_Int ldof_offset = pfes->GetMyDofOffset();
145 
146  Array<HYPRE_Int> glob_J(A_local->NumNonZeroElems());
147  int *J = A_local->GetJ();
148  for (int i = 0; i < glob_J.Size(); i++)
149  {
150  if (J[i] < lvsize)
151  {
152  glob_J[i] = J[i] + ldof_offset;
153  }
154  else
155  {
156  glob_J[i] = face_nbr_glob_ldof[J[i] - lvsize];
157  }
158  }
159 
160  // TODO - construct dA directly in the A format
161  hdA.Reset(
162  new HypreParMatrix(pfes->GetComm(), lvsize, pfes->GlobalVSize(),
163  pfes->GlobalVSize(), A_local->GetI(), glob_J,
164  A_local->GetData(), pfes->GetDofOffsets(),
165  pfes->GetDofOffsets()));
166  // - hdA owns the new HypreParMatrix
167  // - the above constructor copies all input arrays
168  glob_J.DeleteAll();
169  dA.ConvertFrom(hdA);
170  }
171 
172  // TODO - assemble the Dof_TrueDof_Matrix directly in the required format?
173  Ph.ConvertFrom(pfes->Dof_TrueDof_Matrix());
174  // TODO: When Ph.Type() == Operator::ANY_TYPE we want to use the Operator
175  // returned by pfes->GetProlongationMatrix(), however that Operator is a
176  // const Operator, so we cannot store it in OperatorHandle. We need a const
177  // version of class OperatorHandle, e.g. ConstOperatorHandle.
178 
179  A.MakePtAP(dA, Ph);
180 }
181 
183 {
185  ParallelAssemble(Mh, m);
186  Mh.SetOperatorOwner(false);
187  return Mh.As<HypreParMatrix>();
188 }
189 
191 {
192  ParMesh *pmesh = pfes->GetParMesh();
194  Array<int> vdofs1, vdofs2, vdofs_all;
196 
197  int nfaces = pmesh->GetNSharedFaces();
198  for (int i = 0; i < nfaces; i++)
199  {
200  T = pmesh->GetSharedFaceTransformations(i);
201  pfes->GetElementVDofs(T->Elem1No, vdofs1);
202  pfes->GetFaceNbrElementVDofs(T->Elem2No, vdofs2);
203  vdofs1.Copy(vdofs_all);
204  for (int j = 0; j < vdofs2.Size(); j++)
205  {
206  vdofs2[j] += height;
207  }
208  vdofs_all.Append(vdofs2);
209  for (int k = 0; k < fbfi.Size(); k++)
210  {
211  fbfi[k]->AssembleFaceMatrix(*pfes->GetFE(T->Elem1No),
212  *pfes->GetFaceNbrFE(T->Elem2No),
213  *T, elemmat);
214  if (keep_nbr_block)
215  {
216  mat->AddSubMatrix(vdofs_all, vdofs_all, elemmat, skip_zeros);
217  }
218  else
219  {
220  mat->AddSubMatrix(vdofs1, vdofs_all, elemmat, skip_zeros);
221  }
222  }
223  }
224 }
225 
226 void ParBilinearForm::Assemble(int skip_zeros)
227 {
228  if (mat == NULL && fbfi.Size() > 0)
229  {
231  pAllocMat();
232  }
233 
234  BilinearForm::Assemble(skip_zeros);
235 
236  if (fbfi.Size() > 0)
237  {
238  AssembleSharedFaces(skip_zeros);
239  }
240 }
241 
242 void ParBilinearForm
244  HypreParMatrix &A, const HypreParVector &X,
245  HypreParVector &B) const
246 {
247  Array<int> dof_list;
248 
249  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
250 
251  // do the parallel elimination
252  A.EliminateRowsCols(dof_list, X, B);
253 }
254 
257  HypreParMatrix &A) const
258 {
259  Array<int> dof_list;
260 
261  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
262 
263  return A.EliminateRowsCols(dof_list);
264 }
265 
266 void ParBilinearForm::TrueAddMult(const Vector &x, Vector &y, const double a)
267 const
268 {
269  MFEM_VERIFY(fbfi.Size() == 0, "the case of interior face integrators is not"
270  " implemented");
271 
272  if (X.ParFESpace() != pfes)
273  {
274  X.SetSpace(pfes);
275  Y.SetSpace(pfes);
276  }
277 
278  X.Distribute(&x);
279  mat->Mult(X, Y);
280  pfes->Dof_TrueDof_Matrix()->MultTranspose(a, Y, 1.0, y);
281 }
282 
284  const Array<int> &ess_tdof_list, Vector &x, Vector &b,
285  OperatorHandle &A, Vector &X, Vector &B, int copy_interior)
286 {
287  // Finish the matrix assembly and perform BC elimination, storing the
288  // eliminated part of the matrix.
289  FormSystemMatrix(ess_tdof_list, A);
290 
291  const Operator &P = *pfes->GetProlongationMatrix();
292  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
293 
294  // Transform the system and perform the elimination in B, based on the
295  // essential BC values from x. Restrict the BC part of x in X, and set the
296  // non-BC part to zero. Since there is no good initial guess for the Lagrange
297  // multipliers, set X = 0.0 for hybridization.
298  if (static_cond)
299  {
300  // Schur complement reduction to the exposed dofs
301  static_cond->ReduceSystem(x, b, X, B, copy_interior);
302  }
303  else if (hybridization)
304  {
305  // Reduction to the Lagrange multipliers system
306  HypreParVector true_X(pfes), true_B(pfes);
307  P.MultTranspose(b, true_B);
308  R.Mult(x, true_X);
309  p_mat.EliminateBC(p_mat_e, ess_tdof_list, true_X, true_B);
310  R.MultTranspose(true_B, b);
311  hybridization->ReduceRHS(true_B, B);
312  X.SetSize(B.Size());
313  X = 0.0;
314  }
315  else
316  {
317  // Variational restriction with P
318  X.SetSize(pfes->TrueVSize());
319  B.SetSize(X.Size());
320  P.MultTranspose(b, B);
321  R.Mult(x, X);
322  p_mat.EliminateBC(p_mat_e, ess_tdof_list, X, B);
323  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
324  }
325 }
326 
328  OperatorHandle &A)
329 {
330  // Finish the matrix assembly and perform BC elimination, storing the
331  // eliminated part of the matrix.
332  if (static_cond)
333  {
335  {
336  static_cond->SetEssentialTrueDofs(ess_tdof_list);
339  }
341  }
342  else
343  {
344  if (mat)
345  {
346  const int remove_zeros = 0;
347  Finalize(remove_zeros);
348  MFEM_VERIFY(p_mat.Ptr() == NULL && p_mat_e.Ptr() == NULL,
349  "The ParBilinearForm must be updated with Update() before "
350  "re-assembling the ParBilinearForm.");
352  delete mat;
353  mat = NULL;
354  delete mat_e;
355  mat_e = NULL;
356  p_mat_e.EliminateRowsCols(p_mat, ess_tdof_list);
357  }
358  if (hybridization)
359  {
361  }
362  else
363  {
364  A = p_mat;
365  }
366  }
367 }
368 
370  const Vector &X, const Vector &b, Vector &x)
371 {
372  const Operator &P = *pfes->GetProlongationMatrix();
373 
374  if (static_cond)
375  {
376  // Private dofs back solve
377  static_cond->ComputeSolution(b, X, x);
378  }
379  else if (hybridization)
380  {
381  // Primal unknowns recovery
382  HypreParVector true_X(pfes), true_B(pfes);
383  P.MultTranspose(b, true_B);
384  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
385  R.Mult(x, true_X); // get essential b.c. from x
386  hybridization->ComputeSolution(true_B, X, true_X);
387  x.SetSize(P.Height());
388  P.Mult(true_X, x);
389  }
390  else
391  {
392  // Apply conforming prolongation
393  x.SetSize(P.Height());
394  P.Mult(X, x);
395  }
396 }
397 
399 {
400  BilinearForm::Update(nfes);
401 
402  if (nfes)
403  {
404  pfes = dynamic_cast<ParFiniteElementSpace *>(nfes);
405  MFEM_VERIFY(pfes != NULL, "nfes must be a ParFiniteElementSpace!");
406  }
407 
408  p_mat.Clear();
409  p_mat_e.Clear();
410 }
411 
412 
414 {
415  // construct the block-diagonal matrix A
416  HypreParMatrix *A =
422  mat);
423 
426 
427  delete A;
428 
429  return rap;
430 }
431 
433 {
434  // construct the rectangular block-diagonal matrix dA
435  OperatorHandle dA(A.Type());
441  mat);
442 
443  OperatorHandle P_test(A.Type()), P_trial(A.Type());
444 
445  // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
447  P_trial.ConvertFrom(trial_pfes->Dof_TrueDof_Matrix());
448 
449  A.MakeRAP(P_test, dA, P_trial);
450 }
451 
452 /// Compute y += a (P^t A P) x, where x and y are vectors on the true dofs
454  const double a) const
455 {
456  if (X.ParFESpace() != trial_pfes)
457  {
460  }
461 
462  X.Distribute(&x);
463  mat->Mult(X, Y);
465 }
466 
467 
469 {
470  MFEM_ASSERT(mat, "matrix is not assembled");
471  MFEM_ASSERT(mat->Finalized(), "matrix is not finalized");
475  delete RA;
476  return RAP;
477 }
478 
480 const
481 {
482  MFEM_VERIFY(mat->Finalized(), "local matrix needs to be finalized for "
483  "GetParBlocks");
484 
486 
487  blocks.SetSize(range_fes->GetVDim(), domain_fes->GetVDim());
488 
489  RLP->GetBlocks(blocks,
492 
493  delete RLP;
494 }
495 
496 }
497 
498 #endif
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:265
int Size() const
Logical size of the array.
Definition: array.hpp:133
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1368
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:253
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:983
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:626
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition: handle.cpp:199
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:90
int * GetJ()
Definition: table.hpp:114
Array< BilinearFormIntegrator * > fbfi
Set of interior face Integrators to be applied.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1558
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:243
virtual const Operator * GetProlongationMatrix() const
Definition: pfespace.cpp:742
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition: hypre.cpp:1081
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1075
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:171
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:478
HypreParMatrix & GetParallelMatrix()
Return the parallel hybridized matrix.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:2030
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void SetDims(int rows, int nnz)
Definition: table.cpp:142
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:85
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:190
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.
Definition: handle.cpp:302
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
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...
Definition: staticcond.cpp:417
void EliminateReducedTrueDofs(const Array< int > &ess_rtdof_list, Matrix::DiagonalPolicy dpolicy)
Eliminate the given reduced true dofs from the Schur complement matrix S.
Definition: staticcond.cpp:288
OperatorHandle p_mat_e
void LoseData()
Call this if data has been stolen.
Definition: table.hpp:147
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:86
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void AssembleSharedFaces(int skip_zeros=1)
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
Definition: staticcond.cpp:477
const HYPRE_Int * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:332
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:732
void DeleteAll()
Delete whole array.
Definition: array.hpp:709
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
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 ...
Definition: operator.hpp:52
virtual void Update(FiniteElementSpace *nfes=NULL)
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:974
void SetSize(int m, int n)
Definition: array.hpp:308
ParFiniteElementSpace * pfes
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1008
MPI_Comm GetComm() const
Definition: pfespace.hpp:235
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
StaticCondensation * static_cond
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:614
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:224
void Clear()
Definition: table.cpp:373
SparseMatrix * mat
Sparse matrix to be associated with the form.
void ComputeSolution(const Vector &b, const Vector &sol_r, Vector &sol) const
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...
Definition: handle.cpp:251
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:1945
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:267
HypreParMatrix & GetParallelMatrix()
Return the parallel Schur complement matrix.
Definition: staticcond.hpp:159
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:418
virtual void Update(FiniteElementSpace *nfes=NULL)
ParFiniteElementSpace * test_pfes
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:242
double * GetData() const
Return element data, i.e. array A.
Definition: sparsemat.hpp:138
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:1978
Set the diagonal value to one.
Definition: matrix.hpp:35
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
int * GetI() const
Return the array I.
Definition: sparsemat.hpp:134
Dynamic 2D array using row-major layout.
Definition: array.hpp:289
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the &#39;dofs&#39; array to the given &#39;val&#39;.
Definition: vector.cpp:597
bool Finalized() const
Definition: sparsemat.hpp:310
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:245
HypreParMatrix * ParallelAssemble() const
Returns the matrix &quot;assembled&quot; on the true dofs.
SparseMatrix * mat_e
Matrix used to eliminate b.c.
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
ParFiniteElementSpace * domain_fes
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:486
void Finalize()
Finalize the construction of the Schur complement matrix.
Definition: staticcond.cpp:235
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:110
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:131
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
Table * GetFaceToElementTable() const
Definition: mesh.cpp:3880
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:106
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...
Definition: handle.cpp:60
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:322
ParFiniteElementSpace * range_fes
bool HasEliminatedBC() const
Return true if essential boundary conditions have been eliminated from the Schur complement matrix...
Definition: staticcond.hpp:142
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:563
FiniteElementSpace * fes
FE space on which the form lives.
void ParallelEliminateEssentialBC(const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const
Eliminate essential boundary DOFs from a parallel assembled system.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1335
DenseMatrix elemmat
void MakeRectangularBlockDiag(MPI_Comm comm, HYPRE_Int glob_num_rows, HYPRE_Int glob_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel rectangular block-diagonal matrix using the currently set...
Definition: handle.cpp:90
ParFiniteElementSpace * trial_pfes
int GetFaceNbrVSize() const
Definition: pfespace.hpp:327
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:244
OperatorHandle p_mat
Vector data type.
Definition: vector.hpp:48
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
ID for class HypreParMatrix.
Definition: operator.hpp:128
Hybridization * hybridization
void ReduceRHS(const Vector &b, Vector &b_r) const
int * GetI()
Definition: table.hpp:113
const Table & GetElementToDofTable() const
Definition: fespace.hpp:384
void SetEssentialTrueDofs(const Array< int > &ess_tdof_list)
Determine and save internally essential reduced true dofs.
Definition: staticcond.hpp:128
Table * GetFaceToAllElementTable() const
Definition: pmesh.cpp:1856
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
int * GetJ() const
Return the array J.
Definition: sparsemat.hpp:136
Class for parallel meshes.
Definition: pmesh.hpp:32
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition: handle.cpp:122
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
void MakeRAP(OperatorHandle &Rt, OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product R A P, where R = Rt^t.
Definition: handle.cpp:160
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:85