GRINS-0.6.0
test_stokes_poiseuille_flow_parsed_viscosity_parsed_conductivity.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"
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&, // parameters, not needed
47  const std::string&, // sys_name, not needed
48  const std::string&); // unk_name, not needed);
49 
50 libMesh::Gradient
51 exact_derivative( const libMesh::Point& p,
52  const libMesh::Parameters&, // parameters, not needed
53  const std::string&, // sys_name, not needed
54  const std::string&); // unk_name, not needed);
55 
57 {
58 public:
59 
62  { return; };
63 
64  ~ParabolicBCFactory(){return;};
65 
66  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > build_dirichlet( );
67 };
68 
69 int main(int argc, char* argv[])
70 {
71 
72 #ifdef GRINS_USE_GRVY_TIMERS
73  GRVY::GRVY_Timer_Class grvy_timer;
74  grvy_timer.Init("GRINS Timer");
75 #endif
76 
77  // Check command line count.
78  if( argc < 2 )
79  {
80  // TODO: Need more consistent error handling.
81  std::cerr << "Error: Must specify libMesh input file." << std::endl;
82  exit(1); // TODO: something more sophisticated for parallel runs?
83  }
84 
85  // libMesh input file should be first argument
86  std::string libMesh_input_filename = argv[1];
87 
88  // Create our GetPot object.
89  GetPot libMesh_inputfile( libMesh_input_filename );
90 
91 #ifdef GRINS_USE_GRVY_TIMERS
92  grvy_timer.BeginTimer("Initialize Solver");
93 #endif
94 
95  // Initialize libMesh library.
96  libMesh::LibMeshInit libmesh_init(argc, argv);
97 
98  GRINS::SimulationBuilder sim_builder;
99 
100  std::tr1::shared_ptr<ParabolicBCFactory> bc_factory( new ParabolicBCFactory );
101 
102  sim_builder.attach_bc_factory(bc_factory);
103 
104  GRINS::Simulation grins( libMesh_inputfile,
105  sim_builder,
106  libmesh_init.comm() );
107 
108 #ifdef GRINS_USE_GRVY_TIMERS
109  grvy_timer.EndTimer("Initialize Solver");
110 
111  // Attach GRVY timer to solver
112  grins.attach_grvy_timer( &grvy_timer );
113 #endif
114 
115  // Solve
116  grins.run();
117 
118  // Get equation systems to create ExactSolution object
119  std::tr1::shared_ptr<libMesh::EquationSystems> es = grins.get_equation_system();
120 
121  // Create Exact solution object and attach exact solution quantities
122  libMesh::ExactSolution exact_sol(*es);
123 
124  exact_sol.attach_exact_value(&exact_solution);
125  exact_sol.attach_exact_deriv(&exact_derivative);
126 
127  // Compute error and get it in various norms
128  exact_sol.compute_error("GRINS", "u");
129 
130  double l2error = exact_sol.l2_error("GRINS", "u");
131  double h1error = exact_sol.h1_error("GRINS", "u");
132 
133  int return_flag = 0;
134 
135  const double tol = 1.0e-8;
136 
137  if( l2error > tol || h1error > tol )
138  {
139  return_flag = 1;
140 
141  std::cout << "Tolerance exceeded for velocity in Poiseuille test." << std::endl
142  << "l2 error = " << l2error << std::endl
143  << "h1 error = " << h1error << std::endl;
144  }
145 
146  // Compute error and get it in various norms
147  exact_sol.compute_error("GRINS", "p");
148 
149  l2error = exact_sol.l2_error("GRINS", "p");
150  h1error = exact_sol.h1_error("GRINS", "p");
151 
152  if( l2error > tol || h1error > tol )
153  {
154  return_flag = 1;
155 
156  std::cout << "Tolerance exceeded for pressure in Poiseuille test." << std::endl
157  << "l2 error = " << l2error << std::endl
158  << "h1 error = " << h1error << std::endl;
159  }
160 
161  return return_flag;
162 }
163 
164 std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > ParabolicBCFactory::build_dirichlet( )
165 {
166  GRINS::DBCContainer cont;
167  cont.add_var_name( "u" );
168  cont.add_bc_id( 1 );
169  cont.add_bc_id( 3 );
170 
171  std::tr1::shared_ptr<libMesh::FunctionBase<libMesh::Number> > u_func( new GRINS::ParabolicProfile );
172 
173  cont.set_func( u_func );
174 
175  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > mymap;
176 
177  mymap.insert( std::pair<GRINS::PhysicsName, GRINS::DBCContainer >(GRINS::stokes, cont) );
178 
179  return mymap;
180 }
181 
182 libMesh::Number
183 exact_solution( const libMesh::Point& p,
184  const libMesh::Parameters& /*params*/, // parameters, not needed
185  const std::string& /*sys_name*/, // sys_name, not needed
186  const std::string& var ) // unk_name, not needed);
187 {
188  const double x = p(0);
189  const double y = p(1);
190 
191  libMesh::Number f = 0.0;
192  // Hardcoded to velocity in input file.
193  if( var == "u" ) f = 4*y*(1-y);
194  if( var == "p" ) f = 120.0 + (80.0-120.0)/5.0*x;
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 }
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
libMesh::Number exact_solution(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
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.
const PhysicsName stokes
libMesh::Gradient exact_derivative(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
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