https://mooseframework.inl.gov
NodalPatchRecoveryBase.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "NodalPatchRecoveryBase.h"
11 #include "MathUtils.h"
12 
13 #include <Eigen/Dense>
14 
15 // TIMPI includes
16 #include "timpi/communicator.h"
17 #include "timpi/parallel_sync.h"
18 #include "libmesh/parallel_eigen.h"
19 #include <iomanip>
20 
21 // Remove duplicates and sort entries from a vector of element IDs.
22 // This is necessary because user input may contain repeated element IDs,
23 // which is problematic when using these IDs as keys.
24 // In Patch Recovery, we also want to avoid duplicate elements contributing
25 // multiple times to the Ae and be.
26 static std::vector<dof_id_type>
27 removeDuplicateEntries(const std::vector<dof_id_type> & ids)
28 {
29  std::vector<dof_id_type> key = ids;
30  std::sort(key.begin(), key.end());
31  key.erase(std::unique(key.begin(), key.end()), key.end());
32  return key;
33 }
34 
37 {
39 
40  MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH");
42  "patch_polynomial_order",
43  orders,
44  "Polynomial order used in least squares fitting of material property "
45  "over the local patch of elements connected to a given node");
46 
47  params.addRelationshipManager("ElementSideNeighborLayers",
49  [](const InputParameters &, InputParameters & rm_params)
50  { rm_params.set<unsigned short>("layers") = 2; });
51 
52  params.addParamNamesToGroup("patch_polynomial_order", "Advanced");
53 
54  return params;
55 }
56 
58  : ElementUserObject(parameters),
59  _qp(0),
60  _patch_polynomial_order(
61  static_cast<unsigned int>(getParam<MooseEnum>("patch_polynomial_order"))),
62  _multi_index(MathUtils::multiIndex(_mesh.dimension(), _patch_polynomial_order)),
63  _q(_multi_index.size()),
64  _distributed_mesh(_mesh.isDistributedMesh()),
65  _proc_ids(n_processors())
66 {
67  std::iota(_proc_ids.begin(), _proc_ids.end(), 0);
68 }
69 
70 Real
72  const std::vector<dof_id_type> & elem_ids) const
73 {
74  const RealEigenVector coef = getCoefficients(elem_ids); // const version
75  // Compute the fitted nodal value
77  return p.dot(coef);
78 }
79 
80 const RealEigenVector
81 NodalPatchRecoveryBase::getCoefficients(const std::vector<dof_id_type> & elem_ids) const
82 {
83  auto elem_ids_reduced = removeDuplicateEntries(elem_ids);
84 
85  RealEigenVector coef = RealEigenVector::Zero(_q);
86  // Before we go, check if we have enough sample points for solving the least square fitting
87  if (_q_point.size() * elem_ids_reduced.size() < _q)
88  mooseError("There are not enough sample points to recover the nodal value, try reducing the "
89  "polynomial order or using a higher-order quadrature scheme.");
90 
91  // Assemble the least squares problem over the patch
92  RealEigenMatrix A = RealEigenMatrix::Zero(_q, _q);
93  RealEigenVector b = RealEigenVector::Zero(_q);
94  for (auto elem_id : elem_ids_reduced)
95  {
96  const auto elem = _mesh.elemPtr(elem_id);
97  if (elem /*prevent segmentation fault in distributed mesh*/)
98  if (!hasBlocks(elem->subdomain_id()))
99  mooseError("Element with id = ",
100  elem_id,
101  " is not in the block. "
102  "Please use nodalPatchRecovery with elements in the block only.");
103 
104  if (_Ae.find(elem_id) == _Ae.end())
105  mooseError("Missing entry for elem_id = ", elem_id, " in _Ae.");
106  if (_be.find(elem_id) == _be.end())
107  mooseError("Missing entry for elem_id = ", elem_id, " in _be.");
108 
109  A += libmesh_map_find(_Ae, elem_id);
110  b += libmesh_map_find(_be, elem_id);
111  }
112 
113  // Solve the least squares fitting
114  coef = A.completeOrthogonalDecomposition().solve(b);
115 
116  return coef;
117 }
118 
119 const RealEigenVector
120 NodalPatchRecoveryBase::getCachedCoefficients(const std::vector<dof_id_type> & elem_ids)
121 {
122  // Check cache
123  auto key = removeDuplicateEntries(elem_ids);
124 
125  if (key == _cached_elem_ids)
126  return _cached_coef;
127  else
128  {
129  const auto coef = getCoefficients(key); // const version
130 
131  _cached_elem_ids = key; // Update the cached element IDs
132  _cached_coef = coef; // Update the cached coefficients
133 
134  return coef;
135  }
136 }
137 
140 {
141  RealEigenVector p(_q);
143  for (unsigned int r = 0; r < _multi_index.size(); r++)
144  {
145  polynomial = 1.0;
146  mooseAssert(_multi_index[r].size() == _mesh.dimension(), "Wrong multi-index size.");
147  for (unsigned int c = 0; c < _multi_index[r].size(); c++)
148  for (unsigned int p = 0; p < _multi_index[r][c]; p++)
149  polynomial *= q_point(c);
150  p(r) = polynomial;
151  }
152  return p;
153 }
154 
155 void
157 {
158  // Clear cached data to ensure coefficients are correctly recomputed in the next patch recovery
159  // iteration. _Ae and _be must also be reset, as the associated patch elements may differ in the
160  // upcoming iteration.
161 
162  _cached_elem_ids.clear();
163  _cached_coef = RealEigenVector::Zero(_q);
164  _Ae.clear();
165  _be.clear();
166 }
167 
168 void
170 {
171  RealEigenMatrix Ae = RealEigenMatrix::Zero(_q, _q);
172  RealEigenVector be = RealEigenVector::Zero(_q);
173  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
174  {
176  Ae += p * p.transpose();
177  be += computeValue() * p;
178  }
179 
180  dof_id_type elem_id = _current_elem->id();
181 
182  _Ae[elem_id] = Ae;
183  _be[elem_id] = be;
184 }
185 
186 void
188 {
189  const auto & npr = static_cast<const NodalPatchRecoveryBase &>(uo);
190  _Ae.insert(npr._Ae.begin(), npr._Ae.end());
191  _be.insert(npr._be.begin(), npr._be.end());
192 }
193 
194 void
196 {
197  // When calling nodalPatchRecovery, we may need to know _Ae and _be on algebraically ghosted
198  // elements. However, this userobject is only run on local elements, so we need to query those
199  // information from other processors in this finalize() method.
200  sync();
201 }
202 
203 std::unordered_map<processor_id_type, std::vector<dof_id_type>>
204 NodalPatchRecoveryBase::gatherRequestList(const std::vector<dof_id_type> & specific_elems)
205 {
206  std::unordered_map<processor_id_type, std::vector<dof_id_type>> query_ids;
207 
208  typedef std::pair<processor_id_type, dof_id_type> PidElemPair;
209  std::unordered_map<processor_id_type, std::vector<PidElemPair>> push_data;
210 
211  for (const auto & entry : specific_elems)
212  {
213  const auto * elem = _mesh.elemPtr(entry);
214  if (_distributed_mesh)
215  {
216  if (!elem)
217  continue; // Prevent segmentation fault in distributed mesh
218 
219  if (hasBlocks(elem->subdomain_id()))
220  for (processor_id_type pid = 0; pid < n_processors(); ++pid)
221  if (pid != processor_id())
222  push_data[pid].push_back(std::make_pair(elem->processor_id(), elem->id()));
223  }
224  else
225  // Add to query_ids if this element's data is on a different processor
226  addToQuery(elem, query_ids);
227  }
228 
229  if (_distributed_mesh)
230  {
231  auto push_receiver =
232  [&](const processor_id_type, const std::vector<PidElemPair> & received_data)
233  {
234  for (const auto & [pid, id] : received_data)
235  query_ids[pid].push_back(id);
236  };
237 
238  Parallel::push_parallel_vector_data(_mesh.comm(), push_data, push_receiver);
239  }
240 
241  return query_ids;
242 }
243 
244 std::unordered_map<processor_id_type, std::vector<dof_id_type>>
246 {
247  std::unordered_map<processor_id_type, std::vector<dof_id_type>> query_ids;
248 
249  typedef std::pair<processor_id_type, dof_id_type> PidElemPair;
250  std::unordered_map<processor_id_type, std::vector<PidElemPair>> push_data;
251 
252  for (const auto & elem : _fe_problem.getEvaluableElementRange())
253  addToQuery(elem, query_ids);
254 
255  return query_ids;
256 }
257 
258 void
260 {
261  const auto query_ids = gatherRequestList();
262  syncHelper(query_ids);
263 }
264 
265 void
266 NodalPatchRecoveryBase::sync(const std::vector<dof_id_type> & specific_elems)
267 {
268  const auto query_ids = gatherRequestList(specific_elems);
269  syncHelper(query_ids);
270 }
271 
272 void
274  const std::unordered_map<processor_id_type, std::vector<dof_id_type>> & query_ids)
275 {
276  typedef std::pair<RealEigenMatrix, RealEigenVector> AbPair;
277 
278  // Answer queries received from other processors
279  auto gather_data = [this](const processor_id_type /*pid*/,
280  const std::vector<dof_id_type> & elem_ids,
281  std::vector<AbPair> & ab_pairs)
282  {
283  for (const auto & elem_id : elem_ids)
284  ab_pairs.emplace_back(libmesh_map_find(_Ae, elem_id), libmesh_map_find(_be, elem_id));
285  };
286 
287  // Gather answers received from other processors
288  auto act_on_data = [this](const processor_id_type /*pid*/,
289  const std::vector<dof_id_type> & elem_ids,
290  const std::vector<AbPair> & ab_pairs)
291  {
292  for (const auto i : index_range(elem_ids))
293  {
294  const auto elem_id = elem_ids[i];
295  const auto & [Ae, be] = ab_pairs[i];
296  _Ae[elem_id] = Ae;
297  _be[elem_id] = be;
298  }
299  };
300 
301  // the send and receive are called inside the pull_parallel_vector_data function
302  libMesh::Parallel::pull_parallel_vector_data<AbPair>(
303  _communicator, query_ids, gather_data, act_on_data, 0);
304 }
305 
306 void
308  const libMesh::Elem * elem,
309  std::unordered_map<processor_id_type, std::vector<dof_id_type>> & query_ids) const
310 {
311  if (hasBlocks(elem->subdomain_id()) && elem->processor_id() != processor_id())
312  query_ids[elem->processor_id()].push_back(elem->id());
313 }
void syncHelper(const std::unordered_map< processor_id_type, std::vector< dof_id_type >> &query_ids)
Helper function to perform the actual communication of _Ae and _be.
std::vector< int > _proc_ids
The processor IDs vector in the running.
std::map< dof_id_type, RealEigenVector > _be
The element-level b vector.
static InputParameters validParams()
void addToQuery(const libMesh::Elem *elem, std::unordered_map< processor_id_type, std::vector< dof_id_type >> &query_ids) const
Adds an element to the map provided in query_ids if it belongs to a different processor.
void initialize() override
Called before execute() is ever called so that data can be cleared.
const MooseArray< Point > & _q_point
virtual Elem * elemPtr(const dof_id_type i)
Definition: MooseMesh.C:3113
static InputParameters validParams()
std::vector< dof_id_type > _cached_elem_ids
Cache for least-squares coefficients used in nodal patch recovery.
void sync()
Synchronizes local matrices and vectors (_Ae, _be) across processors.
std::map< dof_id_type, RealEigenMatrix > _Ae
The element-level A matrix.
RealEigenVector evaluateBasisFunctions(const Point &q_point) const
Compute the P vector at a given point i.e.
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
static std::vector< dof_id_type > removeDuplicateEntries(const std::vector< dof_id_type > &ids)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const Parallel::Communicator & comm() const
const Parallel::Communicator & _communicator
void addRelationshipManager(const std::string &name, Moose::RelationshipManagerType rm_type, Moose::RelationshipManagerInputParameterCallback input_parameter_callback=nullptr)
Tells MOOSE about a RelationshipManager that this object needs.
virtual Real nodalPatchRecovery(const Point &p, const std::vector< dof_id_type > &elem_ids) const
Solve the least-squares problem.
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
void execute() override
Execute method.
uint8_t processor_id_type
processor_id_type n_processors() const
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
dof_id_type id() const
virtual unsigned int dimension() const
Returns MeshBase::mesh_dimension(), (not MeshBase::spatial_dimension()!) of the underlying libMesh me...
Definition: MooseMesh.C:2928
NodalPatchRecoveryBase(const InputParameters &parameters)
const std::vector< std::vector< unsigned int > > _multi_index
Multi-index table for a polynomial basis.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
bool _distributed_mesh
Whether the mesh is distributed.
R polynomial(const C &c, const T x)
Evaluate a polynomial with the coefficients c at x.
Definition: MathUtils.h:272
virtual Real computeValue()=0
Compute the quantity to recover using nodal patch recovery.
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
void finalize() override
Finalize.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
const QBase *const & _qrule
const Elem *const & _current_elem
The current element pointer (available during execute())
FEProblemBase & _fe_problem
Reference to the FEProblemBase for this user object.
Definition: UserObject.h:211
std::vector< std::vector< unsigned int > > multiIndex(unsigned int dim, unsigned int order)
generate a complete multi index table for given dimension and order i.e.
Definition: MathUtils.C:186
void threadJoin(const UserObject &) override
Must override.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:267
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
const libMesh::ConstElemRange & getEvaluableElementRange()
In general, {evaluable elements} >= {local elements} U {algebraic ghosting elements}.
const RealEigenVector getCoefficients(const std::vector< dof_id_type > &elem_ids) const
Compute coefficients without reading or writing cached values The coefficients returned by this funct...
std::unordered_map< processor_id_type, std::vector< dof_id_type > > gatherRequestList()
Builds a query map of element IDs that require data from other processors.
bool hasBlocks(const SubdomainName &name) const
Test if the supplied block name is valid for this object.
processor_id_type processor_id() const
processor_id_type processor_id() const
const RealEigenVector getCachedCoefficients(const std::vector< dof_id_type > &elem_ids)
Compute coefficients, using cached values if available, and store any newly computed coefficients in ...
const unsigned int _q
Number of basis functions.
void ErrorVector unsigned int
auto index_range(const T &sizable)
Base class for user-specific data.
Definition: UserObject.h:40
uint8_t dof_id_type
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...