GRINS-0.8.0
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
GRINS::NeumannBCFunctionBase< FunctionType, FEShape > Class Template Reference

#include <neumann_bc_function_base.h>

Inheritance diagram for GRINS::NeumannBCFunctionBase< FunctionType, FEShape >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::NeumannBCFunctionBase< FunctionType, FEShape >:
Collaboration graph
[legend]

Public Member Functions

 NeumannBCFunctionBase (VariableIndex var)
 Constructor for function with only one variable. More...
 
 NeumannBCFunctionBase (const std::vector< VariableIndex > &vars)
 Constructor for function with several variables. More...
 
virtual ~NeumannBCFunctionBase ()
 
virtual bool eval_flux (bool compute_jacobian, AssemblyContext &context, libMesh::Real sign, bool is_axisymmetric)
 
- Public Member Functions inherited from GRINS::NeumannBCAbstract
 NeumannBCAbstract ()
 
virtual ~NeumannBCAbstract ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 

Protected Member Functions

FEShape eval_func (AssemblyContext &context, const libMesh::Point &point, libMesh::Real time, unsigned int component, libMesh::FEMFunctionBase< FEShape > &func)
 Helper function to dispatch to FEMFunctionBase API. More...
 
FEShape eval_func (AssemblyContext &, const libMesh::Point &point, libMesh::Real time, unsigned int component, libMesh::FunctionBase< FEShape > &func)
 Helper function to dispatch to FunctionBase API. More...
 

Protected Attributes

std::vector< VariableIndex_vars
 Variable indices for the variables whose Neumann contribution we're computing. More...
 
libMesh::UniquePtr< FunctionType > _func
 Function object for the actual Neumann flux. More...
 

Additional Inherited Members

- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 

Detailed Description

template<typename FunctionType, typename FEShape>
class GRINS::NeumannBCFunctionBase< FunctionType, FEShape >

Definition at line 42 of file neumann_bc_function_base.h.

Constructor & Destructor Documentation

template<typename FunctionType, typename FEShape>
GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::NeumannBCFunctionBase ( VariableIndex  var)
inline

Constructor for function with only one variable.

Definition at line 47 of file neumann_bc_function_base.h.

49  _vars(1,var)
50  {}
std::vector< VariableIndex > _vars
Variable indices for the variables whose Neumann contribution we're computing.
template<typename FunctionType, typename FEShape>
GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::NeumannBCFunctionBase ( const std::vector< VariableIndex > &  vars)
inline

Constructor for function with several variables.

This intended for FunctionTypes that are Composite so we can treat "group" variables, e.g. Displacement, as the same time using a composite function.

Definition at line 56 of file neumann_bc_function_base.h.

58  _vars(vars)
59  {}
std::vector< VariableIndex > _vars
Variable indices for the variables whose Neumann contribution we're computing.
template<typename FunctionType, typename FEShape>
virtual GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::~NeumannBCFunctionBase ( )
inlinevirtual

Definition at line 61 of file neumann_bc_function_base.h.

61 {};

Member Function Documentation

template<typename FunctionType , typename FEShape >
bool GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::eval_flux ( bool  compute_jacobian,
AssemblyContext context,
libMesh::Real  sign,
bool  is_axisymmetric 
)
virtual

Implements GRINS::NeumannBCAbstract.

Definition at line 41 of file neumann_bc_function_base.C.

