12 #include "../tmop.hpp"
14 #include "../gridfunc.hpp"
15 #include "../../general/forall.hpp"
16 #include "../../linalg/kernels.hpp"
30 constexpr
int DIM = 3;
32 const int Q1D = T_Q1D ? T_Q1D : q1d;
37 MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
39 const int Q1D = T_Q1D ? T_Q1D : q1d;
40 MFEM_FOREACH_THREAD(qy,y,Q1D)
42 MFEM_FOREACH_THREAD(qx,x,Q1D)
44 MFEM_FOREACH_THREAD(qz,z,Q1D)
46 kernels::Set(DIM,DIM, 1.0, &W(0,0), &J(0,0,qx,qy,qz,e));
64 constexpr
int DIM = 3;
66 const double detW = w_.
Det();
67 const int D1D = T_D1D ? T_D1D : d1d;
68 const int Q1D = T_Q1D ? T_Q1D : q1d;
76 MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
78 const int D1D = T_D1D ? T_D1D : d1d;
79 const int Q1D = T_Q1D ? T_Q1D : q1d;
81 constexpr
int MQ1 = T_Q1D ? T_Q1D : T_MAX;
82 constexpr
int MD1 = T_D1D ? T_D1D : T_MAX;
84 MFEM_SHARED
double BG[2][MQ1*MD1];
85 MFEM_SHARED
double DDD[3][MD1*MD1*MD1];
86 MFEM_SHARED
double DDQ[6][MD1*MD1*MQ1];
87 MFEM_SHARED
double DQQ[9][MD1*MQ1*MQ1];
88 MFEM_SHARED
double QQQ[9][MQ1*MQ1*MQ1];
90 kernels::internal::LoadX<MD1>(e,D1D,X,DDD);
91 kernels::internal::LoadBG<MD1,MQ1>(D1D,Q1D,
b,g,BG);
93 kernels::internal::GradX<MD1,MQ1>(D1D,Q1D,BG,DDD,DDQ);
94 kernels::internal::GradY<MD1,MQ1>(D1D,Q1D,BG,DDQ,DQQ);
95 kernels::internal::GradZ<MD1,MQ1>(D1D,Q1D,BG,DQQ,QQQ);
97 MFEM_FOREACH_THREAD(qz,z,Q1D)
99 MFEM_FOREACH_THREAD(qy,y,Q1D)
101 MFEM_FOREACH_THREAD(qx,x,Q1D)
104 const double *Wid = &W(0,0);
105 kernels::internal::PullGrad<MQ1>(Q1D,qx,qy,qz,QQQ,Jtr);
106 const double detJ = kernels::Det<3>(Jtr);
107 const double alpha = std::pow(detJ/detW,1./3);
122 MFEM_ASSERT(target_type == IDEAL_SHAPE_UNIT_SIZE || nodes !=
nullptr,
"");
123 const Mesh *mesh = fes.GetMesh();
124 const int NE = mesh->
GetNE();
126 if (NE == 0) {
return true; }
129 "mixed meshes are not supported");
130 MFEM_VERIFY(!fes.IsVariableOrder(),
"variable orders are not supported");
135 const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
138 const int D1D = maps.
ndof;
139 const int Q1D = maps.
nqpt;
140 const int id = (D1D << 4 ) | Q1D;
144 case IDEAL_SHAPE_UNIT_SIZE:
146 MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_UNIT_SIZE_3D_KERNEL,
149 case IDEAL_SHAPE_EQUAL_SIZE:
return false;
150 case IDEAL_SHAPE_GIVEN_SIZE:
152 MFEM_VERIFY(nodes,
"");
154 const Operator *R = fes.GetElementRestriction(ordering);
158 MFEM_ASSERT(nodes->FESpace()->GetVDim() == 3,
"");
159 MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_GIVEN_SIZE_3D_KERNEL,
162 case GIVEN_SHAPE_AND_SIZE:
return false;
163 default:
return false;
Abstract class for all finite elements.
Class for an integration rule - an Array of IntegrationPoint.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
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...
int GetNE() const
Returns number of elements.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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)
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Array< double > B
Basis functions evaluated at quadrature points.
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.
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Lexicographic ordering for tensor-product FiniteElements.
MFEM_HOST_DEVICE void Set(const int height, const int width, const double alpha, const TA *Adata, TB *Bdata)
Compute B = alpha*A, where the matrices A and B are of size height x width with data Adata and Bdata...
double * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
Rank 3 tensor (array of matrices)
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), 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.
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).