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

Classes

class  UnitCubePoissonMGProblem< ELEMENT, MESH >
 

Namespaces

 Global_Parameters
 Namespace for global parameters.
 
 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

HelmholtzSmootherSmoother_Factory_Function_Helper::set_pre_smoother ()
 
HelmholtzSmootherSmoother_Factory_Function_Helper::set_post_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.

601 {
602  //------------------------
603  // Command line arguments
604  //------------------------
605  // Store command line arguments
606  CommandLineArgs::setup(argc,argv);
607 
608  // Number of nodes in one direction of an element. Setting nnode_1d to 2
609  // corresponds to using linear interpolation. Setting nnode_1d=3 corresponds
610  // to quadratic interpolation and setting nnode_1d to 4 corresponds to cubic
611  // interpolation
613  "--nnode_1d",&Global_Parameters::Nnode_1d);
614 
615  // The minimum level of uniform refinement
618 
619  // The additional levels of uniform refinement
622 
623  // The maximum number of adaptive refinements
625  "--n_adapt",&Global_Parameters::N_adaptations);
626 
627  // The additional levels of uniform refinement
630 
631  // The choice of pre-smoother
633  "--presmoother",&Global_Parameters::Pre_smoother_flag);
634 
635  // The choice of post-smoother
637  "--postsmoother",&Global_Parameters::Post_smoother_flag);
638 
639  // The choice of linear solver
641  "--linear_solver",&Global_Parameters::Linear_solver_flag);
642 
643  // The choice to suppress all or some of the MG solver output
646 
647  // Decide whether or not to display convergence information
650 
651  // Parse command line
653 
654  // Document what has been specified on the command line
656 
657  //--------------------------------
658  // Set the documentation directory
659  //--------------------------------
660  // Set output directory
662 
663  //--------------------
664  // Set up the problem
665  //--------------------
666  // Initialise a null pointer to the class Problem
667  Problem* problem_pt=0;
668 
669  // Set the problem pointer depending on the input (defaulted to nnode_1d=2)
671  {
672  // Using linear interpolation
673  typedef RefineableQPoissonElement<3,2> ELEMENT;
674 
675  // Typedef the mesh and template it by the chosen element type
677 
678  // Set the problem pointer
681  }
682  else if (Global_Parameters::Nnode_1d==3)
683  {
684  // Using quadratic interpolation
685  typedef RefineableQPoissonElement<3,3> ELEMENT;
686 
687  // Typedef the mesh and template it by the chosen element type
689 
690  // Set the problem pointer
693  }
694  else if (Global_Parameters::Nnode_1d==4)
695  {
696  // Using cubic interpolation
697  typedef RefineableQPoissonElement<3,4> ELEMENT;
698 
699  // Typedef the mesh and template it by the chosen element type
701 
702  // Set the problem pointer
705  }
706  else
707  {
708  // Throw an error otherwise
709  throw OomphLibError("nnode_1d can only be 2,3 or 4.",
712  }
713 
714  //-------------------
715  // Solve the problem!
716  //-------------------
717  // If the user wishes to use adaptive refinement then we use the Newton solver
718  // with a given argument to indicate how many adaptations can be done
720  {
721  // If the user did not choose to suppress everything
723  {
724  oomph_info << "\n====================Initial Refinement====================\n"
725  << std::endl;
726  }
727 
728  // Refine once at least otherwise the V-cycle only uses SuperLU and converges
729  // instantly
730  problem_pt->refine_uniformly();
731 
732  // Solve the problem
734  }
735  // If the user instead wishes to use uniform refinement
736  else
737  {
738  // Keep refining until the minimum refinement level is reached
739  for (unsigned i=0;i<Global_Parameters::Min_refinement_level;i++)
740  {
741  // If the user did not choose to suppress everything
743  {
744  oomph_info << "\n====================Initial Refinement====================\n"
745  << std::endl;
746  }
747 
748  // Add additional refinement
749  problem_pt->refine_uniformly();
750  }
751 
752  // Solve the problem
753  problem_pt->newton_solve();
754 
755  // Refine and solve until the additional refinements have been completed
756  for (unsigned i=0;i<Global_Parameters::Add_refinement_level;i++)
757  {
758  // If the user did not choose to suppress everything
760  {
761  oomph_info << "===================Additional Refinement==================\n"
762  << std::endl;
763  }
764 
765  // Add additional refinement
766  problem_pt->refine_uniformly();
767 
768  // Solve the problem
769  problem_pt->newton_solve();
770  }
771  } // if (Global_Parameters::Use_adaptation_flag)
772 } // end of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: unit_cube_poisson.cc:156
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
Definition: oomph_definitions.h:222
Definition: problem.h:151
void refine_uniformly(const Vector< unsigned > &nrefine_for_mesh)
Definition: problem.h:2575
void newton_solve()
Use Newton method to solve the problem.
Definition: problem.cc:8783
Definition: refineable_poisson_elements.h:193
Refineable version of simple cubic 3D Brick mesh class.
Definition: simple_cubic_mesh.template.h:169
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
unsigned Output_management_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:88
unsigned Add_refinement_level
The additional levels of uniform refinement.
Definition: two_d_poisson_tanh_flux_bc_validate.cc:56
unsigned Post_smoother_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:74
unsigned Use_adaptation_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:64
DocInfo Doc_info
DocInfo object used for documentation of the solution.
Definition: two_d_poisson_tanh_flux_bc_validate.cc:97
unsigned N_adaptations
The number of adaptations allowed by the Newton solver.
Definition: two_d_poisson_tanh_flux_bc_validate.cc:59
unsigned Nnode_1d
The number of nodes in one direction (default=2)
Definition: two_d_poisson_tanh_flux_bc_validate.cc:50
unsigned Pre_smoother_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:69
unsigned Linear_solver_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:79
unsigned Doc_convergence_flag
Definition: two_d_poisson_tanh_flux_bc_validate.cc:94
unsigned Min_refinement_level
The minimum level of uniform refinement.
Definition: two_d_poisson_tanh_flux_bc_validate.cc:53
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
void specify_command_line_flag(const std::string &command_line_flag, const std::string &doc)
Specify possible argument-free command line flag.
Definition: oomph_utilities.cc:451
void parse_and_assign(int argc, char *argv[], const bool &throw_on_unrecognised_args)
Definition: oomph_utilities.cc:760
void doc_specified_flags()
Document specified command line flags.
Definition: oomph_utilities.cc:610
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Global_Parameters::Add_refinement_level, Global_Parameters::Doc_convergence_flag, Global_Parameters::Doc_info, oomph::CommandLineArgs::doc_specified_flags(), TanhSolnForPoisson::get_source(), i, Global_Parameters::Linear_solver_flag, Global_Parameters::Min_refinement_level, Global_Parameters::N_adaptations, oomph::Problem::newton_solve(), Global_Parameters::Nnode_1d, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, Global_Parameters::Output_management_flag, oomph::CommandLineArgs::parse_and_assign(), Global_Parameters::Post_smoother_flag, Global_Parameters::Pre_smoother_flag, oomph::Problem::refine_uniformly(), oomph::DocInfo::set_directory(), Flag_definition::setup(), oomph::CommandLineArgs::specify_command_line_flag(), and Global_Parameters::Use_adaptation_flag.