26 #include "grins_config.h"
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"
55 class MultiphysicsSystem;
59 const std::string& system_name,
60 const std::string& var,
61 const std::string& norm,
73 { this->_initialized =
true; }
75 virtual libMesh::Number
operator() (
const libMesh::Point&,
const libMesh::Real = 0)
76 { libmesh_not_implemented(); }
78 virtual libMesh::Number
compute_final_value(
const libMesh::DenseVector<libMesh::Number>& u_nu_values ) =0;
81 const libMesh::Real t,
82 libMesh::DenseVector<libMesh::Number>& output)
88 libMesh::Point p_copy(p);
90 p_copy(0) = p_copy(1);
96 p_copy(0) = 1 - p_copy(0);
98 p_copy(0) = 2*p_copy(0);
100 libMesh::DenseVector<libMesh::Number> u_nu_values;
121 {
return u_nu_values(0)/21.995539; }
123 virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
clone()
const
136 {
return u_nu_values(1)/(2.0*21.995539); }
138 virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
clone()
const
170 virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
173 std::vector<std::string>& var_names,
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"));
180 libMesh::UniquePtr<libMesh::CompositeFunction<libMesh::Number> >
183 libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
186 std::vector<GRINS::VariableIndex> vars(1);
187 vars[0] = system.variable_number(var_names[0]);
188 composite_func->attach_subfunction(*bound_func, vars);
190 vars[0] = system.variable_number(var_names[1]);
191 composite_func->attach_subfunction(libMesh::ZeroFunction<libMesh::Number>(), vars);
193 return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >( composite_func.release() );
207 virtual libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
210 std::vector<std::string>& var_names,
213 libmesh_assert_equal_to(var_names.size(), 1);
214 libmesh_assert_equal_to(var_names[0], std::string(
"nu"));
216 libMesh::UniquePtr<libMesh::CompositeFunction<libMesh::Number> >
219 libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >
222 std::vector<GRINS::VariableIndex> vars(1);
223 vars[0] = system.variable_number(var_names[0]);
224 composite_func->attach_subfunction(*bound_func, vars);
226 return libMesh::UniquePtr<libMesh::FunctionBase<libMesh::Number> >( composite_func.release() );
232 int main(
int argc,
char* argv[])
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");
256 libMesh::SerialMesh mesh(libmesh_init.comm());
258 std::string oned_mesh = command_line(
"mesh-1d",
"DIE!");
259 mesh.read(oned_mesh);
262 libMesh::EquationSystems equation_systems (mesh);
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);
271 libMesh::System & turbulent_bc_system = equation_systems.get_system<libMesh::System>(
"flow");
275 turbulent_bc_system.update();
278 equation_systems.print_info();
281 libMesh::UniquePtr<libMesh::MeshFunction> turbulent_bc_values;
283 libMesh::UniquePtr<libMesh::NumericVector<libMesh::Number> > turbulent_bc_soln = libMesh::NumericVector<libMesh::Number>::build(equation_systems.comm());
285 std::vector<libMesh::Number> flow_soln;
287 turbulent_bc_system.update_global_solution(flow_soln);
289 turbulent_bc_soln->init(turbulent_bc_system.solution->size(),
true, libMesh::SERIAL);
291 (*turbulent_bc_soln) = flow_soln;
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);
297 turbulent_bc_values = libMesh::UniquePtr<libMesh::MeshFunction>
298 (
new libMesh::MeshFunction(equation_systems,
300 turbulent_bc_system.get_dof_map(),
301 turbulent_bc_system_variables ));
303 turbulent_bc_values->init();
322 libMesh::ExactSolution exact_sol(*es);
324 libMesh::EquationSystems es_ref( es->get_mesh() );
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);
333 exact_sol.attach_reference_solution( &es_ref );
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++ )
340 vars[v] = command_line(
"vars",
"DIE!", v);
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++ )
348 norms[n] = command_line(
"norms",
"DIE!", n);
349 if( norms[n] != std::string(
"L2") &&
350 norms[n] != std::string(
"H1") )
352 std::cerr <<
"ERROR: Invalid norm input " << norms[n] << std::endl
353 <<
" Valid values are: L2" << std::endl
354 <<
" H1" << std::endl;
361 for(
unsigned int v = 0; v < n_vars; v++ )
363 exact_sol.compute_error(system_name, vars[v]);
368 double tol = command_line(
"tol", 1.0e-10);
371 for(
unsigned int v = 0; v < n_vars; v++ )
373 for(
unsigned int n = 0; n < n_norms; n++ )
375 test_error_norm( exact_sol, system_name, vars[v], norms[n], tol, return_flag );
383 const std::string& system_name,
384 const std::string& var,
385 const std::string& norm,
396 std::cout <<
"==========================================================" << std::endl
397 <<
"Checking variable " << var <<
" using error norm " << norm <<
" with tol " << tol <<
"...";
399 if( norm == std::string(
"L2") )
401 error = exact_sol.l2_error(system_name, var);
403 else if( norm == std::string(
"H1") )
405 error = exact_sol.h1_error(system_name, var);
409 std::cerr <<
"ERROR: Invalid norm " << norm << std::endl;
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;
425 std::cout <<
"PASSED!" << std::endl
426 <<
"==========================================================" << std::endl;
SATurbBCFactoryBase(const std::string &bc_type_name)
TurbBoundFuncBase(libMesh::MeshFunction *turbulent_bc_values)
void init()
Initialize the Simulation objects.
SharedPtr< libMesh::EquationSystems > get_equation_system()
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
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
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
SATurbNuBCFactory(const std::string &bc_type_name)
TurbulentBdyFunctionU(libMesh::MeshFunction *turbulent_bc_values)
DirichletBCFactoryFunctionBase(const std::string &bc_type_name)
const GetPot & get_command_line() const
const std::string & get_multiphysics_system_name() const
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.
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.
Simulation & get_simulation()