13 #include "../../fem/qinterp/dispatch.hpp" 14 #include "../../general/forall.hpp" 15 #include "../../linalg/dtensor.hpp" 27 const int n2 = A.
Size();
28 const int n = sqrt(n2);
33 lu.GetInverseMatrix(n, Ainv.
GetData());
44 for (
int i = 0; i < n; ++i)
46 const double h = gll_pts[i+1] - gll_pts[i];
48 for (
int iq = 0; iq < ir.
Size(); ++iq)
51 const double x = gll_pts[i] + h*ip.
x;
52 const double w = h*ip.
weight;
54 for (
int j = 0; j < n; ++j)
64 const int n = sqrt(B.
Size());
66 for (
int i=0; i<n; ++i)
for (
int j=0; j<n; ++j) { Bt[i+j*n] = B[j+i*n]; }
74 MFEM_VERIFY(fec1,
"Must be L2 finite element space");
76 const int btype = fec1->GetBasisType();
80 if (no_op) {
return; }
86 MFEM_VERIFY(tbe !=
nullptr,
"Must be a tensor element.");
90 const int pp1 =
p + 1;
107 if (no_op) { y = x;
return; }
108 using namespace internal::quadrature_interpolator;
110 TensorValues<QVectorLayout::byVDIM>(ne, 1, dof2quad, x, y);
115 if (no_op) { y = x;
return; }
116 using namespace internal::quadrature_interpolator;
118 TensorValues<QVectorLayout::byVDIM>(ne, 1, dof2quad, x, y);
126 p(fes.GetMaxElementOrder())
130 MFEM_VERIFY(elem_restr != NULL,
"Missing element restriction.");
133 MFEM_VERIFY(rt_fec,
"Must be RT finite element space.");
135 const int cb_type = rt_fec->GetClosedBasisType();
136 const int ob_type = rt_fec->GetOpenBasisType();
140 if (no_op) {
return; }
142 const int pp1 = p + 1;
151 for (
int i = 0; i < pp1; ++i)
153 cbasis.
Eval(cpts2[i],
b);
154 for (
int j = 0; j < pp1; ++j)
156 Bci_1d[i + j*pp1] =
b[j];
167 const double *ChangeOfBasis_RT::GetOpenMap(Mode mode)
const 178 const double *ChangeOfBasis_RT::GetClosedMap(Mode mode)
const 193 const int D1D = p + 1;
194 const int ND = (p+1)*p;
195 const double *BC = GetClosedMap(mode);
196 const double *BO = GetOpenMap(mode);
202 for (
int c = 0; c <
DIM; ++c)
204 const int nx = (c == 0) ? D1D : D1D-1;
205 const int ny = (c == 1) ? D1D : D1D-1;
206 const double *Bx = (c == 0) ? BC : BO;
207 const double *By = (c == 1) ? BC : BO;
209 for (
int i = 0; i < ND; ++i)
211 Y(i + c*ND, e) = 0.0;
213 for (
int iy = 0; iy < ny; ++ iy)
215 double xx[DofQuadLimits::MAX_D1D];
216 for (
int ix = 0; ix < nx; ++ix) { xx[ix] = 0.0; }
217 for (
int jx = 0; jx < nx; ++jx)
219 const double val = X(jx + iy*nx + c*nx*ny, e);
220 for (
int ix = 0; ix < nx; ++ix)
222 xx[ix] += val*Bx[ix + jx*nx];
225 for (
int jy = 0; jy < ny; ++jy)
227 const double b = By[jy + iy*ny];
228 for (
int ix = 0; ix < nx; ++ix)
230 Y(ix + jy*nx + c*nx*ny, e) += xx[ix]*
b;
242 const int D1D = p + 1;
243 const int ND = (p+1)*p*p;
244 const double *BC = GetClosedMap(mode);
245 const double *BO = GetOpenMap(mode);
251 for (
int c = 0; c <
DIM; ++c)
253 const int nx = (c == 0) ? D1D : D1D-1;
254 const int ny = (c == 1) ? D1D : D1D-1;
255 const int nz = (c == 2) ? D1D : D1D-1;
256 const double *Bx = (c == 0) ? BC : BO;
257 const double *By = (c == 1) ? BC : BO;
258 const double *Bz = (c == 2) ? BC : BO;
260 for (
int i = 0; i < ND; ++i)
262 Y(i + c*ND, e) = 0.0;
264 for (
int iz = 0; iz < nz; ++ iz)
266 double xy[DofQuadLimits::MAX_D1D][DofQuadLimits::MAX_D1D];
267 for (
int iy = 0; iy < ny; ++iy)
269 for (
int ix = 0; ix < nx; ++ix)
274 for (
int iy = 0; iy < ny; ++iy)
276 double xx[DofQuadLimits::MAX_D1D];
277 for (
int ix = 0; ix < nx; ++ix) { xx[ix] = 0.0; }
278 for (
int ix = 0; ix < nx; ++ix)
280 const double val = X(ix + iy*nx + iz*nx*ny + c*ND, e);
281 for (
int jx = 0; jx < nx; ++jx)
283 xx[jx] += val*Bx[jx + ix*nx];
286 for (
int jy = 0; jy < ny; ++jy)
288 const double b = By[jy + iy*ny];
289 for (
int jx = 0; jx < nx; ++jx)
291 xy[jy][jx] += xx[jx] *
b;
295 for (
int jz = 0; jz < nz; ++jz)
297 const double b = Bz[jz + iz*nz];
298 for (
int jy = 0; jy < ny; ++jy)
300 for (
int jx = 0; jx < nx; ++jx)
302 Y(jx + jy*nx + jz*nx*ny + c*ND, e) += xy[jy][jx] *
b;
311 void ChangeOfBasis_RT::Mult(
const Vector &x,
Vector &y, Mode mode)
const 313 if (no_op) { y = x;
return; }
332 elem_restr->
Mult(x_l, x_e);
334 if (dim == 2) {
MultRT_2D(x_e, y_e, mode); }
340 if (R) { R->
Mult(y_l, y); }
341 else { MFEM_VERIFY(P == NULL,
"Invalid state."); }
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Abstract class for all finite elements.
const class FiniteElement * FE
The FiniteElement that created and owns this object.
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
void SetSize(int s)
Resize the vector to size s.
T * GetData()
Returns the data.
void MultInverse(const Vector &x, Vector &y) const
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
ChangeOfBasis_RT(FiniteElementSpace &fes)
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void ComputeInverse(const Array< double > &A, Array< double > &Ainv)
Compute the inverse of the matrix A and store the result in Ainv.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
const double * GetPoints(const int p, const int btype)
Get the coordinates of the points of the given BasisType, btype.
int GetMaxElementOrder() const
Return the maximum polynomial order.
const FiniteElementCollection * FEColl() const
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
void SubcellIntegrals(int n, const Poly_1D::Basis &basis, Array< double > &B)
void MultLeftInverse(const Vector &x, Vector &y) const
void MultRT_2D(const Vector &x, Vector &y, Mode mode) const
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
ChangeOfBasis_L2(FiniteElementSpace &fes)
Operator(int s=0)
Construct a square Operator with given size s (default 0).
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
bool IsIdentityProlongation(const Operator *P)
double p(const Vector &x, double t)
Operator that converts FiniteElementSpace L-vectors to E-vectors.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Integrated GLL indicator functions.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Class for integration point with weight.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
void MultRT_3D(const Vector &x, Vector &y, Mode mode) const
Lexicographic ordering for tensor-product FiniteElements.
int Size() const
Return the logical size of the array.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
double u(const Vector &xvec)
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Arbitrary order "L2-conforming" discontinuous finite elements.