13 #include "../../general/forall.hpp" 14 #include "../../fem/pbilinearform.hpp" 24 order(face_fes.GetMaxElementOrder()),
26 edge_fes(face_fes.GetParMesh(), &edge_fec),
36 const int op1 = o + 1;
37 const int nface =
dim*o*o*op1;
42 for (
int c=0; c<
dim; ++c)
44 const int nx = (c == 0) ? op1 : o;
45 const int ny = (c == 1) ? op1 : o;
47 const int c1 = (c+1)%
dim;
48 const int c2 = (c+2)%
dim;
50 const int nx_e1 = (c1 == 0) ? o : op1;
51 const int ny_e1 = (c1 == 1) ? o : op1;
52 const int nx_e2 = (c2 == 0) ? o : op1;
53 const int ny_e2 = (c2 == 1) ? o : op1;
55 for (
int i=0; i<o*o*op1; ++i)
58 const int ix = i1[0] = i2[0] = i%nx;
59 const int iy = i1[1] = i2[1] = (i/nx)%ny;
60 const int iz = i1[2] = i2[2] = i/nx/ny;
62 const int iface = ix + iy*nx + iz*nx*ny + c*o*o*op1;
67 const int ie0 = c1*o*op1*op1 + ix + iy*nx_e1 + iz*nx_e1*ny_e1;
68 const int ie1 = c1*o*op1*op1 + i1[0] + i1[1]*nx_e1 + i1[2]*nx_e1*ny_e1;
70 const int ie2 = c2*o*op1*op1 + ix + iy*nx_e2 + iz*nx_e2*ny_e2;
71 const int ie3 = c2*o*op1*op1 + i2[0] + i2[1]*nx_e2 + i2[2]*nx_e2*ny_e2;
95 const int nnz = 4*nedge_dof;
97 MFEM_FORALL(i, nedge_dof+1, I[i] = 4*i;);
110 MFEM_VERIFY(R_f != NULL && R_e != NULL,
"");
116 const auto offsets_f = R_f->Offsets().Read();
117 const auto indices_f = R_f->Indices().Read();
118 const auto gather_e =
Reshape(R_e->GatherMap().Read(), nedge_per_el, nel_ho);
120 const auto f2e =
Reshape(face2edge.Read(), 4, nface_per_el);
126 auto J = C_local.
WriteJ();
130 MFEM_FORALL(i, nface_dof,
132 const int sj = indices_f[offsets_f[i]];
133 const int j = (sj >= 0) ? sj : -1 - sj;
134 const int sgn_f = (sj >= 0) ? 1 : -1;
135 const int j_loc = j%nface_per_el;
136 const int j_el = j/nface_per_el;
138 for (
int k=0; k<4; ++k)
140 const int je_loc = f2e(k, j_loc);
141 const int sje = gather_e(je_loc, j_el);
142 const int je = (sje >= 0) ? sje : -1 - sje;
143 const int sgn_e = (sje >= 0) ? 1 : -1;
144 const int sgn = (k == 1 || k == 2) ? -1 : 1;
146 V[i*4 + k] = sgn*sgn_f*sgn_e;
163 C_diag.SetOperatorOwner(
false);
198 "The ADS solver is only valid in 3D.");
205 batched_lor.
Assemble(a_ho, ess_tdof_list, A);
206 xyz = lor_ams.StealCoordinateVector();
208 lor_ads.StealCurlMatrix(),
209 lor_ams.StealGradientMatrix(),
210 lor_ams.StealXCoordinate(),
211 lor_ams.StealYCoordinate(),
212 lor_ams.StealZCoordinate());
Create and assemble a low-order refined version of a ParBilinearForm.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
The Auxiliary-space Divergence Solver in hypre.
Efficient batched assembly of LOR discretizations on device.
Pointer to an Operator of a specified type.
int GetNDofs() const
Returns number of degrees of freedom.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Abstract parallel finite element space.
HypreParMatrix * StealCurlMatrix()
Steal ownership of the discrete curl matrix.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Memory< double > & GetMemoryData()
static bool FormIsSupported(BilinearForm &a)
Returns true if the form a supports batched assembly, false otherwise.
Memory< int > & GetMemoryI()
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Memory< int > & GetMemoryJ()
ParFiniteElementSpace edge_fes
The associated Nedelec space.
const OperatorHandle & GetAssembledSystem() const
Returns the assembled LOR system (const version).
void HypreStealOwnership(HypreParMatrix &A_hyp, SparseMatrix &A_diag)
Make A_hyp steal ownership of its diagonal part A_diag.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
ParMesh * GetParMesh() const
int * WriteI(bool on_dev=true)
HypreParMatrix * C
The discrete curl matrix.
ParFiniteElementSpace & face_fes
The RT space.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
static constexpr int dim
Spatial dimension, always 3.
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
int GetNE() const
Returns number of elements in the mesh.
HYPRE_BigInt GlobalTrueVSize() const
void MakeRectangularBlockDiag(MPI_Comm comm, HYPRE_BigInt glob_num_rows, HYPRE_BigInt glob_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, SparseMatrix *diag)
Reset the OperatorHandle to hold a parallel rectangular block-diagonal matrix using the currently set...
const Vector & GetLORVertexCoordinates()
Return the vertices of the LOR mesh in E-vector format.
bool IsIdentityProlongation(const Operator *P)
SolverType & GetSolver()
Access the underlying solver.
Mesh * GetMesh() const
Returns the mesh.
Operator that converts FiniteElementSpace L-vectors to E-vectors.
HYPRE_BigInt GlobalVSize() const
void FormCurlMatrix()
Form the discrete curl matrix (not part of the public API).
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void Assemble(BilinearForm &a, const Array< int > ess_dofs, OperatorHandle &A)
Assemble the given form as a matrix and place the result in A.
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
BatchedLOR_ADS(ParFiniteElementSpace &pfes_ho_, const Vector &X_vert)
Construct the BatchedLOR_AMS object associated with the 3D RT space pfes_ho_.
int height
Dimension of the output / number of rows in the matrix.
int * WriteJ(bool on_dev=true)
HYPRE_BigInt * GetDofOffsets() const
T * HostWrite()
Shortcut for mfem::Write(a.GetMemory(), a.Size(), false).
const int order
Polynomial degree.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
LORSolver(BilinearForm &a_ho, const Array< int > &ess_tdof_list, int ref_type=BasisType::GaussLobatto)
Create a solver of type SolverType, formed using the assembled SparseMatrix of the LOR version of a_h...
Lexicographic ordering for tensor-product FiniteElements.
Operator * Ptr() const
Access the underlying Operator pointer.
ID for class HypreParMatrix.
T * StealPointer(T *&ptr)
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Wrapper for hypre's ParCSR matrix class.
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
HYPRE_BigInt * GetTrueDofOffsets() const
void Form3DFaceToEdge(Array< int > &face2edge)
Form the local elementwise discrete curl matrix.
int width
Dimension of the input / number of columns in the matrix.
double * WriteData(bool on_dev=true)