GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions | Private Attributes
GRINS::IntegratedFunction< Function > Class Template Reference

IntegratedFunction. More...

#include <integrated_function.h>

Inheritance diagram for GRINS::IntegratedFunction< Function >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::IntegratedFunction< Function >:
Collaboration graph
[legend]

Public Member Functions

 IntegratedFunction (unsigned int p_level, SharedPtr< Function > f, RayfireMesh *rayfire, const std::string &qoi_name)
 Constructor. More...
 
 IntegratedFunction (const GetPot &input, unsigned int p_level, SharedPtr< Function > f, const std::string &input_qoi_string, const std::string &qoi_name)
 Constructor. More...
 
 IntegratedFunction (const IntegratedFunction &original)
 Copy Constructor. More...
 
virtual QoIBaseclone () const
 Required to provide clone (deep-copy) for adding QoI object to libMesh objects. More...
 
virtual bool assemble_on_interior () const
 Does the QoI need an element interior assembly loop? More...
 
virtual bool assemble_on_sides () const
 Does the QoI need a domain boundary assembly loop? More...
 
virtual void element_qoi (AssemblyContext &context, const unsigned int qoi_index)
 Compute the qoi value. More...
 
virtual void element_qoi_derivative (AssemblyContext &context, const unsigned int qoi_index)
 Compute the qoi derivative with respect to the solution. More...
 
virtual void init (const GetPot &input, const MultiphysicsSystem &system, unsigned int qoi_num)
 Initializes the rayfire with the mesh from system. More...
 
virtual void reinit (MultiphysicsSystem &system)
 Reinitialize the rayfire. More...
 
const RayfireMeshget_rayfire ()
 
- Public Member Functions inherited from GRINS::QoIBase
 QoIBase (const std::string &qoi_name)
 
virtual ~QoIBase ()
 
virtual void init_context (AssemblyContext &)
 
virtual void side_qoi (AssemblyContext &, const unsigned int)
 Compute the qoi value on the domain boundary. More...
 
virtual void side_qoi_derivative (AssemblyContext &, const unsigned int)
 Compute the qoi derivative with respect to the solution on the domain boundary. More...
 
virtual void parallel_op (const libMesh::Parallel::Communicator &communicator, libMesh::Number &sys_qoi, libMesh::Number &local_qoi)
 Call the parallel operation for this QoI and cache the value. More...
 
virtual void thread_join (libMesh::Number &qoi, const libMesh::Number &other_qoi)
 Call the operation to accumulate this QoI from multiple threads. More...
 
virtual void output_qoi (std::ostream &out) const
 Basic output for computed QoI's. More...
 
libMesh::Number value () const
 Returns the current QoI value. More...
 
const std::string & name () const
 Returns the name of this QoI. More...
 
- 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...
 

Private Member Functions

libMesh::Real qoi_value (Function &f, AssemblyContext &context, const libMesh::Point &xyz)
 Compute the value of a QoI at a QP. More...
 
void qoi_derivative (Function &f, AssemblyContext &context, const libMesh::Point &qp_xyz, const libMesh::Real JxW, const unsigned int qoi_index)
 Compute derivatiuves at QP. More...
 
 IntegratedFunction ()
 User cannot call empty constructor. More...
 
template<>
libMesh::Real qoi_value (FEMFunctionAndDerivativeBase< libMesh::Real > &f, AssemblyContext &context, const libMesh::Point &xyz)
 
template<>
libMesh::Real qoi_value (libMesh::FunctionBase< libMesh::Real > &f, AssemblyContext &, const libMesh::Point &xyz)
 
template<>
void qoi_derivative (FEMFunctionAndDerivativeBase< libMesh::Real > &f, AssemblyContext &context, const libMesh::Point &qp_xyz, const libMesh::Real JxW, const unsigned int qoi_index)
 
template<>
void qoi_derivative (libMesh::FunctionBase< libMesh::Real > &, AssemblyContext &, const libMesh::Point &, const libMesh::Real, const unsigned int)
 

Private Attributes

unsigned int _p_level
 Quadrature order. More...
 
SharedPtr< Function > _f
 Pointer to the template class used for function evaluation. More...
 
libMesh::UniquePtr< RayfireMesh_rayfire
 Pointer to RayfireMesh object. 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...
 
