12#ifndef MFEM_LOR_BATCHED
13#define MFEM_LOR_BATCHED
131static T *GetIntegrator(BilinearForm &
a)
133 Array<BilinearFormIntegrator*> *integs =
a.GetDBFI();
136 for (
auto *i : *integs)
138 if (
auto *ti =
dynamic_cast<T*
>(i))
149template <
typename INTEGRATOR>
152 INTEGRATOR *i = GetIntegrator<INTEGRATOR>(
a);
156 auto *coeff =
const_cast<Coefficient*
>(i->GetCoefficient());
157 if (coeff) { coeff_vector.
Project(*coeff); }
Efficient batched assembly of LOR discretizations on device.
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.
const Vector & GetLORVertexCoordinates()
Return the vertices of the LOR mesh in E-vector 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 for the batched LOR assembly kernels.
FiniteElementSpace & fes_ho
The associated high-order space.
CoefficientVector c2
Coefficient of second integrator.
Vector & X_vert
Mesh coordinate vector.
BatchedLORKernel(FiniteElementSpace &fes_ho_, Vector &X_vert_, Vector &sparse_ij_, Array< int > &sparse_mapping_)
Vector & sparse_ij
Local element sparsity matrix data.
QuadratureSpace qs
Quadrature space for coefficients.
CoefficientVector c1
Coefficient of first integrator.
IntegrationRule ir
Collocated integration rule.
Array< int > & sparse_mapping
Local element sparsity pattern.
Class to represent a coefficient evaluated at quadrature points.
void SetConstant(real_t constant)
Set this vector to the given constant.
void Project(Coefficient &coeff)
Evaluate the given Coefficient at the quadrature points defined by qs.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Class for an integration rule - an Array of IntegrationPoint.
Class used by MFEM to store pointers to host and/or device memory.
int Capacity() const
Return the size of the allocated memory.
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
void Delete()
Delete the owned pointers and reset the Memory object.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Pointer to an Operator of a specified type.
Class representing the storage layout of a QuadratureFunction.
void EnsureCapacity(Memory< T > &mem, int capacity)
Ensure that mem has at least capacity capacity.
CoefficientStorage
Flags that determine what storage optimizations to use in CoefficientVector.
@ COMPRESSED
Enable all above compressions.
void ProjectLORCoefficient(BilinearForm &a, CoefficientVector &coeff_vector)
IntegrationRule GetCollocatedIntRule(FiniteElementSpace &fes)