unit_cube_poisson.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/simple_cubic_mesh.h"

Classes

class  UnitCubePoissonMGProblem< ELEMENT, MESH >
 

Namespaces

 Smoother_Factory_Function_Helper
 Returns a pointer to a smoother of the appropriate type.
 
 TanhSolnForPoisson
 Namespace for exact solution for Poisson equation with "sharp step".
 

Functions

SmootherSmoother_Factory_Function_Helper::set_smoother ()
 
void TanhSolnForPoisson::get_exact_u (const Vector< double > &x, Vector< double > &u)
 Exact solution as a Vector. More...
 
void TanhSolnForPoisson::get_exact_u (const Vector< double > &x, double &u)
 Exact solution as a scalar. More...
 
void TanhSolnForPoisson::get_source (const Vector< double > &x, double &source)
 Source function to make it an exact solution. More...
 
int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int argc  ,
char argv[] 
)

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// Driver for 3D Poisson problem in a unit cube. Solution has a sharp step.

457 {
458  //------------------------------------
459  // Sort out documentation information:
460  //------------------------------------
461  // Create a DocInfo object to document the solution
462  DocInfo doc_info;
463 
464  // Set output directory
465  doc_info.set_directory("RESLT");
466 
467  // Step number
468  doc_info.number()=0;
469 
470  //--------------------
471  // Set up the problem
472  //--------------------
473  // Using linear interpolation
474  typedef RefineableQPoissonElement<3,3> ELEMENT;
475 
476  // Typedef the mesh and template it by the chosen element type
478 
479  // Set the problem
481 
482  //-------------------
483  // Solve the problem:
484  //-------------------
485  // Set the number of times to initially refine the problem. Must be at
486  // least 1 otherwise the multigrid solver uses SuperLU to solve the problem
487  unsigned n_refine=2;
488 
489  // Keep refining until the minimum refinement level is reached
490  for (unsigned i=0;i<n_refine;i++)
491  {
492  // Indicate the problem is about to be refined
493  oomph_info << "\n===================="
494  << "Initial Refinement"
495  << "====================\n"
496  << std::endl;
497 
498  // Refine the problem
499  problem.refine_uniformly();
500  }
501 
502  // Solve the problem
503  problem.newton_solve();
504 
505  // Document the solution
506  problem.doc_solution(doc_info);
507 } // End of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: unit_cube_poisson.cc:156
Definition: oomph_utilities.h:499
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
Definition: refineable_poisson_elements.h:193
Refineable version of simple cubic 3D Brick mesh class.
Definition: simple_cubic_mesh.template.h:169
void get_source(const Vector< double > &x, double &source)
Source function to make it an exact solution.
Definition: extrude_with_macro_element_representation.cc:224
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References TanhSolnForPoisson::get_source(), i, oomph::DocInfo::number(), oomph::oomph_info, problem, and oomph::DocInfo::set_directory().