- Protected Attributes inherited from GRINS::QoIBase
std::string _qoi_name
 
libMesh::Number _qoi_value
 

Detailed Description

template<typename Function>
class GRINS::IntegratedFunction< Function >

IntegratedFunction.

The purpose of this class is to integrate a function along a 1D line through a 2D mesh. This class utilizes the RayfireMesh class to calculate the line and provide 1D elements for evaluating the function.

The template parameter must be FunctionBase or FEMFunctionBase.

Currently, calculating element qoi derivatives is not supported.

Definition at line 53 of file integrated_function.h.

Constructor & Destructor Documentation

template<typename Function>
GRINS::IntegratedFunction< Function >::IntegratedFunction ( unsigned int  p_level,
SharedPtr< Function >  f,
RayfireMesh rayfire,
const std::string &  qoi_name 
)

Constructor.

Parameters
p_levelThe desired Gauss Quadrature level
fA FunctionBase or FEMFunctionBase object for evaluting the QoI
rayfireA RayfireMesh object (will be initialized in init()) Ownership will be taken by an internal libMesh::UniquePtr
qoi_namePassed to the QoIBase

Definition at line 48 of file integrated_function.C.

48  :
49  QoIBase(qoi_name),
50  _p_level(p_level),
51  _f(f),
52  _rayfire(rayfire)
53  {}
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
unsigned int _p_level
Quadrature order.
QoIBase(const std::string &qoi_name)
Definition: qoi_base.C:39
SharedPtr< Function > _f
Pointer to the template class used for function evaluation.
template<typename Function>
GRINS::IntegratedFunction< Function >::IntegratedFunction ( const GetPot &  input,
unsigned int  p_level,
SharedPtr< Function >  f,
const std::string &  input_qoi_string,
const std::string &  qoi_name 
)

Constructor.

Used by the QoIFactory. Passes GetPot through to RayfireMesh for construction

Definition at line 56 of file integrated_function.C.

References GRINS::IntegratedFunction< Function >::_rayfire.

