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

Go to the source code of this file.

Classes

class  AxiConCylBCFactory
 

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 &  sys,
const std::string &  var 
)

Definition at line 222 of file test_axi_ns_con_cyl_flow.C.

Referenced by main().

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

Definition at line 202 of file test_axi_ns_con_cyl_flow.C.

Referenced by main().

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

Definition at line 71 of file test_axi_ns_con_cyl_flow.C.

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

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 }
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 &)
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