https://mooseframework.inl.gov
ProjectedStatefulMaterialNodalPatchRecovery.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 
11 #include "ElementUserObject.h"
12 #include "MaterialBase.h"
13 #include "MathUtils.h"
14 #include "Assembly.h"
15 
16 // TIMPI includes
17 #include "timpi/communicator.h"
18 #include "timpi/parallel_sync.h"
19 #include "libmesh/parallel_eigen.h"
20 
29 
32 {
34  params.addClassDescription(
35  "Prepare patches for use in nodal patch recovery based on a material property for material "
36  "property states projected onto nodal variables.");
37  params.addRelationshipManager("ElementSideNeighborLayers",
39  [](const InputParameters &, InputParameters & rm_params)
40  { rm_params.set<unsigned short>("layers") = 2; });
41 
42  return params;
43 }
44 
45 template <typename T, bool is_ad>
48 {
50  MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH");
51  params.addRequiredParam<MaterialPropertyName>(
52  "property", "Name of the material property to perform nodal patch recovery on");
54  "patch_polynomial_order",
55  orders,
56  "Polynomial order used in least squares fitting of material property "
57  "over the local patch of elements connected to a given node");
58 
59  params.addRelationshipManager("ElementSideNeighborLayers",
61  [](const InputParameters &, InputParameters & rm_params)
62  { rm_params.set<unsigned short>("layers") = 2; });
63 
64  params.addParamNamesToGroup("patch_polynomial_order", "Advanced");
65 
66  return params;
67 }
68 
69 template <typename T, bool is_ad>
73  _qp(0),
74  _n_components(Moose::SerialAccess<T>::size()),
75  _patch_polynomial_order(
76  static_cast<unsigned int>(getParam<MooseEnum>("patch_polynomial_order"))),
77  _multi_index(MathUtils::multiIndex(_mesh.dimension(), _patch_polynomial_order)),
78  _q(_multi_index.size()),
79  _prop(getGenericMaterialProperty<T, is_ad>("property")),
80  _current_subdomain_id(_assembly.currentSubdomainID())
81 {
82 }
83 
84 template <typename T, bool is_ad>
85 void
87 {
88  // get all material classes that provide properties for this object
89  _required_materials = buildRequiredMaterials();
90 }
91 
92 template <typename T, bool is_ad>
93 Real
95  const Point & x, const std::vector<dof_id_type> & elem_ids, std::size_t component) const
96 {
97  // Before we go, check if we have enough sample points for solving the least square fitting
98  if (_q_point.size() * elem_ids.size() < _q)
99  mooseError("There are not enough sample points to recover the nodal value, try reducing the "
100  "polynomial order or using a higher-order quadrature scheme.");
101 
102  // Assemble the least squares problem over the patch
103  RealEigenMatrix A = RealEigenMatrix::Zero(_q, _q);
104  RealEigenVector b = RealEigenVector::Zero(_q);
105  for (const auto elem_id : elem_ids)
106  {
107  const auto abs = libmesh_map_find(_abs, elem_id);
108  A += abs.first;
109  b += abs.second[component];
110  }
111 
112  // Solve the least squares fitting
113  const RealEigenVector coef = A.completeOrthogonalDecomposition().solve(b);
114 
115  // Compute the fitted nodal value
116  const RealEigenVector p = evaluateBasisFunctions(x);
117  return p.dot(coef);
118 }
119 
120 template <typename T, bool is_ad>
123  const Point & q_point) const
124 {
125  RealEigenVector p(_q);
127  for (const auto r : index_range(_multi_index))
128  {
129  polynomial = 1.0;
130  mooseAssert(_multi_index[r].size() == _mesh.dimension(), "Wrong multi-index size.");
131  for (const auto c : index_range(_multi_index[r]))
132  polynomial *= MathUtils::pow(q_point(c), _multi_index[r][c]);
133  p(r) = polynomial;
134  }
135  return p;
136 }
137 
138 template <typename T, bool is_ad>
139 void
141 {
142  _abs.clear();
143 }
144 
145 template <typename T, bool is_ad>
146 void
148 {
149  if (_t_step == 0)
150  for (const auto & mat : _required_materials[_current_subdomain_id])
151  mat->initStatefulProperties(_qrule->size());
152 
153  std::vector<RealEigenVector> bs(_n_components, RealEigenVector::Zero(_q));
154  RealEigenMatrix Ae = RealEigenMatrix::Zero(_q, _q);
155 
156  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
157  {
158  RealEigenVector p = evaluateBasisFunctions(_q_point[_qp]);
159  Ae += p * p.transpose();
160 
161  std::size_t index = 0;
162  for (const auto & v : Moose::serialAccess(_prop[_qp]))
163  bs[index++] += MetaPhysicL::raw_value(v) * p;
164  }
165 
166  const dof_id_type elem_id = _current_elem->id();
167  _abs[elem_id] = {Ae, bs};
168 }
169 
170 template <typename T, bool is_ad>
171 void
173 {
174  const auto & npr =
176  _abs.insert(npr._abs.begin(), npr._abs.end());
177 }
178 
179 template <typename T, bool is_ad>
180 void
182 {
183  // When calling nodalPatchRecovery, we may need to know _Ae and _be on algebraically ghosted
184  // elements. However, this userobject is only run on local elements, so we need to query those
185  // information from other processors in this finalize() method.
186 
187  // Populate algebraically ghosted elements to query
188  std::unordered_map<processor_id_type, std::vector<dof_id_type>> query_ids;
189  const ConstElemRange evaluable_elem_range = _fe_problem.getEvaluableElementRange();
190  for (const auto * const elem : evaluable_elem_range)
191  if (elem->processor_id() != processor_id())
192  query_ids[elem->processor_id()].push_back(elem->id());
193 
194  // Answer queries received from other processors
195  auto gather_data = [this](const processor_id_type /*pid*/,
196  const std::vector<dof_id_type> & elem_ids,
197  std::vector<ElementData> & abs_data)
198  {
199  for (const auto i : index_range(elem_ids))
200  abs_data.emplace_back(libmesh_map_find(_abs, elem_ids[i]));
201  };
202 
203  // Gather answers received from other processors
204  auto act_on_data = [this](const processor_id_type /*pid*/,
205  const std::vector<dof_id_type> & elem_ids,
206  const std::vector<ElementData> & abs_data)
207  {
208  for (const auto i : index_range(elem_ids))
209  _abs[elem_ids[i]] = abs_data[i];
210  };
211 
212  libMesh::Parallel::pull_parallel_vector_data<ElementData>(
213  _communicator, query_ids, gather_data, act_on_data, 0);
214 }
215 
virtual void initialize() override
Called before execute() is ever called so that data can be cleared.
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:42
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:333
static InputParameters validParams()
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real nodalPatchRecovery(const Point &p, const std::vector< dof_id_type > &elem_ids, std::size_t component) const override
Solve the least-squares problem.
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.
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...
RealEigenVector evaluateBasisFunctions(const Point &q_point) const
Compute the P vector at a given point i.e.
uint8_t processor_id_type
virtual void threadJoin(const UserObject &) override
Must override.
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
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
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
SerialAccessRange< T > serialAccess(T &obj)
Definition: SerialAccess.h:151
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("MooseApp", ProjectedStatefulMaterialNodalPatchRecoveryReal)
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 addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
T pow(T x, int e)
Definition: MathUtils.h:91
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...