31#ifndef MFEM_BILININTEG_ELASTICITY_KERNELS_HPP
32#define MFEM_BILININTEG_ELASTICITY_KERNELS_HPP
66void ElasticityAddMultPA(
const int dim,
const int nDofs,
67 const FiniteElementSpace &fespace,
const CoefficientVector &lambda,
68 const CoefficientVector &
mu,
const GeometricFactors &geom,
69 const DofToQuad &maps,
const Vector &x, QuadratureFunction &QVec, Vector &y);
94void ElasticityComponentAddMultPA(
95 const int dim,
const int nDofs,
const FiniteElementSpace &fespace,
96 const CoefficientVector &lambda,
const CoefficientVector &
mu,
97 const GeometricFactors &geom,
const DofToQuad &maps,
const Vector &x,
98 QuadratureFunction &QVec, Vector &y,
const int i_block,
const int j_block);
122void ElasticityAssembleEA(
const int dim,
const int i_block,
const int j_block,
123 const int nDofs,
const IntegrationRule &ir,
124 const CoefficientVector &lambda,
125 const CoefficientVector &
mu,
const GeometricFactors &geom,
126 const DofToQuad &maps, Vector &emat);
138void ElasticityAssembleDiagonalPA(
const int dim,
const int nDofs,
139 const CoefficientVector &lambda,
140 const CoefficientVector &
mu,
const GeometricFactors &geom,
141 const DofToQuad &maps, QuadratureFunction &QVec, Vector &diag);
144template<
int dim,
int i_block = -1,
int j_block = -1>
145void ElasticityAddMultPA_(
const int nDofs,
const FiniteElementSpace &fespace,
146 const CoefficientVector &lambda,
const CoefficientVector &
mu,
147 const GeometricFactors &geom,
const DofToQuad &maps,
const Vector &x,
148 QuadratureFunction &QVec, Vector &y)
150 using future::tensor;
155 static_assert((i_block < 0) == (j_block < 0),
156 "i_block and j_block must both be non-negative or strictly negative.");
157 static constexpr int d =
dim;
158 static constexpr int qLower = (i_block < 0) ? 0 : i_block;
159 static constexpr int qUpper = (i_block < 0) ? d : i_block+1;
160 static constexpr int qSize = qUpper-qLower;
161 static constexpr int aLower = (j_block < 0) ? 0 : j_block;
162 static constexpr int aUpper = (j_block < 0) ? d : j_block+1;
163 static constexpr int aSize = aUpper-aLower;
164 static constexpr bool isComponent = (i_block >= 0);
167 const auto &ir = QVec.GetIntRule(0);
168 const QuadratureInterpolator *E_To_Q_Map = fespace.GetQuadratureInterpolator(
172 E_To_Q_Map->PhysDerivatives(x, QVec);
174 const int numPoints = ir.GetNPoints();
175 const int numEls = fespace.GetNE();
176 const auto lamDev =
Reshape(lambda.Read(), numPoints, numEls);
177 const auto muDev =
Reshape(
mu.Read(), numPoints, numEls);
178 const auto J =
Reshape(geom.J.Read(), numPoints, d, d, numEls);
179 auto Q =
Reshape(QVec.ReadWrite(), numPoints, d, qSize, numEls);
180 const real_t *ipWeights = ir.GetWeights().Read();
184 MFEM_FOREACH_THREAD(
p, x,numPoints)
186 auto invJ =
inv(make_tensor<d, d>(
187 [&](
int i,
int j) {
return J(
p, i, j, e); }));
188 tensor<real_t, aSize, d> gradx;
192 for (
int i = 0; i < d; i++)
194 gradx(0,i) = Q(
p, i, 0, e);
199 for (
int j = 0; j < d; j++)
201 for (
int i = 0; i < d; i++)
203 gradx(i,j) = Q(
p, i, j, e);
209 for (
int i = aLower; i < aUpper; i++)
212 const int iIndex = isComponent ? 0 : i;
213 div += gradx(iIndex,i);
216 for (
int m = 0; m < d; m++)
218 for (
int q = qLower; q < qUpper; q++)
227 for (
int a = 0;
a < d;
a++)
229 contraction += 2*((
a == q)*invJ(m,j_block)
230 + (j_block==q)*invJ(m,
a))*(gradx(0,
a));
235 for (
int a = 0;
a < d;
a++)
237 for (
int b = 0;
b < d;
b++)
239 contraction += ((
a == q)*invJ(m,
b) + (
b == q)*invJ(m,
a))
240 *(gradx(
a,
b) + gradx(
b,
a));
246 const int qIndex = isComponent ? 0 : q;
247 Q(
p,m,qIndex,e) = w*(lamDev(
p, e)*invJ(m,q)*div
248 + 0.5*muDev(
p, e)*contraction);
255 const auto QRead =
Reshape(QVec.Read(), numPoints, d, qSize, numEls);
256 const auto G =
Reshape(maps.G.Read(), numPoints, d, nDofs);
257 auto yDev =
Reshape(y.ReadWrite(), nDofs, qSize, numEls);
260 MFEM_FOREACH_THREAD(i, y, nDofs)
262 MFEM_FOREACH_THREAD(q, x, qSize)
264 const int qIndex = isComponent ? 0 : q;
266 for (
int m = 0; m < d; m++ )
268 for (
int p = 0;
p < numPoints;
p++ )
270 sum += QRead(
p,m,qIndex,e)*G(
p,m,i);
273 yDev(i, qIndex, e) += sum;
281void ElasticityAssembleDiagonalPA_(
const int nDofs,
282 const CoefficientVector &lambda,
283 const CoefficientVector &
mu,
const GeometricFactors &geom,
284 const DofToQuad &maps, QuadratureFunction &QVec, Vector &diag)
286 using future::tensor;
292 const auto &ir = QVec.GetIntRule(0);
293 static constexpr int d =
dim;
294 const int numPoints = ir.GetNPoints();
295 const int numEls = lambda.Size()/numPoints;
296 const auto lamDev =
Reshape(lambda.Read(), numPoints, numEls);
297 const auto muDev =
Reshape(
mu.Read(), numPoints, numEls);
298 const auto J =
Reshape(geom.J.Read(), numPoints, d, d, numEls);
299 auto Q =
Reshape(QVec.ReadWrite(), numPoints, d,d, d, numEls);
300 const real_t *ipWeights = ir.GetWeights().Read();
303 MFEM_FOREACH_THREAD(
p, x,numPoints)
305 auto invJ =
inv(make_tensor<d, d>(
306 [&](
int i,
int j) {
return J(
p, i, j, e); }));
308 for (
int n = 0; n < d; n++)
310 for (
int m = 0; m < d; m++)
312 for (
int q = 0; q < d; q++)
318 for (
int a = 0;
a < d;
a++)
320 for (
int b = 0;
b < d;
b++)
322 contraction += ((
a == q)*invJ(m,
b) + (
b==q)*invJ(m,
a))*((
a == q)
323 *invJ(n,
b) + (
b==q)*invJ(n,
a));
328 Q(
p,m,n,q,e) = w*(lamDev(
p, e)*invJ(m,q)*invJ(n,q)
329 + 0.5*muDev(
p, e)*contraction);
337 const auto QRead =
Reshape(QVec.Read(), numPoints, d, d, d, numEls);
338 auto diagDev =
Reshape(diag.Write(), nDofs, d, numEls);
339 const auto G =
Reshape(maps.G.Read(), numPoints, d, nDofs);
342 MFEM_FOREACH_THREAD(i, y, nDofs)
344 MFEM_FOREACH_THREAD(q, x, d)
347 for (
int n = 0; n < d; n++)
349 for (
int m = 0; m < d; m++)
351 for (
int p = 0;
p < numPoints;
p++ )
353 sum += QRead(
p,m,n,q,e)*G(
p,m,i)*G(
p,n,i);
357 diagDev(i, q, e) = sum;
365void ElasticityAssembleEA_(
const int i_block,
368 const IntegrationRule &ir,
369 const CoefficientVector &lambda,
370 const CoefficientVector &
mu,
371 const GeometricFactors &geom,
372 const DofToQuad &maps,
375 using future::tensor;
381 static constexpr int d =
dim;
382 const int numPoints = ir.GetNPoints();
383 const int numEls = lambda.Size()/numPoints;
384 const auto lamDev =
Reshape(lambda.Read(), numPoints, numEls);
385 const auto muDev =
Reshape(
mu.Read(), numPoints, numEls);
386 const auto J =
Reshape(geom.J.Read(), numPoints, d, d, numEls);
387 const auto G =
Reshape(maps.G.Read(), numPoints, d, nDofs);
388 auto ematDev =
Reshape(emat.Write(), nDofs, nDofs, numEls);
389 const real_t *ipWeights = ir.GetWeights().Read();
392 MFEM_FOREACH_THREAD(JDof, y, nDofs)
394 MFEM_FOREACH_THREAD(IDof, x, nDofs)
397 for (
int p = 0 ;
p < numPoints;
p++)
399 auto invJ =
inv(make_tensor<d, d>(
400 [&](
int i,
int j) {
return J(
p, i, j, e); }));
402 for (
int n = 0; n < d; n++)
404 for (
int m = 0; m < d; m++)
408 for (
int a = 0;
a < d;
a++)
410 for (
int b = 0;
b < d;
b++)
412 contraction += ((
a == i_block)*invJ(m,
b) + (
b==i_block)*invJ(m,
413 a))*((
a == j_block)*invJ(n,
414 b) + (
b==j_block)*invJ(n,
a));
419 sum += w*(lamDev(
p, e)*invJ(m,i_block)*invJ(n,j_block)
420 + 0.5*muDev(
p, e)*contraction)*G(
p,m,IDof)*G(
p,n,JDof);
424 ematDev(IDof, JDof, e) = sum;
MFEM_HOST_DEVICE T det(const tensor< T, 1, 1 > &A)
Returns the determinant of a matrix.
MFEM_HOST_DEVICE constexpr auto make_tensor(lambda_type f) -> tensor< decltype(f())>
Creates a tensor of requested dimension by subsequent calls to a functor Can be thought of as analogo...
MFEM_HOST_DEVICE tensor< T, 1, 1 > inv(const tensor< T, 1, 1 > &A)
Inverts a matrix.
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void forall_2D(int N, int X, int Y, lambda &&body)
real_t p(const Vector &x, real_t t)
Implementation of the tensor class.