https://mooseframework.inl.gov
NodalPatchRecoveryBase.h
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 #pragma once
11 
12 // MOOSE includes
13 #include "ElementUserObject.h"
14 
16 
18 {
19 public:
21 
23 
32  virtual Real nodalPatchRecovery(const Point & p, const std::vector<dof_id_type> & elem_ids) const;
33 
43  const RealEigenVector getCoefficients(const std::vector<dof_id_type> & elem_ids) const;
44 
55  const RealEigenVector getCachedCoefficients(const std::vector<dof_id_type> & elem_ids);
56 
57  void initialize() override;
58  void execute() override;
59  void threadJoin(const UserObject &) override;
60  void finalize() override;
61 
63  const std::vector<std::vector<unsigned int>> & multiIndex() const { return _multi_index; }
64 
65 protected:
67  virtual Real computeValue() = 0;
68 
69  unsigned int _qp;
70 
71 private:
78  std::unordered_map<processor_id_type, std::vector<dof_id_type>>
79  gatherRequestList(const std::vector<dof_id_type> & specific_elems);
80 
86  std::unordered_map<processor_id_type, std::vector<dof_id_type>> gatherRequestList();
87 
100  RealEigenVector evaluateBasisFunctions(const Point & q_point) const;
101 
109  void sync(const std::vector<dof_id_type> & specific_elems);
110 
115  void sync();
116 
120  void
121  syncHelper(const std::unordered_map<processor_id_type, std::vector<dof_id_type>> & query_ids);
122 
129  void
130  addToQuery(const libMesh::Elem * elem,
131  std::unordered_map<processor_id_type, std::vector<dof_id_type>> & query_ids) const;
132 
134  const unsigned int _patch_polynomial_order;
135 
152  const std::vector<std::vector<unsigned int>> _multi_index;
153 
155  const unsigned int _q;
156 
158  std::map<dof_id_type, RealEigenMatrix> _Ae;
159 
161  std::map<dof_id_type, RealEigenVector> _be;
162 
166  mutable std::vector<dof_id_type> _cached_elem_ids;
168 
171 
173  std::vector<int> _proc_ids;
174 
176 };
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.
const unsigned int _patch_polynomial_order
The polynomial order, default is variable order.
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 std::vector< std::vector< unsigned int > > & multiIndex() const
Returns the multi-index table.
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.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
RealEigenVector evaluateBasisFunctions(const Point &q_point) const
Compute the P vector at a given point i.e.
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) const
Solve the least-squares problem.
void execute() override
Execute method.
uint8_t processor_id_type
NodalPatchRecoveryBase(const InputParameters &parameters)
const std::vector< std::vector< unsigned int > > _multi_index
Multi-index table for a polynomial basis.
bool _distributed_mesh
Whether the mesh is distributed.
virtual Real computeValue()=0
Compute the quantity to recover using nodal patch recovery.
void finalize() override
Finalize.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void threadJoin(const UserObject &) override
Must override.
Base class for mesh modifiers modifying element subdomains.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
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.
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.
Base class for user-specific data.
Definition: UserObject.h:40