optimisation/use_coarse_base_meshes/two_d_poisson_adapt.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/simple_rectangular_quadmesh.h"

Classes

class  SimpleRefineableRectangularQuadMesh< ELEMENT >
 
class  RefineablePoissonProblem< ELEMENT >
 

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_source (const Vector< double > &x, double &source)
 Source function to make it an exact solution. More...
 
int main (int argc, char *argv[])
 Driver code for 2D Poisson problem. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char argv[] 
)

Driver code for 2D Poisson problem.

433 {
434 
435 
436  // Number of uniform refinements relative to a 2x2 base mesh
437  unsigned n_refine;
438 
439  // Get number of refinement levels from command line or use default
440  if (argc==1)
441  {
442  n_refine=5;
443  }
444  else if (argc==2)
445  {
446  n_refine=atoi(argv[1]);
447  }
448  else
449  {
450  std::string error_message =
451  "Wrong number of input arguments. The options are: \n";
452  error_message +=
453  "No args: Default number of refinements\n";
454  error_message +=
455  "One arg: Required number of refinements\n";
456 
457  throw OomphLibError(error_message,
460  }
461 
462 
463  // Loop over two versions of refinement: Build a fine base mesh
464  // or refine coarse base mesh uniformly. The former version is
465  // much more expensive because the black-box quadtree forest generator
466  // creates the quadtree forest in the base mesh.
467  for (unsigned i=0;i<2;i++)
468  {
469 
470  // Exponent for number of elements in base mesh
471  int n_power;
472  if (i==0)
473  {
474  cout << std::endl;
475  cout << "/////////////////////////////////////////////////////////"
476  << std::endl;
477  cout << "Building problem with a coarse base mesh" << std::endl;
478  cout << "followed by uniform refinements." << std::endl;
479  cout << "/////////////////////////////////////////////////////////"
480  << std::endl;
481  cout << std::endl;
482  n_power=-n_refine;
483  }
484  else
485  {
486  cout << std::endl;
487  cout << "/////////////////////////////////////////////////////////"
488  << std::endl;
489  cout << "Building problem with a fine base mesh" << std::endl;
490  cout << "/////////////////////////////////////////////////////////"
491  << std::endl;
492  cout << std::endl;
493  n_power=n_refine+1;
494  }
495 
496  //Set up the problem
497  //------------------
498 
499  // Initialise timers
500  clock_t t_start = clock();
501 
502  // Create the problem with 2D nine-node refineable elements from the
503  // RefineableQuadPoissonElement family. Pass pointer to source function.
506 
507  // Finish/doc timing
508  clock_t t_end = clock();
509  double total_time=double(t_end-t_start)/CLOCKS_PER_SEC;
510  cout << std::endl;
511  cout << "======================================================= " << std::endl;
512  cout << "Total time for Problem setup [sec]: " << total_time << std::endl;
513  cout << "======================================================= " << std::endl;
514  cout << std::endl;
515 
516 
517  // Create label for output
518  //------------------------
519  DocInfo doc_info;
520 
521  // Set output directory
522  doc_info.set_directory("RESLT");
523 
524  // Step number
525  doc_info.number()=0;
526 
527  // Set the orientation of the "step" to 45 degrees
529 
530  // Choose a large value for the steepness of the "step"
532 
533 
534  // Solve as linear and nonlinear problem
535  //--------------------------------------
536  problem.set_problem_is_nonlinear();
537 
538  unsigned nstep=2;
539  for (unsigned istep=0;istep<nstep;istep++)
540  {
541 
542  // Initialise timers
543  clock_t t_start = clock();
544 
545  if (problem.is_problem_nonlinear())
546  {
547  cout << std::endl << std::endl;
548  cout << "============================ " << std::endl;
549  cout << "Solving as nonlinear problem " << std::endl;
550  cout << "============================ " << std::endl;
551  cout << std::endl << std::endl;
552  }
553  else
554  {
555  cout << std::endl << std::endl;
556  cout << "============================ " << std::endl;
557  cout << "Solving as linear problem " << std::endl;
558  cout << "============================ " << std::endl;
559  cout << std::endl << std::endl;
560  }
561 
562  // Solve the problem
563  problem.newton_solve();
564 
565  // Finish/doc timing
566  clock_t t_end = clock();
567  double total_time=double(t_end-t_start)/CLOCKS_PER_SEC;
568  cout << "======================================================= "
569  << std::endl;
570  cout << "Total time for Newton solve [sec]: " << total_time << std::endl;
571  cout << "======================================================= "
572  << std::endl;
573 
574  //Output the solution
575  problem.doc_solution(doc_info);
576 
577  //Increment counter for solutions
578  doc_info.number()++;
579 
580  // Next time around solve as linear problem
581  problem.set_problem_is_linear();
582 
583  }
584 
585  // Give user a chance to check the memory usage
586  cout << std::endl << std::endl;
587  cout << "Execution is paused while problem is in core" << std::endl;
588  pause("Have a look at the memory usage now");
589 
590  }
591 
592 } //end of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: optimisation/use_coarse_base_meshes/two_d_poisson_adapt.cc:149
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: oomph_definitions.h:222
double TanPhi
Parameter for angle Phi of "step".
Definition: HypreSolver_test.cc:51
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
double Alpha
Parameter for steepness of step.
Definition: extrude_with_macro_element_representation.cc:185
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
void pause(std::string message)
Pause and display message.
Definition: oomph_utilities.cc:1265
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References TanhSolnForPoisson::Alpha, TanhSolnForPoisson::get_source(), i, oomph::DocInfo::number(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::pause(), problem, oomph::DocInfo::set_directory(), oomph::Global_string_for_annotation::string(), and TanhSolnForPoisson::TanPhi.