MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
lor_batched.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_batched.hpp"
15#include <climits>
16#include "../pbilinearform.hpp"
17
18// Specializations
19#include "lor_h1.hpp"
20#include "lor_nd.hpp"
21#include "lor_rt.hpp"
22
23namespace mfem
24{
25
26template <typename T1, typename T2>
28{
29 Array<BilinearFormIntegrator*> *integs = a.GetDBFI();
30 if (integs == NULL) { return false; }
31 if (integs->Size() == 1)
32 {
33 BilinearFormIntegrator *i0 = (*integs)[0];
34 if (dynamic_cast<T1*>(i0) || dynamic_cast<T2*>(i0)) { return true; }
35 }
36 else if (integs->Size() == 2)
37 {
38 BilinearFormIntegrator *i0 = (*integs)[0];
39 BilinearFormIntegrator *i1 = (*integs)[1];
40 if ((dynamic_cast<T1*>(i0) && dynamic_cast<T2*>(i1)) ||
41 (dynamic_cast<T2*>(i0) && dynamic_cast<T1*>(i1)))
42 {
43 return true;
44 }
45 }
46 return false;
47}
48
50{
51 const FiniteElementCollection *fec = a.FESpace()->FEColl();
52 // TODO: check for maximum supported orders
53
54 // Batched LOR requires all tensor elements
55 if (!UsesTensorBasis(*a.FESpace())) { return false; }
56
57 if (dynamic_cast<const H1_FECollection*>(fec))
58 {
60 }
61 else if (dynamic_cast<const ND_FECollection*>(fec))
62 {
64 }
65 else if (dynamic_cast<const RT_FECollection*>(fec))
66 {
68 }
69 return false;
70}
71
73 Vector &X_vert)
74{
75 Mesh &mesh_ho = *fes_ho.GetMesh();
76 mesh_ho.EnsureNodes();
77
78 // Get nodal points at the LOR vertices
79 const int dim = mesh_ho.Dimension();
80 const int sdim = mesh_ho.SpaceDimension();
81 const int nel_ho = mesh_ho.GetNE();
82 const int order = fes_ho.GetMaxElementOrder();
83 const int nd1d = order + 1;
84 const int ndof_per_el = static_cast<int>(pow(nd1d, dim));
85
86 const GridFunction *nodal_gf = mesh_ho.GetNodes();
87 const FiniteElementSpace *nodal_fes = nodal_gf->FESpace();
88 const Operator *nodal_restriction =
90
91 // Map from nodal L-vector to E-vector
92 Vector nodal_evec(nodal_restriction->Height());
93 nodal_restriction->Mult(*nodal_gf, nodal_evec);
94
96
97 // Map from nodal E-vector to Q-vector at the LOR vertex points
98 X_vert.SetSize(sdim*ndof_per_el*nel_ho);
99 const QuadratureInterpolator *quad_interp =
100 nodal_fes->GetQuadratureInterpolator(ir);
102 quad_interp->Values(nodal_evec, X_vert);
103}
104
105// The following two functions (GetMinElt and GetAndIncrementNnzIndex) are
106// copied from restriction.cpp. Should they be factored out?
107
108// Return the minimal value found in both my_elts and nbr_elts
109static MFEM_HOST_DEVICE int GetMinElt(const int *my_elts, const int n_my_elts,
110 const int *nbr_elts, const int n_nbr_elts)
111{
112 int min_el = INT_MAX;
113 for (int i = 0; i < n_my_elts; i++)
114 {
115 const int e_i = my_elts[i];
116 if (e_i >= min_el) { continue; }
117 for (int j = 0; j < n_nbr_elts; j++)
118 {
119 if (e_i==nbr_elts[j])
120 {
121 min_el = e_i; // we already know e_i < min_el
122 break;
123 }
124 }
125 }
126 return min_el;
127}
128
129// Returns the index where a non-zero entry should be added and increment the
130// number of non-zeros for the row i_L.
131static MFEM_HOST_DEVICE int GetAndIncrementNnzIndex(const int i_L, int* I)
132{
133 int ind = AtomicAdd(I[i_L],1);
134 return ind;
135}
136
138{
139 static constexpr int Max = 16;
140
141 const int nvdof = fes_ho.GetVSize();
142
143 const int ndof_per_el = fes_ho.GetFE(0)->GetDof();
144 const int nel_ho = fes_ho.GetNE();
145 const int nnz_per_row = sparse_mapping.Size()/ndof_per_el;
146
148 const Operator *op = fes_ho.GetElementRestriction(ordering);
149 const ElementRestriction *el_restr =
150 dynamic_cast<const ElementRestriction*>(op);
151 MFEM_VERIFY(el_restr != nullptr, "Bad element restriction");
152
153 const Array<int> &el_dof_lex_ = el_restr->GatherMap();
154 const Array<int> &dof_glob2loc_ = el_restr->Indices();
155 const Array<int> &dof_glob2loc_offsets_ = el_restr->Offsets();
156
157 const auto el_dof_lex = Reshape(el_dof_lex_.Read(), ndof_per_el, nel_ho);
158 const auto dof_glob2loc = dof_glob2loc_.Read();
159 const auto K = dof_glob2loc_offsets_.Read();
160 const auto map = Reshape(sparse_mapping.Read(), nnz_per_row, ndof_per_el);
161
162 auto I = A.WriteI();
163
164 mfem::forall(nvdof + 1, [=] MFEM_HOST_DEVICE (int ii) { I[ii] = 0; });
165 mfem::forall(ndof_per_el*nel_ho, [=] MFEM_HOST_DEVICE (int i)
166 {
167 const int ii_el = i%ndof_per_el;
168 const int iel_ho = i/ndof_per_el;
169 const int sii = el_dof_lex(ii_el, iel_ho);
170 const int ii = (sii >= 0) ? sii : -1 -sii;
171 // Get number and list of elements containing this DOF
172 int i_elts[Max];
173 const int i_offset = K[ii];
174 const int i_next_offset = K[ii+1];
175 const int i_ne = i_next_offset - i_offset;
176 for (int e_i = 0; e_i < i_ne; ++e_i)
177 {
178 const int si_E = dof_glob2loc[i_offset+e_i]; // signed
179 const int i_E = (si_E >= 0) ? si_E : -1 - si_E;
180 i_elts[e_i] = i_E/ndof_per_el;
181 }
182 for (int j = 0; j < nnz_per_row; ++j)
183 {
184 int jj_el = map(j, ii_el);
185 if (jj_el < 0) { continue; }
186 // LDOF index of column
187 const int sjj = el_dof_lex(jj_el, iel_ho); // signed
188 const int jj = (sjj >= 0) ? sjj : -1 - sjj;
189 const int j_offset = K[jj];
190 const int j_next_offset = K[jj+1];
191 const int j_ne = j_next_offset - j_offset;
192 if (i_ne == 1 || j_ne == 1) // no assembly required
193 {
194 AtomicAdd(I[ii], 1);
195 }
196 else // assembly required
197 {
198 int j_elts[Max];
199 for (int e_j = 0; e_j < j_ne; ++e_j)
200 {
201 const int sj_E = dof_glob2loc[j_offset+e_j]; // signed
202 const int j_E = (sj_E >= 0) ? sj_E : -1 - sj_E;
203 const int elt = j_E/ndof_per_el;
204 j_elts[e_j] = elt;
205 }
206 const int min_e = GetMinElt(i_elts, i_ne, j_elts, j_ne);
207 if (iel_ho == min_e) // add the nnz only once
208 {
209 AtomicAdd(I[ii], 1);
210 }
211 }
212 }
213 });
214 // TODO: on device, this is a scan operation
215 // We need to sum the entries of I, we do it on CPU as it is very sequential.
216 auto h_I = A.HostReadWriteI();
217 int sum = 0;
218 for (int i = 0; i < nvdof; i++)
219 {
220 const int nnz = h_I[i];
221 h_I[i] = sum;
222 sum+=nnz;
223 }
224 h_I[nvdof] = sum;
225
226 // Return the number of nnz
227 return h_I[nvdof];
228}
229
231{
232 const int nvdof = fes_ho.GetVSize();
233 const int ndof_per_el = fes_ho.GetFE(0)->GetDof();
234 const int nel_ho = fes_ho.GetNE();
235 const int nnz_per_row = sparse_mapping.Size()/ndof_per_el;
236
238 const Operator *op = fes_ho.GetElementRestriction(ordering);
239 const ElementRestriction *el_restr =
240 dynamic_cast<const ElementRestriction*>(op);
241 MFEM_VERIFY(el_restr != nullptr, "Bad element restriction");
242
243 const Array<int> &el_dof_lex_ = el_restr->GatherMap();
244 const Array<int> &dof_glob2loc_ = el_restr->Indices();
245 const Array<int> &dof_glob2loc_offsets_ = el_restr->Offsets();
246
247 const auto el_dof_lex = Reshape(el_dof_lex_.Read(), ndof_per_el, nel_ho);
248 const auto dof_glob2loc = dof_glob2loc_.Read();
249 const auto K = dof_glob2loc_offsets_.Read();
250
251 const auto V = Reshape(sparse_ij.Read(), nnz_per_row, ndof_per_el, nel_ho);
252 const auto map = Reshape(sparse_mapping.Read(), nnz_per_row, ndof_per_el);
253
254 Array<int> I_(nvdof + 1);
255 const auto I = I_.Write();
256 const auto J = A.WriteJ();
257 auto AV = A.WriteData();
258
259 // Copy A.I into I, use it as a temporary buffer
260 {
261 const auto I2 = A.ReadI();
262 mfem::forall(nvdof + 1, [=] MFEM_HOST_DEVICE (int i) { I[i] = I2[i]; });
263 }
264
265 static constexpr int Max = 16;
266
267 mfem::forall(ndof_per_el*nel_ho, [=] MFEM_HOST_DEVICE (int i)
268 {
269 const int ii_el = i%ndof_per_el;
270 const int iel_ho = i/ndof_per_el;
271 // LDOF index of current row
272 const int sii = el_dof_lex(ii_el, iel_ho); // signed
273 const int ii = (sii >= 0) ? sii : -1 - sii;
274 // Get number and list of elements containing this DOF
275 int i_elts[Max];
276 int i_B[Max];
277 const int i_offset = K[ii];
278 const int i_next_offset = K[ii+1];
279 const int i_ne = i_next_offset - i_offset;
280 for (int e_i = 0; e_i < i_ne; ++e_i)
281 {
282 const int si_E = dof_glob2loc[i_offset+e_i]; // signed
283 const bool plus = si_E >= 0;
284 const int i_E = plus ? si_E : -1 - si_E;
285 i_elts[e_i] = i_E/ndof_per_el;
286 const int i_Bi = i_E % ndof_per_el;
287 i_B[e_i] = plus ? i_Bi : -1 - i_Bi; // encode with sign
288 }
289 for (int j=0; j<nnz_per_row; ++j)
290 {
291 int jj_el = map(j, ii_el);
292 if (jj_el < 0) { continue; }
293 // LDOF index of column
294 const int sjj = el_dof_lex(jj_el, iel_ho); // signed
295 const int jj = (sjj >= 0) ? sjj : -1 - sjj;
296 const int sgn = ((sjj >=0 && sii >= 0) || (sjj < 0 && sii <0)) ? 1 : -1;
297 const int j_offset = K[jj];
298 const int j_next_offset = K[jj+1];
299 const int j_ne = j_next_offset - j_offset;
300 if (i_ne == 1 || j_ne == 1) // no assembly required
301 {
302 const int nnz = GetAndIncrementNnzIndex(ii, I);
303 J[nnz] = jj;
304 AV[nnz] = sgn*V(j, ii_el, iel_ho);
305 }
306 else // assembly required
307 {
308 int j_elts[Max];
309 int j_B[Max];
310 for (int e_j = 0; e_j < j_ne; ++e_j)
311 {
312 const int sj_E = dof_glob2loc[j_offset+e_j]; // signed
313 const bool plus = sj_E >= 0;
314 const int j_E = plus ? sj_E : -1 - sj_E;
315 j_elts[e_j] = j_E/ndof_per_el;
316 const int j_Bj = j_E % ndof_per_el;
317 j_B[e_j] = plus ? j_Bj : -1 - j_Bj; // encode with sign
318 }
319 const int min_e = GetMinElt(i_elts, i_ne, j_elts, j_ne);
320 if (iel_ho == min_e) // add the nnz only once
321 {
322 real_t val = 0.0;
323 for (int k = 0; k < i_ne; k++)
324 {
325 const int iel_ho_2 = i_elts[k];
326 const int sii_el_2 = i_B[k]; // signed
327 const int ii_el_2 = (sii_el_2 >= 0) ? sii_el_2 : -1 -sii_el_2;
328 for (int l = 0; l < j_ne; l++)
329 {
330 const int jel_ho_2 = j_elts[l];
331 if (iel_ho_2 == jel_ho_2)
332 {
333 const int sjj_el_2 = j_B[l]; // signed
334 const int jj_el_2 = (sjj_el_2 >= 0) ? sjj_el_2 : -1 -sjj_el_2;
335 const int sgn_2 = ((sjj_el_2 >=0 && sii_el_2 >= 0)
336 || (sjj_el_2 < 0 && sii_el_2 <0)) ? 1 : -1;
337 int j2 = -1;
338 // find nonzero in matrix of other element
339 for (int m = 0; m < nnz_per_row; ++m)
340 {
341 if (map(m, ii_el_2) == jj_el_2)
342 {
343 j2 = m;
344 break;
345 }
346 }
347 MFEM_ASSERT_KERNEL(j >= 0, "Can't find nonzero");
348 val += sgn_2*V(j2, ii_el_2, iel_ho_2);
349 }
350 }
351 }
352 const int nnz = GetAndIncrementNnzIndex(ii, I);
353 J[nnz] = jj;
354 AV[nnz] = val;
355 }
356 }
357 }
358 });
359}
360
362{
363 const int nvdof = fes_ho.GetVSize();
364
365 // If A contains an existing SparseMatrix, reuse it (and try to reuse its
366 // I, J, A arrays if they are big enough)
367 SparseMatrix *A_mat = A.Is<SparseMatrix>();
368 if (!A_mat)
369 {
370 A_mat = new SparseMatrix;
371 A.Reset(A_mat);
372 }
373
374 A_mat->OverrideSize(nvdof, nvdof);
375
376 A_mat->GetMemoryI().New(nvdof+1, Device::GetDeviceMemoryType());
377 int nnz = FillI(*A_mat);
378
381 FillJAndData(*A_mat);
382}
383
384template <int ORDER, int SDIM, typename LOR_KERNEL>
385static void Assemble_(LOR_KERNEL &kernel, int dim)
386{
387 if (dim == 2) { kernel.template Assemble2D<ORDER,SDIM>(); }
388 else if (dim == 3) { kernel.template Assemble3D<ORDER>(); }
389 else { MFEM_ABORT("Unsupported dimension"); }
390}
391
392template <int ORDER, typename LOR_KERNEL>
393static void Assemble_(LOR_KERNEL &kernel, int dim, int sdim)
394{
395 if (sdim == 2) { Assemble_<ORDER,2>(kernel, dim); }
396 else if (sdim == 3) { Assemble_<ORDER,3>(kernel, dim); }
397 else { MFEM_ABORT("Unsupported space dimension."); }
398}
399
400template <typename LOR_KERNEL>
401static void Assemble_(LOR_KERNEL &kernel, int dim, int sdim, int order)
402{
403 switch (order)
404 {
405 case 1: Assemble_<1>(kernel, dim, sdim); break;
406 case 2: Assemble_<2>(kernel, dim, sdim); break;
407 case 3: Assemble_<3>(kernel, dim, sdim); break;
408 case 4: Assemble_<4>(kernel, dim, sdim); break;
409 case 5: Assemble_<5>(kernel, dim, sdim); break;
410 case 6: Assemble_<6>(kernel, dim, sdim); break;
411 case 7: Assemble_<7>(kernel, dim, sdim); break;
412 case 8: Assemble_<8>(kernel, dim, sdim); break;
413 default: MFEM_ABORT("No kernel order " << order << "!");
414 }
415}
416
417template <typename LOR_KERNEL>
419{
420 LOR_KERNEL kernel(a, fes_ho, X_vert, sparse_ij, sparse_mapping);
421
422 const int dim = fes_ho.GetMesh()->Dimension();
423 const int sdim = fes_ho.GetMesh()->SpaceDimension();
424 const int order = fes_ho.GetMaxElementOrder();
425
426 Assemble_(kernel, dim, sdim, order);
427}
428
430{
431 // Assemble the matrix, depending on what the form is.
432 // This fills in the arrays sparse_ij and sparse_mapping.
434 if (dynamic_cast<const H1_FECollection*>(fec))
435 {
437 {
439 }
440 }
441 else if (dynamic_cast<const ND_FECollection*>(fec))
442 {
444 {
446 }
447 }
448 else if (dynamic_cast<const RT_FECollection*>(fec))
449 {
451 {
453 }
454 }
455
456 return SparseIJToCSR(A);
457}
458
459#ifdef MFEM_USE_MPI
461 BilinearForm &a, const Array<int> &ess_dofs, OperatorHandle &A)
462{
463 // Assemble the system matrix local to this partition
464 OperatorHandle A_local;
465 AssembleWithoutBC(a, A_local);
466
467 ParBilinearForm *pa =
468 dynamic_cast<ParBilinearForm*>(&a);
469
470 pa->ParallelRAP(*A_local.As<SparseMatrix>(), A, true);
471
472 A.As<HypreParMatrix>()->EliminateBC(ess_dofs,
474}
475#endif
476
478 BilinearForm &a, const Array<int> ess_dofs, OperatorHandle &A)
479{
480#ifdef MFEM_USE_MPI
481 if (dynamic_cast<ParFiniteElementSpace*>(&fes_ho))
482 {
483 return ParAssemble(a, ess_dofs, A);
484 }
485#endif
486
488 SparseMatrix *A_mat = A.As<SparseMatrix>();
489
490 A_mat->EliminateBC(ess_dofs,
492}
493
499
501{
503 const Geometry::Type geom = fes.GetMesh()->GetElementGeometry(0);
504 const int nd1d = fes.GetMaxElementOrder() + 1;
505 return irs.Get(geom, 2*nd1d - 3);
506}
507
508} // namespace mfem
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition backends.hpp:92
int Size() const
Return the logical size of the array.
Definition array.hpp:144
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:325
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
void FillJAndData(SparseMatrix &A) const
Fill the J and data arrays of the sparse matrix A.
void SparseIJToCSR(OperatorHandle &A) const
After assembling the "sparse IJ" format, convert it to CSR.
void AssemblyKernel(BilinearForm &a)
Fill in sparse_ij and sparse_mapping using one of the specialized LOR assembly kernel classes.
Array< int > sparse_mapping
The sparsity pattern of the element matrices.
void AssembleWithoutBC(BilinearForm &a, OperatorHandle &A)
Assemble the system without eliminating essential DOFs.
Vector X_vert
LOR vertex coordinates.
Vector sparse_ij
The elementwise LOR matrices in a sparse "ij" format.
static void FormLORVertexCoordinates(FiniteElementSpace &fes_ho, Vector &X_vert)
Compute the vertices of the LOR mesh and place the result in X_vert.
void ParAssemble(BilinearForm &a, const Array< int > &ess_dofs, OperatorHandle &A)
Assemble the system in parallel and place the result in A.
static bool FormIsSupported(BilinearForm &a)
Returns true if the form a supports batched assembly, false otherwise.
FiniteElementSpace & fes_ho
The high-order space.
BatchedLORAssembly(FiniteElementSpace &fes_ho_)
Construct the batched assembly object corresponding to fes_ho_.
void Assemble(BilinearForm &a, const Array< int > ess_dofs, OperatorHandle &A)
Assemble the given form as a matrix and place the result in A.
int FillI(SparseMatrix &A) const
Fill the I array of the sparse matrix A.
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:274
Operator that converts FiniteElementSpace L-vectors to E-vectors.
const Array< int > & GatherMap() const
const Array< int > & Offsets() const
const Array< int > & Indices() const
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition fespace.cpp:1404
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
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition fespace.cpp:1336
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 GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Container class for integration rules.
Definition intrules.hpp:416
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Mesh data type.
Definition mesh.hpp:56
Geometry::Type GetElementGeometry(int i) const
Definition mesh.hpp:1371
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
Definition mesh.cpp:6159
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
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
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 Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
Definition handle.hpp:108
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
@ DIAG_ONE
Set the diagonal value to one.
Definition operator.hpp:50
@ DIAG_KEEP
Keep the diagonal value.
Definition operator.hpp:51
Class for parallel bilinear form.
void ParallelRAP(SparseMatrix &loc_A, OperatorHandle &A, bool steal_loc_A=false)
Compute parallel RAP operator and store it in A as a HypreParMatrix.
Abstract parallel finite element space.
Definition pfespace.hpp:29
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
void SetOutputLayout(QVectorLayout layout) const
Set the desired output Q-vector layout. The default value is QVectorLayout::byNODES.
void Values(const Vector &e_vec, Vector &q_val) const
Interpolate the values of the E-vector e_vec at quadrature points.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
Data type sparse matrix.
Definition sparsemat.hpp:51
void EliminateBC(const Array< int > &ess_dofs, DiagonalPolicy diag_policy)
Eliminate essential (Dirichlet) boundary conditions.
int * WriteJ(bool on_dev=true)
Memory< int > & GetMemoryI()
const int * ReadI(bool on_dev=true) const
int * WriteI(bool on_dev=true)
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Memory< int > & GetMemoryJ()
Memory< real_t > & GetMemoryData()
real_t * WriteData(bool on_dev=true)
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
int dim
Definition ex24.cpp:53
real_t a
Definition lissajous.cpp:41
bool HasIntegrators(BilinearForm &a)
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
float real_t
Definition config.hpp:43
void EliminateBC(const HypreParMatrix &A, const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s....
Definition hypre.cpp:3379
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
IntegrationRule GetCollocatedIntRule(FiniteElementSpace &fes)
void forall(int N, lambda &&body)
Definition forall.hpp:754