GRINS-0.8.0
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
GRINS::AverageNusseltNumber Class Reference

#include <average_nusselt_number.h>

Inheritance diagram for GRINS::AverageNusseltNumber:
Inheritance graph
[legend]
Collaboration diagram for GRINS::AverageNusseltNumber:
Collaboration graph
[legend]

Public Member Functions

 AverageNusseltNumber (const std::string &qoi_name)
 
virtual ~AverageNusseltNumber ()
 
virtual QoIBaseclone () const
 Clone this QoI. 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 side_qoi (AssemblyContext &context, const unsigned int qoi_index)
 Compute the qoi value on the domain boundary. More...
 
virtual void side_qoi_derivative (AssemblyContext &context, const unsigned int qoi_index)
 Compute the qoi derivative with respect to the solution on the domain boundary. More...
 
virtual void init (const GetPot &input, const MultiphysicsSystem &system, unsigned int qoi_num)
 Method to allow QoI to cache any system information needed for QoI calculation, for example, solution variable indices. More...
 
virtual void init_context (AssemblyContext &context)
 
- Public Member Functions inherited from GRINS::QoIBase
 QoIBase (const std::string &qoi_name)
 
virtual ~QoIBase ()
 
virtual void reinit (MultiphysicsSystem &)
 Reinitialize QoI. More...
 
virtual void element_qoi (AssemblyContext &, const unsigned int)
 Compute the qoi value for element interiors. More...
 
virtual void element_qoi_derivative (AssemblyContext &, const unsigned int)
 Compute the qoi derivative with respect to the solution on element interiors. 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...
 

Protected Member Functions

void parse_thermal_conductivity (const GetPot &input)
 

Protected Attributes

libMesh::Real _k
 Thermal conductivity. More...
 
const PrimitiveTempFEVariables_temp_vars
 
std::set< libMesh::boundary_id_type > _bc_ids
 List of boundary ids for which we want to compute this QoI. More...
 
libMesh::Real _scaling
 Scaling constant. More...
 
- Protected Attributes inherited from GRINS::QoIBase
std::string _qoi_name
 
libMesh::Number _qoi_value
 

Private Member Functions

 AverageNusseltNumber ()
 

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

Definition at line 37 of file average_nusselt_number.h.

Constructor & Destructor Documentation

GRINS::AverageNusseltNumber::AverageNusseltNumber ( const std::string &  qoi_name)

Definition at line 44 of file average_nusselt_number.C.

45  : QoIBase(qoi_name)
46  {
47  return;
48  }
QoIBase(const std::string &qoi_name)
Definition: qoi_base.C:39
GRINS::AverageNusseltNumber::~AverageNusseltNumber ( )
virtual

Definition at line 50 of file average_nusselt_number.C.

51  {
52  return;
53  }
GRINS::AverageNusseltNumber::AverageNusseltNumber ( )
private

Referenced by clone().

Member Function Documentation

bool GRINS::AverageNusseltNumber::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 85 of file average_nusselt_number.h.

86  {
87  return false;
88  }
bool GRINS::AverageNusseltNumber::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 91 of file average_nusselt_number.h.

92  {
93  return true;
94  }
QoIBase * GRINS::AverageNusseltNumber::clone ( ) const
virtual

Clone this QoI.

We return a raw pointer, but it is expected for the user to take ownership and delete the object when done.

Implements GRINS::QoIBase.

Definition at line 55 of file average_nusselt_number.C.

References _k, _scaling, AverageNusseltNumber(), and GRINS::ParameterUser::move_parameter().

56  {
57  AverageNusseltNumber *returnval = new AverageNusseltNumber( *this );
58  returnval->move_parameter(_k, returnval->_k);
59  returnval->move_parameter(_scaling, returnval->_scaling);
60  return returnval;
61  }
libMesh::Real _scaling
Scaling constant.
libMesh::Real _k
Thermal conductivity.
void GRINS::AverageNusseltNumber::init ( const GetPot &  ,
const MultiphysicsSystem ,
unsigned int   
)
virtual

Method to allow QoI to cache any system information needed for QoI calculation, for example, solution variable indices.

Reimplemented from GRINS::QoIBase.

Definition at line 64 of file average_nusselt_number.C.

References GRINS::VariablesParsing::QOI, and GRINS::VariablesParsing::temp_variable_name().

