MFEM  v4.4.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-2022, 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_marker,
38  AddIntegratorFn add_integrator,
39  const IntegrationRule *ir)
40 {
41  Array<BilinearFormIntegrator*> *integrators = (a_from.*get_integrators)();
42  Array<Array<int>*> *markers = (a_from.*get_markers)();
43 
44  for (int i=0; i<integrators->Size(); ++i)
45  {
46  if (*markers[i])
47  {
48  (a_to.*add_integrator_marker)((*integrators)[i], *(*markers[i]));
49  }
50  else
51  {
52  (a_to.*add_integrator)((*integrators)[i]);
53  }
54  ir_map[(*integrators)[i]] = ((*integrators)[i])->GetIntegrationRule();
55  if (ir) { ((*integrators)[i])->SetIntegrationRule(*ir); }
56  }
57 }
58 
59 void LORBase::ResetIntegrationRules(GetIntegratorsFn get_integrators)
60 {
61  Array<BilinearFormIntegrator*> *integrators = (a->*get_integrators)();
62  for (int i=0; i<integrators->Size(); ++i)
63  {
64  ((*integrators)[i])->SetIntegrationRule(*ir_map[(*integrators)[i]]);
65  }
66 }
67 
69 {
70  const FiniteElementCollection *fec_ho = fes_ho.FEColl();
71  if (dynamic_cast<const H1_FECollection*>(fec_ho)) { return H1; }
72  else if (dynamic_cast<const ND_FECollection*>(fec_ho)) { return ND; }
73  else if (dynamic_cast<const RT_FECollection*>(fec_ho)) { return RT; }
74  else if (dynamic_cast<const L2_FECollection*>(fec_ho)) { return L2; }
75  else { MFEM_ABORT("Bad LOR space type."); }
76  return INVALID;
77 }
78 
80 {
81  FESpaceType type = GetFESpaceType();
82  return (type == L2 || type == RT) ? 0 : 1;
83 }
84 
86 {
87  FESpaceType type = GetFESpaceType();
88  MFEM_VERIFY(type != H1 && type != L2, "");
89 
90  auto get_dof_map = [](FiniteElementSpace &fes_, int i)
91  {
92  const FiniteElement *fe = fes_.GetFE(i);
93  auto tfe = dynamic_cast<const TensorBasisElement*>(fe);
94  MFEM_ASSERT(tfe != NULL, "");
95  return tfe->GetDofMap();
96  };
97 
98  FiniteElementSpace &fes_lor = *fes;
99  Mesh &mesh_lor = *fes_lor.GetMesh();
100  int dim = mesh_lor.Dimension();
101  const CoarseFineTransformations &cf_tr = mesh_lor.GetRefinementTransforms();
102 
103  using GeomRef = std::pair<Geometry::Type, int>;
104  std::map<GeomRef, int> point_matrices_offsets;
105  perm_.SetSize(fes_lor.GetVSize());
106 
107  Array<int> vdof_ho, vdof_lor;
108  for (int ilor=0; ilor<mesh_lor.GetNE(); ++ilor)
109  {
110  int iho = cf_tr.embeddings[ilor].parent;
111  int p = fes_ho.GetOrder(iho);
112  int lor_index = cf_tr.embeddings[ilor].matrix;
113  // We use the point matrix index to identify the local LOR element index
114  // within the high-order coarse element.
115  //
116  // In variable-order spaces, the point matrices for each order are
117  // concatenated sequentially, so for the given element order, we need to
118  // find the offset that will give us the point matrix index relative to
119  // the current element order only.
120  GeomRef id(mesh_lor.GetElementBaseGeometry(ilor), p);
121  if (point_matrices_offsets.find(id) == point_matrices_offsets.end())
122  {
123  point_matrices_offsets[id] = lor_index;
124  }
125  lor_index -= point_matrices_offsets[id];
126 
127  fes_ho.GetElementVDofs(iho, vdof_ho);
128  fes_lor.GetElementVDofs(ilor, vdof_lor);
129 
130  if (type == L2)
131  {
132  perm_[vdof_lor[0]] = vdof_ho[lor_index];
133  continue;
134  }
135 
136  int p1 = p+1;
137  int ndof_per_dim = (dim == 2) ? p*p1 : type == ND ? p*p1*p1 : p*p*p1;
138 
139  const Array<int> &dofmap_ho = get_dof_map(fes_ho, iho);
140  const Array<int> &dofmap_lor = get_dof_map(fes_lor, ilor);
141 
142  int off_x = lor_index % p;
143  int off_y = (lor_index / p) % p;
144  int off_z = (lor_index / p) / p;
145 
146  auto set_perm = [&](int off_lor, int off_ho, int n1, int n2)
147  {
148  for (int i1=0; i1<2; ++i1)
149  {
150  int m = (dim == 2 || type == RT) ? 1 : 2;
151  for (int i2=0; i2<m; ++i2)
152  {
153  int i;
154  i = dofmap_lor[off_lor + i1 + i2*2];
155  int s1 = i < 0 ? -1 : 1;
156  int idof_lor = vdof_lor[absdof(i)];
157  i = dofmap_ho[off_ho + i1*n1 + i2*n2];
158  int s2 = i < 0 ? -1 : 1;
159  int idof_ho = vdof_ho[absdof(i)];
160  int s3 = idof_lor < 0 ? -1 : 1;
161  int s4 = idof_ho < 0 ? -1 : 1;
162  int s = s1*s2*s3*s4;
163  i = absdof(idof_ho);
164  perm_[absdof(idof_lor)] = s < 0 ? -1-absdof(i) : absdof(i);
165  }
166  }
167  };
168 
169  int offset;
170 
171  if (type == ND)
172  {
173  // x
174  offset = off_x + off_y*p + off_z*p*p1;
175  set_perm(0, offset, p, p*p1);
176  // y
177  offset = ndof_per_dim + off_x + off_y*(p1) + off_z*p1*p;
178  set_perm(dim == 2 ? 2 : 4, offset, 1, p*p1);
179  // z
180  if (dim == 3)
181  {
182  offset = 2*ndof_per_dim + off_x + off_y*p1 + off_z*p1*p1;
183  set_perm(8, offset, 1, p+1);
184  }
185  }
186  else if (type == RT)
187  {
188  // x
189  offset = off_x + off_y*p1 + off_z*p*p1;
190  set_perm(0, offset, 1, 0);
191  // y
192  offset = ndof_per_dim + off_x + off_y*p + off_z*p1*p;
193  set_perm(2, offset, p, 0);
194  // z
195  if (dim == 3)
196  {
197  offset = 2*ndof_per_dim + off_x + off_y*p + off_z*p*p;
198  set_perm(4, offset, p*p, 0);
199  }
200  }
201  }
202 }
203 
205 {
206  FESpaceType type = GetFESpaceType();
207  if (type == H1 || type == L2)
208  {
209  // H1 and L2: no permutation necessary, return identity
211  for (int i=0; i<perm.Size(); ++i) { perm[i] = i; }
212  return;
213  }
214 
215 #ifdef MFEM_USE_MPI
216  ParFiniteElementSpace *pfes_ho
217  = dynamic_cast<ParFiniteElementSpace*>(&fes_ho);
218  ParFiniteElementSpace *pfes_lor = dynamic_cast<ParFiniteElementSpace*>(fes);
219  if (pfes_ho && pfes_lor)
220  {
221  Array<int> l_perm;
223  perm.SetSize(pfes_lor->GetTrueVSize());
224  for (int i=0; i<l_perm.Size(); ++i)
225  {
226  int j = l_perm[i];
227  int s = j < 0 ? -1 : 1;
228  int t_i = pfes_lor->GetLocalTDofNumber(i);
229  int t_j = pfes_ho->GetLocalTDofNumber(absdof(j));
230  // Either t_i and t_j both -1, or both non-negative
231  if ((t_i < 0 && t_j >=0) || (t_j < 0 && t_i >= 0))
232  {
233  MFEM_ABORT("Inconsistent DOF numbering");
234  }
235  if (t_i < 0) { continue; }
236  perm[t_i] = s < 0 ? -1 - t_j : t_j;
237  }
238  }
239  else
240 #endif
241  {
243  }
244 }
245 
247 {
248  if (perm.Size() == 0) { ConstructDofPermutation(); }
249  return perm;
250 }
251 
253 {
254  FESpaceType type = GetFESpaceType();
255  return type == H1 || type == L2;
256 }
257 
259 {
260  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
261  return A;
262 }
263 
265 {
267  AddIntegrators(a_ho, *a, &BilinearForm::GetDBFI,
269  AddIntegrators(a_ho, *a, &BilinearForm::GetFBFI,
271  AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBBFI,
275  AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBFBFI,
279  a->Assemble();
280  a->FormSystemMatrix(ess_dofs, A);
281  ResetIntegrationRules(&BilinearForm::GetDBFI);
282  ResetIntegrationRules(&BilinearForm::GetFBFI);
283  ResetIntegrationRules(&BilinearForm::GetBBFI);
284  ResetIntegrationRules(&BilinearForm::GetBFBFI);
285 }
286 
288 {
289  if (!HasSameDofNumbering())
290  {
291  Array<int> p;
294  }
295  else
296  {
298  }
299 }
300 
301 template <typename FEC>
303 {
304  const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
305  if (fec)
306  {
307  int btype = fec->GetBasisType();
308  if (btype != BasisType::GaussLobatto)
309  {
310  mfem::err << "\nWARNING: Constructing low-order refined "
311  << "discretization with basis type\n"
312  << BasisType::Name(btype) << ". "
313  << "The LOR discretization is only spectrally equivalent\n"
314  << "with Gauss-Lobatto basis.\n" << std::endl;
315  }
316  }
317 }
318 
319 template <typename FEC>
321 {
322  const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
323  if (fec)
324  {
325  int cbtype = fec->GetClosedBasisType();
326  int obtype = fec->GetOpenBasisType();
327  if (cbtype != BasisType::GaussLobatto || obtype != BasisType::IntegratedGLL)
328  {
329  mfem::err << "\nWARNING: Constructing vector low-order refined "
330  << "discretization with basis type \npair ("
331  << BasisType::Name(cbtype) << ", "
332  << BasisType::Name(obtype) << "). "
333  << "The LOR discretization is only spectrally\nequivalent "
334  << "with basis types (Gauss-Lobatto, IntegratedGLL).\n"
335  << std::endl;
336  }
337  }
338 }
339 
341 {
342  CheckScalarBasisType<H1_FECollection>(fes);
343  CheckVectorBasisType<ND_FECollection>(fes);
344  CheckVectorBasisType<RT_FECollection>(fes);
345  // L2 is a bit more complicated, for now don't verify basis type
346 }
347 
349  : irs(0, Quadrature1D::GaussLobatto), fes_ho(fes_ho_)
350 {
351  Mesh &mesh_ = *fes_ho_.GetMesh();
352  int dim = mesh_.Dimension();
353  Array<Geometry::Type> geoms;
354  mesh_.GetGeometries(dim, geoms);
355  if (geoms.Size() == 1 && Geometry::IsTensorProduct(geoms[0]))
356  {
357  ir_el = &irs.Get(geoms[0], 1);
358  ir_face = &irs.Get(Geometry::TensorProductGeometry(dim-1), 1);
359  }
360  else
361  {
362  ir_el = NULL;
363  ir_face = NULL;
364  }
365  a = NULL;
366 }
367 
369 {
370  delete a;
371  delete fes;
372  delete fec;
373  delete mesh;
374 }
375 
377  const Array<int> &ess_tdof_list,
378  int ref_type)
379  : LORDiscretization(*a_ho_.FESpace(), ref_type)
380 {
381  AssembleSystem(a_ho_, ess_tdof_list);
382 }
383 
385  int ref_type) : LORBase(fes_ho)
386 {
387  CheckBasisType(fes_ho);
388 
389  Mesh &mesh_ho = *fes_ho.GetMesh();
390  // For H1, ND and RT spaces, use refinement = element order, for DG spaces,
391  // use refinement = element order + 1 (since LOR is p = 0 in this case).
392  int increment = (GetFESpaceType() == L2) ? 1 : 0;
393  Array<int> refinements(mesh_ho.GetNE());
394  for (int i=0; i<refinements.Size(); ++i)
395  {
396  refinements[i] = fes_ho.GetOrder(i) + increment;
397  }
398  mesh = new Mesh(Mesh::MakeRefined(mesh_ho, refinements, ref_type));
399 
400  fec = fes_ho.FEColl()->Clone(GetLOROrder());
401  fes = new FiniteElementSpace(mesh, fec);
403 
405 }
406 
408  const Array<int> &ess_dofs)
409 {
410  delete a;
411  a = new BilinearForm(&GetFESpace());
412  AssembleSystem_(a_ho, ess_dofs);
413 }
414 
416 {
417  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
418  return *A.As<SparseMatrix>();
419 }
420 
421 #ifdef MFEM_USE_MPI
422 
424  const Array<int> &ess_tdof_list,
425  int ref_type)
426  : ParLORDiscretization(*a_ho_.ParFESpace(), ref_type)
427 {
428  AssembleSystem(a_ho_, ess_tdof_list);
429 }
430 
432  int ref_type) : LORBase(fes_ho)
433 {
434  if (fes_ho.GetMyRank() == 0) { CheckBasisType(fes_ho); }
435  // TODO: support variable-order spaces in parallel
436  MFEM_VERIFY(!fes_ho.IsVariableOrder(),
437  "Cannot construct LOR operators on variable-order spaces");
438 
439  int order = fes_ho.GetMaxElementOrder();
440  if (GetFESpaceType() == L2) { ++order; }
441 
442  ParMesh &mesh_ho = *fes_ho.GetParMesh();
443  ParMesh *pmesh = new ParMesh(ParMesh::MakeRefined(mesh_ho, order, ref_type));
444  mesh = pmesh;
445 
446  fec = fes_ho.FEColl()->Clone(GetLOROrder());
447  ParFiniteElementSpace *pfes = new ParFiniteElementSpace(pmesh, fec);
448  fes = pfes;
450 
452 }
453 
455  const Array<int> &ess_dofs)
456 {
457  delete a;
458  a = new ParBilinearForm(&GetParFESpace());
459  AssembleSystem_(a_ho, ess_dofs);
460 }
461 
463 {
464  MFEM_VERIFY(a != NULL && A.Ptr() != NULL, "No LOR system assembled");
465  return *A.As<HypreParMatrix>();
466 }
467 
469 {
470  return static_cast<ParFiniteElementSpace&>(*fes);
471 }
472 
473 #endif
474 
475 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Integrated GLL indicator functions.
Definition: fe_base.hpp:38
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
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:563
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:376
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system.
Definition: lor.cpp:258
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
Create and assemble a low-order refined version of a ParBilinearForm.
Definition: lor.hpp:156
static Type TensorProductGeometry(int dim)
Definition: geom.hpp:113
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
Definition: lor.cpp:68
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:455
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
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:1366
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:79
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:65
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:3721
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:5888
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
Definition: lor.cpp:246
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:928
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:1062
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition: fe_base.hpp:91
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:100
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition: lor.cpp:462
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1062
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:302
Create and assemble a low-order refined version of a BilinearForm.
Definition: lor.hpp:127
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:289
bool HasSameDofNumbering() const
Definition: lor.cpp:252
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition: lor.cpp:407
ID for class SparseMatrix.
Definition: operator.hpp:258
FiniteElementSpace * fes
Definition: lor.hpp:68
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
Data type sparse matrix.
Definition: sparsemat.hpp:46
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
Mesh * mesh
Definition: lor.hpp:66
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
Definition: lor.hpp:121
void SetupProlongationAndRestriction()
Definition: lor.cpp:287
void CheckBasisType(const FiniteElementSpace &fes)
Definition: lor.cpp:340
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition: fespace.hpp:451
Array< Embedding > embeddings
Fine element positions in their parents.
Definition: ncmesh.hpp:73
Array< int > perm
Definition: lor.hpp:71
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:566
int Dimension() const
Definition: mesh.hpp:999
BilinearForm * a
Definition: lor.hpp:69
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:87
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:88
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:679
int GetOrder(int i) const
Returns the polynomial degree of the i&#39;th finite element.
Definition: fespace.hpp:547
void AssembleSystem(ParBilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition: lor.cpp:454
void ConstructDofPermutation() const
Definition: lor.cpp:204
const CoarseFineTransformations & GetRefinementTransforms()
Definition: mesh.cpp:9849
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:296
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:108
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
Definition: lor.cpp:468
FiniteElementCollection * fec
Definition: lor.hpp:67
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:423
int dim
Definition: ex24.cpp:53
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
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:2783
void ConstructLocalDofPermutation(Array< int > &perm_) const
Definition: lor.cpp:85
Class for parallel bilinear form.
LORBase(FiniteElementSpace &fes_ho_)
Definition: lor.cpp:348
void CheckVectorBasisType(const FiniteElementSpace &fes)
Definition: lor.cpp:320
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
Definition: lor.cpp:415
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:577
ID for class HypreParMatrix.
Definition: operator.hpp:259
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:70
RefCoord s[3]
void AssembleSystem_(BilinearForm &a_ho, const Array< int > &ess_dofs)
Definition: lor.cpp:264
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
OperatorHandle A
Definition: lor.hpp:70
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:132
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.