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 
22 {
24 
25  MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH");
27  "patch_polynomial_order",
28  orders,
29  "Polynomial order used in least squares fitting of material property "
30  "over the local patch of elements connected to a given node");
31 
32  params.addRelationshipManager("ElementSideNeighborLayers",
34  [](const InputParameters &, InputParameters & rm_params)
35  { rm_params.set<unsigned short>("layers") = 2; });
36 
37  params.addParamNamesToGroup("patch_polynomial_order", "Advanced");
38 
39  return params;
40 }
41 
43  : ElementUserObject(parameters),
44  _qp(0),
45  _patch_polynomial_order(
46  static_cast<unsigned int>(getParam<MooseEnum>("patch_polynomial_order"))),
47  _multi_index(MathUtils::multiIndex(_mesh.dimension(), _patch_polynomial_order)),
48  _q(_multi_index.size())
49 {
50 }
51 
52 Real
54  const std::vector<dof_id_type> & elem_ids) const
55 {
56  // Before we go, check if we have enough sample points for solving the least square fitting
57  if (_q_point.size() * elem_ids.size() < _q)
58  mooseError("There are not enough sample points to recover the nodal value, try reducing the "
59  "polynomial order or using a higher-order quadrature scheme.");
60 
61  // Assemble the least squares problem over the patch
62  RealEigenMatrix A = RealEigenMatrix::Zero(_q, _q);
63  RealEigenVector b = RealEigenVector::Zero(_q);
64  for (auto elem_id : elem_ids)
65  {
66  A += libmesh_map_find(_Ae, elem_id);
67  b += libmesh_map_find(_be, elem_id);
68  }
69 
70  // Solve the least squares fitting
71  RealEigenVector coef = A.completeOrthogonalDecomposition().solve(b);
72 
73  // Compute the fitted nodal value
75  return p.dot(coef);
76 }
77 
80 {
83  for (unsigned int r = 0; r < _multi_index.size(); r++)
84  {
85  polynomial = 1.0;
86  mooseAssert(_multi_index[r].size() == _mesh.dimension(), "Wrong multi-index size.");
87  for (unsigned int c = 0; c < _multi_index[r].size(); c++)
88  for (unsigned int p = 0; p < _multi_index[r][c]; p++)
89  polynomial *= q_point(c);
90  p(r) = polynomial;
91  }
92  return p;
93 }
94 
95 void
97 {
98  _Ae.clear();
99  _be.clear();
100 }
101 
102 void
104 {
105  RealEigenMatrix Ae = RealEigenMatrix::Zero(_q, _q);
106  RealEigenVector be = RealEigenVector::Zero(_q);
107  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
108  {
110  Ae += p * p.transpose();
111  be += computeValue() * p;
112  }
113 
114  dof_id_type elem_id = _current_elem->id();
115  _Ae[elem_id] = Ae;
116  _be[elem_id] = be;
117 }
118 
119 void
121 {
122  const auto & npr = static_cast<const NodalPatchRecoveryBase &>(uo);
123  _Ae.insert(npr._Ae.begin(), npr._Ae.end());
124  _be.insert(npr._be.begin(), npr._be.end());
125 }
126 
127 void
129 {
130  // When calling nodalPatchRecovery, we may need to know _Ae and _be on algebraically ghosted
131  // elements. However, this userobject is only run on local elements, so we need to query those
132  // information from other processors in this finalize() method.
133 
134  // Populate algebraically ghosted elements to query
135  std::unordered_map<processor_id_type, std::vector<dof_id_type>> query_ids;
136  const ConstElemRange evaluable_elem_range = _fe_problem.getEvaluableElementRange();
137  for (auto elem : evaluable_elem_range)
138  if (elem->processor_id() != processor_id())
139  query_ids[elem->processor_id()].push_back(elem->id());
140 
141  typedef std::pair<RealEigenMatrix, RealEigenVector> AbPair;
142 
143  // Answer queries received from other processors
144  auto gather_data = [this](const processor_id_type /*pid*/,
145  const std::vector<dof_id_type> & elem_ids,
146  std::vector<AbPair> & ab_pairs)
147  {
148  for (const auto i : index_range(elem_ids))
149  {
150  const auto elem_id = elem_ids[i];
151  ab_pairs.emplace_back(libmesh_map_find(_Ae, elem_id), libmesh_map_find(_be, elem_id));
152  }
153  };
154 
155  // Gather answers received from other processors
156  auto act_on_data = [this](const processor_id_type /*pid*/,
157  const std::vector<dof_id_type> & elem_ids,
158  const std::vector<AbPair> & ab_pairs)
159  {
160  for (const auto i : index_range(elem_ids))
161  {
162  const auto elem_id = elem_ids[i];
163  const auto & [Ae, be] = ab_pairs[i];
164  _Ae[elem_id] = Ae;
165  _be[elem_id] = be;
166  }
167  };
168 
169  libMesh::Parallel::pull_parallel_vector_data<AbPair>(
170  _communicator, query_ids, gather_data, act_on_data, 0);
171 }
std::map< dof_id_type, RealEigenVector > _be
The element-level b vector.
static InputParameters validParams()
virtual void initialize() override
Called before execute() is ever called so that data can be cleared.
const MooseArray< Point > & _q_point
static InputParameters validParams()
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.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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...
virtual void execute() override
Execute method.
uint8_t processor_id_type
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
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
The multi-index table.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
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
virtual void finalize() override
Finalize.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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
virtual 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}.
processor_id_type processor_id() const
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...