MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_pa_jp3.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 "../tmop_tools.hpp"
15 #include "../../general/forall.hpp"
16 #include "../../linalg/kernels.hpp"
17 
18 namespace mfem
19 {
20 
21 MFEM_REGISTER_TMOP_KERNELS(double, MinDetJpr_Kernel_3D,
22  const int NE,
23  const Array<double> &b_,
24  const Array<double> &g_,
25  const Vector &x_,
26  Vector &DetJ,
27  const int d1d,
28  const int q1d)
29 {
30  constexpr int DIM = 3;
31  const int D1D = T_D1D ? T_D1D : d1d;
32  const int Q1D = T_Q1D ? T_Q1D : q1d;
33 
34  const auto b = Reshape(b_.Read(), Q1D, D1D);
35  const auto g = Reshape(g_.Read(), Q1D, D1D);
36  const auto X = Reshape(x_.Read(), D1D, D1D, D1D, DIM, NE);
37 
38  auto E = Reshape(DetJ.Write(), Q1D, Q1D, Q1D, NE);
39 
40  MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
41  {
42  const int D1D = T_D1D ? T_D1D : d1d;
43  const int Q1D = T_Q1D ? T_Q1D : q1d;
44  constexpr int MQ1 = T_Q1D ? T_Q1D : T_MAX;
45  constexpr int MD1 = T_D1D ? T_D1D : T_MAX;
46 
47  MFEM_SHARED double BG[2][MQ1*MD1];
48  MFEM_SHARED double DDD[3][MD1*MD1*MD1];
49  MFEM_SHARED double DDQ[6][MD1*MD1*MQ1];
50  MFEM_SHARED double DQQ[9][MD1*MQ1*MQ1];
51  MFEM_SHARED double QQQ[9][MQ1*MQ1*MQ1];
52 
53  kernels::internal::LoadX<MD1>(e,D1D,X,DDD);
54  kernels::internal::LoadBG<MD1,MQ1>(D1D,Q1D,b,g,BG);
55 
56  kernels::internal::GradX<MD1,MQ1>(D1D,Q1D,BG,DDD,DDQ);
57  kernels::internal::GradY<MD1,MQ1>(D1D,Q1D,BG,DDQ,DQQ);
58  kernels::internal::GradZ<MD1,MQ1>(D1D,Q1D,BG,DQQ,QQQ);
59 
60  MFEM_FOREACH_THREAD(qz,z,Q1D)
61  {
62  MFEM_FOREACH_THREAD(qy,y,Q1D)
63  {
64  MFEM_FOREACH_THREAD(qx,x,Q1D)
65  {
66  double J[9];
67  kernels::internal::PullGrad<MQ1>(Q1D,qx,qy,qz,QQQ,J);
68  E(qx,qy,qz,e) = kernels::Det<3>(J);
69  }
70  }
71  }
72  });
73 
74  return DetJ.Min();
75 }
76 
78  const Vector &X) const
79 {
81  const Operator *R = fes->GetElementRestriction(ordering);
83  XE.UseDevice(true);
84  R->Mult(X, XE);
85 
86  const DofToQuad &maps = fes->GetFE(0)->GetDofToQuad(ir, DofToQuad::TENSOR);
87  const int NE = fes->GetMesh()->GetNE();
88  const int NQ = ir.GetNPoints();
89  const int D1D = maps.ndof;
90  const int Q1D = maps.nqpt;
91  const int id = (D1D << 4 ) | Q1D;
92  const Array<double> &B = maps.B;
93  const Array<double> &G = maps.G;
94 
95  Vector E(NE*NQ);
96  E.UseDevice(true);
97 
98  MFEM_LAUNCH_TMOP_KERNEL(MinDetJpr_Kernel_3D,id,NE,B,G,XE,E);
99 }
100 
101 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
const IntegrationRule & ir
Definition: tmop_tools.hpp:145
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:162
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:174
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
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
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:457
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
double Min() const
Returns the minimal element of the vector.
Definition: vector.cpp:1202
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
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:366
double MinDetJpr_3D(const FiniteElementSpace *, const Vector &) const
Definition: tmop_pa_jp3.cpp:77
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:206
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
Lexicographic ordering for tensor-product FiniteElements.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1259
Vector data type.
Definition: vector.hpp:60
Abstract operator.
Definition: operator.hpp:24
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