MFEM  v4.2.0
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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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 = Memory<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 off-diagonal
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  int Elem2NbrNo = T->Elem2No - pmesh->GetNE();
202  pfes->GetElementVDofs(T->Elem1No, vdofs1);
203  pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
204  vdofs1.Copy(vdofs_all);
205  for (int j = 0; j < vdofs2.Size(); j++)
206  {
207  if (vdofs2[j] >= 0)
208  {
209  vdofs2[j] += height;
210  }
211  else
212  {
213  vdofs2[j] -= height;
214  }
215  }
216  vdofs_all.Append(vdofs2);
217  for (int k = 0; k < fbfi.Size(); k++)
218  {
219  fbfi[k]->AssembleFaceMatrix(*pfes->GetFE(T->Elem1No),
220  *pfes->GetFaceNbrFE(Elem2NbrNo),
221  *T, elemmat);
222  if (keep_nbr_block)
223  {
224  mat->AddSubMatrix(vdofs_all, vdofs_all, elemmat, skip_zeros);
225  }
226  else
227  {
228  mat->AddSubMatrix(vdofs1, vdofs_all, elemmat, skip_zeros);
229  }
230  }
231  }
232 }
233 
234 void ParBilinearForm::Assemble(int skip_zeros)
235 {
236  if (mat == NULL && fbfi.Size() > 0)
237  {
239  pAllocMat();
240  }
241 
242  BilinearForm::Assemble(skip_zeros);
243 
244  if (!ext && fbfi.Size() > 0)
245  {
246  AssembleSharedFaces(skip_zeros);
247  }
248 }
249 
250 void ParBilinearForm
252  HypreParMatrix &A, const HypreParVector &X,
253  HypreParVector &B) const
254 {
255  Array<int> dof_list;
256 
257  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
258 
259  // do the parallel elimination
260  A.EliminateRowsCols(dof_list, X, B);
261 }
262 
265  HypreParMatrix &A) const
266 {
267  Array<int> dof_list;
268 
269  pfes->GetEssentialTrueDofs(bdr_attr_is_ess, dof_list);
270 
271  return A.EliminateRowsCols(dof_list);
272 }
273 
274 void ParBilinearForm::TrueAddMult(const Vector &x, Vector &y, const double a)
275 const
276 {
277  MFEM_VERIFY(fbfi.Size() == 0, "the case of interior face integrators is not"
278  " implemented");
279 
280  if (X.ParFESpace() != pfes)
281  {
282  X.SetSpace(pfes);
283  Y.SetSpace(pfes);
284  }
285 
286  X.Distribute(&x);
287  if (ext)
288  {
289  ext->Mult(X, Y);
290  }
291  else
292  {
293  mat->Mult(X, Y);
294  }
295  pfes->Dof_TrueDof_Matrix()->MultTranspose(a, Y, 1.0, y);
296 }
297 
299  const Array<int> &ess_tdof_list, Vector &x, Vector &b,
300  OperatorHandle &A, Vector &X, Vector &B, int copy_interior)
301 {
302  if (ext)
303  {
304  ext->FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
305  return;
306  }
307 
308  // Finish the matrix assembly and perform BC elimination, storing the
309  // eliminated part of the matrix.
310  FormSystemMatrix(ess_tdof_list, A);
311 
312  const Operator &P = *pfes->GetProlongationMatrix();
313  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
314 
315  // Transform the system and perform the elimination in B, based on the
316  // essential BC values from x. Restrict the BC part of x in X, and set the
317  // non-BC part to zero. Since there is no good initial guess for the Lagrange
318  // multipliers, set X = 0.0 for hybridization.
319  if (static_cond)
320  {
321  // Schur complement reduction to the exposed dofs
322  static_cond->ReduceSystem(x, b, X, B, copy_interior);
323  }
324  else if (hybridization)
325  {
326  // Reduction to the Lagrange multipliers system
327  HypreParVector true_X(pfes), true_B(pfes);
328  P.MultTranspose(b, true_B);
329  R.Mult(x, true_X);
330  p_mat.EliminateBC(p_mat_e, ess_tdof_list, true_X, true_B);
331  R.MultTranspose(true_B, b);
332  hybridization->ReduceRHS(true_B, B);
333  X.SetSize(B.Size());
334  X = 0.0;
335  }
336  else
337  {
338  // Variational restriction with P
339  X.SetSize(pfes->TrueVSize());
340  B.SetSize(X.Size());
341  P.MultTranspose(b, B);
342  R.Mult(x, X);
343  p_mat.EliminateBC(p_mat_e, ess_tdof_list, X, B);
344  if (!copy_interior) { X.SetSubVectorComplement(ess_tdof_list, 0.0); }
345  }
346 }
347 
349  OperatorHandle &A)
350 {
351  if (ext)
352  {
353  ext->FormSystemMatrix(ess_tdof_list, A);
354  return;
355  }
356 
357  // Finish the matrix assembly and perform BC elimination, storing the
358  // eliminated part of the matrix.
359  if (static_cond)
360  {
362  {
363  static_cond->SetEssentialTrueDofs(ess_tdof_list);
366  }
368  }
369  else
370  {
371  if (mat)
372  {
373  const int remove_zeros = 0;
374  Finalize(remove_zeros);
375  MFEM_VERIFY(p_mat.Ptr() == NULL && p_mat_e.Ptr() == NULL,
376  "The ParBilinearForm must be updated with Update() before "
377  "re-assembling the ParBilinearForm.");
379  delete mat;
380  mat = NULL;
381  delete mat_e;
382  mat_e = NULL;
383  p_mat_e.EliminateRowsCols(p_mat, ess_tdof_list);
384  }
385  if (hybridization)
386  {
388  }
389  else
390  {
391  A = p_mat;
392  }
393  }
394 }
395 
397  const Vector &X, const Vector &b, Vector &x)
398 {
399  if (ext)
400  {
401  ext->RecoverFEMSolution(X, b, x);
402  return;
403  }
404 
405  const Operator &P = *pfes->GetProlongationMatrix();
406 
407  if (static_cond)
408  {
409  // Private dofs back solve
410  static_cond->ComputeSolution(b, X, x);
411  }
412  else if (hybridization)
413  {
414  // Primal unknowns recovery
415  HypreParVector true_X(pfes), true_B(pfes);
416  P.MultTranspose(b, true_B);
417  const SparseMatrix &R = *pfes->GetRestrictionMatrix();
418  R.Mult(x, true_X); // get essential b.c. from x
419  hybridization->ComputeSolution(true_B, X, true_X);
420  x.SetSize(P.Height());
421  P.Mult(true_X, x);
422  }
423  else
424  {
425  // Apply conforming prolongation
426  x.SetSize(P.Height());
427  P.Mult(X, x);
428  }
429 }
430 
432 {
433  BilinearForm::Update(nfes);
434 
435  if (nfes)
436  {
437  pfes = dynamic_cast<ParFiniteElementSpace *>(nfes);
438  MFEM_VERIFY(pfes != NULL, "nfes must be a ParFiniteElementSpace!");
439  }
440 
441  p_mat.Clear();
442  p_mat_e.Clear();
443 }
444 
445 
447 {
448  // construct the block-diagonal matrix A
449  HypreParMatrix *A =
455  mat);
456 
459 
460  delete A;
461 
462  return rap;
463 }
464 
466 {
467  // construct the rectangular block-diagonal matrix dA
468  OperatorHandle dA(A.Type());
474  mat);
475 
476  OperatorHandle P_test(A.Type()), P_trial(A.Type());
477 
478  // TODO - construct the Dof_TrueDof_Matrix directly in the required format.
480  P_trial.ConvertFrom(trial_pfes->Dof_TrueDof_Matrix());
481 
482  A.MakeRAP(P_test, dA, P_trial);
483 }
484 
485 /// Compute y += a (P^t A P) x, where x and y are vectors on the true dofs
487  const double a) const
488 {
489  if (X.ParFESpace() != trial_pfes)
490  {
493  }
494 
495  X.Distribute(&x);
496  mat->Mult(X, Y);
498 }
499 
501  const Array<int>
502  &trial_tdof_list,
503  const Array<int> &test_tdof_list,
504  OperatorHandle &A)
505 {
506  if (ext)
507  {
508  ext->FormRectangularSystemOperator(trial_tdof_list, test_tdof_list, A);
509  return;
510  }
511 
512  if (mat)
513  {
514  Finalize();
516  delete mat;
517  mat = NULL;
518  delete mat_e;
519  mat_e = NULL;
520  HypreParMatrix *temp =
521  p_mat.As<HypreParMatrix>()->EliminateCols(trial_tdof_list);
522  p_mat.As<HypreParMatrix>()->EliminateRows(test_tdof_list);
523  p_mat_e.Reset(temp, true);
524  }
525 
526  A = p_mat;
527 }
528 
530  const Array<int>
531  &trial_tdof_list,
532  const Array<int> &test_tdof_list, Vector &x,
533  Vector &b, OperatorHandle &A, Vector &X,
534  Vector &B)
535 {
536  if (ext)
537  {
538  ext->FormRectangularLinearSystem(trial_tdof_list, test_tdof_list,
539  x, b, A, X, B);
540  return;
541  }
542 
543  FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list, A);
544 
545  const Operator *test_P = test_pfes->GetProlongationMatrix();
546  const SparseMatrix *trial_R = trial_pfes->GetRestrictionMatrix();
547 
550  test_P->MultTranspose(b, B);
551  trial_R->Mult(x, X);
552 
553  p_mat_e.As<HypreParMatrix>()->Mult(-1.0, X, 1.0, B);
554  B.SetSubVector(test_tdof_list, 0.0);
555 }
556 
558 {
559  MFEM_ASSERT(mat, "Matrix is not assembled");
560  MFEM_ASSERT(mat->Finalized(), "Matrix is not finalized");
564  delete RA;
565  return RAP;
566 }
567 
569 const
570 {
571  MFEM_VERIFY(mat->Finalized(), "Local matrix needs to be finalized for "
572  "GetParBlocks");
573 
575 
576  blocks.SetSize(range_fes->GetVDim(), domain_fes->GetVDim());
577 
578  RLP->GetBlocks(blocks,
581 
582  delete RLP;
583 }
584 
585 }
586 
587 #endif
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:412
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1368
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1359
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
void ConvertFrom(OperatorHandle &A)
Convert the given OperatorHandle A to the currently set type id.
Definition: handle.cpp:200
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:96
int * GetJ()
Definition: table.hpp:114
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
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:1640
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
ParGridFunction X
Auxiliary objects used in TrueAddMult().
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:247
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
BilinearFormExtension * ext
Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA)...
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:459
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:1045
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:173
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:382
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
HypreParMatrix & GetParallelMatrix()
Return the parallel hybridized matrix.
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:2508
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:178
void SetDims(int rows, int nnz)
Definition: table.cpp:144
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:91
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:831
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
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:340
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
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:418
int * GetI()
Return the array I.
Definition: sparsemat.hpp:173
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
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:154
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:87
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AssembleSharedFaces(int skip_zeros=1)
void ComputeSolution(const Vector &b, const Vector &sc_sol, Vector &sol) const
Definition: staticcond.cpp:478
const HYPRE_Int * GetFaceNbrGlobalDofMap()
Definition: pfespace.hpp:349
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:885
void DeleteAll()
Delete the whole array.
Definition: array.hpp:821
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:92
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)=0
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:183
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:1165
void SetSize(int m, int n)
Definition: array.hpp:352
ParFiniteElementSpace * pfes
Points to the same object as fes.
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1199
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
virtual void FormRectangularSystemOperator(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)=0
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
StaticCondensation * static_cond
Owned.
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:726
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
void Clear()
Definition: table.cpp:381
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
double b
Definition: lissajous.cpp:42
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:252
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:2399
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
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:414
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
ParFiniteElementSpace * test_pfes
Points to the same object as test_fes.
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
SparseMatrix * mat
Owned.
Set the diagonal value to one.
Definition: operator.hpp:48
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2436
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Definition: vector.cpp:656
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:466
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:249
HypreParMatrix * ParallelAssemble() const
Returns the matrix &quot;assembled&quot; on the true dofs.
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
ParFiniteElementSpace * domain_fes
Points to the same object as trial_fes.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
MixedBilinearFormExtension * ext
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
void Finalize()
Finalize the construction of the Schur complement matrix.
Definition: staticcond.cpp:235
ParGridFunction X
Auxiliary objects used in TrueAddMult().
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:116
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition: operator.cpp:85
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:27
Table * GetFaceToElementTable() const
Definition: mesh.cpp:4903
double a
Definition: lissajous.cpp:41
virtual 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:112
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:339
SparseMatrix * mat_e
Owned.
ParFiniteElementSpace * range_fes
Points to the same object as test_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:721
FiniteElementSpace * fes
FE space on which the form lives. Not owned.
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 in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Form the parallel linear system A X = B, corresponding to this mixed bilinear form and the linear for...
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:91
ParFiniteElementSpace * trial_pfes
Points to the same object as trial_fes.
int GetFaceNbrVSize() const
Definition: pfespace.hpp:344
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:248
OperatorHandle p_mat
Vector data type.
Definition: vector.hpp:51
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)=0
ID for class HypreParMatrix.
Definition: operator.hpp:241
Hybridization * hybridization
Owned.
void ReduceRHS(const Vector &b, Vector &b_r) const
int * GetI()
Definition: table.hpp:113
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:543
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:2310
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
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.
Class for parallel meshes.
Definition: pmesh.hpp:32
OperatorHandle p_mat
Matrix and eliminated matrix.
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
void MakePtAP(OperatorHandle &A, OperatorHandle &P)
Reset the OperatorHandle to hold the product P^t A P.
Definition: handle.cpp:123
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
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:161
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:137
virtual void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)=0
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108