two_d_poisson_tanh_flux_bc.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/rectangular_quadmesh.h"

Classes

class  FluxPoissonMGProblem< ELEMENT, MESH >
 

Namespaces

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

Functions

void TanhSolnForPoisson::get_exact_u (const Vector< double > &x, Vector< double > &u)
 Exact solution as a Vector. More...
 
void TanhSolnForPoisson::source_function (const Vector< double > &x, double &source)
 Source function required to make the solution above an exact solution. More...
 
void TanhSolnForPoisson::prescribed_flux_on_fixed_x_boundary (const Vector< double > &x, double &flux)
 Flux required by the exact solution on a boundary on which x is fixed. More...
 
SmootherSmoother_Factory_Function_Helper::set_smoother ()
 
int main (int argc, char **argv)
 

Function Documentation

◆ main()

int main ( int argc  ,
char **  argv 
)

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// Demonstrate how to solve 2D Poisson problem with flux boundary conditions.

563 {
564  //------------------------------------
565  // Sort out documentation information:
566  //------------------------------------
567  // Create a DocInfo object to document the solution
568  DocInfo doc_info;
569 
570  // Set output directory
571  doc_info.set_directory("RESLT");
572 
573  // Step number
574  doc_info.number()=0;
575 
576  //--------------------
577  // Set up the problem:
578  //--------------------
579  // Use 2D quad elements with quadratic interpolation
580  typedef RefineableQPoissonElement<2,3> ELEMENT;
581 
582  // Use a rectangular quad mesh templated by the chosen element type
584 
585  // Set up the problem with the chosen element and mesh type
588 
589  //-------------------
590  // Solve the problem:
591  //-------------------
592  // Set the number of times to initially refine the problem. Must be at
593  // least 1 otherwise the multigrid solver uses SuperLU to solve the problem
594  unsigned n_refine=1;
595 
596  // Keep refining until the minimum refinement level is reached
597  for (unsigned i=0;i<n_refine;i++)
598  {
599  // Indicate the problem is about to be refined
600  oomph_info << "\n===================="
601  << "Initial Refinement"
602  << "====================\n"
603  << std::endl;
604 
605  // Refine the problem
606  problem.refine_uniformly();
607  }
608 
609  // Set the orientation of the "step" to 45 degrees
611 
612  // Initial value for the steepness of the "step"
614 
615  // Do a couple of solutions for different forcing functions
616  //---------------------------------------------------------
617  // Choose the number of choices of Alpha to use
618  unsigned nstep=4;
619 
620  // Loop over the different choices of Alpha
621  for (unsigned istep=0;istep<nstep;istep++)
622  {
623  // Increase the steepness of the step:
625 
626  // Tell the user what problem we're solving
627  std::cout << "\n\nSolving for TanhSolnForPoisson::Alpha="
628  << TanhSolnForPoisson::Alpha << "\n" << std::endl;
629 
630  // Solve the problem
631  problem.newton_solve();
632 
633  // Document the solution
634  problem.doc_solution(doc_info);
635  }
636 } // End of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: two_d_poisson_tanh_flux_bc.cc:115
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
Definition: rectangular_quadmesh.template.h:326
double TanPhi
Parameter for angle Phi of "step".
Definition: HypreSolver_test.cc:51
void source_function(const Vector< double > &x, double &source)
Source function required to make the solution above an exact solution.
Definition: HypreSolver_test.cc:60
double Alpha
Parameter for steepness of step.
Definition: extrude_with_macro_element_representation.cc:185
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References TanhSolnForPoisson::Alpha, i, oomph::DocInfo::number(), oomph::oomph_info, problem, oomph::DocInfo::set_directory(), TanhSolnForPoisson::source_function(), and TanhSolnForPoisson::TanPhi.