MFEM
v4.5.2
Finite element discretization library
|
Efficient batched assembly of LOR discretizations on device. More...
#include <lor_batched.hpp>
Public Member Functions | |
BatchedLORAssembly (FiniteElementSpace &fes_ho_) | |
Construct the batched assembly object corresponding to fes_ho_. More... | |
void | Assemble (BilinearForm &a, const Array< int > ess_dofs, OperatorHandle &A) |
Assemble the given form as a matrix and place the result in A. More... | |
const Vector & | GetLORVertexCoordinates () |
Return the vertices of the LOR mesh in E-vector format. More... | |
GPU kernel functions | |
These functions should be considered protected, but they contain MFEM_FORALL kernels, and so they must be public (this is a compiler limitation). | |
int | FillI (SparseMatrix &A) const |
Fill the I array of the sparse matrix A. More... | |
void | FillJAndData (SparseMatrix &A) const |
Fill the J and data arrays of the sparse matrix A. More... | |
void | ParAssemble (BilinearForm &a, const Array< int > &ess_dofs, OperatorHandle &A) |
Assemble the system in parallel and place the result in A. More... | |
Static Public Member Functions | |
static bool | FormIsSupported (BilinearForm &a) |
Returns true if the form a supports batched assembly, false otherwise. More... | |
static void | FormLORVertexCoordinates (FiniteElementSpace &fes_ho, Vector &X_vert) |
Compute the vertices of the LOR mesh and place the result in X_vert. More... | |
Protected Member Functions | |
void | SparseIJToCSR (OperatorHandle &A) const |
After assembling the "sparse IJ" format, convert it to CSR. More... | |
void | AssembleWithoutBC (BilinearForm &a, OperatorHandle &A) |
Assemble the system without eliminating essential DOFs. More... | |
template<typename LOR_KERNEL > | |
void | AssemblyKernel (BilinearForm &a) |
Fill in sparse_ij and sparse_mapping using one of the specialized LOR assembly kernel classes. More... | |
Protected Attributes | |
FiniteElementSpace & | fes_ho |
The high-order space. More... | |
Vector | X_vert |
LOR vertex coordinates. More... | |
Vector | sparse_ij |
The elementwise LOR matrices in a sparse "ij" format. More... | |
Array< int > | sparse_mapping |
The sparsity pattern of the element matrices. More... | |
Efficient batched assembly of LOR discretizations on device.
This class should typically be used by the user-facing classes LORDiscretization and ParLORDiscretization. Only certain bilinear forms are supported, currently:
Whether a form is supported can be checked with the static member function BatchedLORAssembly::FormIsSupported.
Definition at line 33 of file lor_batched.hpp.
mfem::BatchedLORAssembly::BatchedLORAssembly | ( | FiniteElementSpace & | fes_ho_ | ) |
Construct the batched assembly object corresponding to fes_ho_.
Definition at line 488 of file lor_batched.cpp.
void mfem::BatchedLORAssembly::Assemble | ( | BilinearForm & | a, |
const Array< int > | ess_dofs, | ||
OperatorHandle & | A | ||
) |
Assemble the given form as a matrix and place the result in A.
In serial, the result will be a SparseMatrix. In parallel, the result will be a HypreParMatrix.
Definition at line 471 of file lor_batched.cpp.
|
protected |
Assemble the system without eliminating essential DOFs.
Definition at line 423 of file lor_batched.cpp.
|
protected |
Fill in sparse_ij and sparse_mapping using one of the specialized LOR assembly kernel classes.
Definition at line 384 of file lor_batched.cpp.
int mfem::BatchedLORAssembly::FillI | ( | SparseMatrix & | A | ) | const |
Fill the I array of the sparse matrix A.
Definition at line 136 of file lor_batched.cpp.
void mfem::BatchedLORAssembly::FillJAndData | ( | SparseMatrix & | A | ) | const |
Fill the J and data arrays of the sparse matrix A.
Definition at line 229 of file lor_batched.cpp.
|
static |
Returns true if the form a supports batched assembly, false otherwise.
Definition at line 49 of file lor_batched.cpp.
|
static |
Compute the vertices of the LOR mesh and place the result in X_vert.
Definition at line 72 of file lor_batched.cpp.
|
inline |
Return the vertices of the LOR mesh in E-vector format.
Definition at line 74 of file lor_batched.hpp.
void mfem::BatchedLORAssembly::ParAssemble | ( | BilinearForm & | a, |
const Array< int > & | ess_dofs, | ||
OperatorHandle & | A | ||
) |
Assemble the system in parallel and place the result in A.
Definition at line 454 of file lor_batched.cpp.
|
protected |
After assembling the "sparse IJ" format, convert it to CSR.
Definition at line 360 of file lor_batched.cpp.
|
protected |
The high-order space.
Definition at line 36 of file lor_batched.hpp.
|
protected |
The elementwise LOR matrices in a sparse "ij" format.
This is interpreted to have shape (nnz_per_row, ndof_per_el, nel_ho). For index (i, j, k), this represents row j of the kth element matrix. The column index is given by sparse_mapping(i, j).
Definition at line 45 of file lor_batched.hpp.
|
protected |
The sparsity pattern of the element matrices.
This array should be interpreted as having shape (nnz_per_row, ndof_per_el). For local DOF index j, sparse_mapping(i, j) is the column index of the ith nonzero in the jth row. If the index is negative, that entry should be skipped (there is no corresponding nonzero).
Definition at line 54 of file lor_batched.hpp.
|
protected |
LOR vertex coordinates.
Definition at line 38 of file lor_batched.hpp.