www.mooseframework.org
NodalPatchRecovery.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "NodalPatchRecovery.h"
11 #include "SwapBackSentinel.h"
12 
14 
17 {
19  MooseEnum orders("CONSTANT FIRST SECOND THIRD FOURTH");
20  params.addParam<MooseEnum>("patch_polynomial_order",
21  orders,
22  "Polynomial order used in least squares fitting of material property "
23  "over the local patch of elements connected to a given node");
24  params.addParamNamesToGroup("patch_polynomial_order", "Advanced");
25 
26  // TODO make this class work with relationship manager
27  // params.registerRelationshipManagers("ElementSideNeighborLayers");
28  // params.addPrivateParam<unsigned int short>("element_side_neighbor_layers", 2);
29 
30  return params;
31 }
32 
34  : AuxKernel(parameters),
35  _patch_polynomial_order(
36  isParamValid("patch_polynomial_order")
37  ? static_cast<unsigned int>(getParam<MooseEnum>("patch_polynomial_order"))
38  : static_cast<unsigned int>(_var.order())),
39  _fe_problem(*parameters.get<FEProblemBase *>("_fe_problem_base")),
40  _multi_index(multiIndex(_mesh.dimension(), _patch_polynomial_order))
41 {
42  // if the patch polynomial order is lower than the variable order,
43  // it is very likely that the patch recovery is not used at its max accuracy
44  if (_patch_polynomial_order < static_cast<unsigned int>(_var.order()))
45  mooseWarning("Specified 'patch_polynomial_order' is lower than the AuxVariable's order");
46 
47  // TODO remove the manual ghosting once relationship manager is working correctly
48  // no need to ghost if this aux is elemental
49  if (isNodal())
50  {
51  MeshBase & meshhelper = _mesh.getMesh();
52  meshhelper.allow_renumbering(false);
53  for (const auto & elem :
54  as_range(meshhelper.semilocal_elements_begin(), meshhelper.semilocal_elements_end()))
55  _fe_problem.addGhostedElem(elem->id());
56  }
57 }
58 
59 std::vector<std::vector<unsigned int>>
60 NodalPatchRecovery::nChooseK(unsigned int N, unsigned int K)
61 {
62  std::vector<std::vector<unsigned int>> n_choose_k;
63  std::vector<unsigned int> row;
64  std::string bitmask(K, 1); // K leading 1's
65  bitmask.resize(N, 0); // N-K trailing 0's
66 
67  do
68  {
69  row.clear();
70  row.push_back(0);
71  for (unsigned int i = 0; i < N; ++i) // [0..N-1] integers
72  if (bitmask[i])
73  row.push_back(i + 1);
74  row.push_back(N + 1);
75  n_choose_k.push_back(row);
76  } while (std::prev_permutation(bitmask.begin(), bitmask.end()));
77 
78  return n_choose_k;
79 }
80 
81 std::vector<std::vector<unsigned int>>
82 NodalPatchRecovery::multiIndex(unsigned int dim, unsigned int order)
83 {
84  // first row all zero
85  std::vector<std::vector<unsigned int>> multi_index;
86  std::vector<std::vector<unsigned int>> n_choose_k;
87  std::vector<unsigned int> row(dim, 0);
88  multi_index.push_back(row);
89 
90  if (dim == 1)
91  for (unsigned int q = 1; q <= order; q++)
92  {
93  row[0] = q;
94  multi_index.push_back(row);
95  }
96  else
97  for (unsigned int q = 1; q <= order; q++)
98  {
99  n_choose_k = nChooseK(dim + q - 1, dim - 1);
100  for (unsigned int r = 0; r < n_choose_k.size(); r++)
101  {
102  row.clear();
103  for (unsigned int c = 1; c < n_choose_k[0].size(); c++)
104  row.push_back(n_choose_k[r][c] - n_choose_k[r][c - 1] - 1);
105  multi_index.push_back(row);
106  }
107  }
108 
109  return multi_index;
110 }
111 
112 void
114 {
115  _P.resize(_multi_index.size());
116  Real polynomial;
117  for (unsigned int r = 0; r < _multi_index.size(); r++)
118  {
119  polynomial = 1.0;
120  for (unsigned int c = 0; c < _multi_index[0].size(); c++)
121  for (unsigned int p = 0; p < _multi_index[r][c]; p++)
122  polynomial *= q_point(c);
123  _P(r) = polynomial;
124  }
125 }
126 
127 void
129 {
130  DenseMatrix<Number> PxP;
131  PxP.resize(_multi_index.size(), _multi_index.size());
132  PxP.outer_product(_P, _P);
133  _A.add(1.0, PxP);
134 }
135 
136 void
138 {
139  _B.add(val, _P);
140 }
141 
142 void
144 {
145  // clean and resize P, A, B, coef
146  _P.resize(_multi_index.size());
147  _A.resize(_multi_index.size(), _multi_index.size());
148  _B.resize(_multi_index.size());
149  _coef.resize(_multi_index.size());
150 
151  // activate dependent material properties
152  std::set<unsigned int> needed_mat_props;
153  const std::set<unsigned int> & mp_deps = getMatPropDependencies();
154  needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
155  _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid);
156 }
157 
158 void
160 {
161  // Fall back on standard procedure for element variables
162  if (!isNodal())
163  {
165  return;
166  }
167 
168  // Limit current use of NodalPatchRecovery to a single processor
169  if (_communicator.size() > 1)
170  mooseError("The nodal patch recovery option, which calculates the Zienkiewicz-Zhu patch "
171  "recovery for nodal variables (family = LAGRANGE), is not currently implemented for "
172  "parallel runs. Run in serial if you must use the nodal patch capability");
173 
174  // Use Zienkiewicz-Zhu patch recovery for nodal variables
175  reinitPatch();
176 
177  // get node-to-conneted-elem map
178  const std::map<dof_id_type, std::vector<dof_id_type>> & node_to_elem_map = _mesh.nodeToElemMap();
179  auto node_to_elem_pair = node_to_elem_map.find(_current_node->id());
180  mooseAssert(node_to_elem_pair != node_to_elem_map.end(), "Missing entry in node to elem map");
181  std::vector<dof_id_type> elem_ids = node_to_elem_pair->second;
182 
183  // consider the case for corner node
184  if (elem_ids.size() == 1)
185  {
186  dof_id_type elem_id = elem_ids[0];
187  for (auto & n : _mesh.elemPtr(elem_id)->node_ref_range())
188  {
189  node_to_elem_pair = node_to_elem_map.find(n.id());
190  std::vector<dof_id_type> elem_ids_candidate = node_to_elem_pair->second;
191  if (elem_ids_candidate.size() > elem_ids.size())
192  elem_ids = elem_ids_candidate;
193  }
194  }
195 
196  // before we go, check if we have enough sample points for solving the least square fitting
197  if (_q_point.size() * elem_ids.size() < _multi_index.size())
198  mooseError("There are not enough sample points to recover the nodal value, try reducing the "
199  "polynomial order or using a higher-order quadrature scheme.");
200 
201  // general treatment for side nodes and internal nodes
202  for (auto elem_id : elem_ids)
203  {
204  const Elem * elem = _mesh.elemPtr(elem_id);
205 
206  if (!elem->is_semilocal(_mesh.processor_id()))
207  mooseError("skipped non local elem!");
208 
209  _fe_problem.prepare(elem, _tid);
210  _fe_problem.reinitElem(elem, _tid);
211 
212  // Set up Sentinel class so that, even if reinitMaterials() throws, we
213  // still remember to swap back during stack unwinding.
215  _fe_problem.reinitMaterials(elem->subdomain_id(), _tid);
216 
217  for (_qp = 0; _qp < _q_point.size(); _qp++)
218  {
222  }
223  }
224 
226 
227  // solve for the coefficients of least square fitting
228  // as we are on the physical coordinates, A is not always SPD
229  _A.svd_solve(_B, _coef);
230 
231  // compute fitted nodal value
233  Real nodal_value = _P.dot(_coef);
234 
235  // set nodal value
237  dof_id_type dof = _var.nodalDofIndex();
238  {
239  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
240  _solution.set(dof, nodal_value);
241  }
242  _var.setNodalValue(nodal_value);
243 }
MooseMesh::elemPtr
virtual Elem * elemPtr(const dof_id_type i)
Definition: MooseMesh.C:2289
FEProblemBase::reinitMaterials
virtual void reinitMaterials(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
Definition: FEProblemBase.C:2933
MooseMesh::nodeToElemMap
const std::map< dof_id_type, std::vector< dof_id_type > > & nodeToElemMap()
If not already created, creates a map from every node to all elements to which they are connected.
Definition: MooseMesh.C:709
MooseObject::mooseError
void mooseError(Args &&... args) const
Definition: MooseObject.h:141
NodalPatchRecovery::_A
DenseMatrix< Number > _A
Definition: NodalPatchRecovery.h:143
AuxKernelTempl::_tid
THREAD_ID _tid
Thread ID.
Definition: AuxKernel.h:169
SwapBackSentinel.h
defineLegacyParams
defineLegacyParams(NodalPatchRecovery)
FEProblemBase::clearActiveMaterialProperties
virtual void clearActiveMaterialProperties(THREAD_ID tid) override
Clear the active material properties.
Definition: FEProblemBase.C:4383
MooseEnum
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
InputParameters::addParam
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object.
Definition: InputParameters.h:1198
FEProblemBase::swapBackMaterials
virtual void swapBackMaterials(THREAD_ID tid)
Definition: FEProblemBase.C:3072
NodalPatchRecovery::_fe_problem
FEProblemBase & _fe_problem
Definition: NodalPatchRecovery.h:140
AuxKernelTempl
Base class for creating new auxiliary kernels and auxiliary boundary conditions.
Definition: AuxKernel.h:35
NodalPatchRecovery::computePVector
virtual void computePVector(Point q_point)
compute the P vector at given point i.e.
Definition: NodalPatchRecovery.C:113
NodalPatchRecovery::compute
virtual void compute() override
solve the equation Ac = B where c is the coefficients vector from least square fitting nodal value is...
Definition: NodalPatchRecovery.C:159
NodalPatchRecovery::accumulateBVector
virtual void accumulateBVector(Real val)
calculate the patch load vector as \sum_1^n P^Tval where n is the number of quadrature points in the ...
Definition: NodalPatchRecovery.C:137
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
AuxKernelTempl::compute
virtual void compute()
Computes the value and stores it in the solution vector.
Definition: AuxKernel.C:274
AuxKernelTempl::_q_point
const MooseArray< Point > & _q_point
Dimension of the problem being solved.
Definition: AuxKernel.h:200
MooseVariableBase::order
Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int.
Definition: MooseVariableBase.C:123
AuxKernelTempl::_mesh
MooseMesh & _mesh
Mesh this kernel is active on.
Definition: AuxKernel.h:195
NodalPatchRecovery::_B
DenseVector< Number > _B
Definition: NodalPatchRecovery.h:144
MooseMesh::getMesh
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:2599
MooseVariableFE::setNodalValue
void setNodalValue(const OutputType &value, unsigned int idx=0)
Set nodal value.
Definition: MooseVariableFE.C:678
FEProblemBase::reinitElem
virtual void reinitElem(const Elem *elem, THREAD_ID tid) override
Definition: FEProblemBase.C:1620
NodalPatchRecovery::_coef
DenseVector< Number > _coef
Definition: NodalPatchRecovery.h:145
NodalPatchRecovery::nChooseK
virtual std::vector< std::vector< unsigned int > > nChooseK(unsigned int N, unsigned int K)
Find out how many different ways you can choose K items from N items set without repetition and witho...
Definition: NodalPatchRecovery.C:60
AuxKernelTempl::validParams
static InputParameters validParams()
Definition: AuxKernel.C:28
NodalPatchRecovery::reinitPatch
virtual void reinitPatch()
reserve space for A, B, P, and prepare required material properties
Definition: NodalPatchRecovery.C:143
AuxKernelTempl::isNodal
bool isNodal() const
Nodal or elemental kernel?
Definition: AuxKernel.h:89
FEProblemBase::reinitNode
virtual void reinitNode(const Node *node, THREAD_ID tid) override
Definition: FEProblemBase.C:1670
AuxKernelTempl::_current_node
const Node *const & _current_node
Current node (valid only for nodal kernels)
Definition: AuxKernel.h:218
InputParameters::addParamNamesToGroup
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...
Definition: InputParameters.C:590
FEProblemBase::prepare
virtual void prepare(const Elem *elem, THREAD_ID tid) override
Definition: FEProblemBase.C:1227
AuxKernelTempl::computeValue
virtual ComputeValueType computeValue()=0
Compute and return the value of the aux variable.
NodalPatchRecovery
Definition: NodalPatchRecovery.h:30
FEProblemBase::addGhostedElem
virtual void addGhostedElem(dof_id_type elem_id) override
Will make sure that all dofs connected to elem_id are ghosted to this processor.
Definition: FEProblemBase.C:1534
NodalPatchRecovery.h
NodalPatchRecovery::_P
DenseVector< Number > _P
Definition: NodalPatchRecovery.h:142
SwapBackSentinel
The "SwapBackSentinel" class's destructor guarantees that FEProblemBase::swapBackMaterials{Face,...
Definition: SwapBackSentinel.h:32
FEProblemBase::setActiveMaterialProperties
virtual void setActiveMaterialProperties(const std::set< unsigned int > &mat_prop_ids, THREAD_ID tid) override
Record and set the material properties required by the current computing thread.
Definition: FEProblemBase.C:4373
NodalPatchRecovery::_patch_polynomial_order
const unsigned int _patch_polynomial_order
polynomial order, default is variable order
Definition: NodalPatchRecovery.h:138
NodalPatchRecovery::validParams
static InputParameters validParams()
Definition: NodalPatchRecovery.C:16
AuxKernelTempl::_qp
unsigned int _qp
Quadrature point index.
Definition: AuxKernel.h:227
MaterialPropertyInterface::getMatPropDependencies
const std::set< unsigned int > & getMatPropDependencies() const
Retrieve the set of material properties that this object depends on.
Definition: MaterialPropertyInterface.h:187
N
PetscInt N
Definition: PetscDMMoose.C:1504
AuxKernelTempl::_var
MooseVariableFE< ComputeValueType > & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: AuxKernel.h:172
MooseArray::size
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
FEProblemBase
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
Definition: FEProblemBase.h:139
AuxKernelTempl::_solution
NumericVector< Number > & _solution
reference to the solution vector of auxiliary system
Definition: AuxKernel.h:224
NodalPatchRecovery::multiIndex
virtual 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: NodalPatchRecovery.C:82
MooseObject::mooseWarning
void mooseWarning(Args &&... args) const
Definition: MooseObject.h:150
NodalPatchRecovery::NodalPatchRecovery
NodalPatchRecovery(const InputParameters &parameters)
Definition: NodalPatchRecovery.C:33
NodalPatchRecovery::_multi_index
std::vector< std::vector< unsigned int > > _multi_index
Definition: NodalPatchRecovery.h:141
NodalPatchRecovery::accumulateAMatrix
virtual void accumulateAMatrix()
Calculate the patch stiffness matrix as \sum_1^n P^TP where n is the number of quadrature points in t...
Definition: NodalPatchRecovery.C:128
MooseVariableFE::nodalDofIndex
const dof_id_type & nodalDofIndex() const override
Definition: MooseVariableFE.h:164
n
PetscInt n
Definition: PetscDMMoose.C:1504