MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
lor.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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 "lor_batched.hpp"
14#include "../restriction.hpp"
15#include "../pbilinearform.hpp"
17
18namespace mfem
19{
20
21void LORBase::AddIntegrators(BilinearForm &a_from,
22 BilinearForm &a_to,
23 GetIntegratorsFn get_integrators,
24 AddIntegratorFn add_integrator,
25 const IntegrationRule *ir)
26{
27 Array<BilinearFormIntegrator*> *integrators = (a_from.*get_integrators)();
28 for (int i=0; i<integrators->Size(); ++i)
29 {
30 BilinearFormIntegrator *integrator = (*integrators)[i];
31 (a_to.*add_integrator)(integrator);
32 ir_map[integrator] = integrator->GetIntegrationRule();
33 if (ir) { integrator->SetIntegrationRule(*ir); }
34 }
35}
36
37void LORBase::AddIntegratorsAndMarkers(BilinearForm &a_from,
38 BilinearForm &a_to,
39 GetIntegratorsFn get_integrators,
40 GetMarkersFn get_markers,
41 AddIntegratorMarkersFn add_integrator_marker,
42 AddIntegratorFn add_integrator,
43 const IntegrationRule *ir)
44{
45 Array<BilinearFormIntegrator*> *integrators = (a_from.*get_integrators)();
46 Array<Array<int>*> *markers = (a_from.*get_markers)();
47
48 for (int i=0; i<integrators->Size(); ++i)
49 {
50 BilinearFormIntegrator *integrator = (*integrators)[i];
51 if (*markers[i])
52 {
53 (a_to.*add_integrator_marker)(integrator, *(*markers[i]));
54 }
55 else
56 {
57 (a_to.*add_integrator)(integrator);
58 }
59 ir_map[integrator] = integrator->GetIntegrationRule();
60 if (ir) { integrator->SetIntegrationRule(*ir); }
61 }
62}
63
64void LORBase::ResetIntegrationRules(GetIntegratorsFn get_integrators)
65{
66 Array<BilinearFormIntegrator*> *integrators = (a->*get_integrators)();
67 for (int i=0; i<integrators->Size(); ++i)
68 {
69 ((*integrators)[i])->SetIntRule(ir_map[(*integrators)[i]]);
70 }
71}
72
74{
75 const FiniteElementCollection *fec_ho = fes_ho.FEColl();
76 if (dynamic_cast<const H1_FECollection*>(fec_ho)) { return H1; }
77 else if (dynamic_cast<const ND_FECollection*>(fec_ho)) { return ND; }
78 else if (dynamic_cast<const RT_FECollection*>(fec_ho)) { return RT; }
79 else if (dynamic_cast<const L2_FECollection*>(fec_ho)) { return L2; }
80 else { MFEM_ABORT("Bad LOR space type."); }
81 return INVALID;
82}
83
85{
87 return (type == L2 || type == RT) ? 0 : 1;
88}
89
91{
93 MFEM_VERIFY(type != H1 && type != L2, "");
94
95 auto get_dof_map = [](FiniteElementSpace &fes_, int i)
96 {
97 const FiniteElement *fe = fes_.GetFE(i);
98 auto tfe = dynamic_cast<const TensorBasisElement*>(fe);
99 MFEM_ASSERT(tfe != NULL, "");
100 return tfe->GetDofMap();
101 };
102
103 FiniteElementSpace &fes_lor = GetFESpace();
104 Mesh &mesh_lor = *fes_lor.GetMesh();
105 int dim = mesh_lor.Dimension();
106 const CoarseFineTransformations &cf_tr = mesh_lor.GetRefinementTransforms();
107
108 using GeomRef = std::pair<Geometry::Type, int>;
109 std::map<GeomRef, int> point_matrices_offsets;
110 perm_.SetSize(fes_lor.GetVSize());
111
112 Array<int> vdof_ho, vdof_lor;
113 for (int ilor=0; ilor<mesh_lor.GetNE(); ++ilor)
114 {
115 int iho = cf_tr.embeddings[ilor].parent;
116 int p = fes_ho.GetOrder(iho);
117 int lor_index = cf_tr.embeddings[ilor].matrix;
118 // We use the point matrix index to identify the local LOR element index
119 // within the high-order coarse element.
120 //
121 // In variable-order spaces, the point matrices for each order are
122 // concatenated sequentially, so for the given element order, we need to
123 // find the offset that will give us the point matrix index relative to
124 // the current element order only.
125 GeomRef id(mesh_lor.GetElementBaseGeometry(ilor), p);
126 if (point_matrices_offsets.find(id) == point_matrices_offsets.end())
127 {
128 point_matrices_offsets[id] = lor_index;
129 }
130 lor_index -= point_matrices_offsets[id];
131
132 fes_ho.GetElementVDofs(iho, vdof_ho);
133 fes_lor.GetElementVDofs(ilor, vdof_lor);
134
135 if (type == L2)
136 {
137 perm_[vdof_lor[0]] = vdof_ho[lor_index];
138 continue;
139 }
140
141 int p1 = p+1;
142 int ndof_per_dim = (dim == 2) ? p*p1 : type == ND ? p*p1*p1 : p*p*p1;
143
144 const Array<int> &dofmap_ho = get_dof_map(fes_ho, iho);
145 const Array<int> &dofmap_lor = get_dof_map(fes_lor, ilor);
146
147 int off_x = lor_index % p;
148 int off_y = (lor_index / p) % p;
149 int off_z = (lor_index / p) / p;
150
151 auto set_perm = [&](int off_lor, int off_ho, int n1, int n2)
152 {
153 for (int i1=0; i1<2; ++i1)
154 {
155 int m = (dim == 2 || type == RT) ? 1 : 2;
156 for (int i2=0; i2<m; ++i2)
157 {
158 int i;
159 i = dofmap_lor[off_lor + i1 + i2*2];
160 int s1 = i < 0 ? -1 : 1;
161 int idof_lor = vdof_lor[absdof(i)];
162 i = dofmap_ho[off_ho + i1*n1 + i2*n2];
163 int s2 = i < 0 ? -1 : 1;
164 int idof_ho = vdof_ho[absdof(i)];
165 int s3 = idof_lor < 0 ? -1 : 1;
166 int s4 = idof_ho < 0 ? -1 : 1;
167 int s = s1*s2*s3*s4;
168 i = absdof(idof_ho);
169 perm_[absdof(idof_lor)] = s < 0 ? -1-absdof(i) : absdof(i);
170 }
171 }
172 };
173
174 int offset;
175
176 if (type == ND)
177 {
178 // x
179 offset = off_x + off_y*p + off_z*p*p1;
180 set_perm(0, offset, p, p*p1);
181 // y
182 offset = ndof_per_dim + off_x + off_y*(p1) + off_z*p1*p;
183 set_perm(dim == 2 ? 2 : 4, offset, 1, p*p1);
184 // z
185 if (dim == 3)
186 {
187 offset = 2*ndof_per_dim + off_x + off_y*p1 + off_z*p1*p1;
188 set_perm(8, offset, 1, p+1);
189 }
190 }
191 else if (type == RT)
192 {
193 // x
194 offset = off_x + off_y*p1 + off_z*p*p1;
195 set_perm(0, offset, 1, 0);
196 // y
197 offset = ndof_per_dim + off_x + off_y*p + off_z*p1*p;
198 set_perm(2, offset, p, 0);
199 // z
200 if (dim == 3)
201 {
202 offset = 2*ndof_per_dim + off_x + off_y*p + off_z*p*p;
203 set_perm(4, offset, p*p, 0);
204 }
205 }
206 }
207}
208
210{
212 if (type == H1 || type == L2)
213 {
214 // H1 and L2: no permutation necessary, return identity
216 for (int i=0; i<perm.Size(); ++i) { perm[i] = i; }
217 return;
218 }
219
220#ifdef MFEM_USE_MPI
221 ParFiniteElementSpace *pfes_ho
222 = dynamic_cast<ParFiniteElementSpace*>(&fes_ho);
223 ParFiniteElementSpace *pfes_lor
224 = dynamic_cast<ParFiniteElementSpace*>(&GetFESpace());
225 if (pfes_ho && pfes_lor)
226 {
227 Array<int> l_perm;
229 perm.SetSize(pfes_lor->GetTrueVSize());
230 for (int i=0; i<l_perm.Size(); ++i)
231 {
232 int j = l_perm[i];
233 int s = j < 0 ? -1 : 1;
234 int t_i = pfes_lor->GetLocalTDofNumber(i);
235 int t_j = pfes_ho->GetLocalTDofNumber(absdof(j));
236 // Either t_i and t_j both -1, or both non-negative
237 if ((t_i < 0 && t_j >=0) || (t_j < 0 && t_i >= 0))
238 {
239 MFEM_ABORT("Inconsistent DOF numbering");
240 }
241 if (t_i < 0) { continue; }
242 perm[t_i] = s < 0 ? -1 - t_j : t_j;
243 }
244 }
245 else
246#endif
247 {
249 }
250}
251
253{
254 if (perm.Size() == 0) { ConstructDofPermutation(); }
255 return perm;
256}
257
259{
261 return type == H1 || type == L2;
262}
263
265{
266 MFEM_VERIFY(A.Ptr() != NULL, "No LOR system assembled");
267 return A;
268}
269
271{
272 MFEM_VERIFY(A.Ptr() != NULL, "No LOR system assembled");
273 return A;
274}
275
289
290template <typename FEC>
292{
293 const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
294 if (fec)
295 {
296 int btype = fec->GetBasisType();
297 if (btype != BasisType::GaussLobatto)
298 {
299 mfem::err << "\nWARNING: Constructing low-order refined "
300 << "discretization with basis type\n"
301 << BasisType::Name(btype) << ". "
302 << "The LOR discretization is only spectrally equivalent\n"
303 << "with Gauss-Lobatto basis.\n" << std::endl;
304 }
305 }
306}
307
308template <typename FEC>
310{
311 const FEC *fec = dynamic_cast<const FEC*>(fes.FEColl());
312 if (fec)
313 {
314 int cbtype = fec->GetClosedBasisType();
315 int obtype = fec->GetOpenBasisType();
316 if (cbtype != BasisType::GaussLobatto || obtype != BasisType::IntegratedGLL)
317 {
318 mfem::err << "\nWARNING: Constructing vector low-order refined "
319 << "discretization with basis type \npair ("
320 << BasisType::Name(cbtype) << ", "
321 << BasisType::Name(obtype) << "). "
322 << "The LOR discretization is only spectrally\nequivalent "
323 << "with basis types (Gauss-Lobatto, IntegratedGLL).\n"
324 << std::endl;
325 }
326 }
327}
328
330{
334 // L2 is a bit more complicated, for now don't verify basis type
335}
336
337LORBase::LORBase(FiniteElementSpace &fes_ho_, int ref_type_)
338 : irs(0, Quadrature1D::GaussLobatto), ref_type(ref_type_), fes_ho(fes_ho_)
339{
340 Mesh &mesh_ = *fes_ho_.GetMesh();
341 int dim = mesh_.Dimension();
343 mesh_.GetGeometries(dim, geoms);
344 if (geoms.Size() == 1 && Geometry::IsTensorProduct(geoms[0]))
345 {
346 ir_el = &irs.Get(geoms[0], 1);
347 ir_face = &irs.Get(Geometry::TensorProductGeometry(dim-1), 1);
348 }
349 else
350 {
351 ir_el = NULL;
352 ir_face = NULL;
353 }
354 a = NULL;
355}
356
358{
359 // In the case of "batched assembly", the creation of the LOR mesh and
360 // space can be completely omitted (for efficiency). In this case, the
361 // fes object is NULL, and we need to create it when requested.
362 if (fes == NULL) { const_cast<LORBase*>(this)->FormLORSpace(); }
363 return *fes;
364}
365
367{
368 A.Clear();
369 delete a;
371 {
372 // Skip forming the space
373 a = nullptr;
374 if (batched_lor == nullptr)
375 {
377 }
378 batched_lor->Assemble(a_ho, ess_dofs, A);
379 }
380 else
381 {
382 LegacyAssembleSystem(a_ho, ess_dofs);
383 }
384}
385
387 const Array<int> &ess_dofs)
388{
389 // TODO: use AssemblyLevel::FULL here instead of AssemblyLevel::LEGACY.
390 // This is waiting for parallel assembly + BCs with AssemblyLevel::FULL.
391 // In that case, maybe "LegacyAssembleSystem" is not a very clear name.
392
393 // If the space is not formed already, it will be constructed lazily in
394 // GetFESpace.
395 FiniteElementSpace &fes_lor = GetFESpace();
396#ifdef MFEM_USE_MPI
397 if (auto *pfes = dynamic_cast<ParFiniteElementSpace*>(&fes_lor))
398 {
399 a = new ParBilinearForm(pfes);
400 }
401 else
402#endif
403 {
404 a = new BilinearForm(&fes_lor);
405 }
406
408 AddIntegrators(a_ho, *a, &BilinearForm::GetDBFI,
410 AddIntegrators(a_ho, *a, &BilinearForm::GetFBFI,
412 AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBBFI,
416 AddIntegratorsAndMarkers(a_ho, *a, &BilinearForm::GetBFBFI,
420
421 a->Assemble();
422 a->FormSystemMatrix(ess_dofs, A);
423
424 ResetIntegrationRules(&BilinearForm::GetDBFI);
425 ResetIntegrationRules(&BilinearForm::GetFBFI);
426 ResetIntegrationRules(&BilinearForm::GetBBFI);
427 ResetIntegrationRules(&BilinearForm::GetBFBFI);
428}
429
431{
432 delete batched_lor;
433 delete a;
434 delete fes;
435 delete fec;
436 delete mesh;
437}
438
440 const Array<int> &ess_tdof_list,
441 int ref_type_)
442 : LORBase(*a_ho_.FESpace(), ref_type_)
443{
446 AssembleSystem(a_ho_, ess_tdof_list);
447}
448
450 int ref_type_) : LORBase(fes_ho, ref_type_)
451{
454}
455
457{
458 Mesh &mesh_ho = *fes_ho.GetMesh();
459 // For H1, ND and RT spaces, use refinement = element order, for DG spaces,
460 // use refinement = element order + 1 (since LOR is p = 0 in this case).
461 int increment = (GetFESpaceType() == L2) ? 1 : 0;
462 Array<int> refinements(mesh_ho.GetNE());
463 for (int i=0; i<refinements.Size(); ++i)
464 {
465 refinements[i] = fes_ho.GetOrder(i) + increment;
466 }
467 mesh = new Mesh(Mesh::MakeRefined(mesh_ho, refinements, ref_type));
468
470 const int vdim = fes_ho.GetVDim();
471 const Ordering::Type ordering = fes_ho.GetOrdering();
472 fes = new FiniteElementSpace(mesh, fec, vdim, ordering);
474}
475
477{
478 MFEM_VERIFY(A.Ptr() != nullptr, "No LOR system assembled");
479 return *A.As<SparseMatrix>();
480}
481
482#ifdef MFEM_USE_MPI
483
485 const Array<int> &ess_tdof_list,
486 int ref_type_) : LORBase(*a_ho_.ParFESpace(), ref_type_)
487{
488 ParFiniteElementSpace *pfes_ho = a_ho_.ParFESpace();
489 if (pfes_ho->GetMyRank() == 0) { CheckBasisType(fes_ho); }
491 AssembleSystem(a_ho_, ess_tdof_list);
492}
493
495 ParFiniteElementSpace &fes_ho, int ref_type_) : LORBase(fes_ho, ref_type_)
496{
497 if (fes_ho.GetMyRank() == 0) { CheckBasisType(fes_ho); }
499}
500
502{
503 ParFiniteElementSpace &pfes_ho = static_cast<ParFiniteElementSpace&>(fes_ho);
504 // TODO: support variable-order spaces in parallel
505 MFEM_VERIFY(!pfes_ho.IsVariableOrder(),
506 "Cannot construct LOR operators on variable-order spaces");
507
508 int order = pfes_ho.GetMaxElementOrder();
509 if (GetFESpaceType() == L2) { ++order; }
510
511 ParMesh &mesh_ho = *pfes_ho.GetParMesh();
512 ParMesh *pmesh = new ParMesh(ParMesh::MakeRefined(mesh_ho, order, ref_type));
513 mesh = pmesh;
514
515 fec = pfes_ho.FEColl()->Clone(GetLOROrder());
516 const int vdim = fes_ho.GetVDim();
517 const Ordering::Type ordering = fes_ho.GetOrdering();
518 fes = new ParFiniteElementSpace(pmesh, fec, vdim, ordering);
520}
521
523{
524 MFEM_VERIFY(A.Ptr() != nullptr, "No LOR system assembled");
525 return *A.As<HypreParMatrix>();
526}
527
532
533#endif // MFEM_USE_MPI
534
535} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ IntegratedGLL
Integrated GLL indicator functions.
Definition fe_base.hpp:39
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition fe_base.hpp:92
Efficient batched assembly of LOR discretizations on device.
static bool FormIsSupported(BilinearForm &a)
Returns true if the form a supports batched assembly, false otherwise.
void Assemble(BilinearForm &a, const Array< int > ess_dofs, OperatorHandle &A)
Assemble the given form as a matrix and place the result in A.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
Array< Array< int > * > * GetBBFI_Marker()
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the i...
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
void UseExternalIntegrators()
Indicate that integrators are not owned by the BilinearForm.
Array< Array< int > * > * GetBFBFI_Marker()
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the in...
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition fe_coll.cpp:371
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:581
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes, const Array< int > *perm)
Copies the prolongation and restriction matrices from fes.
Definition fespace.cpp:98
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
Definition fespace.hpp:696
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition fespace.hpp:577
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
Abstract class for all finite elements.
Definition fe_base.hpp:239
static bool IsTensorProduct(Type geom)
Definition geom.hpp:108
static Type TensorProductGeometry(int dim)
Definition geom.hpp:113
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition lor.hpp:23
LORBase(FiniteElementSpace &fes_ho_, int ref_type_)
Construct the LORBase object for the given FE space and refinement type.
Definition lor.cpp:337
int GetLOROrder() const
Returns the order of the LOR space. 1 for H1 or ND, 0 for L2 or RT.
Definition lor.cpp:84
void ConstructLocalDofPermutation(Array< int > &perm_) const
Definition lor.cpp:90
bool HasSameDofNumbering() const
Definition lor.cpp:258
void ConstructDofPermutation() const
Definition lor.cpp:209
FiniteElementSpace * fes
Definition lor.hpp:69
void AssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho.
Definition lor.cpp:366
FESpaceType GetFESpaceType() const
Returns the type of finite element space: H1, ND, RT or L2.
Definition lor.cpp:73
FiniteElementCollection * fec
Definition lor.hpp:68
void SetupProlongationAndRestriction()
Definition lor.cpp:276
class BatchedLORAssembly * batched_lor
Definition lor.hpp:71
FiniteElementSpace & fes_ho
Definition lor.hpp:66
OperatorHandle A
Definition lor.hpp:72
BilinearForm * a
Definition lor.hpp:70
virtual void FormLORSpace()=0
Construct the LOR space (overridden for serial and parallel versions).
Array< int > perm
Definition lor.hpp:73
int ref_type
Definition lor.hpp:65
const Array< int > & GetDofPermutation() const
Returns the permutation that maps LOR DOFs to high-order DOFs.
Definition lor.cpp:252
FiniteElementSpace & GetFESpace() const
Returns the low-order refined finite element space.
Definition lor.cpp:357
virtual ~LORBase()
Definition lor.cpp:430
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system (const version).
Definition lor.cpp:270
Mesh * mesh
Definition lor.hpp:67
void LegacyAssembleSystem(BilinearForm &a_ho, const Array< int > &ess_dofs)
Assembles the LOR system corresponding to a_ho using the legacy method.
Definition lor.cpp:386
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:439
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
Definition lor.cpp:456
SparseMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a SparseMatrix.
Definition lor.cpp:476
Mesh data type.
Definition mesh.hpp:56
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:6972
const CoarseFineTransformations & GetRefinementTransforms() const
Definition mesh.cpp:11082
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
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:4290
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1385
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:465
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition handle.hpp:132
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Definition handle.hpp:124
@ MFEM_SPARSEMAT
ID for class SparseMatrix.
Definition operator.hpp:286
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
Type
Ordering methods:
Definition fespace.hpp:34
Class for parallel bilinear form.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
Abstract parallel finite element space.
Definition pfespace.hpp:29
int GetLocalTDofNumber(int ldof) const
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
ParFiniteElementSpace & GetParFESpace() const
Return the LOR ParFiniteElementSpace.
Definition lor.cpp:528
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:484
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition lor.cpp:522
void FormLORSpace() override
Construct the LOR space (overridden for serial and parallel versions).
Definition lor.cpp:501
Class for parallel meshes.
Definition pmesh.hpp:34
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:1375
A class container for 1D quadrature type constants.
Definition intrules.hpp:394
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
Data type sparse matrix.
Definition sparsemat.hpp:51
int dim
Definition ex24.cpp:53
void CheckScalarBasisType(const FiniteElementSpace &fes)
Definition lor.cpp:291
void CheckBasisType(const FiniteElementSpace &fes)
Definition lor.cpp:329
void CheckVectorBasisType(const FiniteElementSpace &fes)
Definition lor.cpp:309
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
real_t p(const Vector &x, real_t t)
RefCoord s[3]
Defines the coarse-fine transformations of all fine elements.
Definition ncmesh.hpp:72
Array< Embedding > embeddings
Fine element positions in their parents.
Definition ncmesh.hpp:74