25 #include "grins_config.h"
34 #include "libmesh/exact_solution.h"
35 #include "libmesh/parsed_function.h"
36 #include "libmesh/composite_function.h"
39 const std::string& system_name,
40 const std::string& var,
41 const std::string& norm,
42 const libMesh::Real exact_error,
45 int main(
int argc,
char* argv[])
51 if( !command_line.have_variable(
"test_data") )
53 std::cerr <<
"ERROR: Must specify solution test data on command line with test_data=<file>."
58 if( !command_line.have_variable(
"vars") )
60 std::cerr <<
"ERROR: Must specify variables on command line with vars='var1 var2'"
65 if( !command_line.have_variable(
"norms") )
67 std::cerr <<
"ERROR: Must specify error norms on command line with norms='L2 H1'"
72 if( !command_line.have_variable(
"tol") )
74 std::cerr <<
"ERROR: Must specify test tolerance on command line with tol=<tol>"
80 unsigned int n_vars = command_line.vector_variable_size(
"vars");
81 std::vector<std::string> vars(n_vars);
82 for(
unsigned int v = 0; v < n_vars; v++ )
84 vars[v] = command_line(
"vars",
"DIE!", v);
88 unsigned int n_norms = command_line.vector_variable_size(
"norms");
89 std::vector<std::string> norms(n_norms);
90 for(
unsigned int n = 0; n < n_norms; n++ )
92 norms[n] = command_line(
"norms",
"DIE!", n);
93 if( norms[n] != std::string(
"L2") )
95 std::cerr <<
"ERROR: Invalid norm input " << norms[n] << std::endl
96 <<
" Valid values are: L2" << std::endl;
101 std::vector<std::vector<libMesh::Real> > error_values;
102 error_values.resize(n_vars);
104 for(
unsigned int v = 0; v < n_vars; v++ )
106 std::string var = vars[v];
108 error_values[v].resize(n_norms);
110 for(
unsigned int n = 0; n < n_norms; n++ )
112 std::string norm_cli = var+
"_"+norms[n]+
"_error";
114 if( !command_line.have_variable(norm_cli) )
116 std::cerr <<
"ERROR: Must specify "+norms[n]+
" errors on command line with "+
117 var+
"_"+norms[n]+
"_error='<error>'" << std::endl;
121 error_values[v][n] = command_line( norm_cli, -1.0 );
123 if( error_values[v][n] < 0.0 )
125 std::cerr <<
"ERROR: Invalid value of error --- norms should be positive!" << std::endl
126 <<
" Found error value: " << error_values[v][n] << std::endl;
136 GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = mesh_builder.
build(inputfile,libmesh_init.comm());
138 libMesh::EquationSystems es(*mesh);
141 std::string test_data = command_line(
"test_data",
"DIE!" );
144 std::ifstream i(test_data.c_str());
147 std::cerr <<
"Error: Could not read from test_data file "
148 << test_data << std::endl;
154 libMesh::EquationSystems::READ_HEADER |
155 libMesh::EquationSystems::READ_DATA |
156 libMesh::EquationSystems::READ_ADDITIONAL_DATA);
158 std::string system_name =
"GRINS-TEST";
159 libMesh::System& system = es.get_system(system_name);
162 libMesh::ExactSolution exact_sol(es);
171 for(
unsigned int v = 0; v < n_vars; v++ )
173 std::string var = vars[v];
174 std::string soln_cli = var+
"_exact_soln";
176 if( !command_line.have_variable(soln_cli) )
178 std::cerr <<
"ERROR: Must specify the exact solution for the variable"+var+
"on" << std::endl
179 <<
" the command line with "+soln_cli+
"=<expression>"
184 unsigned int var_index = system.variable_number(var);
186 std::string parsed_soln = command_line(soln_cli,
"NULL");
188 std::vector<unsigned int> index_map(1,var_index);
194 exact_sol.attach_exact_value(0, &exact_sols);
198 for(
unsigned int v = 0; v < n_vars; v++ )
200 exact_sol.compute_error(system_name, vars[v]);
203 double tol = command_line(
"tol", 1.0e-10);
209 for(
unsigned int v = 0; v < n_vars; v++ )
211 for(
unsigned int n = 0; n < n_norms; n++ )
213 test_flag =
test_error_norm( exact_sol, system_name, vars[v], norms[n], error_values[v][n], tol );
224 const std::string& system_name,
225 const std::string& var,
226 const std::string& norm,
227 const libMesh::Real exact_error,
234 std::cout <<
"==========================================================" << std::endl
235 <<
"Checking variable " << var <<
" using error norm " << norm <<
" with tol " << tol <<
"...";
237 if( norm == std::string(
"L2") )
239 error = exact_sol.l2_error(system_name, var);
241 else if( norm == std::string(
"H1") )
243 error = exact_sol.h1_error(system_name, var);
247 std::cout <<
"ERROR: Invalid norm " << norm << std::endl;
251 if( std::abs(error-exact_error) > tol )
255 std::cout <<
"Tolerance exceeded for generic regression test!" << std::endl
256 <<
"tolerance = " << tol << std::endl
257 <<
"norm of error = " << error << std::endl
258 <<
"exact error = " << exact_error << std::endl
259 <<
"norm type = " << norm << std::endl
260 <<
"var = " << var << std::endl;
264 std::cout <<
"PASSED!" << std::endl
265 <<
"==========================================================" << std::endl;
const GetPot & get_input_file() const
const libMesh::LibMeshInit & get_libmesh_init() const
const GetPot & get_command_line() const
int main(int argc, char *argv[])
SharedPtr< libMesh::UnstructuredMesh > build(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
Builds the libMesh::Mesh according to input options.
Class to encapsulate initializing and running GRINS Simulation.
int test_error_norm(libMesh::ExactSolution &exact_sol, const std::string &system_name, const std::string &var, const std::string &norm, const libMesh::Real exact_error, const double tol)