67  {
68  this->parse_thermal_conductivity(input);
69 
70  this->set_parameter
71  ( _scaling, input, "QoI/NusseltNumber/scaling", 1.0 );
72 
73  if( this->_k < 0.0 )
74  {
75  std::cerr << "Error: thermal conductivity for AverageNusseltNumber must be positive." << std::endl
76  << "Found k = " << _k << std::endl;
77  libmesh_error();
78  }
79 
80  // Read boundary ids for which we want to compute
81  int num_bcs = input.vector_variable_size("QoI/NusseltNumber/bc_ids");
82 
83  if( num_bcs <= 0 )
84  {
85  std::cerr << "Error: Must specify at least one boundary id to compute"
86  << " average Nusselt number." << std::endl
87  << "Found: " << num_bcs << std::endl;
88  libmesh_error();
89  }
90 
91  for( int i = 0; i < num_bcs; i++ )
92  {
93  _bc_ids.insert( input("QoI/NusseltNumber/bc_ids", -1, i ) );
94  }
95 
96  _temp_vars = &GRINSPrivate::VariableWarehouse::get_variable_subclass<PrimitiveTempFEVariables>(VariablesParsing::temp_variable_name(input,std::string("NusseltNumber"),VariablesParsing::QOI));
97  }
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.
std::set< libMesh::boundary_id_type > _bc_ids
List of boundary ids for which we want to compute this QoI.
const PrimitiveTempFEVariables * _temp_vars
libMesh::Real _scaling
Scaling constant.
static std::string temp_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
libMesh::Real _k
Thermal conductivity.
void parse_thermal_conductivity(const GetPot &input)
void GRINS::AverageNusseltNumber::init_context ( AssemblyContext context)
virtual

Reimplemented from GRINS::QoIBase.

Definition at line 99 of file average_nusselt_number.C.

References _temp_vars, and GRINS::PrimitiveTempFEVariables::T().

100  {
101  libMesh::FEBase* T_fe;
102 
103  context.get_side_fe<libMesh::Real>(this->_temp_vars->T(), T_fe);
104 
105  T_fe->get_dphi();
106  T_fe->get_JxW();
107 
108  return;
109  }
const PrimitiveTempFEVariables * _temp_vars
VariableIndex T() const
void GRINS::AverageNusseltNumber::parse_thermal_conductivity ( const GetPot &  input)
protected

Definition at line 201 of file average_nusselt_number.C.

References _k, GRINS::MaterialsParsing::dep_input_warning(), GRINS::MaterialsParsing::duplicate_input_test(), and GRINS::ParameterUser::set_parameter().

202  {
203  std::string material = input("QoI/NusseltNumber/material", "NoMaterial!");
204 
206  "Materials/"+material+"/ThermalConductivity/value",
207  "QoI/NusseltNumber/thermal_conductivity");
208 
209  // Parse the old version
210  if( input.have_variable("QoI/NusseltNumber/thermal_conductivity") )
211  {
212  MaterialsParsing::dep_input_warning( "QoI/NusseltNumber/thermal_conductivity",
213  "ThermalConductivity/value" );
214 
215  this->set_parameter
216  ( _k, input, "QoI/NusseltNumber/thermal_conductivity", -1.0 );
217  }
218  // Parse new version
219  else if( input.have_variable("Materials/"+material+"/ThermalConductivity/value") )
220  {
221  // This is currently only valid for constant thermal conductivity models
222  if( input("Materials/"+material+"/ThermalConductivity/model", "DIE!")
223  != std::string("constant") )
224  {
225  libmesh_error_msg("ERROR: Only constant ThermalConductivity model supported in NusseltNumber!");
226  }
227 
228  this->set_parameter
229  ( _k, input, "Materials/"+material+"/ThermalConductivity/value", -1.0 );
230  }
231  else
232  libmesh_error_msg("ERROR: Could not find valid thermal conducitivity value!");
233 
234  // thermal conducivity should be positive
235  if( _k <= 0.0 )
236  libmesh_error_msg("ERROR: Detected non-positive thermal conductivity!");
237  }
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.
libMesh::Real _k
Thermal conductivity.
static void dep_input_warning(const std::string &old_option, const std::string &property)
Helper function for parsing/maintaing backward compatibility.
static void duplicate_input_test(const GetPot &input, const std::string &option1, const std::string &option2)
Helper function for parsing/maintaing backward compatibility.
void GRINS::AverageNusseltNumber::side_qoi ( AssemblyContext ,
const unsigned int   
)
virtual

Compute the qoi value on the domain boundary.

Override this method if your QoI is defined on the domain boundary

Reimplemented from GRINS::QoIBase.

Definition at line 111 of file average_nusselt_number.C.

References _bc_ids, _k, _scaling, _temp_vars, and GRINS::PrimitiveTempFEVariables::T().

