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

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