45  {
46  libmesh_assert(!_vars.empty());
47 
48  unsigned int n_vars = this->_vars.size();
49 
50  // The number of local degrees of freedom in each variable.
51  // We're assuming that each var in our set is the same FEType
52  const unsigned int n_var_dofs = context.get_dof_indices(this->_vars[0]).size();
53 #ifndef NDEBUG
54  for( unsigned int v = 1; v < n_vars; v++ )
55  libmesh_assert_equal_to( n_var_dofs, context.get_dof_indices(this->_vars[v]).size() );
56 #endif
57 
58  libMesh::FEGenericBase<FEShape>* side_fe = NULL;
59  context.get_side_fe( this->_vars[0], side_fe );
60 
61  // Element Jacobian * quadrature weight for side integration.
62  const std::vector<libMesh::Real> &JxW_side = side_fe->get_JxW();
63 
64  // The var shape functions at side quadrature points.
65  const std::vector<std::vector<FEShape> >& var_phi_side =
66  side_fe->get_phi();
67 
68  // Physical location of the quadrature points
69  const std::vector<libMesh::Point>& var_qpoint =
70  side_fe->get_xyz();
71 
72  std::vector<libMesh::DenseSubVector<libMesh::Number>*> F_vars(_vars.size(),NULL);
73 
74  for( unsigned int v = 0; v < n_vars; v++ )
75  F_vars[v] = &(context.get_elem_residual(this->_vars[v])); // residual
76 
77  unsigned int n_qpoints = context.get_side_qrule().n_points();
78 
79  FEShape value = 0.0;
80 
81  for (unsigned int qp=0; qp != n_qpoints; qp++)
82  {
83  libMesh::Real jac = JxW_side[qp];
84 
85  if( is_axisymmetric )
86  {
87  const libMesh::Number r = var_qpoint[qp](0);
88  jac *= r;
89  }
90 
91  for( unsigned int v = 0; v < n_vars; v++ )
92  {
93  value = this->eval_func( context, var_qpoint[qp], context.get_time(), this->_vars[v], *_func );
94 
95  for (unsigned int i=0; i != n_var_dofs; i++)
96  (*F_vars[v])(i) += sign*(value*var_phi_side[i][qp])*jac;
97  }
98  }
99 
100  if( ParsedFunctionTraits<FunctionType>::is_fem_function && compute_jacobian )
101  libmesh_error_msg("ERROR: Analytical Jacobians of FEMFunctionBase objects currently not supported!");
102 
103  return compute_jacobian;
104  }
std::vector< VariableIndex > _vars
Variable indices for the variables whose Neumann contribution we're computing.
FEShape eval_func(AssemblyContext &context, const libMesh::Point &point, libMesh::Real time, unsigned int component, libMesh::FEMFunctionBase< FEShape > &func)
Helper function to dispatch to FEMFunctionBase API.
libMesh::UniquePtr< FunctionType > _func
Function object for the actual Neumann flux.
template<typename FunctionType, typename FEShape>
FEShape GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::eval_func ( AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real  time,
unsigned int  component,
libMesh::FEMFunctionBase< FEShape > &  func 
)
inlineprotected

Helper function to dispatch to FEMFunctionBase API.

Definition at line 71 of file neumann_bc_function_base.h.

74  {
75  return func.component(context,component,point,time);
76  }
template<typename FunctionType, typename FEShape>
FEShape GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::eval_func ( AssemblyContext ,
const libMesh::Point &  point,
libMesh::Real  time,
unsigned int  component,
libMesh::FunctionBase< FEShape > &  func 
)
inlineprotected

Helper function to dispatch to FunctionBase API.

Definition at line 79 of file neumann_bc_function_base.h.

82  {
83  return func.component(component,point,time);
84  }

Member Data Documentation

template<typename FunctionType, typename FEShape>
libMesh::UniquePtr<FunctionType> GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::_func
protected

Function object for the actual Neumann flux.

Subclasses should initialize this function appropriately at construction time.

Definition at line 95 of file neumann_bc_function_base.h.

template<typename FunctionType, typename FEShape>
std::vector<VariableIndex> GRINS::NeumannBCFunctionBase< FunctionType, FEShape >::_vars
protected

Variable indices for the variables whose Neumann contribution we're computing.

We're assuming all _vars use the same FunctionType object. This allows us to use Composite type functions and treat "group" variables, like Displacement, at the same time.

Definition at line 90 of file neumann_bc_function_base.h.


The documentation for this class was generated from the following files:

Generated on Tue Dec 19 2017 12:47:31 for GRINS-0.8.0 by  doxygen 1.8.9.1