mpi/distribution/hp_adaptive_poisson/two_d_poisson_hp_adapt.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/simple_rectangular_quadmesh.h"

Classes

class  SimpleRefineableRectangularQuadMesh< ELEMENT >
 
class  RefineableTwoDPoissonProblem< 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_exact_gradient (const Vector< double > &x, Vector< double > &dudx)
 Exact gradient 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 hp-adaptive solution to the 2D Poisson problem. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char **  argv 
)

Driver code for hp-adaptive solution to the 2D Poisson problem.

338 {
339 
340 #ifdef OOMPH_HAS_MPI
341 
342  // Initialise MPI
343  MPI_Helpers::init(argc,argv);
344 
345 #endif
346 
347  // Store command line arguments
348  CommandLineArgs::setup(argc,argv);
349 
350  // Create label for output
351  //------------------------
352  DocInfo doc_info;
353 
354  // Set output directory
355  doc_info.set_directory("RESLT");
356 
357  //Set up the problem
358  //------------------
359 
360  // Create the problem with 2D hp-refineable elements from the
361  // PRefineableQPoissonElement family. Pass pointer to source function.
364 
365  //Are there command-line arguments?
366  if(CommandLineArgs::Argc==1)
367  {
368 
369 #ifdef OOMPH_HAS_MPI
370 
371  // Provide storage for each element's partition number
372  const unsigned n_element=problem.mesh_pt()->nelement();
373  Vector<unsigned> out_element_partition(n_element);
374 
375  // Distribute the problem
376  bool report_stats=true;
377  out_element_partition=problem.distribute(report_stats);
378 
379  // Write partition to disk
380  std::ofstream output_file;
381  char filename[100];
382  sprintf(filename,"out_hp_adaptive_poisson_partition.dat");
383  output_file.open(filename);
384  for (unsigned e=0;e<n_element;e++)
385  {
386  output_file << out_element_partition[e] << std::endl;
387  }
388 
389  // Check halo schemes (optional)
390  problem.check_halo_schemes(doc_info);
391 
392 #endif
393 
394  }
395  else // Validation run - read in partition from file
396  {
397 
398 #ifdef OOMPH_HAS_MPI
399 
400  DocInfo mesh_doc_info;
401  mesh_doc_info.set_directory("RESLT_MESH");
402  mesh_doc_info.number()=0;
403  std::ifstream input_file;
404  char filename[100];
405 
406  // Get the partition to be used from file
407  const unsigned n_element=problem.mesh_pt()->nelement();
408  Vector<unsigned> element_partition(n_element);
409  sprintf(filename,"hp_adaptive_poisson_partition.dat");
410  input_file.open(filename);
411  std::string input_string;
412  for (unsigned e=0;e<n_element;e++)
413  {
414  getline(input_file,input_string,'\n');
415  element_partition[e]=atoi(input_string.c_str());
416  }
417 
418  // Now perform the distribution
419  bool report_stats=true;
420  problem.distribute(element_partition,mesh_doc_info,report_stats);
421 
422 #endif
423 
424  } // end of mesh distribution
425 
426 
427  // Check if we're ready to go:
428  //----------------------------
429  cout << "\n\n\nProblem self-test ";
430  if (problem.self_test()==0)
431  {
432  cout << "passed: Problem can be solved." << std::endl;
433  }
434  else
435  {
436  throw OomphLibError("Self test failed",
439  }
440 
441  // Set the orientation of the "step" to 45 degrees
443 
444  // Choose a large value for the steepness of the "step"
446 
447  // Initially uniformly refine the mesh
448  for (unsigned p=0; p<0; p++)
449  {
450  cout << "p-refining:" << endl;
451  problem.p_refine_uniformly();
452  }
453  for (unsigned h=0; h<0; h++)
454  {
455  cout << "h-refining:" << endl;
456  problem.refine_uniformly();
457  }
458 
459  problem.newton_solve();
460 
461  //Output the solution
462  doc_info.number()=1;
463  problem.doc_solution(doc_info);
464 
465  problem.p_adapt();
466  problem.newton_solve();
467  doc_info.number()=2;
468  problem.doc_solution(doc_info);
469 
470  problem.adapt();
471  problem.newton_solve();
472  doc_info.number()=3;
473  problem.doc_solution(doc_info);
474 
475  problem.p_adapt();
476  problem.newton_solve();
477  doc_info.number()=4;
478  problem.doc_solution(doc_info);
479 
480  problem.adapt();
481  problem.newton_solve();
482  doc_info.number()=5;
483  problem.doc_solution(doc_info);
484 
485  problem.p_adapt();
486  problem.newton_solve();
487  doc_info.number()=6;
488  problem.doc_solution(doc_info);
489 
490  problem.adapt();
491  problem.newton_solve();
492  doc_info.number()=7;
493  problem.doc_solution(doc_info);
494 
495  problem.p_adapt();
496  problem.newton_solve();
497  doc_info.number()=8;
498  problem.doc_solution(doc_info);
499 
500  problem.adapt();
501  problem.newton_solve();
502  doc_info.number()=9;
503  problem.doc_solution(doc_info);
504 
505 
506  // Count hanging nodes
507  unsigned hang_no=0;
508  for (unsigned n=0; n<problem.mesh_pt()->nnode(); n++)
509  {
510  if (problem.mesh_pt()->node_pt(n)->is_hanging()) hang_no++;
511  }
512  cout << hang_no << " hanging nodes in mesh." << endl;
513 
514  //Output the solution
515  doc_info.number()=0;
516  problem.doc_solution(doc_info);
517 
518 
519 // Finalise MPI
520 #ifdef OOMPH_HAS_MPI
521 
522  MPI_Helpers::finalize();
523 
524 #endif
525 
526 } //end of main
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
float * p
Definition: Tutorial_Map_using.cpp:9
Definition: mpi/distribution/hp_adaptive_poisson/two_d_poisson_hp_adapt.cc:127
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
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
string filename
Definition: MergeRestartFiles.py:39
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
int Argc
Number of arguments + 1.
Definition: oomph_utilities.cc:407
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#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, oomph::CommandLineArgs::Argc, e(), MergeRestartFiles::filename, TanhSolnForPoisson::get_source(), n, oomph::DocInfo::number(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, p, problem, oomph::DocInfo::set_directory(), Flag_definition::setup(), oomph::Global_string_for_annotation::string(), and TanhSolnForPoisson::TanPhi.