MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
tmop_pa_da3.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
17namespace mfem
18{
19
21 const int NE,
22 const int ncomp,
23 const int sizeidx,
24 const real_t input_min_size,
25 const DenseMatrix &w_,
26 const Array<real_t> &b_,
27 const Vector &x_,
28 const Vector &nc_reduce,
29 DenseTensor &j_,
30 const int d1d,
31 const int q1d)
32{
33 MFEM_VERIFY(ncomp==1,"");
34 constexpr int DIM = 3;
35 const int D1D = T_D1D ? T_D1D : d1d;
36 const int Q1D = T_Q1D ? T_Q1D : q1d;
37 MFEM_VERIFY(D1D <= Q1D, "");
38
39 const auto b = Reshape(b_.Read(), Q1D, D1D);
40 const auto W = Reshape(w_.Read(), DIM,DIM);
41 const auto X = Reshape(x_.Read(), D1D, D1D, D1D, ncomp, NE);
42 auto J = Reshape(j_.Write(), DIM,DIM, Q1D,Q1D,Q1D, NE);
43
44 const real_t infinity = std::numeric_limits<real_t>::infinity();
45 MFEM_VERIFY(sizeidx == 0,"");
46 MFEM_VERIFY(MFEM_CUDA_BLOCKS==256,"");
47
48 const real_t *nc_red = nc_reduce.Read();
49
50 mfem::forall_3D(NE, Q1D, Q1D, Q1D, [=] MFEM_HOST_DEVICE (int e)
51 {
52 const int D1D = T_D1D ? T_D1D : d1d;
53 const int Q1D = T_Q1D ? T_Q1D : q1d;
54 constexpr int MQ1 = T_Q1D ? T_Q1D : T_MAX;
55 constexpr int MD1 = T_D1D ? T_D1D : T_MAX;
56 constexpr int MDQ = (MQ1 > MD1) ? MQ1 : MD1;
57
58 MFEM_SHARED real_t sB[MQ1*MD1];
59 MFEM_SHARED real_t sm0[MDQ*MDQ*MDQ];
60 MFEM_SHARED real_t sm1[MDQ*MDQ*MDQ];
61
62 kernels::internal::LoadB<MD1,MQ1>(D1D,Q1D,b,sB);
63
64 ConstDeviceMatrix B(sB, D1D, Q1D);
65 DeviceCube DDD(sm0, MD1,MD1,MD1);
66 DeviceCube DDQ(sm1, MD1,MD1,MQ1);
67 DeviceCube DQQ(sm0, MD1,MQ1,MQ1);
68 DeviceCube QQQ(sm1, MQ1,MQ1,MQ1);
69
70 kernels::internal::LoadX(e,D1D,sizeidx,X,DDD);
71
72 real_t min;
73 MFEM_SHARED real_t min_size[MFEM_CUDA_BLOCKS];
74 DeviceTensor<3,real_t> M((real_t*)(min_size),D1D,D1D,D1D);
75 const DeviceTensor<3,const real_t> D((real_t*)(DDD+sizeidx),D1D,D1D,D1D);
76 MFEM_FOREACH_THREAD(t,x,MFEM_CUDA_BLOCKS) { min_size[t] = infinity; }
77 MFEM_SYNC_THREAD;
78 MFEM_FOREACH_THREAD(dz,z,D1D)
79 {
80 MFEM_FOREACH_THREAD(dy,y,D1D)
81 {
82 MFEM_FOREACH_THREAD(dx,x,D1D)
83 {
84 M(dx,dy,dz) = D(dx,dy,dz);
85 }
86 }
87 }
88 MFEM_SYNC_THREAD;
89 for (int wrk = MFEM_CUDA_BLOCKS >> 1; wrk > 0; wrk >>= 1)
90 {
91 MFEM_FOREACH_THREAD(t,x,MFEM_CUDA_BLOCKS)
92 {
93 if (t < wrk && MFEM_THREAD_ID(y)==0 && MFEM_THREAD_ID(z)==0)
94 {
95 min_size[t] = fmin(min_size[t], min_size[t+wrk]);
96 }
97 }
98 MFEM_SYNC_THREAD;
99 }
100 min = min_size[0];
101 if (input_min_size > 0.) { min = input_min_size; }
102 kernels::internal::EvalX(D1D,Q1D,B,DDD,DDQ);
103 kernels::internal::EvalY(D1D,Q1D,B,DDQ,DQQ);
104 kernels::internal::EvalZ(D1D,Q1D,B,DQQ,QQQ);
105 MFEM_FOREACH_THREAD(qx,x,Q1D)
106 {
107 MFEM_FOREACH_THREAD(qy,y,Q1D)
108 {
109 MFEM_FOREACH_THREAD(qz,z,Q1D)
110 {
111 real_t T;
112 kernels::internal::PullEval(qx,qy,qz,QQQ,T);
113 const real_t shape_par_vals = T;
114 const real_t size = fmax(shape_par_vals, min) / nc_red[e];
115 const real_t alpha = std::pow(size, 1.0/DIM);
116 for (int i = 0; i < DIM; i++)
117 {
118 for (int j = 0; j < DIM; j++)
119 {
120 J(i,j,qx,qy,qz,e) = alpha * W(i,j);
121 }
122 }
123 }
124 }
125 }
126 });
127}
128
129// PA.Jtr Size = (dim, dim, PA.ne*PA.nq);
131 const IntegrationRule &ir,
132 const Vector &xe,
133 DenseTensor &Jtr) const
134{
135 MFEM_VERIFY(target_type == IDEAL_SHAPE_GIVEN_SIZE ||
137
138 MFEM_VERIFY(tspec_fesv, "No target specifications have been set.");
139 const FiniteElementSpace *fes = tspec_fesv;
140
141 // Cases that are not implemented below
142 if (skewidx != -1 ||
143 aspectratioidx != -1 ||
144 orientationidx != -1 ||
145 fes->GetMesh()->Dimension() != 3 ||
146 sizeidx == -1)
147 {
148 return ComputeAllElementTargets_Fallback(pa_fes, ir, xe, Jtr);
149 }
150
151 const Mesh *mesh = fes->GetMesh();
152 const int NE = mesh->GetNE();
153 // Quick return for empty processors:
154 if (NE == 0) { return; }
155 const int dim = mesh->Dimension();
156 MFEM_VERIFY(mesh->GetNumGeometries(dim) <= 1,
157 "mixed meshes are not supported");
158 MFEM_VERIFY(!fes->IsVariableOrder(), "variable orders are not supported");
159 const FiniteElement &fe = *fes->GetFE(0);
162 const DofToQuad &maps = fe.GetDofToQuad(ir, mode);
163 const Array<real_t> &B = maps.B;
164 const int D1D = maps.ndof;
165 const int Q1D = maps.nqpt;
166 const real_t input_min_size = lim_min_size;
167
168 Vector nc_size_red(NE, Device::GetDeviceMemoryType());
169 nc_size_red.HostWrite();
170 NCMesh *ncmesh = tspec_fesv->GetMesh()->ncmesh;
171 for (int e = 0; e < NE; e++)
172 {
173 nc_size_red(e) = (ncmesh) ? ncmesh->GetElementSizeReduction(e) : 1.0;
174 }
175
176 Vector tspec_e;
178 const Operator *R = fes->GetElementRestriction(ordering);
179 MFEM_VERIFY(R,"");
180 MFEM_VERIFY(R->Height() == NE*ncomp*D1D*D1D*D1D,"");
182 tspec_e.UseDevice(true);
183 tspec.UseDevice(true);
184 R->Mult(tspec, tspec_e);
185 const int id = (D1D << 4 ) | Q1D;
186 MFEM_LAUNCH_TMOP_KERNEL(DatcSize,id,NE,ncomp,sizeidx,input_min_size,W,B,
187 tspec_e, nc_size_red, Jtr);
188}
189
190} // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:317
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:466
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:1104
real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(GetMemory(), TotalSize(), on_dev).
Definition: densemat.hpp:1245
A basic generic Tensor class, appropriate for use on the GPU.
Definition: dtensor.hpp:82
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:274
FiniteElementSpace * tspec_fesv
Definition: tmop.hpp:1542
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.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition: fe_base.hpp:137
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition: fe_base.hpp:150
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:161
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:189
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:174
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:178
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:220
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:581
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition: fespace.cpp:3168
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1336
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:559
Abstract class for all finite elements.
Definition: fe_base.hpp:239
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:326
const DenseMatrix & GetGeomToPerfGeomJac(int GeomType) const
Definition: geom.hpp:98
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:100
Mesh data type.
Definition: mesh.hpp:56
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1160
NCMesh * ncmesh
Optional nonconforming mesh extension.
Definition: mesh.hpp:291
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:6961
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
Definition: ncmesh.hpp:123
int GetElementSizeReduction(int i) const
Definition: ncmesh.cpp:5229
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:1360
void ComputeAllElementTargets_Fallback(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition: tmop.cpp:1562
Vector data type.
Definition: vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:474
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:136
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:538
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:486
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_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)
Definition: tmop_pa_da3.cpp:20
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:45
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:775
float real_t
Definition: config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:75
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
RefCoord t[3]