GRINS-0.7.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-2016 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"
33 
34 // libMesh
35 #include "libmesh/getpot.h"
36 #include "libmesh/fem_system.h"
37 #include "libmesh/quadrature.h"
38 
39 namespace GRINS
40 {
41  AverageNusseltNumber::AverageNusseltNumber( const std::string& qoi_name )
42  : QoIBase(qoi_name)
43  {
44  return;
45  }
46 
48  {
49  return;
50  }
51 
53  {
54  AverageNusseltNumber *returnval = new AverageNusseltNumber( *this );
55  returnval->move_parameter(_k, returnval->_k);
56  returnval->move_parameter(_scaling, returnval->_scaling);
57  return returnval;
58  }
59 
61  (const GetPot& input,
62  const MultiphysicsSystem& system,
63  unsigned int /*qoi_num*/ )
64  {
65  this->parse_thermal_conductivity(input);
66 
67  this->set_parameter
68  ( _scaling, input, "QoI/NusseltNumber/scaling", 1.0 );
69 
70  if( this->_k < 0.0 )
71  {
72  std::cerr << "Error: thermal conductivity for AverageNusseltNumber must be positive." << std::endl
73  << "Found k = " << _k << std::endl;
74  libmesh_error();
75  }
76 
77  // Read boundary ids for which we want to compute
78  int num_bcs = input.vector_variable_size("QoI/NusseltNumber/bc_ids");
79 
80  if( num_bcs <= 0 )
81  {
82  std::cerr << "Error: Must specify at least one boundary id to compute"
83  << " average Nusselt number." << std::endl
84  << "Found: " << num_bcs << std::endl;
85  libmesh_error();
86  }
87 
88  for( int i = 0; i < num_bcs; i++ )
89  {
90  _bc_ids.insert( input("QoI/NusseltNumber/bc_ids", -1, i ) );
91  }
92 
93  // Grab temperature variable index
94  std::string T_var_name = input( "Physics/VariableNames/Temperature",
96 
97  this->_T_var = system.variable_number(T_var_name);
98 
99  return;
100  }
101 
103  {
104  libMesh::FEBase* T_fe;
105 
106  context.get_side_fe<libMesh::Real>(this->_T_var, T_fe);
107 
108  T_fe->get_dphi();
109  T_fe->get_JxW();
110 
111  return;
112  }
113 
115  const unsigned int qoi_index )
116  {
117  bool on_correct_side = false;
118 
119  for (std::set<libMesh::boundary_id_type>::const_iterator id =
120  _bc_ids.begin(); id != _bc_ids.end(); id++ )
121  if( context.has_side_boundary_id( (*id) ) )
122  {
123  on_correct_side = true;
124  break;
125  }
126 
127  if (!on_correct_side)
128  return;
129 
130  libMesh::FEBase* side_fe;
131  context.get_side_fe<libMesh::Real>(this->_T_var, side_fe);
132 
133  const std::vector<libMesh::Real> &JxW = side_fe->get_JxW();
134 
135  const std::vector<libMesh::Point>& normals = side_fe->get_normals();
136 
137  unsigned int n_qpoints = context.get_side_qrule().n_points();
138 
139  libMesh::Number& qoi = context.get_qois()[qoi_index];
140 
141  // Loop over quadrature points
142 
143  for (unsigned int qp = 0; qp != n_qpoints; qp++)
144  {
145  // Get the solution value at the quadrature point
146  libMesh::Gradient grad_T = 0.0;
147  context.side_gradient(this->_T_var, qp, grad_T);
148 
149  // Update the elemental increment dR for each qp
150  qoi += (this->_scaling)*(this->_k)*(grad_T*normals[qp])*JxW[qp];
151 
152  } // quadrature loop
153  }
154 
156  const unsigned int qoi_index )
157  {
158 
159  for( std::set<libMesh::boundary_id_type>::const_iterator id = _bc_ids.begin();
160  id != _bc_ids.end(); id++ )
161  {
162  if( context.has_side_boundary_id( (*id) ) )
163  {
164  libMesh::FEBase* T_side_fe;
165  context.get_side_fe<libMesh::Real>(this->_T_var, T_side_fe);
166 
167  const std::vector<libMesh::Real> &JxW = T_side_fe->get_JxW();
168 
169  const std::vector<libMesh::Point>& normals = T_side_fe->get_normals();
170 
171  unsigned int n_qpoints = context.get_side_qrule().n_points();
172 
173  const unsigned int n_T_dofs = context.get_dof_indices(_T_var).size();
174 
175  const std::vector<std::vector<libMesh::Gradient> >& T_gradphi = T_side_fe->get_dphi();
176 
177  libMesh::DenseSubVector<libMesh::Number>& dQ_dT =
178  context.get_qoi_derivatives(qoi_index, _T_var);
179 
180  // Loop over quadrature points
181  for (unsigned int qp = 0; qp != n_qpoints; qp++)
182  {
183  // Get the solution value at the quadrature point
184  libMesh::Gradient grad_T = 0.0;
185  context.side_gradient(this->_T_var, qp, grad_T);
186 
187  // Update the elemental increment dR for each qp
188  //qoi += (this->_scaling)*(this->_k)*(grad_T*normals[qp])*JxW[qp];
189 
190  for( unsigned int i = 0; i != n_T_dofs; i++ )
191  {
192  dQ_dT(i) += _scaling*_k*T_gradphi[i][qp]*normals[qp]*JxW[qp];
193  }
194 
195  } // quadrature loop
196 
197  } // end check on boundary id
198 
199  }
200 
201  return;
202  }
203 
205  {
206  std::string material = input("QoI/NusseltNumber/material", "NoMaterial!");
207 
209  "Materials/"+material+"/ThermalConductivity/value",
210  "QoI/NusseltNumber/thermal_conductivity");
211 
212  // Parse the old version
213  if( input.have_variable("QoI/NusseltNumber/thermal_conductivity") )
214  {
215  MaterialsParsing::dep_input_warning( "QoI/NusseltNumber/thermal_conductivity",
216  "ThermalConductivity/value" );
217 
218  this->set_parameter
219  ( _k, input, "QoI/NusseltNumber/thermal_conductivity", -1.0 );
220  }
221  // Parse new version
222  else if( input.have_variable("Materials/"+material+"/ThermalConductivity/value") )
223  {
224  // This is currently only valid for constant thermal conductivity models
225  if( input("Materials/"+material+"/ThermalConductivity/model", "DIE!")
226  != std::string("constant") )
227  {
228  libmesh_error_msg("ERROR: Only constant ThermalConductivity model supported in NusseltNumber!");
229  }
230 
231  this->set_parameter
232  ( _k, input, "Materials/"+material+"/ThermalConductivity/value", -1.0 );
233  }
234  else
235  libmesh_error_msg("ERROR: Could not find valid thermal conducitivity value!");
236 
237  // thermal conducivity should be positive
238  if( _k <= 0.0 )
239  libmesh_error_msg("ERROR: Detected non-positive thermal conductivity!");
240  }
241 
242 } //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.
virtual void init_context(AssemblyContext &context)
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.
const std::string T_var_name_default
temperature
virtual QoIBase * clone() const
Clone this QoI.
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.
VariableIndex _T_var
Temperature variable index.
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 Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10