56  :
57  QoIBase(qoi_name),
58  _p_level(p_level),
59  _f(f)
60  {
61  _rayfire.reset(new RayfireMesh(input,input_qoi_string));
62  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
unsigned int _p_level
Quadrature order.
QoIBase(const std::string &qoi_name)
Definition: qoi_base.C:39
SharedPtr< Function > _f
Pointer to the template class used for function evaluation.
template<typename Function>
GRINS::IntegratedFunction< Function >::IntegratedFunction ( const IntegratedFunction< Function > &  original)

Copy Constructor.

Required to deep-copy the UniquePtr RayfireMesh object

Definition at line 65 of file integrated_function.C.

References GRINS::IntegratedFunction< Function >::_rayfire.

65  :
66  QoIBase(original),
67  _p_level(original._p_level),
68  _f(original._f)
69  {
70  // call RayfireMesh copy constructor
71  this->_rayfire.reset(new RayfireMesh( *((original._rayfire).get())) );
72  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
unsigned int _p_level
Quadrature order.
QoIBase(const std::string &qoi_name)
Definition: qoi_base.C:39
SharedPtr< Function > _f
Pointer to the template class used for function evaluation.
template<typename Function>
GRINS::IntegratedFunction< Function >::IntegratedFunction ( )
private

User cannot call empty constructor.

Member Function Documentation

template<typename Function >
bool GRINS::IntegratedFunction< Function >::assemble_on_interior ( ) const
inlinevirtual

Does the QoI need an element interior assembly loop?

This is pure virtual to force to user to specify.

Implements GRINS::QoIBase.

Definition at line 132 of file integrated_function.h.

133  {
134  return true;
135  }
template<typename Function >
bool GRINS::IntegratedFunction< Function >::assemble_on_sides ( ) const
inlinevirtual

Does the QoI need a domain boundary assembly loop?

This is pure virtual to force to user to specify.

Implements GRINS::QoIBase.

Definition at line 139 of file integrated_function.h.

140  {
141  return false;
142  }
template<typename Function >
QoIBase * GRINS::IntegratedFunction< Function >::clone ( ) const
virtual

Required to provide clone (deep-copy) for adding QoI object to libMesh objects.

Implements GRINS::QoIBase.

Reimplemented in GRINS::SpectroscopicAbsorption.

Definition at line 75 of file integrated_function.C.

76  {
77  return new IntegratedFunction<Function>( *this );
78  }
template<typename Function >
void GRINS::IntegratedFunction< Function >::element_qoi ( AssemblyContext context,
const unsigned int  qoi_index 
)
virtual

Compute the qoi value.

Reimplemented from GRINS::QoIBase.

Definition at line 96 of file integrated_function.C.

98  {
99  const libMesh::Elem & original_elem = context.get_elem();
100  const libMesh::Elem * rayfire_elem = _rayfire->map_to_rayfire_elem(original_elem.id());
101 
102  // rayfire_elem will be NULL if the main_elem
103  // is not in the rayfire
104  if (rayfire_elem)
105  {
106  // create and init the quadrature base on the rayfire elem
107  libMesh::QGauss qbase(rayfire_elem->dim(),libMesh::Order(_p_level));
108  qbase.init(rayfire_elem->type(),libMesh::Order(_p_level));
109 
110  // need the QP coordinates and JxW
111  libMesh::UniquePtr< libMesh::FEGenericBase<libMesh::Real> > fe = libMesh::FEGenericBase<libMesh::Real>::build(rayfire_elem->dim(), libMesh::FEType(libMesh::FIRST, libMesh::LAGRANGE) );
112 
113  fe->attach_quadrature_rule( &qbase );
114  const std::vector<libMesh::Real> & JxW = fe->get_JxW();
115  const std::vector<libMesh::Point> & xyz = fe->get_xyz();
116 
117  fe->reinit(rayfire_elem);
118 
119  const unsigned int n_qpoints = fe->n_quadrature_points();
120 
121  libMesh::Number & qoi = context.get_qois()[qoi_index];
122 
123  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
124  qoi += this->qoi_value((*_f),context,xyz[qp])*JxW[qp];
125 
126  }
127  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
libMesh::Real qoi_value(Function &f, AssemblyContext &context, const libMesh::Point &xyz)
Compute the value of a QoI at a QP.
unsigned int _p_level
Quadrature order.
SharedPtr< Function > _f
Pointer to the template class used for function evaluation.
template<typename Function >
void GRINS::IntegratedFunction< Function >::element_qoi_derivative ( AssemblyContext context,
const unsigned int  qoi_index 
)
virtual

Compute the qoi derivative with respect to the solution.

Currently not implemented

Reimplemented from GRINS::QoIBase.

Definition at line 130 of file integrated_function.C.

132  {
133  const libMesh::Elem& original_elem = context.get_elem();
134  const libMesh::Elem* rayfire_elem = _rayfire->map_to_rayfire_elem(original_elem.id());
135 
136  // rayfire_elem will be NULL if the main_elem
137  // is not in the rayfire
138  if (rayfire_elem)
139  {
140  // create and init the quadrature base on the rayfire elem
141  libMesh::QGauss qbase(rayfire_elem->dim(),libMesh::Order(_p_level));
142  qbase.init(rayfire_elem->type(),libMesh::Order(_p_level));
143 
144  // need the QP coordinates and JxW
145  libMesh::UniquePtr< libMesh::FEGenericBase<libMesh::Real> > fe = libMesh::FEGenericBase<libMesh::Real>::build(rayfire_elem->dim(), libMesh::FEType(libMesh::FIRST, libMesh::LAGRANGE) );
146 
147  fe->attach_quadrature_rule( &qbase );
148  const std::vector<libMesh::Real> & JxW = fe->get_JxW();
149  const std::vector<libMesh::Point> & xyz = fe->get_xyz();
150 
151  fe->reinit(rayfire_elem);
152 
153  const unsigned int n_qpoints = fe->n_quadrature_points();
154 
155  // FIXME the inverse_map() is ugly AF
156  for (unsigned int qp = 0; qp != n_qpoints; qp++)
157  this->qoi_derivative((*_f),context,xyz[qp],JxW[qp],qoi_index);
158 
159  }
160  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
unsigned int _p_level
Quadrature order.
void qoi_derivative(Function &f, AssemblyContext &context, const libMesh::Point &qp_xyz, const libMesh::Real JxW, const unsigned int qoi_index)
Compute derivatiuves at QP.
SharedPtr< Function > _f
Pointer to the template class used for function evaluation.
template<typename Function>
const RayfireMesh& GRINS::IntegratedFunction< Function >::get_rayfire ( )
inline

Definition at line 103 of file integrated_function.h.

104  {
105  return *(_rayfire.get());
106  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
template<typename Function >
void GRINS::IntegratedFunction< Function >::init ( const GetPot &  input,
const MultiphysicsSystem system,
unsigned int  qoi_num 
)
virtual

Initializes the rayfire with the mesh from system.

Reimplemented from GRINS::QoIBase.

Definition at line 82 of file integrated_function.C.

85  {
86  _rayfire->init(system.get_mesh());
87  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.
template<typename Function>
void GRINS::IntegratedFunction< Function >::qoi_derivative ( Function &  f,
AssemblyContext context,
const libMesh::Point &  qp_xyz,
const libMesh::Real  JxW,
const unsigned int  qoi_index 
)
private

Compute derivatiuves at QP.

template<>
void GRINS::IntegratedFunction< FEMFunctionAndDerivativeBase< libMesh::Real > >::qoi_derivative ( FEMFunctionAndDerivativeBase< libMesh::Real > &  f,
AssemblyContext context,
const libMesh::Point &  qp_xyz,
const libMesh::Real  JxW,
const unsigned int  qoi_index 
)
private

Definition at line 177 of file integrated_function.C.

References GRINS::FEMFunctionAndDerivativeBase< Output >::derivatives().

179  {
180  f.derivatives(context,qp_xyz,JxW,qoi_index);
181  }
virtual void derivatives(libMesh::FEMContext &context, const libMesh::Point &qp_xyz, const libMesh::Real &JxW, const unsigned int qoi_index, const libMesh::Real time=0.)=0
Function Derivative Evaluation.
template<>
void GRINS::IntegratedFunction< libMesh::FunctionBase< libMesh::Real > >::qoi_derivative ( libMesh::FunctionBase< libMesh::Real > &  ,
AssemblyContext ,
const libMesh::Point &  ,
const libMesh::Real  ,
const unsigned  int 
)
private

Definition at line 184 of file integrated_function.C.

186  {
187  // derivatives are always zero for FunctionBase
188  }
template<typename Function>
libMesh::Real GRINS::IntegratedFunction< Function >::qoi_value ( Function &  f,
AssemblyContext context,
const libMesh::Point &  xyz 
)
private

Compute the value of a QoI at a QP.

template<>
libMesh::Real GRINS::IntegratedFunction< FEMFunctionAndDerivativeBase< libMesh::Real > >::qoi_value ( FEMFunctionAndDerivativeBase< libMesh::Real > &  f,
AssemblyContext context,
const libMesh::Point &  xyz 
)
private

Definition at line 164 of file integrated_function.C.

165  {
166  return f(context,xyz);
167  }
template<>
libMesh::Real GRINS::IntegratedFunction< libMesh::FunctionBase< libMesh::Real > >::qoi_value ( libMesh::FunctionBase< libMesh::Real > &  f,
AssemblyContext ,
const libMesh::Point &  xyz 
)
private

Definition at line 170 of file integrated_function.C.

171  {
172  return f(xyz);
173  }
template<typename Function >
void GRINS::IntegratedFunction< Function >::reinit ( MultiphysicsSystem system)
virtual

Reinitialize the rayfire.

Reimplemented from GRINS::QoIBase.

Definition at line 90 of file integrated_function.C.

91  {
92  _rayfire->reinit(system.get_mesh());
93  }
libMesh::UniquePtr< RayfireMesh > _rayfire
Pointer to RayfireMesh object.

Member Data Documentation

template<typename Function>
SharedPtr<Function> GRINS::IntegratedFunction< Function >::_f
private

Pointer to the template class used for function evaluation.

Definition at line 113 of file integrated_function.h.

template<typename Function>
unsigned int GRINS::IntegratedFunction< Function >::_p_level
private

Quadrature order.

Definition at line 110 of file integrated_function.h.

template<typename Function>
libMesh::UniquePtr<RayfireMesh> GRINS::IntegratedFunction< Function >::_rayfire
private

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