https://mooseframework.inl.gov
NodalPatchRecovery.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 #include "AuxKernel.h"
13 #include "FEProblem.h"
14 
26 {
27 public:
29 
31  virtual ~NodalPatchRecovery(){};
40  template <typename T>
41  const MaterialProperty<T> & getMaterialProperty(const std::string & name);
42  template <typename T, bool is_ad>
44 
53  template <typename T>
54  const MaterialProperty<T> & getMaterialPropertyOld(const std::string & name);
63  template <typename T>
64  const MaterialProperty<T> & getMaterialPropertyOlder(const std::string & name);
65 
66 protected:
75  virtual std::vector<std::vector<unsigned int>> nChooseK(unsigned int N, unsigned int K);
76 
94  virtual std::vector<std::vector<unsigned int>> multiIndex(unsigned int dim, unsigned int order);
95 
108  virtual void computePVector(Point q_point);
109 
114  virtual void accumulateAMatrix();
115 
122  virtual void accumulateBVector(Real val);
123 
125  virtual void reinitPatch();
126 
132  virtual void compute() override;
133 
134 private:
136  const unsigned int _patch_polynomial_order;
137 
139  std::vector<std::vector<unsigned int>> _multi_index;
144 };
145 
146 template <typename T>
147 const MaterialProperty<T> &
149 {
150  return MaterialPropertyInterface::getMaterialProperty<T>(name);
151 }
152 
153 template <typename T, bool is_ad>
156 {
157  return MaterialPropertyInterface::getGenericMaterialProperty<T, is_ad>(name);
158 }
159 
160 template <typename T>
161 const MaterialProperty<T> &
163 {
164  return MaterialPropertyInterface::getMaterialPropertyOld<T>(name);
165 }
166 
167 template <typename T>
168 const MaterialProperty<T> &
170 {
171  return MaterialPropertyInterface::getMaterialPropertyOlder<T>(name);
172 }
NodalPatchRecovery(const InputParameters &parameters)
const MaterialProperty< T > & getMaterialProperty(const std::string &name)
This function overrides the one implemented in AuxKernel.C to suppress warnings when retrieving mater...
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:127
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:159
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.
FEProblemBase & _fe_problem
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.
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:99
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
DenseVector< Number > _P
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...
DenseVector< Number > _B
static InputParameters validParams()
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
DenseMatrix< Number > _A
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