MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tmop_pa.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 "../linearform.hpp"
14 #include "../pgridfunc.hpp"
15 #include "../tmop_tools.hpp"
16 #include "../quadinterpolator.hpp"
17 #include "../../general/forall.hpp"
18 #include "../../linalg/kernels.hpp"
19 
20 namespace mfem
21 {
22 
24  const FiniteElementSpace &fes)
25 {
26  MFEM_VERIFY(PA.enabled, "PA extension setup has not been done!");
27  MFEM_VERIFY(PA.fes == &fes, "");
28  // TODO: we need a more robust way to check that the 'fes' used when
29  // AssemblePA() was called has not been modified or completely destroyed and
30  // a new object created at the same address.
31 
32  if (PA.Jtr_needs_update || targetC->UsesPhysicalCoordinates())
33  {
35  PA.Jtr_debug_grad = true;
36  }
37 
38  if (PA.dim == 2)
39  {
41  if (lim_coeff) { AssembleGradPA_C0_2D(xe); }
42  }
43 
44  if (PA.dim == 3)
45  {
47  if (lim_coeff) { AssembleGradPA_C0_3D(xe); }
48  }
49 }
50 
52 {
53  const MemoryType mt = (pa_mt == MemoryType::DEFAULT) ?
55  // Return immediately if limiting is not enabled
56  if (lim_coeff == nullptr) { return; }
57  MFEM_VERIFY(lim_nodes0, "internal error");
58 
59  MFEM_VERIFY(PA.enabled, "AssemblePA_Limiting but PA is not enabled!");
60  MFEM_VERIFY(lim_func, "No TMOP_LimiterFunction specification!")
61  MFEM_VERIFY(dynamic_cast<TMOP_QuadraticLimiter*>(lim_func),
62  "Only TMOP_QuadraticLimiter is supported");
63 
64  const FiniteElementSpace *fes = PA.fes;
65  const int NE = PA.ne;
66  if (NE == 0) { return; } // Quick return for empty processors
67  const IntegrationRule &ir = *PA.ir;
68 
70 
71  // H0 for lim_coeff, (dim x dim) Q-vector
72  PA.H0.UseDevice(true);
73  PA.H0.SetSize(PA.dim * PA.dim * PA.nq * NE, mt);
74 
75  // lim_coeff -> PA.C0 (Q-vector)
76  PA.C0.UseDevice(true);
77  if (ConstantCoefficient* cQ =
78  dynamic_cast<ConstantCoefficient*>(lim_coeff))
79  {
80  PA.C0.SetSize(1, Device::GetMemoryType());
81  PA.C0.HostWrite();
82  PA.C0(0) = cQ->constant;
83  }
84  else
85  {
86  PA.C0.SetSize(PA.nq * PA.ne, Device::GetMemoryType());
87  auto C0 = Reshape(PA.C0.HostWrite(), PA.nq, PA.ne);
88  for (int e = 0; e < NE; ++e)
89  {
91  for (int q = 0; q < ir.GetNPoints(); ++q)
92  {
93  C0(q,e) = lim_coeff->Eval(T, ir.IntPoint(q));
94  }
95  }
96  }
97 
98  // lim_nodes0 -> PA.X0 (E-vector)
99  MFEM_VERIFY(lim_nodes0->FESpace() == fes, "");
100  const Operator *n0_R = fes->GetElementRestriction(ordering);
101  PA.X0.SetSize(n0_R->Height(), Device::GetMemoryType());
102  PA.X0.UseDevice(true);
103  n0_R->Mult(*lim_nodes0, PA.X0);
104 
105  // Limiting distances: lim_dist -> PA.LD (E-vector)
106  // TODO: remove the hack for the case lim_dist == NULL.
107  const FiniteElementSpace *limfes = (lim_dist) ? lim_dist->FESpace() : fes;
108  const FiniteElement &lim_fe = *limfes->GetFE(0);
109  PA.maps_lim = &lim_fe.GetDofToQuad(ir, DofToQuad::TENSOR);
110  PA.LD.SetSize(NE*lim_fe.GetDof(), Device::GetMemoryType());
111  PA.LD.UseDevice(true);
112  if (lim_dist)
113  {
114  const Operator *ld_R = limfes->GetElementRestriction(ordering);
115  ld_R->Mult(*lim_dist, PA.LD);
116  }
117  else
118  {
119  PA.LD = 1.0;
120  }
121 }
122 
124  const IntegrationRule &ir,
125  const Vector &xe,
126  DenseTensor &Jtr) const
127 {
128  MFEM_VERIFY(Jtr.SizeI() == Jtr.SizeJ() && Jtr.SizeI() > 1, "");
129  const int dim = Jtr.SizeI();
130  bool done = false;
131  if (dim == 2) { done = ComputeAllElementTargets<2>(fes, ir, xe, Jtr); }
132  if (dim == 3) { done = ComputeAllElementTargets<3>(fes, ir, xe, Jtr); }
133 
134  if (!done) { ComputeAllElementTargets_Fallback(fes, ir, xe, Jtr); }
135 }
136 
138  const IntegrationRule &ir,
139  const Vector &xe,
140  DenseTensor &Jtr) const
141 {
142  ComputeAllElementTargets_Fallback(fes, ir, xe, Jtr);
143 }
144 
145 
146 // Code paths leading to ComputeElementTargets:
147 // - GetElementEnergy(elfun) which is done through GetGridFunctionEnergyPA(x)
148 // - AssembleElementVectorExact(elfun)
149 // - AssembleElementGradExact(elfun)
150 // - EnableNormalization(x) -> ComputeNormalizationEnergies(x)
151 // - (AssembleElementVectorFD(elfun))
152 // - (AssembleElementGradFD(elfun))
153 // ============================================================================
154 // - TargetConstructor():
155 // - IDEAL_SHAPE_UNIT_SIZE: Wideal
156 // - IDEAL_SHAPE_EQUAL_SIZE: α * Wideal
157 // - IDEAL_SHAPE_GIVEN_SIZE: β * Wideal
158 // - GIVEN_SHAPE_AND_SIZE: β * Wideal
159 // - AnalyticAdaptTC(elfun):
160 // - GIVEN_FULL: matrix_tspec->Eval(Jtr(elfun))
161 // - DiscreteAdaptTC():
162 // - IDEAL_SHAPE_GIVEN_SIZE: size^{1.0/dim} * Jtr(i) (size)
163 // - GIVEN_SHAPE_AND_SIZE: Jtr(i) *= D_rho (ratio)
164 // Jtr(i) *= Q_phi (skew)
165 // Jtr(i) *= R_theta (orientation)
167 {
168  PA.Jtr_needs_update = false;
169  PA.Jtr_debug_grad = false;
170  const FiniteElementSpace *fes = PA.fes;
171  if (PA.ne == 0) { return; } // Quick return for empty processors
172  const IntegrationRule &ir = *PA.ir;
173 
174  // Compute PA.Jtr for all elements
175  targetC->ComputeAllElementTargets(*fes, ir, xe, PA.Jtr);
176 }
177 
179 {
180  const MemoryType mt = (pa_mt == MemoryType::DEFAULT) ?
182  PA.enabled = true;
183  PA.fes = &fes;
184  Mesh *mesh = fes.GetMesh();
185  const int ne = PA.ne = mesh->GetNE();
186  if (ne == 0) { return; } // Quick return for empty processors
187  const int dim = PA.dim = mesh->Dimension();
188  MFEM_VERIFY(PA.dim == 2 || PA.dim == 3, "Not yet implemented!");
189  MFEM_VERIFY(mesh->GetNumGeometries(dim) <= 1,
190  "mixed meshes are not supported");
191  MFEM_VERIFY(!fes.IsVariableOrder(), "variable orders are not supported");
192  const FiniteElement &fe = *fes.GetFE(0);
193  PA.ir = &EnergyIntegrationRule(fe);
194  const IntegrationRule &ir = *PA.ir;
195  MFEM_VERIFY(fes.GetOrdering() == Ordering::byNODES,
196  "PA Only supports Ordering::byNODES!");
197 
198  const int nq = PA.nq = ir.GetNPoints();
199  const DofToQuad::Mode mode = DofToQuad::TENSOR;
200  PA.maps = &fe.GetDofToQuad(ir, mode);
202 
203  // Energy vector, scalar Q-vector
204  PA.E.UseDevice(true);
205  PA.E.SetSize(ne*nq, Device::GetDeviceMemoryType());
206 
207  // H for Grad, (dim x dim) Q-vector
208  PA.H.UseDevice(true);
209  PA.H.SetSize(dim*dim * dim*dim * nq*ne, mt);
210 
211  // Scalar Q-vector of '1', used to compute sums via dot product
212  PA.O.SetSize(ne*nq, Device::GetDeviceMemoryType());
213  PA.O = 1.0;
214 
215  // Setup ref->target Jacobians, PA.Jtr, (dim x dim) Q-vector, DenseTensor
216  PA.Jtr.SetSize(dim, dim, PA.ne*PA.nq, mt);
217  PA.Jtr_needs_update = true;
218  PA.Jtr_debug_grad = false;
219 
220  // Limiting: lim_coeff -> PA.C0, lim_nodes0 -> PA.X0, lim_dist -> PA.LD, PA.H0
221  if (lim_coeff) { AssemblePA_Limiting(); }
222 }
223 
225 {
226  // This method must be called after AssembleGradPA().
227 
228  MFEM_VERIFY(PA.Jtr_needs_update == false, "");
229 
231  {
232  MFEM_VERIFY(PA.Jtr_debug_grad == true, "AssembleGradPA() was not called"
233  " or Jtr was overwritten by another method!");
234  }
235 
236  if (PA.dim == 2)
237  {
239  if (lim_coeff) { AssembleDiagonalPA_C0_2D(de); }
240  }
241 
242  if (PA.dim == 3)
243  {
245  if (lim_coeff) { AssembleDiagonalPA_C0_3D(de); }
246  }
247 }
248 
249 void TMOP_Integrator::AddMultPA(const Vector &xe, Vector &ye) const
250 {
251  // This method must be called after AssemblePA().
252 
253  if (PA.Jtr_needs_update || targetC->UsesPhysicalCoordinates())
254  {
256  }
257 
258  if (PA.dim == 2)
259  {
260  AddMultPA_2D(xe,ye);
261  if (lim_coeff) { AddMultPA_C0_2D(xe,ye); }
262  }
263 
264  if (PA.dim == 3)
265  {
266  AddMultPA_3D(xe,ye);
267  if (lim_coeff) { AddMultPA_C0_3D(xe,ye); }
268  }
269 }
270 
271 void TMOP_Integrator::AddMultGradPA(const Vector &re, Vector &ce) const
272 {
273  // This method must be called after AssembleGradPA().
274 
275  MFEM_VERIFY(PA.Jtr_needs_update == false, "");
276 
278  {
279  MFEM_VERIFY(PA.Jtr_debug_grad == true, "AssembleGradPA() was not called or "
280  "Jtr was overwritten by another method!");
281  }
282 
283  if (PA.dim == 2)
284  {
285  AddMultGradPA_2D(re,ce);
286  if (lim_coeff) { AddMultGradPA_C0_2D(re,ce); }
287  }
288 
289  if (PA.dim == 3)
290  {
291  AddMultGradPA_3D(re,ce);
292  if (lim_coeff) { AddMultGradPA_C0_3D(re,ce); }
293  }
294 }
295 
297 {
298  // This method must be called after AssemblePA().
299 
300  double energy = 0.0;
301 
302  if (PA.Jtr_needs_update || targetC->UsesPhysicalCoordinates())
303  {
305  }
306 
307  if (PA.dim == 2)
308  {
309  energy = GetLocalStateEnergyPA_2D(xe);
310  if (lim_coeff) { energy += GetLocalStateEnergyPA_C0_2D(xe); }
311  }
312 
313  if (PA.dim == 3)
314  {
315  energy = GetLocalStateEnergyPA_3D(xe);
316  if (lim_coeff) { energy += GetLocalStateEnergyPA_C0_3D(xe); }
317  }
318 
319  return energy;
320 }
321 
322 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void AddMultGradPA(const Vector &, Vector &) const
Method for partially assembled gradient action.
Definition: tmop_pa.cpp:271
void AssembleGradPA_2D(const Vector &) const
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:599
void AddMultGradPA_C0_3D(const Vector &, Vector &) const
void AddMultPA_C0_2D(const Vector &, Vector &) const
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
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
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition: mesh.cpp:840
double GetLocalStateEnergyPA_3D(const Vector &) const
Definition: tmop_pa_w3.cpp:155
const GridFunction * lim_dist
Definition: tmop.hpp:1590
struct mfem::TMOP_Integrator::@23 PA
virtual double GetLocalStateEnergyPA(const Vector &) const
Compute the local (to the MPI rank) energy with partial assembly.
Definition: tmop_pa.cpp:296
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5908
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...
Definition: tmop_pa.cpp:137
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).
const IntegrationRule * ir
Definition: tmop.hpp:1676
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
virtual void AssembleGradPA(const Vector &, const FiniteElementSpace &)
Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at th...
Definition: tmop_pa.cpp:23
void AddMultGradPA_C0_2D(const Vector &, Vector &) const
void AssembleDiagonalPA_C0_2D(Vector &) const
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
void AssembleDiagonalPA_3D(Vector &) const
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:277
void AssembleDiagonalPA_2D(Vector &) const
const TargetConstructor * targetC
Definition: tmop.hpp:1573
const GridFunction * lim_nodes0
Definition: tmop.hpp:1587
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
Definition: tmop_pa.cpp:178
double GetLocalStateEnergyPA_2D(const Vector &) const
Definition: tmop_pa_w2.cpp:145
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:668
int Dimension() const
Definition: mesh.hpp:1006
void AssemblePA_Limiting()
Definition: tmop_pa.cpp:51
int SizeI() const
Definition: densemat.hpp:999
Coefficient * lim_coeff
Definition: tmop.hpp:1588
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:647
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
bool UsesPhysicalCoordinates() const
Return true if the methods ComputeElementTargets(), ComputeAllElementTargets(), and ComputeElementTar...
Definition: tmop.hpp:1259
double GetLocalStateEnergyPA_C0_2D(const Vector &) const
void AddMultGradPA_2D(const Vector &, Vector &) const
void AddMultPA_3D(const Vector &, Vector &) const
Definition: tmop_pa_p3.cpp:192
void AddMultPA_2D(const Vector &, Vector &) const
Definition: tmop_pa_p2.cpp:167
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
int SizeJ() const
Definition: densemat.hpp:1000
void AssembleGradPA_C0_2D(const Vector &) const
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition: fe_base.hpp:149
void AssembleDiagonalPA_C0_3D(Vector &) const
const FiniteElementSpace * fes
Definition: tmop.hpp:1675
const IntegrationRule & EnergyIntegrationRule(const FiniteElement &el) const
Definition: tmop.hpp:1737
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
int dim
Definition: ex24.cpp:53
void AssembleGradPA_3D(const Vector &) const
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 AddMultGradPA_3D(const Vector &, Vector &) const
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.
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1259
virtual void AssembleGradDiagonalPA(Vector &) const
Method for computing the diagonal of the gradient with partial assembly.
Definition: tmop_pa.cpp:224
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void ComputeAllElementTargets(const Vector &xe=Vector()) const
Definition: tmop_pa.cpp:166
Vector data type.
Definition: vector.hpp:60
double GetLocalStateEnergyPA_C0_3D(const Vector &) const
void AddMultPA_C0_3D(const Vector &, Vector &) const
Abstract operator.
Definition: operator.hpp:24
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:953
TMOP_LimiterFunction * lim_func
Definition: tmop.hpp:1592
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
Definition: tmop_pa.cpp:249
void AssembleGradPA_C0_3D(const Vector &) const
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
bool ComputeAllElementTargets(const FiniteElementSpace &fes, const IntegrationRule &ir, const Vector &xe, DenseTensor &Jtr) const
Definition: tmop_pa.cpp:123