GRINS-0.8.0
test_turbulent_channel.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 #include "grins_config.h"
27 
28 #include <iostream>
29 
30 // GRINS
31 #include "grins/runner.h"
32 #include "grins/mesh_builder.h"
33 #include "grins/multiphysics_sys.h"
35 #include "grins/var_typedefs.h"
36 
37 //libMesh
38 #include "libmesh/dirichlet_boundaries.h"
39 #include "libmesh/dof_map.h"
40 #include "libmesh/fe_base.h"
41 #include "libmesh/fe_interface.h"
42 #include "libmesh/mesh_function.h"
43 #include "libmesh/serial_mesh.h"
44 #include "libmesh/edge_edge2.h"
45 #include "libmesh/mesh_generation.h"
46 #include "libmesh/linear_implicit_system.h"
47 #include "libmesh/gmv_io.h"
48 #include "libmesh/exact_solution.h"
49 #include "libmesh/zero_function.h"
50 #include "libmesh/composite_function.h"
51 #include "libmesh/zero_function.h"
52 
53 namespace GRINS
54 {
55  class MultiphysicsSystem;
56 }
57 
58 void test_error_norm( libMesh::ExactSolution& exact_sol,
59  const std::string& system_name,
60  const std::string& var,
61  const std::string& norm,
62  const double tol,
63  int& return_flag );
64 
65 namespace GRINSTesting
66 {
67 
68  class TurbBoundFuncBase : public libMesh::FunctionBase<libMesh::Number>
69  {
70  public:
71  TurbBoundFuncBase (libMesh::MeshFunction* turbulent_bc_values)
72  : _turbulent_bc_values(turbulent_bc_values)
73  { this->_initialized = true; }
74 
75  virtual libMesh::Number operator() (const libMesh::Point&, const libMesh::Real = 0)
76  { libmesh_not_implemented(); }
77 
78  virtual libMesh::Number compute_final_value( const libMesh::DenseVector<libMesh::Number>& u_nu_values ) =0;
79 
80  virtual void operator() (const libMesh::Point& p,
81  const libMesh::Real t,
82  libMesh::DenseVector<libMesh::Number>& output)
83  {
84  output.resize(1);
85  output.zero();
86 
87  // Since the turbulent_bc_values object has a solution from a 1-d problem, we have to zero out the y coordinate of p
88  libMesh::Point p_copy(p);
89  // Also, the 1-d solution provided is on the domain [0, 1] on the x axis and we need to map this to the corresponding point on the y axis
90  p_copy(0) = p_copy(1);
91  p_copy(1)= 0.0;
92  // Also, the 1-d solution provided is actually a symmetry solution, so we have to make the following map
93  // x_GRINS < 0.5 => x_meshfunction = 2*x_GRINS , x_GRINS >= 0.5 => x_GRINS = 1 - x_GRINS, x_meshfunction = 2*x_GRINS
94  if(p_copy(0) > 0.5)
95  {
96  p_copy(0) = 1 - p_copy(0);
97  }
98  p_copy(0) = 2*p_copy(0);
99 
100  libMesh::DenseVector<libMesh::Number> u_nu_values;
101  _turbulent_bc_values->operator()(p_copy, t, u_nu_values);
102 
103  output(0) = this->compute_final_value(u_nu_values);
104  }
105 
106  protected:
107 
108  libMesh::MeshFunction* _turbulent_bc_values;
109 
110  };
111 
112  // Class to construct the Dirichlet boundary object and operator for the inlet u velocity and nu profiles
114  {
115  public:
116  TurbulentBdyFunctionU (libMesh::MeshFunction* turbulent_bc_values)
117  : TurbBoundFuncBase(turbulent_bc_values)
118  {}
119 
120  virtual libMesh::Number compute_final_value( const libMesh::DenseVector<libMesh::Number>& u_nu_values )
121  { return u_nu_values(0)/21.995539; }
122 
123  virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> > clone() const
124  { return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> > (new TurbulentBdyFunctionU(_turbulent_bc_values)); }
125  };
126 
127  // Class to construct the Dirichlet boundary object and operator for the inlet u velocity and nu profiles
129  {
130  public:
131  TurbulentBdyFunctionNu (libMesh::MeshFunction* turbulent_bc_values)
132  : TurbBoundFuncBase(turbulent_bc_values)
133  {}
134 
135  virtual libMesh::Number compute_final_value( const libMesh::DenseVector<libMesh::Number>& u_nu_values )
136  { return u_nu_values(1)/(2.0*21.995539); }
137 
138  virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> > clone() const
139  { return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> > (new TurbulentBdyFunctionNu(_turbulent_bc_values)); }
140 
141  };
142 
143  class SATurbBCFactoryBase : public GRINS::DirichletBCFactoryFunctionBase<libMesh::FunctionBase<libMesh::Number> >
144  {
145  public:
146  SATurbBCFactoryBase( const std::string& bc_type_name )
147  : DirichletBCFactoryFunctionBase<libMesh::FunctionBase<libMesh::Number> >(bc_type_name)
148  {}
149 
150  static void set_turb_bc_values( libMesh::MeshFunction* turbulent_bc_values )
151  { _turbulent_bc_values = turbulent_bc_values; }
152 
153  protected:
154 
155  static libMesh::MeshFunction* _turbulent_bc_values;
156  };
157 
158  libMesh::MeshFunction* SATurbBCFactoryBase::_turbulent_bc_values = NULL;
159 
161  {
162  public:
163 
164  SATurbUBCFactory( const std::string& bc_type_name )
165  : SATurbBCFactoryBase(bc_type_name)
166  {}
167 
168  protected:
169 
170  virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
171  build_func( const GetPot& /*input*/,
173  std::vector<std::string>& var_names,
174  const std::string& /*section*/ )
175  {
176  libmesh_assert_equal_to(var_names.size(), 2);
177  libmesh_assert_equal_to(var_names[0], std::string("u"));
178  libmesh_assert_equal_to(var_names[1], std::string("v"));
179 
180  libMesh::UniquePtr<libMesh::CompositeFunction<libMesh::Number> >
181  composite_func( new libMesh::CompositeFunction<libMesh::Number> );
182 
183  libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
185 
186  std::vector<GRINS::VariableIndex> vars(1);
187  vars[0] = system.variable_number(var_names[0]);
188  composite_func->attach_subfunction(*bound_func, vars);
189 
190  vars[0] = system.variable_number(var_names[1]);
191  composite_func->attach_subfunction(libMesh::ZeroFunction<libMesh::Number>(), vars);
192 
193  return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >( composite_func.release() );
194  }
195  };
196 
198  {
199  public:
200 
201  SATurbNuBCFactory( const std::string& bc_type_name )
202  : SATurbBCFactoryBase(bc_type_name)
203  {}
204 
205  protected:
206 
207  virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
208  build_func( const GetPot& /*input*/,
210  std::vector<std::string>& var_names,
211  const std::string& /*section*/ )
212  {
213  libmesh_assert_equal_to(var_names.size(), 1);
214  libmesh_assert_equal_to(var_names[0], std::string("nu"));
215 
216  libMesh::UniquePtr<libMesh::CompositeFunction<libMesh::Number> >
217  composite_func( new libMesh::CompositeFunction<libMesh::Number> );
218 
219  libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
220  bound_func( new TurbulentBdyFunctionNu(_turbulent_bc_values) );
221 
222  std::vector<GRINS::VariableIndex> vars(1);
223  vars[0] = system.variable_number(var_names[0]);
224  composite_func->attach_subfunction(*bound_func, vars);
225 
226  return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >( composite_func.release() );
227  }
228  };
229 
230 } // end namespace GRINSTesting
231 
232 int main(int argc, char* argv[])
233 {
234  // Factories needed for run
235  GRINSTesting::SATurbUBCFactory grins_factory_testing_turb_u_bc("testing_turb_u");
236  GRINSTesting::SATurbNuBCFactory grins_factory_testing_turb_nu_bc("testing_turb_nu");
237 
238  // Need only
239  GRINS::Runner grins(argc,argv);
240 
241  const GetPot & command_line = grins.get_command_line();
242 
243  const GetPot & inputfile = grins.get_input_file();
244 
245  // Don't flag our command-line-specific variables as UFOs later
246  inputfile.have_variable("mesh-1d");
247  inputfile.have_variable("data-1d");
248  inputfile.have_variable("soln-data");
249  inputfile.have_variable("vars");
250  inputfile.have_variable("norms");
251  inputfile.have_variable("tol");
252 
253  const libMesh::LibMeshInit & libmesh_init = grins.get_libmesh_init();
254 
255  // Build a 1-d turbulent_bc_system to get the bc data from files
256  libMesh::SerialMesh mesh(libmesh_init.comm());
257 
258  std::string oned_mesh = command_line("mesh-1d", "DIE!");
259  mesh.read(oned_mesh);
260 
261  // And an EquationSystems to run on it
262  libMesh::EquationSystems equation_systems (mesh);
263 
264  std::string oned_data = command_line("data-1d", "DIE!");
265  equation_systems.read(oned_data, libMesh::XdrMODE::READ,
266  libMesh::EquationSystems::READ_HEADER |
267  libMesh::EquationSystems::READ_DATA |
268  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
269 
270  // Get a reference to the system so that we can call update() on it
271  libMesh::System & turbulent_bc_system = equation_systems.get_system<libMesh::System>("flow");
272 
273  // We need to call update to put system in a consistent state
274  // with the solution that was read in
275  turbulent_bc_system.update();
276 
277  // Print information about the system to the screen.
278  equation_systems.print_info();
279 
280  // Prepare a global solution and a MeshFunction of the Turbulent system
281  libMesh::UniquePtr<libMesh::MeshFunction> turbulent_bc_values;
282 
283  libMesh::UniquePtr<libMesh::NumericVector<libMesh::Number> > turbulent_bc_soln = libMesh::NumericVector<libMesh::Number>::build(equation_systems.comm());
284 
285  std::vector<libMesh::Number> flow_soln;
286 
287  turbulent_bc_system.update_global_solution(flow_soln);
288 
289  turbulent_bc_soln->init(turbulent_bc_system.solution->size(), true, libMesh::SERIAL);
290 
291  (*turbulent_bc_soln) = flow_soln;
292 
293  std::vector<unsigned int>turbulent_bc_system_variables;
294  turbulent_bc_system_variables.push_back(0);
295  turbulent_bc_system_variables.push_back(1);
296 
297  turbulent_bc_values = libMesh::UniquePtr<libMesh::MeshFunction>
298  (new libMesh::MeshFunction(equation_systems,
299  *turbulent_bc_soln,
300  turbulent_bc_system.get_dof_map(),
301  turbulent_bc_system_variables ));
302 
303  turbulent_bc_values->init();
304 
305  GRINSTesting::SATurbBCFactoryBase::set_turb_bc_values( turbulent_bc_values.get() );
306 
307  // Initialize
308  grins.init();
309 
310  // Solve
311  grins.run();
312 
313  GRINS::Simulation & sim = grins.get_simulation();
314 
315  // Get equation systems to create ExactSolution object
316  GRINS::SharedPtr<libMesh::EquationSystems> es = sim.get_equation_system();
317 
318  // Create Exact solution object and attach exact solution quantities
319  //libMesh::ExactSolution exact_sol(*es);
320 
321  // Create Exact solution object and attach exact solution quantities
322  libMesh::ExactSolution exact_sol(*es);
323 
324  libMesh::EquationSystems es_ref( es->get_mesh() );
325 
326  // Filename of file where comparison solution is stashed
327  std::string solution_file = command_line("soln-data", "DIE!");
328  es_ref.read( solution_file, libMesh::XdrMODE::DECODE,
329  libMesh::EquationSystems::READ_HEADER |
330  libMesh::EquationSystems::READ_DATA |
331  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
332 
333  exact_sol.attach_reference_solution( &es_ref );
334 
335  // Now grab the variables for which we want to compare
336  unsigned int n_vars = command_line.vector_variable_size("vars");
337  std::vector<std::string> vars(n_vars);
338  for( unsigned int v = 0; v < n_vars; v++ )
339  {
340  vars[v] = command_line("vars", "DIE!", v);
341  }
342 
343  // Now grab the norms to compute for each variable error
344  unsigned int n_norms = command_line.vector_variable_size("norms");
345  std::vector<std::string> norms(n_norms);
346  for( unsigned int n = 0; n < n_norms; n++ )
347  {
348  norms[n] = command_line("norms", "DIE!", n);
349  if( norms[n] != std::string("L2") &&
350  norms[n] != std::string("H1") )
351  {
352  std::cerr << "ERROR: Invalid norm input " << norms[n] << std::endl
353  << " Valid values are: L2" << std::endl
354  << " H1" << std::endl;
355  }
356  }
357 
358  const std::string& system_name = sim.get_multiphysics_system_name();
359 
360  // Now compute error for each variable
361  for( unsigned int v = 0; v < n_vars; v++ )
362  {
363  exact_sol.compute_error(system_name, vars[v]);
364  }
365 
366  int return_flag = 0;
367 
368  double tol = command_line("tol", 1.0e-10);
369 
370  // Now test error for each variable, for each norm
371  for( unsigned int v = 0; v < n_vars; v++ )
372  {
373  for( unsigned int n = 0; n < n_norms; n++ )
374  {
375  test_error_norm( exact_sol, system_name, vars[v], norms[n], tol, return_flag );
376  }
377  }
378 
379  return return_flag;
380 }
381 
382 void test_error_norm( libMesh::ExactSolution& exact_sol,
383  const std::string& system_name,
384  const std::string& var,
385  const std::string& norm,
386  const double tol,
387  int& return_flag )
388 {
389  // We don't set return_flag unless we are setting it 1
390  // since this function gets called multiple times and we don't
391  // want to overwrite a previous "fail" (return_flag = 1) with
392  // a "pass" (return_flag = 0)
393 
394  double error = 0.0;
395 
396  std::cout << "==========================================================" << std::endl
397  << "Checking variable " << var << " using error norm " << norm << " with tol " << tol << "...";
398 
399  if( norm == std::string("L2") )
400  {
401  error = exact_sol.l2_error(system_name, var);
402  }
403  else if( norm == std::string("H1") )
404  {
405  error = exact_sol.h1_error(system_name, var);
406  }
407  else
408  {
409  std::cerr << "ERROR: Invalid norm " << norm << std::endl;
410  exit(1);
411  }
412 
413  if( error > tol )
414  {
415  return_flag = 1;
416 
417  std::cerr << "Tolerance exceeded for generic regression test!" << std::endl
418  << "tolerance = " << tol << std::endl
419  << "norm of error = " << error << std::endl
420  << "norm type = " << norm << std::endl
421  << "var = " << var << std::endl;
422  }
423  else
424  {
425  std::cout << "PASSED!" << std::endl
426  << "==========================================================" << std::endl;
427  }
428 
429  return;
430 }
SATurbBCFactoryBase(const std::string &bc_type_name)
TurbBoundFuncBase(libMesh::MeshFunction *turbulent_bc_values)
void init()
Initialize the Simulation objects.
Definition: runner.C:53
SharedPtr< libMesh::EquationSystems > get_equation_system()
Definition: simulation.C:393
virtual libMesh::UniquePtr< libMesh::FunctionBase< libMesh::Number > > build_func(const GetPot &, GRINS::MultiphysicsSystem &system, std::vector< std::string > &var_names, const std::string &)
Builds the FunctionBase object for boundary condition.
virtual libMesh::UniquePtr< libMesh::FunctionBase< libMesh::Number > > clone() const
libMesh::MeshFunction * _turbulent_bc_values
const GetPot & get_input_file() const
Definition: runner.h:60
virtual libMesh::UniquePtr< libMesh::FunctionBase< libMesh::Number > > build_func(const GetPot &, GRINS::MultiphysicsSystem &system, std::vector< std::string > &var_names, const std::string &)
Builds the FunctionBase object for boundary condition.
const libMesh::LibMeshInit & get_libmesh_init() const
Definition: runner.h:74
virtual libMesh::Number operator()(const libMesh::Point &, const libMesh::Real=0)
virtual libMesh::UniquePtr< libMesh::FunctionBase< libMesh::Number > > clone() const
virtual libMesh::Number compute_final_value(const libMesh::DenseVector< libMesh::Number > &u_nu_values)
static libMesh::MeshFunction * _turbulent_bc_values
GRINS namespace.
SATurbNuBCFactory(const std::string &bc_type_name)
TurbulentBdyFunctionU(libMesh::MeshFunction *turbulent_bc_values)
const GetPot & get_command_line() const
Definition: runner.h:63
const std::string & get_multiphysics_system_name() const
Definition: simulation.h:185
void test_error_norm(libMesh::ExactSolution &exact_sol, const std::string &system_name, const std::string &var, const std::string &norm, const double tol, int &return_flag)
virtual libMesh::Number compute_final_value(const libMesh::DenseVector< libMesh::Number > &u_nu_values)=0
Interface with libMesh for solving Multiphysics problems.
static void set_turb_bc_values(libMesh::MeshFunction *turbulent_bc_values)
TurbulentBdyFunctionNu(libMesh::MeshFunction *turbulent_bc_values)
Class to encapsulate initializing and running GRINS Simulation.
Definition: runner.h:47
SATurbUBCFactory(const std::string &bc_type_name)
int main(int argc, char *argv[])
virtual libMesh::Number compute_final_value(const libMesh::DenseVector< libMesh::Number > &u_nu_values)
void run()
Runs the simulation that was setup at construction time.
Definition: runner.C:104
Simulation & get_simulation()
Definition: runner.h:66

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