GRINS-0.7.0
generic_solution_regression.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 #include "grins_config.h"
26 
27 // GRINS
28 #include "grins/mesh_builder.h"
29 #include "grins/simulation.h"
31 #include "grins/multiphysics_sys.h"
32 
33 //libMesh
34 #include "libmesh/exact_solution.h"
35 
36 void test_error_norm( libMesh::ExactSolution& exact_sol,
37  const std::string& system_name,
38  const std::string& var,
39  const std::string& norm,
40  const double tol,
41  int& return_flag );
42 
43 int main(int argc, char* argv[])
44 {
45  GetPot command_line(argc,argv);
46 
47  if( !command_line.have_variable("input") )
48  {
49  std::cerr << "ERROR: Must specify input file on command line with input=<file>." << std::endl;
50  exit(1);
51  }
52 
53  if( !command_line.have_variable("soln-data") )
54  {
55  std::cerr << "ERROR: Must specify solution data on command line with soln-data=<file>." << std::endl;
56  exit(1);
57  }
58 
59  if( !command_line.have_variable("vars") )
60  {
61  std::cerr << "ERROR: Must specify variables on command line with vars='var1 var2'" << 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'" << std::endl;
68  exit(1);
69  }
70 
71  if( !command_line.have_variable("tol") )
72  {
73  std::cerr << "ERROR: Must specify test tolerance on command line with tol=<tol>" << std::endl;
74  exit(1);
75  }
76 
77  // libMesh input file should be first argument
78  std::string libMesh_input_filename = command_line("input", "DIE!");
79 
80  {
81  std::ifstream i(libMesh_input_filename.c_str());
82  if (!i)
83  {
84  std::cerr << "Error: Could not read from libMesh input file "
85  << libMesh_input_filename << std::endl;
86  exit(1);
87  }
88  }
89 
90  // Create our GetPot object.
91  GetPot libMesh_inputfile( libMesh_input_filename );
92 
93  // Initialize libMesh library.
94  libMesh::LibMeshInit libmesh_init(argc, argv);
95 
96  GRINS::SimulationBuilder sim_builder;
97 
98  GRINS::Simulation grins( libMesh_inputfile,
99  command_line,
100  sim_builder,
101  libmesh_init.comm() );
102 
103  // Do solve here
104  grins.run();
105 
106  // Get equation systems to create ExactSolution object
107  GRINS::SharedPtr<libMesh::EquationSystems> es = grins.get_equation_system();
108 
109  // Create Exact solution object and attach exact solution quantities
110  libMesh::ExactSolution exact_sol(*es);
111 
112  libMesh::EquationSystems es_ref( es->get_mesh() );
113 
114  // Filename of file where comparison solution is stashed
115  std::string solution_file = command_line("soln-data", "DIE!");
116  es_ref.read( solution_file );
117 
118  exact_sol.attach_reference_solution( &es_ref );
119 
120  // Now grab the variables for which we want to compare
121  unsigned int n_vars = command_line.vector_variable_size("vars");
122  std::vector<std::string> vars(n_vars);
123  for( unsigned int v = 0; v < n_vars; v++ )
124  {
125  vars[v] = command_line("vars", "DIE!", v);
126  }
127 
128  // Now grab the norms to compute for each variable error
129  unsigned int n_norms = command_line.vector_variable_size("norms");
130  std::vector<std::string> norms(n_norms);
131  for( unsigned int n = 0; n < n_norms; n++ )
132  {
133  norms[n] = command_line("norms", "DIE!", n);
134  if( norms[n] != std::string("L2") &&
135  norms[n] != std::string("H1") )
136  {
137  std::cerr << "ERROR: Invalid norm input " << norms[n] << std::endl
138  << " Valid values are: L2" << std::endl
139  << " H1" << std::endl;
140  }
141  }
142 
143  const std::string& system_name = grins.get_multiphysics_system_name();
144 
145  // Now compute error for each variable
146  for( unsigned int v = 0; v < n_vars; v++ )
147  {
148  exact_sol.compute_error(system_name, vars[v]);
149  }
150 
151  int return_flag = 0;
152 
153  double tol = command_line("tol", 1.0e-10);
154 
155  // Now test error for each variable, for each norm
156  for( unsigned int v = 0; v < n_vars; v++ )
157  {
158  for( unsigned int n = 0; n < n_norms; n++ )
159  {
160  test_error_norm( exact_sol, system_name, vars[v], norms[n], tol, return_flag );
161  }
162  }
163 
164  // Now grab any "gold" QoI values
165  unsigned int n_qoi_vals = command_line.vector_variable_size("qois");
166  for( unsigned int n = 0; n < n_qoi_vals; n++ )
167  {
168 
169  std::cout << "==========================================================" << std::endl
170  << "Checking qoi " << n << " with tol " << tol << "...";
171 
172  libMesh::Number gold_qoi = command_line("qois", libMesh::Number(0), n);
173 
174  libMesh::Number computed_qoi =
175  grins.get_multiphysics_system()->qoi[n];
176 
177  double error = computed_qoi - gold_qoi;
178 
179  if (std::abs(error) > tol)
180  {
181  std::cerr << "Tolerance exceeded for generic regression test!" << std::endl
182  << "tolerance = " << tol << std::endl
183  << "error = " << error << std::endl
184  << "qoi index = " << n << std::endl;
185  return_flag = 1;
186  }
187  else
188  std::cout << "PASSED!" << std::endl
189  << "==========================================================" << std::endl;
190 
191  }
192 
193 
194  return return_flag;
195 }
196 
197 void test_error_norm( libMesh::ExactSolution& exact_sol,
198  const std::string& system_name,
199  const std::string& var,
200  const std::string& norm,
201  const double tol,
202  int& return_flag )
203 {
204  // We don't set return_flag unless we are setting it 1
205  // since this function gets called multiple times and we don't
206  // want to overwrite a previous "fail" (return_flag = 1) with
207  // a "pass" (return_flag = 0)
208 
209  double error = 0.0;
210 
211  std::cout << "==========================================================" << std::endl
212  << "Checking variable " << var << " using error norm " << norm << " with tol " << tol << "...";
213 
214  if( norm == std::string("L2") )
215  {
216  error = exact_sol.l2_error(system_name, var);
217  }
218  else if( norm == std::string("H1") )
219  {
220  error = exact_sol.h1_error(system_name, var);
221  }
222  else
223  {
224  std::cerr << "ERROR: Invalid norm " << norm << std::endl;
225  exit(1);
226  }
227 
228  if( error > tol )
229  {
230  return_flag = 1;
231 
232  std::cerr << "Tolerance exceeded for generic regression test!" << std::endl
233  << "tolerance = " << tol << std::endl
234  << "norm of error = " << error << std::endl
235  << "norm type = " << norm << std::endl
236  << "var = " << var << std::endl;
237  }
238  else
239  {
240  std::cout << "PASSED!" << std::endl
241  << "==========================================================" << std::endl;
242  }
243 
244  return;
245 }
int main(int argc, char *argv[])
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)

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10