13 #include "../../fem/quadinterpolator.hpp"
14 #include "../../general/forall.hpp"
16 #include "../pbilinearform.hpp"
26 template <
typename T1,
typename T2>
30 if (integs == NULL) {
return false; }
31 if (integs->
Size() == 1)
34 if (dynamic_cast<T1*>(i0) ||
dynamic_cast<T2*
>(i0)) {
return true; }
36 else if (integs->
Size() == 2)
40 if ((dynamic_cast<T1*>(i0) && dynamic_cast<T2*>(i1)) ||
41 (
dynamic_cast<T2*
>(i0) && dynamic_cast<T1*>(i1)))
57 if (dynamic_cast<const H1_FECollection*>(fec))
59 if (HasIntegrators<DiffusionIntegrator, MassIntegrator>(a)) {
return true; }
61 else if (dynamic_cast<const ND_FECollection*>(fec))
63 if (HasIntegrators<CurlCurlIntegrator, VectorFEMassIntegrator>(a)) {
return true; }
65 else if (dynamic_cast<const RT_FECollection*>(fec))
67 if (HasIntegrators<DivDivIntegrator, VectorFEMassIntegrator>(a)) {
return true; }
80 const int nel_ho = mesh_ho.
GetNE();
82 const int nd1d = order + 1;
83 const int ndof_per_el =
static_cast<int>(pow(nd1d, dim));
92 nodal_restriction->
Mult(*nodal_gf, nodal_evec);
97 X_vert.
SetSize(dim*ndof_per_el*nel_ho);
101 quad_interp->
Values(nodal_evec, X_vert);
108 static MFEM_HOST_DEVICE
int GetMinElt(
const int *my_elts,
const int n_my_elts,
109 const int *nbr_elts,
const int n_nbr_elts)
111 int min_el = INT_MAX;
112 for (
int i = 0; i < n_my_elts; i++)
114 const int e_i = my_elts[i];
115 if (e_i >= min_el) {
continue; }
116 for (
int j = 0; j < n_nbr_elts; j++)
118 if (e_i==nbr_elts[j])
130 static MFEM_HOST_DEVICE
int GetAndIncrementNnzIndex(
const int i_L,
int* I)
138 static constexpr
int Max = 16;
150 MFEM_VERIFY(el_restr !=
nullptr,
"Bad element restriction");
156 const auto el_dof_lex =
Reshape(el_dof_lex_.
Read(), ndof_per_el, nel_ho);
157 const auto dof_glob2loc = dof_glob2loc_.
Read();
158 const auto K = dof_glob2loc_offsets_.
Read();
163 MFEM_FORALL(ii, nvdof + 1, I[ii] = 0;);
164 MFEM_FORALL(i, ndof_per_el*nel_ho,
166 const int ii_el = i%ndof_per_el;
167 const int iel_ho = i/ndof_per_el;
168 const int sii = el_dof_lex(ii_el, iel_ho);
169 const int ii = (sii >= 0) ? sii : -1 -sii;
172 const int i_offset = K[ii];
173 const int i_next_offset = K[ii+1];
174 const int i_ne = i_next_offset - i_offset;
175 for (
int e_i = 0; e_i < i_ne; ++e_i)
177 const int si_E = dof_glob2loc[i_offset+e_i];
178 const int i_E = (si_E >= 0) ? si_E : -1 - si_E;
179 i_elts[e_i] = i_E/ndof_per_el;
181 for (
int j = 0; j < nnz_per_row; ++j)
183 int jj_el = map(j, ii_el);
184 if (jj_el < 0) {
continue; }
186 const int sjj = el_dof_lex(jj_el, iel_ho);
187 const int jj = (sjj >= 0) ? sjj : -1 - sjj;
188 const int j_offset = K[jj];
189 const int j_next_offset = K[jj+1];
190 const int j_ne = j_next_offset - j_offset;
191 if (i_ne == 1 || j_ne == 1)
198 for (
int e_j = 0; e_j < j_ne; ++e_j)
200 const int sj_E = dof_glob2loc[j_offset+e_j];
201 const int j_E = (sj_E >= 0) ? sj_E : -1 - sj_E;
202 const int elt = j_E/ndof_per_el;
205 const int min_e = GetMinElt(i_elts, i_ne, j_elts, j_ne);
217 for (
int i = 0; i < nvdof; i++)
219 const int nnz = h_I[i];
240 MFEM_VERIFY(el_restr !=
nullptr,
"Bad element restriction");
246 const auto el_dof_lex =
Reshape(el_dof_lex_.
Read(), ndof_per_el, nel_ho);
247 const auto dof_glob2loc = dof_glob2loc_.
Read();
248 const auto K = dof_glob2loc_offsets_.
Read();
254 const auto I = I_.Write();
255 const auto J = A.
WriteJ();
260 const auto I2 = A.
ReadI();
261 MFEM_FORALL(i, nvdof + 1, I[i] = I2[i];);
264 static constexpr
int Max = 16;
266 MFEM_FORALL(i, ndof_per_el*nel_ho,
268 const int ii_el = i%ndof_per_el;
269 const int iel_ho = i/ndof_per_el;
271 const int sii = el_dof_lex(ii_el, iel_ho);
272 const int ii = (sii >= 0) ? sii : -1 - sii;
276 const int i_offset = K[ii];
277 const int i_next_offset = K[ii+1];
278 const int i_ne = i_next_offset - i_offset;
279 for (
int e_i = 0; e_i < i_ne; ++e_i)
281 const int si_E = dof_glob2loc[i_offset+e_i];
282 const bool plus = si_E >= 0;
283 const int i_E = plus ? si_E : -1 - si_E;
284 i_elts[e_i] = i_E/ndof_per_el;
285 const int i_Bi = i_E % ndof_per_el;
286 i_B[e_i] = plus ? i_Bi : -1 - i_Bi;
288 for (
int j=0; j<nnz_per_row; ++j)
290 int jj_el = map(j, ii_el);
291 if (jj_el < 0) {
continue; }
293 const int sjj = el_dof_lex(jj_el, iel_ho);
294 const int jj = (sjj >= 0) ? sjj : -1 - sjj;
295 const int sgn = ((sjj >=0 && sii >= 0) || (sjj < 0 && sii <0)) ? 1 : -1;
296 const int j_offset = K[jj];
297 const int j_next_offset = K[jj+1];
298 const int j_ne = j_next_offset - j_offset;
299 if (i_ne == 1 || j_ne == 1)
301 const int nnz = GetAndIncrementNnzIndex(ii, I);
303 AV[nnz] = sgn*V(j, ii_el, iel_ho);
309 for (
int e_j = 0; e_j < j_ne; ++e_j)
311 const int sj_E = dof_glob2loc[j_offset+e_j];
312 const bool plus = sj_E >= 0;
313 const int j_E = plus ? sj_E : -1 - sj_E;
314 j_elts[e_j] = j_E/ndof_per_el;
315 const int j_Bj = j_E % ndof_per_el;
316 j_B[e_j] = plus ? j_Bj : -1 - j_Bj;
318 const int min_e = GetMinElt(i_elts, i_ne, j_elts, j_ne);
322 for (
int k = 0; k < i_ne; k++)
324 const int iel_ho_2 = i_elts[k];
325 const int sii_el_2 = i_B[k];
326 const int ii_el_2 = (sii_el_2 >= 0) ? sii_el_2 : -1 -sii_el_2;
327 for (
int l = 0; l < j_ne; l++)
329 const int jel_ho_2 = j_elts[l];
330 if (iel_ho_2 == jel_ho_2)
332 const int sjj_el_2 = j_B[l];
333 const int jj_el_2 = (sjj_el_2 >= 0) ? sjj_el_2 : -1 -sjj_el_2;
334 const int sgn_2 = ((sjj_el_2 >=0 && sii_el_2 >= 0)
335 || (sjj_el_2 < 0 && sii_el_2 <0)) ? 1 : -1;
338 for (
int m = 0; m < nnz_per_row; ++m)
340 if (map(m, ii_el_2) == jj_el_2)
346 MFEM_ASSERT_KERNEL(j >= 0,
"Can't find nonzero");
347 val += sgn_2*V(j2, ii_el_2, iel_ho_2);
351 const int nnz = GetAndIncrementNnzIndex(ii, I);
376 int nnz =
FillI(*A_mat);
383 template <
typename LOR_KERNEL>
395 case 1: kernel.template Assemble2D<1>();
break;
396 case 2: kernel.template Assemble2D<2>();
break;
397 case 3: kernel.template Assemble2D<3>();
break;
398 case 4: kernel.template Assemble2D<4>();
break;
399 case 5: kernel.template Assemble2D<5>();
break;
400 case 6: kernel.template Assemble2D<6>();
break;
401 case 7: kernel.template Assemble2D<7>();
break;
402 case 8: kernel.template Assemble2D<8>();
break;
403 default: MFEM_ABORT(
"No kernel order " << order <<
"!");
410 case 1: kernel.template Assemble3D<1>();
break;
411 case 2: kernel.template Assemble3D<2>();
break;
412 case 3: kernel.template Assemble3D<3>();
break;
413 case 4: kernel.template Assemble3D<4>();
break;
414 case 5: kernel.template Assemble3D<5>();
break;
415 case 6: kernel.template Assemble3D<6>();
break;
416 case 7: kernel.template Assemble3D<7>();
break;
417 case 8: kernel.template Assemble3D<8>();
break;
418 default: MFEM_ABORT(
"No kernel order " << order <<
"!");
428 if (dynamic_cast<const H1_FECollection*>(fec))
430 if (HasIntegrators<DiffusionIntegrator, MassIntegrator>(a))
432 AssemblyKernel<BatchedLOR_H1>(
a);
435 else if (dynamic_cast<const ND_FECollection*>(fec))
437 if (HasIntegrators<CurlCurlIntegrator, VectorFEMassIntegrator>(a))
439 AssemblyKernel<BatchedLOR_ND>(
a);
442 else if (dynamic_cast<const RT_FECollection*>(fec))
444 if (HasIntegrators<DivDivIntegrator, VectorFEMassIntegrator>(a))
446 AssemblyKernel<BatchedLOR_RT>(
a);
467 Operator::DiagonalPolicy::DIAG_ONE);
475 if (dynamic_cast<ParFiniteElementSpace*>(&
fes_ho))
485 Operator::DiagonalPolicy::DIAG_KEEP);
499 return irs.
Get(geom, 2*nd1d - 3);
int Size() const
Return the logical size of the array.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
void ParAssemble(BilinearForm &a, const Array< int > &ess_dofs, OperatorHandle &A)
Assemble the system in parallel and place the result in A.
FiniteElementSpace & fes_ho
The high-order space.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
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. B.
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetSize(int s)
Resize the vector to size s.
const Array< int > & Indices() const
Pointer to an Operator of a specified type.
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Container class for integration rules.
int GetNE() const
Returns number of elements.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Memory< double > & GetMemoryData()
static bool FormIsSupported(BilinearForm &a)
Returns true if the form a supports batched assembly, false otherwise.
Memory< int > & GetMemoryI()
static void FormLORVertexCoordinates(FiniteElementSpace &fes_ho, Vector &X_vert)
Compute the vertices of the LOR mesh and place the result in X_vert.
Memory< int > & GetMemoryJ()
bool HasIntegrators(BilinearForm &a)
int GetNE() const
Returns number of elements in the mesh.
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Mesh * GetMesh() const
Returns the mesh.
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Vector sparse_ij
The elementwise LOR matrices in a sparse "ij" format.
Vector X_vert
LOR vertex coordinates.
Array< int > sparse_mapping
The sparsity pattern of the element matrices.
int GetMaxElementOrder() const
Return the maximum polynomial order.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
int * WriteI(bool on_dev=true)
FiniteElementSpace * FESpace()
const Array< int > & Offsets() const
int FillI(SparseMatrix &A) const
Fill the I array of the sparse matrix A.
void FillJAndData(SparseMatrix &A) const
Fill the J and data arrays of the sparse matrix A.
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
void Assemble(BilinearForm &a, const Array< int > ess_dofs, OperatorHandle &A)
Assemble the given form as a matrix and place the result in A.
const Array< int > & GatherMap() const
Geometry::Type GetElementGeometry(int i) const
const int * ReadI(bool on_dev=true) const
int * WriteJ(bool on_dev=true)
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Lexicographic ordering for tensor-product FiniteElements.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
void AssembleWithoutBC(BilinearForm &a, OperatorHandle &A)
Assemble the system without eliminating essential DOFs.
void SparseIJToCSR(OperatorHandle &A) const
After assembling the "sparse IJ" format, convert it to CSR.
const FiniteElementCollection * FEColl() const
void GetNodes(Vector &node_coord) const
BatchedLORAssembly(FiniteElementSpace &fes_ho_)
Construct the batched assembly object corresponding to fes_ho_.
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.
Wrapper for hypre's ParCSR matrix class.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
IntegrationRule GetCollocatedIntRule(FiniteElementSpace &fes)
void EliminateBC(const Array< int > &ess_dofs, DiagonalPolicy diag_policy)
Eliminate essential (Dirichlet) boundary conditions.
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
double * WriteData(bool on_dev=true)
void AssemblyKernel(BilinearForm &a)
Fill in sparse_ij and sparse_mapping using one of the specialized LOR assembly kernel classes...