GRINS-0.6.0
test_ns_poiseuille_flow.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"
36 
37 //libMesh
38 #include "libmesh/exact_solution.h"
39 
40 // GRVY
41 #ifdef GRINS_HAVE_GRVY
42 #include "grvy.h"
43 #endif
44 
45 libMesh::Number
46 exact_solution( const libMesh::Point& p,
47  const libMesh::Parameters&, // parameters, not needed
48  const std::string&, // sys_name, not needed
49  const std::string&); // unk_name, not needed);
50 
51 libMesh::Gradient
52 exact_derivative( const libMesh::Point& p,
53  const libMesh::Parameters&, // parameters, not needed
54  const std::string&, // sys_name, not needed
55  const std::string&); // unk_name, not needed);
56 
58 {
59 public:
60 
63  { return; };
64 
65  ~ParabolicBCFactory(){return;};
66 
67  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > build_dirichlet( );
68 };
69 
70 int main(int argc, char* argv[])
71 {
72 
73 #ifdef GRINS_USE_GRVY_TIMERS
74  GRVY::GRVY_Timer_Class grvy_timer;
75  grvy_timer.Init("GRINS Timer");
76 #endif
77 
78  // Check command line count.
79  if( argc < 2 )
80  {
81  // TODO: Need more consistent error handling.
82  std::cerr << "Error: Must specify libMesh input file." << std::endl;
83  exit(1); // TODO: something more sophisticated for parallel runs?
84  }
85 
86  // libMesh input file should be first argument
87  std::string libMesh_input_filename = argv[1];
88 
89  // Create our GetPot object.
90  GetPot libMesh_inputfile( libMesh_input_filename );
91 
92 #ifdef GRINS_USE_GRVY_TIMERS
93  grvy_timer.BeginTimer("Initialize Solver");
94 #endif
95 
96  // Initialize libMesh library.
97  libMesh::LibMeshInit libmesh_init(argc, argv);
98 
99  GRINS::SimulationBuilder sim_builder;
100 
101  std::tr1::shared_ptr<ParabolicBCFactory> bc_factory( new ParabolicBCFactory );
102 
103  sim_builder.attach_bc_factory(bc_factory);
104 
105  GRINS::Simulation grins( libMesh_inputfile,
106  sim_builder,
107  libmesh_init.comm() );
108 
109 #ifdef GRINS_USE_GRVY_TIMERS
110  grvy_timer.EndTimer("Initialize Solver");
111 
112  // Attach GRVY timer to solver
113  grins.attach_grvy_timer( &grvy_timer );
114 #endif
115 
116  // Solve
117  grins.run();
118 
119  // Get equation systems to create ExactSolution object
120  std::tr1::shared_ptr<libMesh::EquationSystems> es = grins.get_equation_system();
121 
122  // Create Exact solution object and attach exact solution quantities
123  libMesh::ExactSolution exact_sol(*es);
124 
125  exact_sol.attach_exact_value(&exact_solution);
126  exact_sol.attach_exact_deriv(&exact_derivative);
127 
128  // Compute error and get it in various norms
129  exact_sol.compute_error("GRINS", "u");
130 
131  double l2error = exact_sol.l2_error("GRINS", "u");
132  double h1error = exact_sol.h1_error("GRINS", "u");
133 
134  int return_flag = 0;
135 
136  if( l2error > 1.0e-9 || h1error > 1.0e-9 )
137  {
138  return_flag = 1;
139 
140  std::cout << "Tolerance exceeded for velocity in Poiseuille test." << std::endl
141  << "l2 error = " << l2error << std::endl
142  << "h1 error = " << h1error << std::endl;
143  }
144 
145  // Compute error and get it in various norms
146  exact_sol.compute_error("GRINS", "p");
147 
148  l2error = exact_sol.l2_error("GRINS", "p");
149  h1error = exact_sol.h1_error("GRINS", "p");
150 
151  if( l2error > 2.0e-9 || h1error > 2.0e-9 )
152  {
153  return_flag = 1;
154 
155  std::cout << "Tolerance exceeded for pressure in Poiseuille test." << std::endl
156  << "l2 error = " << l2error << std::endl
157  << "h1 error = " << h1error << std::endl;
158  }
159 
160  return return_flag;
161 }
162 
163 std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > ParabolicBCFactory::build_dirichlet( )
164 {
165  GRINS::DBCContainer cont;
166  cont.add_var_name( "u" );
167  cont.add_bc_id( 1 );
168  cont.add_bc_id( 3 );
169 
170  std::tr1::shared_ptr<libMesh::FunctionBase<libMesh::Number> > u_func( new GRINS::ParabolicProfile );
171 
172  cont.set_func( u_func );
173 
174  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > mymap;
175 
176  mymap.insert( std::pair<GRINS::PhysicsName, GRINS::DBCContainer >(GRINS::incompressible_navier_stokes, cont) );
177 
178  return mymap;
179 }
180 
181 libMesh::Number
182 exact_solution( const libMesh::Point& p,
183  const libMesh::Parameters& /*params*/, // parameters, not needed
184  const std::string& /*sys_name*/, // sys_name, not needed
185  const std::string& var ) // unk_name, not needed);
186 {
187  const double x = p(0);
188  const double y = p(1);
189 
190  libMesh::Number f = 0;
191  // Hardcoded to velocity in input file.
192  if( var == "u" ) f = 4*y*(1-y);
193  else if( var == "p" ) f = 120.0 + (80.0-120.0)/5.0*x;
194  else libmesh_assert(false);
195 
196  return f;
197 }
198 
199 libMesh::Gradient
200 exact_derivative( const libMesh::Point& p,
201  const libMesh::Parameters& /*params*/, // parameters, not needed
202  const std::string& /*sys_name*/, // sys_name, not needed
203  const std::string& var) // unk_name, not needed);
204 {
205  const double y = p(1);
206 
207  libMesh::Gradient g;
208 
209  // Hardcoded to velocity in input file.
210  if( var == "u" )
211  {
212  g(0) = 0.0;
213  g(1) = 4*(1-y) - 4*y;
214 
215 #if LIBMESH_DIM > 2
216  g(2) = 0.0;
217 #endif
218  }
219 
220  if( var == "p" )
221  {
222  g(0) = (80.0-120.0)/5.0;
223  g(1) = 0.0;
224 
225 #if LIBMESH_DIM > 2
226  g(2) = 0.0;
227 #endif
228  }
229  return g;
230 }
const PhysicsName incompressible_navier_stokes
Object for constructing boundary condition function objects.
Definition: bc_factory.h:50
Simple helper class to setup general Dirichlet boundary conditions.
Definition: dbc_container.h:49
Parabolic profile.
void add_var_name(const GRINS::VariableName &var)
Add variables that are constrained by the Dirichlet bc.
Definition: dbc_container.C:46
GRINS namespace.
libMesh::Gradient exact_derivative(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
libMesh::Number exact_solution(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
int main(int argc, char *argv[])
void add_bc_id(const GRINS::BoundaryID bc_id)
Add boundary id's for which this Dirichlet bc is to be applied.
Definition: dbc_container.C:52
void set_func(std::tr1::shared_ptr< libMesh::FunctionBase< libMesh::Number > > func)
Add the Dirichlet bc functor.
Definition: dbc_container.C:58
std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > build_dirichlet()
Builds all required libMesh::DirichletBoundary objects and adds them to DofMap.
void attach_bc_factory(std::tr1::shared_ptr< BoundaryConditionsFactory > bc_factory)

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