GRINS-0.8.0
generic_exact_solution_testing_app.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 #include "grins_config.h"
26 
27 // GRINS
28 #include "grins/runner.h"
29 #include "grins/grins_enums.h"
30 #include "grins/mesh_builder.h"
31 #include "grins/multiphysics_sys.h"
32 
33 //libMesh
34 #include "libmesh/exact_solution.h"
35 #include "libmesh/parsed_function.h"
36 #include "libmesh/composite_function.h"
37 
38 int test_error_norm( libMesh::ExactSolution& exact_sol,
39  const std::string& system_name,
40  const std::string& var,
41  const std::string& norm,
42  const libMesh::Real exact_error,
43  const double tol );
44 
45 int main(int argc, char* argv[])
46 {
47  GRINS::Runner grins(argc,argv);
48 
49  const GetPot & command_line = grins.get_command_line();
50 
51  if( !command_line.have_variable("test_data") )
52  {
53  std::cerr << "ERROR: Must specify solution test data on command line with test_data=<file>."
54  << std::endl;
55  exit(1);
56  }
57 
58  if( !command_line.have_variable("vars") )
59  {
60  std::cerr << "ERROR: Must specify variables on command line with vars='var1 var2'"
61  << std::endl;
62  exit(1);
63  }
64 
65  if( !command_line.have_variable("norms") )
66  {
67  std::cerr << "ERROR: Must specify error norms on command line with norms='L2 H1'"
68  << std::endl;
69  exit(1);
70  }
71 
72  if( !command_line.have_variable("tol") )
73  {
74  std::cerr << "ERROR: Must specify test tolerance on command line with tol=<tol>"
75  << std::endl;
76  exit(1);
77  }
78 
79  // Grab the variables for which we want to compare the errors
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++ )
83  {
84  vars[v] = command_line("vars", "DIE!", v);
85  }
86 
87  // Now grab the norms to compute for each variable error
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++ )
91  {
92  norms[n] = command_line("norms", "DIE!", n);
93  if( norms[n] != std::string("L2") )
94  {
95  std::cerr << "ERROR: Invalid norm input " << norms[n] << std::endl
96  << " Valid values are: L2" << std::endl;
97  }
98  }
99 
100  // Now that we have the norms, we need to grab the error values for each var, for each norm from the CLI
101  std::vector<std::vector<libMesh::Real> > error_values;
102  error_values.resize(n_vars);
103 
104  for( unsigned int v = 0; v < n_vars; v++ )
105  {
106  std::string var = vars[v];
107 
108  error_values[v].resize(n_norms);
109 
110  for( unsigned int n = 0; n < n_norms; n++ )
111  {
112  std::string norm_cli = var+"_"+norms[n]+"_error";
113 
114  if( !command_line.have_variable(norm_cli) )
115  {
116  std::cerr << "ERROR: Must specify "+norms[n]+" errors on command line with "+
117  var+"_"+norms[n]+"_error='<error>'" << std::endl;
118  exit(1);
119  }
120 
121  error_values[v][n] = command_line( norm_cli, -1.0 );
122 
123  if( error_values[v][n] < 0.0 )
124  {
125  std::cerr << "ERROR: Invalid value of error --- norms should be positive!" << std::endl
126  << " Found error value: " << error_values[v][n] << std::endl;
127  exit(1);
128  }
129  }
130  }
131 
132  const GetPot & inputfile = grins.get_input_file();
133  const libMesh::LibMeshInit & libmesh_init = grins.get_libmesh_init();
134 
135  GRINS::MeshBuilder mesh_builder;
136  GRINS::SharedPtr<libMesh::UnstructuredMesh> mesh = mesh_builder.build(inputfile,libmesh_init.comm());
137 
138  libMesh::EquationSystems es(*mesh);
139 
140  // This needs to match the read counter-part in GRINS::Simulation
141  std::string test_data = command_line("test_data", "DIE!" );
142 
143  {
144  std::ifstream i(test_data.c_str());
145  if (!i)
146  {
147  std::cerr << "Error: Could not read from test_data file "
148  << test_data << std::endl;
149  exit(1);
150  }
151  }
152 
153  es.read(test_data,
154  libMesh::EquationSystems::READ_HEADER |
155  libMesh::EquationSystems::READ_DATA |
156  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
157 
158  std::string system_name = "GRINS-TEST";
159  libMesh::System& system = es.get_system(system_name);
160 
161  // Create Exact solution object and attach exact solution quantities
162  libMesh::ExactSolution exact_sol(es);
163 
164  // Grab the exact solution expression for each variable
165  // We'll directly construct the ParsedFunction for each of the variable exact solutions provided
166  // They are all packed together in the CompositeFunction
167  // We need to do this here because we need the variable index from the equation system
168  {
170 
171  for( unsigned int v = 0; v < n_vars; v++ )
172  {
173  std::string var = vars[v];
174  std::string soln_cli = var+"_exact_soln";
175 
176  if( !command_line.have_variable(soln_cli) )
177  {
178  std::cerr << "ERROR: Must specify the exact solution for the variable"+var+"on" << std::endl
179  << " the command line with "+soln_cli+"=<expression>"
180  << std::endl;
181  exit(1);
182  }
183 
184  unsigned int var_index = system.variable_number(var);
185 
186  std::string parsed_soln = command_line(soln_cli,"NULL");
187 
188  std::vector<unsigned int> index_map(1,var_index);
189 
190  exact_sols.attach_subfunction( libMesh::ParsedFunction<libMesh::Real>(parsed_soln), index_map );
191  }
192 
193  // This is assuming there's only 1 system
194  exact_sol.attach_exact_value(0, &exact_sols);
195  }
196 
197  // Now compute error for each variable
198  for( unsigned int v = 0; v < n_vars; v++ )
199  {
200  exact_sol.compute_error(system_name, vars[v]);
201  }
202 
203  double tol = command_line("tol", 1.0e-10);
204 
205  int return_flag = 0;
206  int test_flag = 0;
207 
208  // Now test error for each variable, for each norm
209  for( unsigned int v = 0; v < n_vars; v++ )
210  {
211  for( unsigned int n = 0; n < n_norms; n++ )
212  {
213  test_flag = test_error_norm( exact_sol, system_name, vars[v], norms[n], error_values[v][n], tol );
214 
215  if( test_flag == 1 )
216  return_flag = 1;
217  }
218  }
219 
220  return return_flag;
221 }
222 
223 int test_error_norm( libMesh::ExactSolution& exact_sol,
224  const std::string& system_name,
225  const std::string& var,
226  const std::string& norm,
227  const libMesh::Real exact_error,
228  const double tol )
229 {
230  int return_flag = 0;
231 
232  double error = 0.0;
233 
234  std::cout << "==========================================================" << std::endl
235  << "Checking variable " << var << " using error norm " << norm << " with tol " << tol << "...";
236 
237  if( norm == std::string("L2") )
238  {
239  error = exact_sol.l2_error(system_name, var);
240  }
241  else if( norm == std::string("H1") )
242  {
243  error = exact_sol.h1_error(system_name, var);
244  }
245  else
246  {
247  std::cout << "ERROR: Invalid norm " << norm << std::endl;
248  exit(1);
249  }
250 
251  if( std::abs(error-exact_error) > tol )
252  {
253  return_flag = 1;
254 
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;
261  }
262  else
263  {
264  std::cout << "PASSED!" << std::endl
265  << "==========================================================" << std::endl;
266  }
267 
268  return return_flag;
269 }
const GetPot & get_input_file() const
Definition: runner.h:60
const libMesh::LibMeshInit & get_libmesh_init() const
Definition: runner.h:74
const GetPot & get_command_line() const
Definition: runner.h:63
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.
Definition: mesh_builder.C:46
Class to encapsulate initializing and running GRINS Simulation.
Definition: runner.h:47
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)

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