12 #ifndef MFEM_ELASTICITY_KERNELS_HPP
13 #define MFEM_ELASTICITY_KERNELS_HPP
18 using mfem::internal::tensor;
19 using mfem::internal::make_tensor;
24 namespace ElasticityKernels
48 template <
int d1d,
int q1d,
typename material_type>
static inline
49 void Apply3D(
const int ne,
58 static constexpr
int dim = 3;
61 const tensor<double, q1d, d1d> &B =
62 make_tensor<q1d, d1d>([&](
int i,
int j) {
return B_[i + q1d*j]; });
64 const tensor<double, q1d, d1d> &G =
65 make_tensor<q1d, d1d>([&](
int i,
int j) {
return G_[i + q1d*j]; });
67 const auto qweights =
Reshape(W_.
Read(), q1d, q1d, q1d);
69 const auto detJ =
Reshape(detJ_.
Read(), q1d, q1d, q1d, ne);
73 MFEM_FORALL_3D(e, ne, q1d, q1d, q1d,
76 MFEM_SHARED tensor<double, 2, 3, q1d, q1d, q1d> smem;
78 MFEM_SHARED_3D_BLOCK_TENSOR(invJ_sigma_detJw,
double, q1d, q1d, q1d, dim, dim);
80 MFEM_SHARED_3D_BLOCK_TENSOR(dudxi,
double, q1d, q1d, q1d, dim, dim);
82 const auto U_el =
Reshape(&U(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
83 KernelHelpers::CalcGrad(B, G, smem, U_el, dudxi);
85 MFEM_FOREACH_THREAD(qx, x, q1d)
87 MFEM_FOREACH_THREAD(qy, y, q1d)
89 MFEM_FOREACH_THREAD(qz, z, q1d)
91 auto invJqp = inv(make_tensor<dim, dim>(
92 [&](
int i,
int j) {
return J(qx, qy, qz, i, j, e); }));
94 auto dudx = dudxi(qz, qy, qx) * invJqp;
98 auto sigma = material.stress(dudx);
100 invJ_sigma_detJw(qx, qy, qz) =
101 invJqp *
sigma * detJ(qx, qy, qz, e) * qweights(qx, qy, qz);
107 auto F =
Reshape(&force(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
108 KernelHelpers::CalcGradTSum(B, G, smem, invJ_sigma_detJw, F);
141 template <
int d1d,
int q1d,
typename material_type>
static inline
142 void ApplyGradient3D(
const int ne,
146 const Vector &U_,
const material_type &material,
147 const bool use_cache_,
const bool recompute_cache_,
150 static constexpr
int dim = 3;
153 const tensor<double, q1d, d1d> &B =
154 make_tensor<q1d, d1d>([&](
int i,
int j) {
return B_[i + q1d*j]; });
156 const tensor<double, q1d, d1d> &G =
157 make_tensor<q1d, d1d>([&](
int i,
int j) {
return G_[i + q1d*j]; });
159 const auto qweights =
Reshape(W_.
Read(), q1d, q1d, q1d);
161 const auto detJ =
Reshape(detJ_.
Read(), q1d, q1d, q1d, ne);
169 MFEM_FORALL_3D(e, ne, q1d, q1d, q1d,
172 MFEM_SHARED tensor<double, 2, 3, q1d, q1d, q1d> smem;
174 MFEM_SHARED tensor<double, q1d, q1d, q1d, dim, dim> invJ_dsigma_detJw;
176 MFEM_SHARED_3D_BLOCK_TENSOR( dudxi,
double, q1d, q1d, q1d, dim, dim);
177 MFEM_SHARED_3D_BLOCK_TENSOR(ddudxi,
double, q1d, q1d, q1d, dim, dim);
179 const auto U_el =
Reshape(&U(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
180 KernelHelpers::CalcGrad(B, G, smem, U_el, dudxi);
182 const auto dU_el =
Reshape(&dU(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
183 KernelHelpers::CalcGrad(B, G, smem, dU_el, ddudxi);
185 MFEM_FOREACH_THREAD(qx, x, q1d)
187 MFEM_FOREACH_THREAD(qy, y, q1d)
189 MFEM_FOREACH_THREAD(qz, z, q1d)
191 auto invJqp = inv(make_tensor<dim, dim>(
192 [&](
int i,
int j) {
return J(qx, qy, qz, i, j, e); }));
194 auto dudx = dudxi(qz, qy, qx) * invJqp;
195 auto ddudx = ddudxi(qz, qy, qx) * invJqp;
200 tensor<double, dim, dim, dim, dim> C;
202 auto C_cache = make_tensor<dim, dim, dim, dim>(
203 [&](
int i,
int j,
int k,
int l) {
return dsigma_cache(e, qx, qy, qz, i, j, k, l); });
205 if (recompute_cache_)
207 C = material.gradient(dudx);
208 for (
int i = 0; i <
dim; i++)
210 for (
int j = 0; j <
dim; j++)
212 for (
int k = 0; k <
dim; k++)
214 for (
int l = 0; l <
dim; l++)
216 dsigma_cache(e, qx, qy, qz, i, j, k, l) = C(i, j, k, l);
227 invJ_dsigma_detJw(qx, qy, qz) =
228 invJqp * ddot(C, ddudx) * detJ(qx, qy, qz, e) * qweights(qx, qy, qz);
232 auto dsigma = material.action_of_gradient(dudx, ddudx);
233 invJ_dsigma_detJw(qx, qy, qz) =
234 invJqp * dsigma * detJ(qx, qy, qz, e) * qweights(qx, qy, qz);
240 auto F =
Reshape(&force(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
241 KernelHelpers::CalcGradTSum(B, G, smem, invJ_dsigma_detJw, F);
268 template <
int d1d,
int q1d,
typename material_type>
static inline
269 void AssembleGradientDiagonal3D(
const int ne,
277 const material_type &material)
279 static constexpr
int dim = 3;
282 const tensor<double, q1d, d1d> &B =
283 make_tensor<q1d, d1d>([&](
int i,
int j) {
return B_[i + q1d*j]; });
285 const tensor<double, q1d, d1d> &G =
286 make_tensor<q1d, d1d>([&](
int i,
int j) {
return G_[i + q1d*j]; });
288 const auto qweights =
Reshape(W_.
Read(), q1d, q1d, q1d);
290 const auto detJ =
Reshape(detJ_.
Read(), q1d, q1d, q1d, ne);
296 MFEM_FORALL_3D(e, ne, q1d, q1d, q1d,
299 MFEM_SHARED tensor<double, 2, 3, q1d, q1d, q1d> smem;
302 MFEM_SHARED_3D_BLOCK_TENSOR(dudxi,
double, q1d, q1d, q1d, dim, dim);
303 MFEM_SHARED_3D_BLOCK_TENSOR(Ke_diag,
double, d1d, d1d, d1d, dim, dim);
305 const auto U_el =
Reshape(&U(0, 0, 0, 0, e), d1d, d1d, d1d, dim);
306 KernelHelpers::CalcGrad(B, G, smem, U_el, dudxi);
308 MFEM_FOREACH_THREAD(qx, x, q1d)
310 MFEM_FOREACH_THREAD(qy, y, q1d)
312 MFEM_FOREACH_THREAD(qz, z, q1d)
314 const auto invJqp = inv(make_tensor<dim, dim>([&](
int i,
int j)
316 return J(qx, qy, qz, i, j, e);
319 const auto dudx = dudxi(qz, qy, qx) * invJqp;
321 const auto dsigma_ddudx = material.gradient(dudx);
323 const double JxW = detJ(qx, qy, qz, e) * qweights(qx, qy, qz);
324 const auto dphidx = KernelHelpers::GradAllShapeFunctions(qx, qy, qz, B, G, invJqp);
326 for (
int dx = 0; dx < d1d; dx++)
328 for (
int dy = 0; dy < d1d; dy++)
330 for (
int dz = 0; dz < d1d; dz++)
334 const auto phi_i = dphidx(dx, dy, dz);
335 const auto val = JxW * ( phi_i * dsigma_ddudx * phi_i);
336 for (
int l = 0; l <
dim; l++)
338 for (
int m = 0; m <
dim; m++)
340 AtomicAdd(Ke_diag(dx, dy, dz, l, m), val[l][m]);
351 MFEM_FOREACH_THREAD(i, x, d1d)
353 MFEM_FOREACH_THREAD(j, y, d1d)
355 MFEM_FOREACH_THREAD(k, z, d1d)
357 for (
int l = 0; l <
dim; l++)
359 for (
int m = 0; m <
dim; m++)
361 Ke_diag_m(i, j, k, l, e, m) = Ke_diag(i, j, k, l, m);
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
void CheckMemoryRestriction(int d1d, int q1d)
Runtime check for memory restrictions that are determined at compile time.
int material(Vector &x, Vector &xmin, Vector &xmax)
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
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.
double sigma(const Vector &x)