GRINS-0.8.0
average_nusselt_number.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
28 
29 // GRINS
30 #include "grins/multiphysics_sys.h"
31 #include "grins/assembly_context.h"
35 #include "grins/single_variable.h"
36 
37 // libMesh
38 #include "libmesh/getpot.h"
39 #include "libmesh/fem_system.h"
40 #include "libmesh/quadrature.h"
41 
42 namespace GRINS
43 {
44  AverageNusseltNumber::AverageNusseltNumber( const std::string& qoi_name )
45  : QoIBase(qoi_name)
46  {
47  return;
48  }
49 
51  {
52  return;
53  }
54 
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  }
62 
64  (const GetPot& input,
65  const MultiphysicsSystem& /*system*/,
66  unsigned int /*qoi_num*/ )
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  }
98 
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  }
110 
112  const unsigned int qoi_index )
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  }
151 
153  const unsigned int qoi_index )
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  }
200 
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  }
238 
239 } //namespace GRINS
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
virtual void init_context(AssemblyContext &context)
VariableIndex T() const
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...
GRINS namespace.
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)

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