MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lor.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 "lor.hpp"
13 #include "pbilinearform.hpp"
14 
15 namespace mfem
16 {
17 
18 void LORBase::AddIntegrators(BilinearForm &a_from,
19  BilinearForm &a_to,
20  GetIntegratorsFn get_integrators,
21  AddIntegratorFn add_integrator,
22  const IntegrationRule *ir)
23 {
24  Array<BilinearFormIntegrator*> *integrators = (a_from.*get_integrators)();
25  for (int i=0; i<integrators->Size(); ++i)
26  {
27  (a_to.*add_integrator)((*integrators)[i]);
28  ir_map[(*integrators)[i]] = ((*integrators)[i])->GetIntegrationRule();
29  if (ir) { ((*integrators)[i])->SetIntegrationRule(*ir); }
30  }
31 }
32 
33 void LORBase::AddIntegratorsAndMarkers(BilinearForm &a_from,
34  BilinearForm &a_to,
35  GetIntegratorsFn get_integrators,
36  GetMarkersFn get_markers,
37  AddIntegratorMarkersFn add_integrator,
38  const IntegrationRule *ir)
39 {
40  Array<BilinearFormIntegrator*> *integrators = (a_from.*get_integrators)();
41  Array<Array<int>*> *markers = (a_from.*get_markers)();
42 
43  for (int i=0; i<integrators->Size(); ++i)
44  {
45  (a_to.*add_integrator)((*integrators)[i], *(*markers[i]));
46  ir_map[(*integrators)[i]] = ((*integrators)[i])->GetIntegrationRule();
47  if (ir) { ((*integrators)[i])->SetIntegrationRule(*ir); }
48  }
49 }
50 
51 void LORBase::ResetIntegrationRules(GetIntegratorsFn get_integrators)
52 {
53  Array<BilinearFormIntegrator*> *integrators = (a->*get_integrators)();
54  for (int i=0; i<integrators->Size(); ++i)
55  {
56  ((*integrators)[i])->SetIntegrationRule(*ir_map[(*integrators)[i]]);
57  }
58 }
59 
61 {
63  if (dynamic_cast<const H1_FECollection*>(fec)) { return H1; }
64  else if (dynamic_cast<const ND_FECollection*>(fec)) { return ND; }
65  else if (dynamic_cast<const RT_FECollection*>(fec)) { return RT; }
66  else if (dynamic_cast<const L2_FECollection*>(fec)) { return L2; }
67  else { MFEM_ABORT("Bad LOR space type."); }
68  return INVALID;
69 }
70 
72 {
73  FESpaceType type = GetFESpaceType();
74  return (type == L2 || type == RT) ? 0 : 1;
75 }
76 
78 {
79  FESpaceType type = GetFESpaceType();
80  MFEM_VERIFY(type != H1 && type != L2, "");
81 
82  auto get_dof_map = [](FiniteElementSpace &fes, int i)
83  {
84  const FiniteElement *fe = fes.GetFE(i);
85  auto tfe = dynamic_cast<const TensorBasisElement*>(fe);
86  MFEM_ASSERT(tfe != NULL, "");
87  return tfe->GetDofMap();
88  };
89 
90  FiniteElementSpace &fes_lor = *fes;
91  Mesh &mesh_lor = *fes_lor.GetMesh();
92  int dim = mesh_lor.Dimension();
93  const CoarseFineTransformations &cf_tr = mesh_lor.GetRefinementTransforms();
94 
95  perm_.SetSize(fes_lor.GetVSize());
96 
97  Array<int> vdof_ho, vdof_lor;
98  for (int ilor=0; ilor<mesh_lor.GetNE(); ++ilor)
99  {
100  int iho = cf_tr.embeddings[ilor].parent;
101  int lor_index = cf_tr.embeddings[ilor].matrix;
102 
103  fes_ho.GetElementVDofs(iho, vdof_ho);
104  fes_lor.GetElementVDofs(ilor, vdof_lor);
105 
106  if (type == L2)
107  {
108  perm_[vdof_lor[0]] = vdof_ho[lor_index];
109  continue;
110  }
111 
112  int p = fes_ho.GetOrder(iho);
113  int p1 = p+1;
114  int ndof_per_dim = (dim == 2) ? p*p1 : type == ND ? p*p1*p1 : p*p*p1;
115 
116  const Array<int> &dofmap_ho = get_dof_map(fes_ho, iho);
117  const Array<int> &dofmap_lor = get_dof_map(fes_lor, ilor);
118 
119  int off_x = lor_index % p;
120  int off_y = (lor_index / p) % p;
121  int off_z = (lor_index / p) / p;
122 
123  auto set_perm = [&](int off_lor, int off_ho, int n1, int n2)
124  {
125  for (int i1=0; i1<2; ++i1)
126  {
127  int m = (dim == 2 || type == RT) ? 1 : 2;
128  for (int i2=0; i2<m; ++i2)
129  {
130  int i;
131  i = dofmap_lor[off_lor + i1 + i2*2];
132  int s1 = i < 0 ? -1 : 1;
133  int idof_lor = vdof_lor[absdof(i)];
134  i = dofmap_ho[off_ho + i1*n1 + i2*n2];
135  int s2 = i < 0 ? -1 : 1;
136  int idof_ho = vdof_ho[absdof(i)];
137  int s3 = idof_lor < 0 ? -1 : 1;
138  int s4 = idof_ho < 0 ? -1 : 1;
139  int s = s1*s2*s3*s4;
140  i = absdof(idof_ho);
141  perm_[absdof(idof_lor)] = s < 0 ? -1-absdof(i) : absdof(i);
142  }
143  }
144  };
145 
146  int offset;
147 
148  if (type == ND)
149  {
150  // x
151  offset = off_x + off_y*p + off_z*p*p1;
152  set_perm(0, offset, p, p*p1);
153  // y
154  offset = ndof_per_dim + off_x + off_y*(p1) + off_z*p1*p;
155  set_perm(dim == 2 ? 2 : 4, offset, 1, p*p1);
156  // z
157  if (dim == 3)
158  {
159  offset = 2*ndof_per_dim + off_x + off_y*p1 + off_z*p1*p1;
160  set_perm(8, offset, 1, p+1);
161  }
162  }
163  else if (type == RT)
164  {
165  // x
166  offset = off_x + off_y*p1 + off_z*p*p1;
167  set_perm(0, offset, 1, 0);
168  // y
169  offset = ndof_per_dim + off_x + off_y*p + off_z*p1*p;
170  set_perm(2, offset, p, 0);
171  // z
172  if (dim == 3)
173  {
174  offset = 2*ndof_per_dim + off_x + off_y*p + off_z*p*p;
175  set_perm(4, offset, p*p, 0);
176  }
177  }
178  }
179 }
180 
182 {
183  FESpaceType type = GetFESpaceType();
184  if (type == H1 || type == L2 || nonconforming)
185  {
186  // H1 and L2: no permutation necessary, return identity
188  for (int i=0; i<perm.Size(); ++i) { perm[i] = i; }
189  return;
190  }
191 
192 #ifdef MFEM_USE_MPI
193  ParFiniteElementSpace *pfes_ho
194  = dynamic_cast<ParFiniteElementSpace*>(&fes_ho);
195  ParFiniteElementSpace *pfes_lor = dynamic_cast<ParFiniteElementSpace*>(fes);
196  if (pfes_ho && pfes_lor)
197  {
198  Array<int> l_perm;
200  perm.SetSize(pfes_lor->GetTrueVSize());
201  for (int i=0; i<l_perm.Size(); ++i)
202  {
203  int j = l_perm[i];
204  int s = j < 0 ? -1 : 1;
205  int t_i = pfes_lor->GetLocalTDofNumber(i);
206  int t_j = pfes_ho->GetLocalTDofNumber(absdof(j));
207  // Either t_i and t_j both -1, or both non-negative
208  if ((t_i < 0 && t_j >=0) || (t_j < 0 && t_i >= 0))
209  {
210  MFEM_ABORT("Inconsistent DOF numbering");
211  }
212  if (t_i < 0) { continue; }
213  perm[t_i] = s < 0 ? -1 - t_j : t_j;
214  }
215  }
216  else
217 #endif
218  {
220  }
221 }
222 
224 {
225  if (perm.Size() == 0) { ConstructDofPermutation(); }
226  return perm;
227 }
228 
230 {
231  FESpaceType type = GetFESpaceType();
232  return (type == H1 || type == L2 || nonconforming) ? false : true;
233 }
234 
236 {
237  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
238  return A;
239 }
240 
241 void LORBase::AssembleSystem(BilinearForm &a_ho, const Array<int> &ess_dofs)
242 {
244  AddIntegrators(a_ho, *a, &BilinearForm::GetDBFI,
246  AddIntegrators(a_ho, *a, &BilinearForm::GetFBFI,
248  AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBBFI,
251  AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBFBFI,
254  a->Assemble();
256  {
257  const Array<int> &p = GetDofPermutation();
258  // Form inverse permutation: given high-order dof i, pi[i] is corresp. LO
259  Array<int> pi(p.Size());
260  for (int i=0; i<p.Size(); ++i)
261  {
262  pi[absdof(p[i])] = i;
263  }
264  Array<int> ess_dofs_perm(ess_dofs.Size());
265  for (int i=0; i<ess_dofs.Size(); ++i)
266  {
267  ess_dofs_perm[i] = pi[ess_dofs[i]];
268  }
269  a->FormSystemMatrix(ess_dofs_perm, A);
270  }
271  else
272  {
273  a->FormSystemMatrix(ess_dofs, A);
274  }
275  ResetIntegrationRules(&BilinearForm::GetDBFI);
276  ResetIntegrationRules(&BilinearForm::GetFBFI);
277  ResetIntegrationRules(&BilinearForm::GetBBFI);
278  ResetIntegrationRules(&BilinearForm::GetBFBFI);
279 }
280 
282 {
284  {
285  Array<int> p;
288  }
289  else
290  {
292  }
293  nonconforming = true;
294 }
295 
296 template <typename FEC>
298 {
299  const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
300  if (fec)
301  {
302  int btype = fec->GetBasisType();
303  if (btype != BasisType::GaussLobatto)
304  {
305  mfem::err << "\nWARNING: Constructing low-order refined "
306  << "discretization with basis type\n"
307  << BasisType::Name(btype) << ". "
308  << "The LOR discretization is only spectrally equivalent\n"
309  << "with Gauss-Lobatto basis.\n" << std::endl;
310  }
311  }
312 }
313 
314 template <typename FEC>
316 {
317  const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
318  if (fec)
319  {
320  int cbtype = fec->GetClosedBasisType();
321  int obtype = fec->GetOpenBasisType();
322  if (cbtype != BasisType::GaussLobatto || obtype != BasisType::IntegratedGLL)
323  {
324  mfem::err << "\nWARNING: Constructing vector low-order refined "
325  << "discretization with basis type \npair ("
326  << BasisType::Name(cbtype) << ", "
327  << BasisType::Name(obtype) << "). "
328  << "The LOR discretization is only spectrally\nequivalent "
329  << "with basis types (Gauss-Lobatto, IntegratedGLL).\n"
330  << std::endl;
331  }
332  }
333 }
334 
336 {
337  CheckScalarBasisType<H1_FECollection>(fes);
338  CheckVectorBasisType<ND_FECollection>(fes);
339  CheckVectorBasisType<RT_FECollection>(fes);
340  // L2 is a bit more complicated, for now don't verify basis type
341 }
342 
344  : irs(0, Quadrature1D::GaussLobatto), fes_ho(fes_ho_)
345 {
346  Mesh &mesh_ = *fes_ho_.GetMesh();
347  int dim = mesh_.Dimension();
348  Array<Geometry::Type> geoms;
349  mesh_.GetGeometries(dim, geoms);
350  if (geoms.Size() == 1 && Geometry::IsTensorProduct(geoms[0]))
351  {
352  ir_el = &irs.Get(geoms[0], 1);
353  ir_face = &irs.Get(Geometry::TensorProductGeometry(dim-1), 1);
354  }
355  else
356  {
357  ir_el = NULL;
358  ir_face = NULL;
359  }
360  a = NULL;
361 }
362 
364 {
365  delete a;
366  delete fes;
367  delete fec;
368  delete mesh;
369 }
370 
372  const Array<int> &ess_tdof_list,
373  int ref_type)
374  : LORDiscretization(*a_ho_.FESpace(), ref_type)
375 {
376  a = new BilinearForm(fes);
377  AssembleSystem(a_ho_, ess_tdof_list);
378 }
379 
381  int ref_type) : LORBase(fes_ho)
382 {
383  CheckBasisType(fes_ho);
384 
385  // TODO: support variable-order spaces
386  MFEM_VERIFY(!fes_ho.IsVariableOrder(),
387  "Cannot construct LOR operators on variable-order spaces");
388 
389  int order = fes_ho.GetMaxElementOrder();
390  if (GetFESpaceType() == L2) { ++order; }
391 
392  Mesh &mesh_ho = *fes_ho.GetMesh();
393  mesh = new Mesh(Mesh::MakeRefined(mesh_ho, order, ref_type));
394 
395  fec = fes_ho.FEColl()->Clone(GetLOROrder());
396  fes = new FiniteElementSpace(mesh, fec);
397  if (fes_ho.Nonconforming()) { SetupNonconforming(); }
398 
400 }
401 
403 {
404  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
405  return *A.As<SparseMatrix>();
406 }
407 
408 #ifdef MFEM_USE_MPI
409 
411  const Array<int> &ess_tdof_list,
412  int ref_type)
413  : ParLORDiscretization(*a_ho_.ParFESpace(), ref_type)
414 {
415  a = new ParBilinearForm(static_cast<ParFiniteElementSpace*>(fes));
416  AssembleSystem(a_ho_, ess_tdof_list);
417 }
418 
420  int ref_type) : LORBase(fes_ho)
421 {
422  if (fes_ho.GetMyRank() == 0) { CheckBasisType(fes_ho); }
423  // TODO: support variable-order spaces
424  MFEM_VERIFY(!fes_ho.IsVariableOrder(),
425  "Cannot construct LOR operators on variable-order spaces");
426 
427  int order = fes_ho.GetMaxElementOrder();
428  if (GetFESpaceType() == L2) { ++order; }
429 
430  ParMesh &mesh_ho = *fes_ho.GetParMesh();
431  ParMesh *pmesh = new ParMesh(ParMesh::MakeRefined(mesh_ho, order, ref_type));
432  mesh = pmesh;
433 
434  fec = fes_ho.FEColl()->Clone(GetLOROrder());
435  ParFiniteElementSpace *pfes = new ParFiniteElementSpace(pmesh, fec);
436  fes = pfes;
437  if (fes_ho.Nonconforming()) { SetupNonconforming(); }
438 
440 }
441 
443 {
444  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
445  return *A.As<HypreParMatrix>();
446 }
447 
449 {
450  return static_cast<ParFiniteElementSpace&>(*fes);
451 }
452 
453 #endif
454 
455 } // namespace mfem
Abstract class for all finite elements.
Definition: fe.hpp:243
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition: lor.hpp:22
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
LORDiscretization(BilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Construct the low-order refined version of a_ho using the given list of essential DOFs...
Definition: lor.cpp:371
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system.
Definition: lor.cpp:235
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:99
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
Create and assemble a low-order refined version of a ParBilinearForm.
Definition: lor.hpp:154
static Type TensorProductGeometry(int dim)
Definition: geom.hpp:112
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
Definition: lor.cpp:60
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:432
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
ParMesh * GetParMesh() const
Definition: pfespace.hpp:267
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
Definition: pmesh.cpp:1361
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
Definition: lor.cpp:71
Integrated GLL indicator functions.
Definition: fe.hpp:41
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:262
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
FiniteElementSpace & fes_ho
Definition: lor.hpp:64
static Mesh MakeRefined(Mesh &orig_mesh, int ref_factor, int ref_type)
Create a refined (by any factor) version of orig_mesh.
Definition: mesh.cpp:3307
void GetGeometries(int dim, Array< Geometry::Type > &el_geoms) const
Return all element geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5439
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
Definition: lor.cpp:223
Array< Array< int > * > * GetBFBFI_Marker()
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the in...
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
Abstract parallel finite element space.
Definition: pfespace.hpp:28
A class container for 1D quadrature type constants.
Definition: intrules.hpp:289
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:819
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe.hpp:94
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes, const Array< int > *perm)
Copies the prolongation and restriction matrices from fes.
Definition: fespace.cpp:96
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition: lor.cpp:442
Array< Array< int > * > * GetBBFI_Marker()
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the i...
void CheckScalarBasisType(const FiniteElementSpace &fes)
Definition: lor.cpp:297
Create and assemble a low-order refined version of a BilinearForm.
Definition: lor.hpp:128
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:279
ID for class SparseMatrix.
Definition: operator.hpp:255
FiniteElementSpace * fes
Definition: lor.hpp:67
Data type sparse matrix.
Definition: sparsemat.hpp:41
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
Mesh * mesh
Definition: lor.hpp:65
void CheckBasisType(const FiniteElementSpace &fes)
Definition: lor.cpp:335
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition: fespace.hpp:428
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:65
Array< int > perm
Definition: lor.hpp:70
void SetupNonconforming()
Definition: lor.cpp:281
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:543
int Dimension() const
Definition: mesh.hpp:911
BilinearForm * a
Definition: lor.hpp:68
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
void UseExternalIntegrators()
Indicate that integrators are not owned by the BilinearForm.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
int GetOrder(int i) const
Returns the polynomial degree of the i&#39;th finite element.
Definition: fespace.hpp:524
void ConstructDofPermutation() const
Definition: lor.cpp:181
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9255
bool Nonconforming() const
Definition: pfespace.hpp:398
bool Nonconforming() const
Definition: fespace.hpp:417
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:293
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
static bool IsTensorProduct(Type geom)
Definition: geom.hpp:107
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
Definition: lor.cpp:448
FiniteElementCollection * fec
Definition: lor.hpp:66
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
ParLORDiscretization(ParBilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Construct the low-order refined version of a_ho using the given list of essential DOFs...
Definition: lor.cpp:410
int dim
Definition: ex24.cpp:53
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system.
Definition: lor.cpp:241
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
virtual 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:2388
void ConstructLocalDofPermutation(Array< int > &perm_) const
Definition: lor.cpp:77
Class for parallel bilinear form.
LORBase(FiniteElementSpace &fes_ho_)
Definition: lor.cpp:343
void CheckVectorBasisType(const FiniteElementSpace &fes)
Definition: lor.cpp:315
bool RequiresDofPermutation() const
Definition: lor.cpp:229
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
Definition: lor.cpp:402
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:554
ID for class HypreParMatrix.
Definition: operator.hpp:256
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:60
RefCoord s[3]
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
OperatorHandle A
Definition: lor.hpp:69
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:127
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
bool nonconforming
Definition: lor.hpp:71