42 template <
typename T,
bool is_ad>
75 virtual std::vector<std::vector<unsigned int>>
nChooseK(
unsigned int N,
unsigned int K);
94 virtual std::vector<std::vector<unsigned int>>
multiIndex(
unsigned int dim,
unsigned int order);
132 virtual void compute()
override;
146 template <
typename T>
150 return MaterialPropertyInterface::getMaterialProperty<T>(
name);
153 template <
typename T,
bool is_ad>
157 return MaterialPropertyInterface::getGenericMaterialProperty<T, is_ad>(
name);
160 template <
typename T>
164 return MaterialPropertyInterface::getMaterialPropertyOld<T>(
name);
167 template <
typename T>
171 return MaterialPropertyInterface::getMaterialPropertyOlder<T>(
name);
NodalPatchRecovery(const InputParameters ¶meters)
const MaterialProperty< T > & getMaterialProperty(const std::string &name)
This function overrides the one implemented in AuxKernel.C to suppress warnings when retrieving mater...
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty(const std::string &name)
This class uses patch recovery technique to recover material properties to smooth nodal fields...
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
virtual const std::string & name() const
Get the name of the class.
FEProblemBase & _fe_problem
virtual ~NodalPatchRecovery()
const MaterialProperty< T > & getMaterialPropertyOlder(const std::string &name)
This function overrides the one implemented in AuxKernel.C to suppress warnings when retrieving mater...
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.
virtual void accumulateBVector(Real val)
calculate the patch load vector as ^n P^Tval where n is the number of quadrature points in the elemen...
typename GenericMaterialPropertyStruct< T, is_ad >::type GenericMaterialProperty
virtual void compute() override
solve the equation Ac = B where c is the coefficients vector from least square fitting nodal value is...
const MaterialProperty< T > & getMaterialPropertyOld(const std::string &name)
This function overrides the one implemented in AuxKernel.C to suppress warnings when retrieving mater...
static InputParameters validParams()
const InputParameters & parameters() const
Get the parameters of the object.
const unsigned int _patch_polynomial_order
polynomial order, default is variable order
virtual void accumulateAMatrix()
Calculate the patch stiffness matrix as ^n P^TP where n is the number of quadrature points in the ele...
virtual void reinitPatch()
reserve space for A, B, P, and prepare required material properties
std::vector< std::vector< unsigned int > > _multi_index
virtual void computePVector(Point q_point)
compute the P vector at given point i.e.
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...
DenseVector< Number > _coef