MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_pa_da3.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 "../../general/forall.hpp"
15 #include "../../linalg/kernels.hpp"
16 
17 namespace mfem
18 {
19 
21  const int NE,
22  const int ncomp,
23  const int sizeidx,
24  const DenseMatrix &w_,
25  const Array<double> &b_,
26  const Vector &x_,
27  DenseTensor &j_,
28  const int d1d,
29  const int q1d)
30 {
31  MFEM_VERIFY(ncomp==1,"");
32  constexpr int DIM = 3;
33  const int D1D = T_D1D ? T_D1D : d1d;
34  const int Q1D = T_Q1D ? T_Q1D : q1d;
35  MFEM_VERIFY(D1D <= Q1D, "");
36 
37  const auto b = Reshape(b_.Read(), Q1D, D1D);
38  const auto W = Reshape(w_.Read(), DIM,DIM);
39  const auto X = Reshape(x_.Read(), D1D, D1D, D1D, ncomp, NE);
40  auto J = Reshape(j_.Write(), DIM,DIM, Q1D,Q1D,Q1D, NE);
41 
43  MFEM_VERIFY(sizeidx == 0,"");
44  MFEM_VERIFY(MFEM_CUDA_BLOCKS==256,"");
45 
46  MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
47  {
48  const int D1D = T_D1D ? T_D1D : d1d;
49  const int Q1D = T_Q1D ? T_Q1D : q1d;
50  constexpr int MQ1 = T_Q1D ? T_Q1D : T_MAX;
51  constexpr int MD1 = T_D1D ? T_D1D : T_MAX;
52  constexpr int MDQ = (MQ1 > MD1) ? MQ1 : MD1;
53 
54  MFEM_SHARED double sB[MQ1*MD1];
55  MFEM_SHARED double sm0[MDQ*MDQ*MDQ];
56  MFEM_SHARED double sm1[MDQ*MDQ*MDQ];
57 
58  kernels::internal::LoadB<MD1,MQ1>(D1D,Q1D,b,sB);
59 
60  ConstDeviceMatrix B(sB, D1D, Q1D);
61  DeviceCube DDD(sm0, MD1,MD1,MD1);
62  DeviceCube DDQ(sm1, MD1,MD1,MQ1);
63  DeviceCube DQQ(sm0, MD1,MQ1,MQ1);
64  DeviceCube QQQ(sm1, MQ1,MQ1,MQ1);
65 
66  kernels::internal::LoadX(e,D1D,sizeidx,X,DDD);
67 
68  double min;
69  MFEM_SHARED double min_size[MFEM_CUDA_BLOCKS];
70  DeviceTensor<3,double> M((double*)(min_size),D1D,D1D,D1D);
71  const DeviceTensor<3,const double> D((double*)(DDD+sizeidx),D1D,D1D,D1D);
72  MFEM_FOREACH_THREAD(t,x,MFEM_CUDA_BLOCKS) { min_size[t] = infinity; }
73  MFEM_SYNC_THREAD;
74  MFEM_FOREACH_THREAD(dz,z,D1D)
75  {
76  MFEM_FOREACH_THREAD(dy,y,D1D)
77  {
78  MFEM_FOREACH_THREAD(dx,x,D1D)
79  {
80  M(dx,dy,dz) = D(dx,dy,dz);
81  }
82  }
83  }
84  MFEM_SYNC_THREAD;
85  for (int wrk = MFEM_CUDA_BLOCKS >> 1; wrk > 0; wrk >>= 1)
86  {
87  MFEM_FOREACH_THREAD(t,x,MFEM_CUDA_BLOCKS)
88  {
89  if (t < wrk && MFEM_THREAD_ID(y)==0 && MFEM_THREAD_ID(z)==0)
90  {
91  min_size[t] = fmin(min_size[t], min_size[t+wrk]);
92  }
93  }
94  MFEM_SYNC_THREAD;
95  }
96  min = min_size[0];
97 
98  kernels::internal::EvalX(D1D,Q1D,B,DDD,DDQ);
99  kernels::internal::EvalY(D1D,Q1D,B,DDQ,DQQ);
100  kernels::internal::EvalZ(D1D,Q1D,B,DQQ,QQQ);
101  MFEM_FOREACH_THREAD(qx,x,Q1D)
102  {
103  MFEM_FOREACH_THREAD(qy,y,Q1D)
104  {
105  MFEM_FOREACH_THREAD(qz,z,Q1D)
106  {
107  double T;
108  kernels::internal::PullEval(qx,qy,qz,QQQ,T);
109  const double shape_par_vals = T;
110  const double size = fmax(shape_par_vals, min);
111  const double alpha = std::pow(size, 1.0/DIM);
112  for (int i = 0; i < DIM; i++)
113  {
114  for (int j = 0; j < DIM; j++)
115  {
116  J(i,j,qx,qy,qz,e) = alpha * W(i,j);
117  }
118  }
119  }
120  }
121  }
122  });
123 }
124 
125 // PA.Jtr Size = (dim, dim, PA.ne*PA.nq);
127  const IntegrationRule &ir,
128  const Vector &xe,
129  DenseTensor &Jtr) const
130 {
131  MFEM_VERIFY(target_type == IDEAL_SHAPE_GIVEN_SIZE ||
133 
134  MFEM_VERIFY(tspec_fesv, "No target specifications have been set.");
135  const FiniteElementSpace *fes = tspec_fesv;
136 
137  // Cases that are not implemented below
138  if (skewidx != -1 ||
139  aspectratioidx != -1 ||
140  orientationidx != -1 ||
141  fes->GetMesh()->Dimension() != 3 ||
142  sizeidx == -1)
143  {
144  return ComputeAllElementTargets_Fallback(pa_fes, ir, xe, Jtr);
145  }
146 
147  const Mesh *mesh = fes->GetMesh();
148  const int NE = mesh->GetNE();
149  // Quick return for empty processors:
150  if (NE == 0) { return; }
151  const int dim = mesh->Dimension();
152  MFEM_VERIFY(mesh->GetNumGeometries(dim) <= 1,
153  "mixed meshes are not supported");
154  MFEM_VERIFY(!fes->IsVariableOrder(), "variable orders are not supported");
155  const FiniteElement &fe = *fes->GetFE(0);
156  const DenseMatrix &W = Geometries.GetGeomToPerfGeomJac(fe.GetGeomType());
157  const DofToQuad::Mode mode = DofToQuad::TENSOR;
158  const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
159  const Array<double> &B = maps.B;
160  const int D1D = maps.ndof;
161  const int Q1D = maps.nqpt;
162 
163  Vector tspec_e;
165  const Operator *R = fes->GetElementRestriction(ordering);
166  MFEM_VERIFY(R,"");
167  MFEM_VERIFY(R->Height() == NE*ncomp*D1D*D1D*D1D,"");
168  tspec_e.SetSize(R->Height(), Device::GetDeviceMemoryType());
169  tspec_e.UseDevice(true);
170  tspec.UseDevice(true);
171  R->Mult(tspec, tspec_e);
172  const int id = (D1D << 4 ) | Q1D;
173  MFEM_LAUNCH_TMOP_KERNEL(DatcSize,id,NE,ncomp,sizeidx,W,B,tspec_e,Jtr);
174 }
175 
176 } // 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
const TargetType target_type
Definition: tmop.hpp:1193
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:463
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
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
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
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
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
FiniteElementSpace * tspec_fesv
Definition: tmop.hpp:1375
int Dimension() const
Definition: mesh.hpp:1006
A basic generic Tensor class, appropriate for use on the GPU.
Definition: dtensor.hpp:81
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
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition: tmop.cpp:1372
Lexicographic ordering for tensor-product FiniteElements.
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1259
virtual void ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Computes reference-to-target transformation Jacobians for all quadrature points in all elements...
RefCoord t[3]
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