38 #include "libmesh/getpot.h"
39 #include "libmesh/fem_system.h"
40 #include "libmesh/quadrature.h"
68 this->parse_thermal_conductivity(input);
71 ( _scaling, input,
"QoI/NusseltNumber/scaling", 1.0 );
75 std::cerr <<
"Error: thermal conductivity for AverageNusseltNumber must be positive." << std::endl
76 <<
"Found k = " << _k << std::endl;
81 int num_bcs = input.vector_variable_size(
"QoI/NusseltNumber/bc_ids");
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;
91 for(
int i = 0; i < num_bcs; i++ )
93 _bc_ids.insert( input(
"QoI/NusseltNumber/bc_ids", -1, i ) );
101 libMesh::FEBase* T_fe;
103 context.get_side_fe<libMesh::Real>(this->
_temp_vars->
T(), T_fe);
112 const unsigned int qoi_index )
114 bool on_correct_side =
false;
116 for (std::set<libMesh::boundary_id_type>::const_iterator
id =
118 if( context.has_side_boundary_id( (*
id) ) )
120 on_correct_side =
true;
124 if (!on_correct_side)
127 libMesh::FEBase* side_fe;
128 context.get_side_fe<libMesh::Real>(this->
_temp_vars->
T(), side_fe);
130 const std::vector<libMesh::Real> &JxW = side_fe->get_JxW();
132 const std::vector<libMesh::Point>& normals = side_fe->get_normals();
134 unsigned int n_qpoints = context.get_side_qrule().n_points();
136 libMesh::Number& qoi = context.get_qois()[qoi_index];
140 for (
unsigned int qp = 0; qp != n_qpoints; qp++)
143 libMesh::Gradient grad_T = 0.0;
144 context.side_gradient(this->
_temp_vars->
T(), qp, grad_T);
147 qoi += (this->
_scaling)*(this->
_k)*(grad_T*normals[qp])*JxW[qp];
153 const unsigned int qoi_index )
156 for( std::set<libMesh::boundary_id_type>::const_iterator
id =
_bc_ids.begin();
159 if( context.has_side_boundary_id( (*
id) ) )
161 libMesh::FEBase* T_side_fe;
162 context.get_side_fe<libMesh::Real>(this->
_temp_vars->
T(), T_side_fe);
164 const std::vector<libMesh::Real> &JxW = T_side_fe->get_JxW();
166 const std::vector<libMesh::Point>& normals = T_side_fe->get_normals();
168 unsigned int n_qpoints = context.get_side_qrule().n_points();
170 const unsigned int n_T_dofs = context.get_dof_indices(this->
_temp_vars->
T()).size();
172 const std::vector<std::vector<libMesh::Gradient> >& T_gradphi = T_side_fe->get_dphi();
174 libMesh::DenseSubVector<libMesh::Number>& dQ_dT =
175 context.get_qoi_derivatives(qoi_index, this->
_temp_vars->
T());
178 for (
unsigned int qp = 0; qp != n_qpoints; qp++)
181 libMesh::Gradient grad_T = 0.0;
182 context.side_gradient(this->
_temp_vars->
T(), qp, grad_T);
187 for(
unsigned int i = 0; i != n_T_dofs; i++ )
189 dQ_dT(i) +=
_scaling*
_k*T_gradphi[i][qp]*normals[qp]*JxW[qp];
203 std::string material = input(
"QoI/NusseltNumber/material",
"NoMaterial!");
206 "Materials/"+material+
"/ThermalConductivity/value",
207 "QoI/NusseltNumber/thermal_conductivity");
210 if( input.have_variable(
"QoI/NusseltNumber/thermal_conductivity") )
213 "ThermalConductivity/value" );
216 (
_k, input,
"QoI/NusseltNumber/thermal_conductivity", -1.0 );
219 else if( input.have_variable(
"Materials/"+material+
"/ThermalConductivity/value") )
222 if( input(
"Materials/"+material+
"/ThermalConductivity/model",
"DIE!")
223 != std::string(
"constant") )
225 libmesh_error_msg(
"ERROR: Only constant ThermalConductivity model supported in NusseltNumber!");
229 (
_k, input,
"Materials/"+material+
"/ThermalConductivity/value", -1.0 );
232 libmesh_error_msg(
"ERROR: Could not find valid thermal conducitivity value!");
236 libmesh_error_msg(
"ERROR: Detected non-positive thermal conductivity!");
virtual void set_parameter(libMesh::Number ¶m_variable, const GetPot &input, const std::string ¶m_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
virtual void init_context(AssemblyContext &context)
virtual ~AverageNusseltNumber()
libMesh::Real _scaling
Scaling constant.
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...
virtual QoIBase * clone() const
Clone this QoI.
static std::string temp_variable_name(const GetPot &input, const std::string &subsection_name, const SECTION_TYPE section_type)
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.
libMesh::Real _k
Thermal conductivity.
virtual void side_qoi(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi value on the domain boundary.
Interface with libMesh for solving Multiphysics problems.
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.
virtual void move_parameter(const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
When cloning an object, we need to update parameter pointers.
void parse_thermal_conductivity(const GetPot &input)