12 #include "../tmop.hpp" 14 #include "../../general/forall.hpp" 15 #include "../../linalg/kernels.hpp" 31 MFEM_VERIFY(ncomp==1,
"");
32 constexpr
int DIM = 3;
33 const int D1D = T_D1D ? T_D1D : d1d;
34 const int Q1D = T_Q1D ? T_Q1D : q1d;
35 MFEM_VERIFY(D1D <= Q1D,
"");
39 const auto X =
Reshape(x_.
Read(), D1D, D1D, D1D, ncomp, NE);
43 MFEM_VERIFY(sizeidx == 0,
"");
44 MFEM_VERIFY(MFEM_CUDA_BLOCKS==256,
"");
46 MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
48 const int D1D = T_D1D ? T_D1D : d1d;
49 const int Q1D = T_Q1D ? T_Q1D : q1d;
50 constexpr
int MQ1 = T_Q1D ? T_Q1D : T_MAX;
51 constexpr
int MD1 = T_D1D ? T_D1D : T_MAX;
52 constexpr
int MDQ = (MQ1 > MD1) ? MQ1 : MD1;
54 MFEM_SHARED
double sB[MQ1*MD1];
55 MFEM_SHARED
double sm0[MDQ*MDQ*MDQ];
56 MFEM_SHARED
double sm1[MDQ*MDQ*MDQ];
58 kernels::internal::LoadB<MD1,MQ1>(D1D,Q1D,
b,sB);
66 kernels::internal::LoadX(e,D1D,sizeidx,X,DDD);
69 MFEM_SHARED
double min_size[MFEM_CUDA_BLOCKS];
72 MFEM_FOREACH_THREAD(
t,x,MFEM_CUDA_BLOCKS) { min_size[
t] =
infinity; }
74 MFEM_FOREACH_THREAD(dz,z,D1D)
76 MFEM_FOREACH_THREAD(dy,y,D1D)
78 MFEM_FOREACH_THREAD(dx,x,D1D)
80 M(dx,dy,dz) = D(dx,dy,dz);
85 for (
int wrk = MFEM_CUDA_BLOCKS >> 1; wrk > 0; wrk >>= 1)
87 MFEM_FOREACH_THREAD(
t,x,MFEM_CUDA_BLOCKS)
89 if (
t < wrk && MFEM_THREAD_ID(y)==0 && MFEM_THREAD_ID(z)==0)
91 min_size[
t] = fmin(min_size[
t], min_size[
t+wrk]);
98 kernels::internal::EvalX(D1D,Q1D,B,DDD,DDQ);
99 kernels::internal::EvalY(D1D,Q1D,B,DDQ,DQQ);
100 kernels::internal::EvalZ(D1D,Q1D,B,DQQ,QQQ);
101 MFEM_FOREACH_THREAD(qx,x,Q1D)
103 MFEM_FOREACH_THREAD(qy,y,Q1D)
105 MFEM_FOREACH_THREAD(qz,z,Q1D)
108 kernels::internal::PullEval(qx,qy,qz,QQQ,T);
109 const double shape_par_vals = T;
110 const double size = fmax(shape_par_vals, min);
111 const double alpha = std::pow(size, 1.0/
DIM);
112 for (
int i = 0; i <
DIM; i++)
114 for (
int j = 0; j <
DIM; j++)
116 J(i,j,qx,qy,qz,e) =
alpha * W(i,j);
134 MFEM_VERIFY(
tspec_fesv,
"No target specifications have been set.");
148 const int NE = mesh->
GetNE();
150 if (NE == 0) {
return; }
153 "mixed meshes are not supported");
154 MFEM_VERIFY(!fes->
IsVariableOrder(),
"variable orders are not supported");
158 const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
160 const int D1D = maps.
ndof;
161 const int Q1D = maps.
nqpt;
167 MFEM_VERIFY(R->
Height() == NE*
ncomp*D1D*D1D*D1D,
"");
172 const int id = (D1D << 4 ) | Q1D;
173 MFEM_LAUNCH_TMOP_KERNEL(DatcSize,
id,NE,
ncomp,
sizeidx,W,B,tspec_e,Jtr);
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Abstract class for all finite elements.
virtual void ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Computes reference-to-target transformation Jacobians for all quadrature points in all elements...
Class for an integration rule - an Array of IntegrationPoint.
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
const TargetType target_type
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
void SetSize(int s)
Resize the vector to size s.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Data type dense matrix using column-major storage.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual void Mult(const Vector &x, Vector &y) const =0
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...
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
MFEM_REGISTER_TMOP_KERNELS(void, DatcSize, const int NE, const int ncomp, const int sizeidx, const DenseMatrix &w_, const Array< double > &b_, const Vector &x_, DenseTensor &j_, const int d1d, const int q1d)
FiniteElementSpace * tspec_fesv
A basic generic Tensor class, appropriate for use on the GPU.
Mesh * GetMesh() const
Returns the mesh.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Array< double > B
Basis functions evaluated at quadrature points.
int GetNE() const
Returns number of elements.
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Lexicographic ordering for tensor-product FiniteElements.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
double * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Rank 3 tensor (array of matrices)
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const