GRINS-0.6.0
test_axi_ns_con_cyl_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 #include "libmesh/zero_function.h"
40 
41 // GRVY
42 #ifdef GRINS_HAVE_GRVY
43 #include "grvy.h"
44 #endif
45 
46 libMesh::Number
47 exact_solution( const libMesh::Point& p,
48  const libMesh::Parameters&, // parameters, not needed
49  const std::string&, // sys_name, not needed
50  const std::string&); // unk_name, not needed);
51 
52 libMesh::Gradient
53 exact_derivative( const libMesh::Point& p,
54  const libMesh::Parameters&, // parameters, not needed
55  const std::string&, // sys_name, not needed
56  const std::string&); // unk_name, not needed);
57 
59 {
60 public:
61 
64  { return; };
65 
66  ~AxiConCylBCFactory(){return;};
67 
68  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > build_dirichlet( );
69 };
70 
71 int main(int argc, char* argv[])
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  // This is a tough problem to get to converge without PETSc
100  libmesh_example_requires
101  (libMesh::default_solver_package() != libMesh::LASPACK_SOLVERS,
102  "--enable-petsc");
103 
104  GRINS::SimulationBuilder sim_builder;
105 
106  std::tr1::shared_ptr<AxiConCylBCFactory> bc_factory( new AxiConCylBCFactory );
107 
108  sim_builder.attach_bc_factory(bc_factory);
109 
110  GRINS::Simulation grins( libMesh_inputfile,
111  sim_builder,
112  libmesh_init.comm() );
113 
114 #ifdef GRINS_USE_GRVY_TIMERS
115  grvy_timer.EndTimer("Initialize Solver");
116 
117  // Attach GRVY timer to solver
118  grins.attach_grvy_timer( &grvy_timer );
119 #endif
120 
121  // Do solve here
122  grins.run();
123 
124  // Get equation systems to create ExactSolution object
125  std::tr1::shared_ptr<libMesh::EquationSystems> es = grins.get_equation_system();
126 
127  // Create Exact solution object and attach exact solution quantities
128  libMesh::ExactSolution exact_sol(*es);
129 
130  exact_sol.attach_exact_value(&exact_solution);
131  exact_sol.attach_exact_deriv(&exact_derivative);
132 
133  // Compute error and get it in various norms
134  exact_sol.compute_error("GRINS", "z_vel");
135 
136  double l2error = exact_sol.l2_error("GRINS", "z_vel");
137  double h1error = exact_sol.h1_error("GRINS", "z_vel");
138 
139  int return_flag = 0;
140 
141  if( l2error > 2.0e-9 || h1error > 4.0e-7 )
142  {
143  return_flag = 1;
144 
145  std::cout << "Tolerance exceeded for velocity in Poiseuille test." << std::endl
146  << "l2 error = " << l2error << std::endl
147  << "h1 error = " << h1error << std::endl;
148  }
149 
150  /*
151  // Compute error and get it in various norms
152  exact_sol.compute_error("GRINS", "p");
153 
154  l2error = exact_sol.l2_error("GRINS", "p");
155  h1error = exact_sol.h1_error("GRINS", "p");
156 
157  if( l2error > 1.0e-11 || h1error > 1.0e-11 )
158  {
159  return_flag = 1;
160 
161  std::cout << "Tolerance exceeded for pressure in Poiseuille test." << std::endl
162  << "l2 error = " << l2error << std::endl
163  << "h1 error = " << h1error << std::endl;
164  }
165  */
166 
167  return return_flag;
168 }
169 
170 std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > AxiConCylBCFactory::build_dirichlet( )
171 {
172  GRINS::DBCContainer cont;
173  cont.add_var_name( "z_vel" );
174  cont.add_bc_id( 0 );
175  cont.add_bc_id( 2 );
176 
177  std::tr1::shared_ptr<libMesh::FunctionBase<libMesh::Number> > vel_func( new GRINS::ConcentricCylinderProfile );
178 
179  cont.set_func( vel_func );
180 
181 
182  GRINS::DBCContainer cont2;
183  cont2.add_var_name( "z_vel" );
184  cont2.add_bc_id( 0 );
185  cont2.add_bc_id( 2 );
186 
187  std::tr1::shared_ptr<libMesh::FunctionBase<libMesh::Number> >
188  vel_func2( new libMesh::ZeroFunction<libMesh::Number> );
189 
190  cont2.set_func( vel_func2 );
191 
192  std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > mymap;
193 
194  mymap.insert( std::pair<GRINS::PhysicsName, GRINS::DBCContainer >(GRINS::incompressible_navier_stokes, cont) );
195 
196  mymap.insert( std::pair<GRINS::PhysicsName, GRINS::DBCContainer >(GRINS::incompressible_navier_stokes, cont2) );
197 
198  return mymap;
199 }
200 
201 libMesh::Number
202 exact_solution( const libMesh::Point& p,
203  const libMesh::Parameters& /*params*/, // parameters, not needed
204  const std::string& /*sys_name*/, // sys_name, not needed
205  const std::string& var ) // unk_name, not needed);
206 {
207  const double r = p(0);
208 
209  const double r0 = 1.0;
210  const double r1 = 2.0;
211  const double u0 = 2.0;
212 
213  libMesh::Number f = 0;
214  // Hardcoded to velocity in input file.
215  if( var == "z_vel" ) f = u0*std::log( r1/r )/std::log( r1/r0 );
216  else libmesh_assert(false);
217 
218  return f;
219 }
220 
221 libMesh::Gradient
222 exact_derivative( const libMesh::Point& p,
223  const libMesh::Parameters& /*params*/, // parameters, not needed
224  const std::string& /*sys_name*/, // sys_name, not needed
225  const std::string& var) // unk_name, not needed);
226 {
227  const double r = p(0);
228 
229  const double r0 = 1.0;
230  const double r1 = 2.0;
231  const double u0 = 2.0;
232 
233  libMesh::Gradient g;
234 
235  // Hardcoded to velocity in input file.
236  if( var == "z_vel" )
237  {
238  g(0) = -u0/std::log(r1/r0)*r/r1*r1/(r*r);
239  g(1) = 0.0;
240  }
241 
242  return g;
243 }
const PhysicsName incompressible_navier_stokes
Object for constructing boundary condition function objects.
Definition: bc_factory.h:50
libMesh::Number exact_solution(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
Simple helper class to setup general Dirichlet boundary conditions.
Definition: dbc_container.h:49
libMesh::Gradient exact_derivative(const libMesh::Point &p, const libMesh::Parameters &, const std::string &, const std::string &)
std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > build_dirichlet()
Builds all required libMesh::DirichletBoundary objects and adds them to DofMap.
void add_var_name(const GRINS::VariableName &var)
Add variables that are constrained by the Dirichlet bc.
Definition: dbc_container.C:46
GRINS namespace.
Profile for flow between axially moving concentric cylinders.
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
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