GRINS-0.8.0
low_mach_cavity_benchmark_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-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 // C++
27 #include <iostream>
28 
29 // GRINS
30 #include "grins/runner.h"
31 
32 // Function for getting initial temperature field
33 libMesh::Real
34 initial_values( const libMesh::Point& p, const libMesh::Parameters &params,
35  const std::string& system_name, const std::string& unknown_name );
36 
37 int main(int argc, char* argv[])
38 {
39  GRINS::Runner grins(argc,argv);
40 
41  // This is a tough problem to get to converge without PETSc
42  libmesh_example_requires
43  (libMesh::default_solver_package() != libMesh::LASPACK_SOLVERS,
44  "--enable-petsc");
45 
46  grins.init();
47 
48  const GetPot & inputfile = grins.get_input_file();
49 
50  GRINS::Simulation & sim = grins.get_simulation();
51 
52  //FIXME: We need to move this to within the Simulation object somehow...
53  std::string restart_file = inputfile( "restart-options/restart_file", "none" );
54 
55  if( restart_file == "none" )
56  {
57  // Asssign initial temperature value
58  std::string system_name = inputfile( "screen-options/system_name", "GRINS" );
59  GRINS::SharedPtr<libMesh::EquationSystems> es = sim.get_equation_system();
60  const libMesh::System& system = es->get_system(system_name);
61 
62  libMesh::Parameters &params = es->parameters;
63  libMesh::Real T_init = inputfile("Materials/TestMaterial/ReferenceTemperature/value", 0.0);
64  libMesh::Real p0_init = inputfile("Materials/TestMaterial/ThermodynamicPressure/value", 0.0);
65 
66  libMesh::Real& dummy_T = params.set<libMesh::Real>("T_init");
67  dummy_T = T_init;
68 
69  libMesh::Real& dummy_p0 = params.set<libMesh::Real>("p0_init");
70  dummy_p0 = p0_init;
71 
72  system.project_solution( initial_values, NULL, params );
73  }
74 
75  grins.run();
76 
77  libMesh::Real qoi = sim.get_qoi_value(0);
78 
79  // Note that this is a *really* coarse mesh. This is just for testing
80  // and not even close to the real QoI for this problem.
81 
82  // Erroneous value from libMesh 0.9.2.2
83  // const libMesh::Real exact_qoi = 4.8158910676325055;
84 
85  // Value after libMesh 7acb6fc9 bugfix
86  const libMesh::Real exact_qoi = 4.8654229502012685;
87 
88  const libMesh::Real tol = 1.0e-9;
89 
90  int return_flag = 0;
91 
92  libMesh::Real rel_error = std::fabs( (qoi-exact_qoi)/exact_qoi );
93 
94  if( rel_error > tol )
95  {
96  // Skip this test until we know what changed
97  // return_flag = 1;
98  return_flag = 77;
99 
100  std::cerr << std::setprecision(16)
101  << std::scientific
102  << "Error: QoI value mismatch." << std::endl
103  << "Computed qoi = " << qoi << std::endl
104  << "Exact qoi = " << exact_qoi << std::endl
105  << "Relative error = " << rel_error << std::endl;
106  }
107 
108  return return_flag;
109 }
110 
111 libMesh::Real
112 initial_values( const libMesh::Point&, const libMesh::Parameters &params,
113  const std::string& , const std::string& unknown_name )
114 {
115  libMesh::Real value = 0.0;
116 
117  if( unknown_name == "T" )
118  value = params.get<libMesh::Real>("T_init");
119 
120  else if( unknown_name == "p0" )
121  value = params.get<libMesh::Real>("p0_init");
122 
123  else
124  value = 0.0;
125 
126  return value;
127 }
void init()
Initialize the Simulation objects.
Definition: runner.C:53
SharedPtr< libMesh::EquationSystems > get_equation_system()
Definition: simulation.C:393
int main(int argc, char *argv[])
const GetPot & get_input_file() const
Definition: runner.h:60
libMesh::Number get_qoi_value(unsigned int qoi_index) const
Definition: simulation.C:398
libMesh::Real initial_values(const libMesh::Point &p, const libMesh::Parameters &params, const std::string &system_name, const std::string &unknown_name)
Class to encapsulate initializing and running GRINS Simulation.
Definition: runner.h:47
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