eighth_sphere_poisson_hp_adapt.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/eighth_sphere_mesh.h"

Classes

class  EighthSpherePoissonProblem< ELEMENT >
 Poisson problem in refineable eighth of a sphere mesh. More...
 

Namespaces

 TanhSolnForPoisson
 Namespace for exact solution for Poisson equation with "sharp step".
 

Functions

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 eighth of a sphere. Solution has a sharp step. If there are any command line arguments, we regard this as a validation run and perform only a single adaptation.

351 {
352 
353  // Store command line arguments
354  CommandLineArgs::setup(argc,argv);
355 
356  // Set up the problem with 27-node brick elements, pass pointer to
357  // source function
360 
361  // Setup labels for output
362  DocInfo doc_info;
363 
364  // Output directory
365  doc_info.set_directory("RESLT");
366 
367  // Step number
368  doc_info.number()=0;
369 
370  // Check if we're ready to go
371  std::cout << "Self test: " << problem.self_test() << std::endl;
372 
373  // Solve the problem
374  problem.newton_solve();
375 
376  //Output solution
377  problem.doc_solution(doc_info);
378 
379  //Increment counter for solutions
380  doc_info.number()++;
381 
382  // Perform initial mesh refinement and solve
383  problem.refine_uniformly();
384  problem.refine_uniformly();
385  problem.p_refine_uniformly();
386 
387  //std::cout << "Self test: " << problem.self_test() << std::endl;
388  problem.newton_solve();
389  problem.doc_solution(doc_info);
390  doc_info.number()++;
391 
392  // Adapt and solve
393  problem.adapt();
394  //std::cout << "Self test: " << problem.self_test() << std::endl;
395  problem.newton_solve();
396  problem.doc_solution(doc_info);
397  doc_info.number()++;
398 
399  problem.p_adapt();
400  //std::cout << "Self test: " << problem.self_test() << std::endl;
401  problem.newton_solve();
402  problem.doc_solution(doc_info);
403  doc_info.number()++;
404 
405 } // end of main
Poisson problem in refineable eighth of a sphere mesh.
Definition: eighth_sphere_poisson.cc:132
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
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
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
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References TanhSolnForPoisson::get_source(), oomph::DocInfo::number(), problem, oomph::DocInfo::set_directory(), and Flag_definition::setup().