GRINS-0.6.0
Classes | Functions
test_stokes_poiseuille_flow.C File Reference
#include "grins_config.h"
#include <iostream>
#include "grins/simulation.h"
#include "grins/simulation_builder.h"
#include "grins/multiphysics_sys.h"
#include "grins/parabolic_profile.h"
#include "libmesh/exact_solution.h"
Include dependency graph for test_stokes_poiseuille_flow.C:

Go to the source code of this file.

Classes

class  ParabolicBCFactory
 

Functions

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[])
 

Function Documentation

libMesh::Gradient exact_derivative ( const libMesh::Point &  p,
const libMesh::Parameters &  ,
const std::string &  ,
const std::string &  var 
)

Definition at line 200 of file test_stokes_poiseuille_flow.C.

Referenced by main().

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 }
libMesh::Number exact_solution ( const libMesh::Point &  p,
const libMesh::Parameters &  ,
const std::string &  ,
const std::string &  var 
)

Definition at line 183 of file test_stokes_poiseuille_flow.C.

Referenced by main().

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 }
int main ( int  argc,
char *  argv[] 
)

Definition at line 69 of file test_stokes_poiseuille_flow.C.

References GRINS::SimulationBuilder::attach_bc_factory(), exact_derivative(), exact_solution(), and GRINS::Simulation::run().

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  const double errortol = 1.0e-8;
134 
135  int return_flag = 0;
136 
137  if( l2error > errortol || h1error > errortol )
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 > errortol || h1error > errortol )
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 }
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 &)
void attach_bc_factory(std::tr1::shared_ptr< BoundaryConditionsFactory > bc_factory)

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