MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
tmop_pa_tc3.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "../tmop.hpp"
13#include "tmop_pa.hpp"
14#include "../gridfunc.hpp"
17
18using namespace mfem;
19
20namespace mfem
21{
22
23MFEM_REGISTER_TMOP_KERNELS(bool, TC_IDEAL_SHAPE_UNIT_SIZE_3D_KERNEL,
24 const int NE,
25 const DenseMatrix &w_,
26 DenseTensor &j_,
27 const int d1d,
28 const int q1d)
29{
30 constexpr int DIM = 3;
31
32 const int Q1D = T_Q1D ? T_Q1D : q1d;
33
34 const auto W = Reshape(w_.Read(), DIM,DIM);
35 auto J = Reshape(j_.Write(), DIM,DIM, Q1D,Q1D,Q1D, NE);
36
37 mfem::forall_3D(NE, Q1D, Q1D, Q1D, [=] MFEM_HOST_DEVICE (int e)
38 {
39 const int Q1D = T_Q1D ? T_Q1D : q1d;
40 MFEM_FOREACH_THREAD(qy,y,Q1D)
41 {
42 MFEM_FOREACH_THREAD(qx,x,Q1D)
43 {
44 MFEM_FOREACH_THREAD(qz,z,Q1D)
45 {
46 kernels::Set(DIM,DIM, 1.0, &W(0,0), &J(0,0,qx,qy,qz,e));
47 }
48 }
49 }
50 });
51 return true;
52}
53
54MFEM_REGISTER_TMOP_KERNELS(bool, TC_IDEAL_SHAPE_GIVEN_SIZE_3D_KERNEL,
55 const int NE,
56 const Array<real_t> &b_,
57 const Array<real_t> &g_,
58 const DenseMatrix &w_,
59 const Vector &x_,
60 DenseTensor &j_,
61 const int d1d,
62 const int q1d)
63{
64 constexpr int DIM = 3;
65
66 const real_t detW = w_.Det();
67 const int D1D = T_D1D ? T_D1D : d1d;
68 const int Q1D = T_Q1D ? T_Q1D : q1d;
69
70 const auto b = Reshape(b_.Read(), Q1D, D1D);
71 const auto g = Reshape(g_.Read(), Q1D, D1D);
72 const auto W = Reshape(w_.Read(), DIM,DIM);
73 const auto X = Reshape(x_.Read(), D1D, D1D, D1D, DIM, NE);
74 auto J = Reshape(j_.Write(), DIM,DIM, Q1D,Q1D,Q1D, NE);
75
76 mfem::forall_3D(NE, Q1D, Q1D, Q1D, [=] MFEM_HOST_DEVICE (int e)
77 {
78 const int D1D = T_D1D ? T_D1D : d1d;
79 const int Q1D = T_Q1D ? T_Q1D : q1d;
80
81 constexpr int MQ1 = T_Q1D ? T_Q1D : T_MAX;
82 constexpr int MD1 = T_D1D ? T_D1D : T_MAX;
83
84 MFEM_SHARED real_t BG[2][MQ1*MD1];
85 MFEM_SHARED real_t DDD[3][MD1*MD1*MD1];
86 MFEM_SHARED real_t DDQ[6][MD1*MD1*MQ1];
87 MFEM_SHARED real_t DQQ[9][MD1*MQ1*MQ1];
88 MFEM_SHARED real_t QQQ[9][MQ1*MQ1*MQ1];
89
90 kernels::internal::LoadX<MD1>(e,D1D,X,DDD);
91 kernels::internal::LoadBG<MD1,MQ1>(D1D,Q1D,b,g,BG);
92
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);
96
97 MFEM_FOREACH_THREAD(qz,z,Q1D)
98 {
99 MFEM_FOREACH_THREAD(qy,y,Q1D)
100 {
101 MFEM_FOREACH_THREAD(qx,x,Q1D)
102 {
103 real_t Jtr[9];
104 const real_t *Wid = &W(0,0);
105 kernels::internal::PullGrad<MQ1>(Q1D,qx,qy,qz,QQQ,Jtr);
106 const real_t detJ = kernels::Det<3>(Jtr);
107 const real_t alpha = std::pow(detJ/detW,1./3);
108 kernels::Set(DIM,DIM,alpha,Wid,&J(0,0,qx,qy,qz,e));
109 }
110 }
111 }
112 });
113 return true;
114}
115
116template<> bool
118 const IntegrationRule &ir,
119 const Vector &,
120 DenseTensor &Jtr) const
121{
122 MFEM_ASSERT(target_type == IDEAL_SHAPE_UNIT_SIZE || nodes != nullptr, "");
123 const Mesh *mesh = fes.GetMesh();
124 const int NE = mesh->GetNE();
125 const int dim = mesh->Dimension();
126 MFEM_VERIFY(mesh->GetNumGeometries(dim) <= 1,
127 "mixed meshes are not supported");
128 MFEM_VERIFY(!fes.IsVariableOrder(), "variable orders are not supported");
129 const FiniteElement &fe = *fes.GetTypicalFE();
130 MFEM_VERIFY(fe.GetGeomType() == Geometry::CUBE, "");
133 const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
134 const Array<real_t> &B = maps.B;
135 const Array<real_t> &G = maps.G;
136 const int D1D = maps.ndof;
137 const int Q1D = maps.nqpt;
138 const int id = (D1D << 4 ) | Q1D;
139
140 switch (target_type)
141 {
142 case IDEAL_SHAPE_UNIT_SIZE: // Jtr(i) = Wideal;
143 {
144 MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_UNIT_SIZE_3D_KERNEL,
145 id,NE,W,Jtr);
146 }
147 case IDEAL_SHAPE_EQUAL_SIZE: return false;
149 {
150 MFEM_VERIFY(nodes, "");
152 const Operator *R = fes.GetElementRestriction(ordering);
154 X.UseDevice(true);
155 R->Mult(*nodes, X);
156 MFEM_ASSERT(nodes->FESpace()->GetVDim() == 3, "");
157 MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_GIVEN_SIZE_3D_KERNEL,
158 id,NE,B,G,W,X,Jtr);
159 }
160 case GIVEN_SHAPE_AND_SIZE: return false;
161 default: return false;
162 }
163 return false;
164}
165
166} // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:337
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
Definition densemat.hpp:474
real_t Det() const
Definition densemat.cpp:516
Rank 3 tensor (array of matrices)
real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:274
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:141
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition fe_base.hpp:214
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition fe_base.hpp:154
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:165
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:193
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:178
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:182
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:709
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition fespace.cpp:1510
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:841
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
Abstract class for all finite elements.
Definition fe_base.hpp:244
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:365
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:331
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition geom.hpp:102
FiniteElementSpace * FESpace()
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Mesh data type.
Definition mesh.hpp:64
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:7243
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
const TargetType target_type
Definition tmop.hpp:1507
const GridFunction * nodes
Definition tmop.hpp:1504
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:494
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:144
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
constexpr int DIM
MFEM_HOST_DEVICE void Set(const int height, const int width, const real_t 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.
Definition kernels.hpp:326
MFEM_HOST_DEVICE T Det(const T *data)
Compute the determinant of a square matrix of size dim with given data.
Definition kernels.hpp:237
MFEM_REGISTER_TMOP_KERNELS(void, DatcSize, const int NE, const int ncomp, const int sizeidx, const real_t input_min_size, const DenseMatrix &w_, const Array< real_t > &b_, const Vector &x_, const Vector &nc_reduce, DenseTensor &j_, const int d1d, const int q1d)
Geometry Geometries
Definition fe.cpp:49
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
void forall_3D(int N, int X, int Y, int Z, lambda &&body)
Definition forall.hpp:774
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:83