113  {
114  bool on_correct_side = false;
115 
116  for (std::set<libMesh::boundary_id_type>::const_iterator id =
117  _bc_ids.begin(); id != _bc_ids.end(); id++ )
118  if( context.has_side_boundary_id( (*id) ) )
119  {
120  on_correct_side = true;
121  break;
122  }
123 
124  if (!on_correct_side)
125  return;
126 
127  libMesh::FEBase* side_fe;
128  context.get_side_fe<libMesh::Real>(this->_temp_vars->T(), side_fe);
129 
130  const std::vector<libMesh::Real> &JxW = side_fe->get_JxW();
131 
132  const std::vector<libMesh::Point>& normals = side_fe->get_normals();
133 
134  unsigned int n_qpoints = context.get_side_qrule().n_points();
135 
136  libMesh::Number& qoi = context.get_qois()[qoi_index];
137 
138  // Loop over quadrature points
139 
140  for (unsigned int qp = 0; qp != n_qpoints; qp++)
141  {
142  // Get the solution value at the quadrature point
143  libMesh::Gradient grad_T = 0.0;
144  context.side_gradient(this->_temp_vars->T(), qp, grad_T);
145 
146  // Update the elemental increment dR for each qp
147  qoi += (this->_scaling)*(this->_k)*(grad_T*normals[qp])*JxW[qp];
148 
149  } // quadrature loop
150  }
std::set< libMesh::boundary_id_type > _bc_ids
List of boundary ids for which we want to compute this QoI.
const PrimitiveTempFEVariables * _temp_vars
VariableIndex T() const
libMesh::Real _scaling
Scaling constant.
libMesh::Real _k
Thermal conductivity.
void GRINS::AverageNusseltNumber::side_qoi_derivative ( AssemblyContext ,
const unsigned int   
)
virtual

Compute the qoi derivative with respect to the solution on the domain boundary.

Override this method if your QoI is defined on the domain boundary

Reimplemented from GRINS::QoIBase.

Definition at line 152 of file average_nusselt_number.C.

References _bc_ids, _k, _scaling, _temp_vars, and GRINS::PrimitiveTempFEVariables::T().

154  {
155 
156  for( std::set<libMesh::boundary_id_type>::const_iterator id = _bc_ids.begin();
157  id != _bc_ids.end(); id++ )
158  {
159  if( context.has_side_boundary_id( (*id) ) )
160  {
161  libMesh::FEBase* T_side_fe;
162  context.get_side_fe<libMesh::Real>(this->_temp_vars->T(), T_side_fe);
163 
164  const std::vector<libMesh::Real> &JxW = T_side_fe->get_JxW();
165 
166  const std::vector<libMesh::Point>& normals = T_side_fe->get_normals();
167 
168  unsigned int n_qpoints = context.get_side_qrule().n_points();
169 
170  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars->T()).size();
171 
172  const std::vector<std::vector<libMesh::Gradient> >& T_gradphi = T_side_fe->get_dphi();
173 
174  libMesh::DenseSubVector<libMesh::Number>& dQ_dT =
175  context.get_qoi_derivatives(qoi_index, this->_temp_vars->T());
176 
177  // Loop over quadrature points
178  for (unsigned int qp = 0; qp != n_qpoints; qp++)
179  {
180  // Get the solution value at the quadrature point
181  libMesh::Gradient grad_T = 0.0;
182  context.side_gradient(this->_temp_vars->T(), qp, grad_T);
183 
184  // Update the elemental increment dR for each qp
185  //qoi += (this->_scaling)*(this->_k)*(grad_T*normals[qp])*JxW[qp];
186 
187  for( unsigned int i = 0; i != n_T_dofs; i++ )
188  {
189  dQ_dT(i) += _scaling*_k*T_gradphi[i][qp]*normals[qp]*JxW[qp];
190  }
191 
192  } // quadrature loop
193 
194  } // end check on boundary id
195 
196  }
197 
198  return;
199  }
std::set< libMesh::boundary_id_type > _bc_ids
List of boundary ids for which we want to compute this QoI.
const PrimitiveTempFEVariables * _temp_vars
VariableIndex T() const
libMesh::Real _scaling
Scaling constant.
libMesh::Real _k
Thermal conductivity.

Member Data Documentation

std::set<libMesh::boundary_id_type> GRINS::AverageNusseltNumber::_bc_ids
protected

List of boundary ids for which we want to compute this QoI.

Definition at line 73 of file average_nusselt_number.h.

Referenced by side_qoi(), and side_qoi_derivative().

libMesh::Real GRINS::AverageNusseltNumber::_k
protected

Thermal conductivity.

Definition at line 68 of file average_nusselt_number.h.

Referenced by clone(), parse_thermal_conductivity(), side_qoi(), and side_qoi_derivative().

libMesh::Real GRINS::AverageNusseltNumber::_scaling
protected

Scaling constant.

Definition at line 76 of file average_nusselt_number.h.

Referenced by clone(), side_qoi(), and side_qoi_derivative().

const PrimitiveTempFEVariables* GRINS::AverageNusseltNumber::_temp_vars
protected

Definition at line 70 of file average_nusselt_number.h.

Referenced by init_context(), side_qoi(), and side_qoi_derivative().


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

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