GRINS-0.6.0
test_ns_couette_flow_2d_y.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2015 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 #include "grins_config.h"
27 
28 #include <iostream>
29 
30 // GRINS
31 #include "grins/mesh_builder.h"
32 #include "grins/simulation.h"
34 #include "grins/multiphysics_sys.h"
35 
36 //libMesh
37 #include "libmesh/exact_solution.h"
38 
39 // GRVY
40 #ifdef GRINS_HAVE_GRVY
41 #include "grvy.h"
42 #endif
43 
44 libMesh::Number
45 exact_solution( const libMesh::Point& p,
46  const libMesh::Parameters& params, // parameters, not needed
47  const std::string& sys, // sys_name, not needed
48  const std::string& var ); // unk_name, not needed);
49 
50 libMesh::Gradient
51 exact_derivative( const libMesh::Point& p,
52  const libMesh::Parameters& params, // parameters, not needed
53  const std::string& sys, // sys_name, not needed
54  const std::string& var ); // unk_name, not needed);
55 
56 int main(int argc, char* argv[])
57 {
58 
59 #ifdef GRINS_USE_GRVY_TIMERS
60  GRVY::GRVY_Timer_Class grvy_timer;
61  grvy_timer.Init("GRINS Timer");
62 #endif
63 
64  // Check command line count.
65  if( argc < 2 )
66  {
67  // TODO: Need more consistent error handling.
68  std::cerr << "Error: Must specify libMesh input file." << std::endl;
69  exit(1); // TODO: something more sophisticated for parallel runs?
70  }
71 
72  // libMesh input file should be first argument
73  std::string libMesh_input_filename = argv[1];
74 
75  // Create our GetPot object.
76  GetPot libMesh_inputfile( libMesh_input_filename );
77 
78 #ifdef GRINS_USE_GRVY_TIMERS
79  grvy_timer.BeginTimer("Initialize Solver");
80 #endif
81 
82  // Initialize libMesh library.
83  libMesh::LibMeshInit libmesh_init(argc, argv);
84 
85  GRINS::SimulationBuilder sim_builder;
86 
87  GRINS::Simulation grins( libMesh_inputfile,
88  sim_builder,
89  libmesh_init.comm() );
90 
91 #ifdef GRINS_USE_GRVY_TIMERS
92  grvy_timer.EndTimer("Initialize Solver");
93 
94  // Attach GRVY timer to solver
95  grins.attach_grvy_timer( &grvy_timer );
96 #endif
97 
98  grins.run();
99 
100 #ifdef GRINS_USE_GRVY_TIMERS
101  grvy_timer.Finalize();
102 #endif
103 
104  std::tr1::shared_ptr<libMesh::EquationSystems> es = grins.get_equation_system();
105 
106  // Create Exact solution object and attach exact solution quantities
107  libMesh::ExactSolution exact_sol(*es);
108 
109  exact_sol.attach_exact_value(&exact_solution);
110  exact_sol.attach_exact_deriv(&exact_derivative);
111 
112  // Compute error and get it in various norms
113  exact_sol.compute_error("GRINS", "v");
114 
115  double l2error = exact_sol.l2_error("GRINS", "v");
116  double h1error = exact_sol.h1_error("GRINS", "v");
117 
118  int return_flag = 0;
119 
120  if( l2error > 4.0e-9 || h1error > 4.0e-9 )
121  {
122  return_flag = 1;
123 
124  std::cout << "Tolerance exceeded for Couette flow test." << std::endl
125  << "l2 error = " << l2error << std::endl
126  << "h1 error = " << h1error << std::endl;
127  }
128 
129  return return_flag;
130 }
131 
132 libMesh::Number
133 exact_solution( const libMesh::Point& p,
134  const libMesh::Parameters& /*params*/, // parameters, not needed
135  const std::string& sys, // sys_name, not needed
136  const std::string& var ) // unk_name, not needed);
137 {
138  const double x = p(0);
139 
140  if( sys != "GRINS" || var != "v" )
141  std::cout << "sys = " << sys << ", var = " << var << std::endl;
142 
143  // Hardcoded to velocity in input file.
144  libMesh::Number f = 10.0*x;
145 
146  return f;
147 }
148 
149 libMesh::Gradient
150 exact_derivative( const libMesh::Point& /*p*/,
151  const libMesh::Parameters& /*params*/, // parameters, not needed
152  const std::string& sys, // sys_name, not needed
153  const std::string& var ) // unk_name, not needed);
154 {
155  libMesh::Gradient g;
156 
157  if( sys != "GRINS" || var != "v" )
158  std::cout << "sys = " << sys << ", var = " << var << std::endl;
159 
160  // Hardcoded to velocity in input file.
161  g(0) = 10.0;
162  g(1) = 0.0;
163 
164 #if LIBMESH_DIM > 2
165  g(2) = 0.0;
166 #endif
167 
168  return g;
169 }
libMesh::Number exact_solution(const libMesh::Point &p, const libMesh::Parameters &params, const std::string &sys, const std::string &var)
int main(int argc, char *argv[])
libMesh::Gradient exact_derivative(const libMesh::Point &p, const libMesh::Parameters &params, const std::string &sys, const std::string &var)

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1