MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_pa_tc3.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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"
15 #include "../../general/forall.hpp"
16 #include "../../linalg/kernels.hpp"
17 
18 using namespace mfem;
19 
20 namespace mfem
21 {
22 
23 MFEM_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(e, NE, Q1D, Q1D, Q1D,
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 
54 MFEM_REGISTER_TMOP_KERNELS(bool, TC_IDEAL_SHAPE_GIVEN_SIZE_3D_KERNEL,
55  const int NE,
56  const Array<double> &b_,
57  const Array<double> &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 double 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(e, NE, Q1D, Q1D, Q1D,
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 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];
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  double Jtr[9];
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);
108  kernels::Set(DIM,DIM,alpha,Wid,&J(0,0,qx,qy,qz,e));
109  }
110  }
111  }
112  });
113  return true;
114 }
115 
116 template<> bool
117 TargetConstructor::ComputeAllElementTargets<3>(const FiniteElementSpace &fes,
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  // Quick return for empty processors:
126  if (NE == 0) { return true; }
127  const int dim = mesh->Dimension();
128  MFEM_VERIFY(mesh->GetNumGeometries(dim) <= 1,
129  "mixed meshes are not supported");
130  MFEM_VERIFY(!fes.IsVariableOrder(), "variable orders are not supported");
131  const FiniteElement &fe = *fes.GetFE(0);
132  MFEM_VERIFY(fe.GetGeomType() == Geometry::CUBE, "");
134  const DofToQuad::Mode mode = DofToQuad::TENSOR;
135  const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
136  const Array<double> &B = maps.B;
137  const Array<double> &G = maps.G;
138  const int D1D = maps.ndof;
139  const int Q1D = maps.nqpt;
140  const int id = (D1D << 4 ) | Q1D;
141 
142  switch (target_type)
143  {
144  case IDEAL_SHAPE_UNIT_SIZE: // Jtr(i) = Wideal;
145  {
146  MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_UNIT_SIZE_3D_KERNEL,
147  id,NE,W,Jtr);
148  }
149  case IDEAL_SHAPE_EQUAL_SIZE: return false;
150  case IDEAL_SHAPE_GIVEN_SIZE:
151  {
152  MFEM_VERIFY(nodes, "");
154  const Operator *R = fes.GetElementRestriction(ordering);
156  X.UseDevice(true);
157  R->Mult(*nodes, X);
158  MFEM_ASSERT(nodes->FESpace()->GetVDim() == 3, "");
159  MFEM_LAUNCH_TMOP_KERNEL(TC_IDEAL_SHAPE_GIVEN_SIZE_3D_KERNEL,
160  id,NE,B,G,W,X,Jtr);
161  }
162  case GIVEN_SHAPE_AND_SIZE: return false;
163  default: return false;
164  }
165  return false;
166 }
167 
168 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:162
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5908
double Det() const
Definition: densemat.cpp:449
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:174
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:170
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
constexpr int DIM
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:118
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:98
Geometry Geometries
Definition: fe.cpp:49
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
double b
Definition: lissajous.cpp:42
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
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)
Definition: tmop_pa_da3.cpp:20
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
int Dimension() const
Definition: mesh.hpp:1006
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:185
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition: fe_base.hpp:149
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
int dim
Definition: ex24.cpp:53
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:206
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...
Definition: kernels.hpp:326
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
double * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
Definition: densemat.hpp:1095
Abstract operator.
Definition: operator.hpp:24
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:953
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
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
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read( GetMemory(), TotalSize(), on_dev).
Definition: densemat.hpp:443