Go to the documentation of this file.
19 MooseEnum orders(
"CONSTANT FIRST SECOND THIRD FOURTH");
22 "Polynomial order used in least squares fitting of material property "
23 "over the local patch of elements connected to a given node");
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))
45 mooseWarning(
"Specified 'patch_polynomial_order' is lower than the AuxVariable's order");
52 meshhelper.allow_renumbering(
false);
53 for (
const auto & elem :
54 as_range(meshhelper.semilocal_elements_begin(), meshhelper.semilocal_elements_end()))
59 std::vector<std::vector<unsigned int>>
62 std::vector<std::vector<unsigned int>> n_choose_k;
63 std::vector<unsigned int> row;
64 std::string bitmask(K, 1);
71 for (
unsigned int i = 0; i <
N; ++i)
75 n_choose_k.push_back(row);
76 }
while (std::prev_permutation(bitmask.begin(), bitmask.end()));
81 std::vector<std::vector<unsigned int>>
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);
91 for (
unsigned int q = 1; q <= order; q++)
94 multi_index.push_back(row);
97 for (
unsigned int q = 1; q <= order; q++)
99 n_choose_k =
nChooseK(dim + q - 1, dim - 1);
100 for (
unsigned int r = 0; r < n_choose_k.size(); r++)
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);
120 for (
unsigned int c = 0; c <
_multi_index[0].size(); c++)
122 polynomial *= q_point(c);
130 DenseMatrix<Number> PxP;
132 PxP.outer_product(
_P,
_P);
152 std::set<unsigned int> needed_mat_props;
154 needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
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");
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;
184 if (elem_ids.size() == 1)
186 dof_id_type elem_id = elem_ids[0];
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;
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.");
202 for (
auto elem_id : elem_ids)
206 if (!elem->is_semilocal(
_mesh.processor_id()))
233 Real nodal_value =
_P.dot(
_coef);
239 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
virtual Elem * elemPtr(const dof_id_type i)
virtual void reinitMaterials(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
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.
void mooseError(Args &&... args) const
defineLegacyParams(NodalPatchRecovery)
virtual void clearActiveMaterialProperties(THREAD_ID tid) override
Clear the active material properties.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
virtual void swapBackMaterials(THREAD_ID tid)
FEProblemBase & _fe_problem
Base class for creating new auxiliary kernels and auxiliary boundary conditions.
virtual void computePVector(Point q_point)
compute the P vector at given point i.e.
virtual void compute() override
solve the equation Ac = B where c is the coefficients vector from least square fitting nodal value is...
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 ...
virtual void compute()
Computes the value and stores it in the solution vector.
const MooseArray< Point > & _q_point
Dimension of the problem being solved.
Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int.
MooseMesh & _mesh
Mesh this kernel is active on.
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
void setNodalValue(const OutputType &value, unsigned int idx=0)
Set nodal value.
virtual void reinitElem(const Elem *elem, THREAD_ID tid) override
DenseVector< Number > _coef
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...
static InputParameters validParams()
virtual void reinitPatch()
reserve space for A, B, P, and prepare required material properties
bool isNodal() const
Nodal or elemental kernel?
virtual void reinitNode(const Node *node, THREAD_ID tid) override
const Node *const & _current_node
Current node (valid only for nodal kernels)
virtual void prepare(const Elem *elem, THREAD_ID tid) override
virtual ComputeValueType computeValue()=0
Compute and return the value of the aux variable.
virtual void addGhostedElem(dof_id_type elem_id) override
Will make sure that all dofs connected to elem_id are ghosted to this processor.
The "SwapBackSentinel" class's destructor guarantees that FEProblemBase::swapBackMaterials{Face,...
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.
const unsigned int _patch_polynomial_order
polynomial order, default is variable order
static InputParameters validParams()
unsigned int _qp
Quadrature point index.
const std::set< unsigned int > & getMatPropDependencies() const
Retrieve the set of material properties that this object depends on.
MooseVariableFE< ComputeValueType > & _var
This is a regular kernel so we cast to a regular MooseVariable.
unsigned int size() const
The number of elements that can currently be stored in the array.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
NumericVector< Number > & _solution
reference to the solution vector of auxiliary system
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.
void mooseWarning(Args &&... args) const
NodalPatchRecovery(const InputParameters ¶meters)
std::vector< std::vector< unsigned int > > _multi_index
virtual void accumulateAMatrix()
Calculate the patch stiffness matrix as \sum_1^n P^TP where n is the number of quadrature points in t...
const dof_id_type & nodalDofIndex() const override