structured_cubic_point_source.cc File Reference
#include <fenv.h>
#include "math.h"
#include <complex>
#include "generic.h"
#include "pml_helmholtz.h"
#include "meshes/simple_cubic_mesh.h"

Classes

class  GlobalParameters::TestPMLMapping
 
class  PMLStructuredCubicHelmholtz< ELEMENT >
 Problem class. More...
 

Namespaces

 GlobalParameters
 Global parameters.
 
 Smoother_Factory_Function_Helper
 Returns a pointer to a smoother of the appropriate type.
 

Functions

void GlobalParameters::update_parameters ()
 Update parameters. More...
 
void GlobalParameters::get_simple_exact_u (const Vector< double > &x, Vector< double > &u)
 Get the exact solution, u, at the spatial position, x. More...
 
bool GlobalParameters::is_in_pinned_region (const Vector< double > &x)
 
HelmholtzSmootherSmoother_Factory_Function_Helper::set_pre_smoother ()
 
HelmholtzSmootherSmoother_Factory_Function_Helper::set_post_smoother ()
 
int main (int argc, char **argv)
 Solve 3D Helmholtz problem for a point source in a unit cube. More...
 

Variables

unsigned GlobalParameters::Nnode_1d =2
 The number of nodes in one direction (default=2) More...
 
unsigned GlobalParameters::Min_refinement_level =1
 The minimum level of uniform refinement. More...
 
unsigned GlobalParameters::Add_refinement_level =0
 The additional levels of uniform refinement. More...
 
unsigned GlobalParameters::N_adaptations =1
 The number of adaptations allowed by the Newton solver. More...
 
unsigned GlobalParameters::Use_adaptation_flag =0
 
unsigned GlobalParameters::Pre_smoother_flag =0
 
unsigned GlobalParameters::Post_smoother_flag =0
 
unsigned GlobalParameters::Linear_solver_flag =1
 
unsigned GlobalParameters::Output_management_flag =0
 
unsigned GlobalParameters::Doc_convergence_flag =0
 
std::ostream * GlobalParameters::Stream_pt
 
double GlobalParameters::Lx =1.0
 
double GlobalParameters::Ly =1.0
 
double GlobalParameters::Lz =1.0
 
unsigned GlobalParameters::Nx =7
 Number of elements in each direction (used by SimpleCubicMesh) More...
 
unsigned GlobalParameters::Ny =7
 
unsigned GlobalParameters::Nz =7
 
double GlobalParameters::Element_width =Lx/double(Nx)
 The element width. More...
 
double GlobalParameters::Pml_thickness =Element_width
 Length of cube in each direction. More...
 
double GlobalParameters::Pi =MathematicalConstants::Pi
 Store the value of Pi. More...
 
double GlobalParameters::Alpha_shift =0.0
 
double GlobalParameters::Wavenumber =sqrt(K_squared)
 Wavenumber (also known as k),k=omega/c. More...
 
double GlobalParameters::Centre =Lx/2.0
 The x and y coordinate of the centre of the cube. More...
 
FiniteElement::SteadyExactSolutionFctPt GlobalParameters::simple_exact_u_pt =&get_simple_exact_u
 
TestPMLMapping * GlobalParameters::Test_pml_mapping_pt =new TestPMLMapping
 Set the new PML mapping. More...
 
unsigned GlobalParameters::Enable_test_pml_mapping_flag =0
 
double GlobalParameters::Eps =1.0e-12
 The tolerance for a point relative to the bounding inner square. More...
 
double Smoother_Factory_Function_Helper::Omega =0.4
 The value of the damping factor for the damped Jacobi smoother. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char **  argv 
)

Solve 3D Helmholtz problem for a point source in a unit cube.

904 {
905  //------------------------
906  // Command line arguments
907  //------------------------
908  // Store command line arguments
909  CommandLineArgs::setup(argc,argv);
910 
911  // Choose the number of nodes in one direction of an element;
912  // Values of nnode_1d:
913  // 2: Bilinear interpolation
914  // 3: Biquadratic interpolation
915  // 4: Bicubic interpolation
917  "--nnode_1d",&GlobalParameters::Nnode_1d);
918 
919  // Choose the minimum level of uniform refinement
922 
923  // Choose the additional levels of uniform refinement
926 
927  // Choose the maximum number of adaptive refinements
929  "--n_adapt",&GlobalParameters::N_adaptations);
930 
931  // Choose how many additional levels of uniform refinement to use
934 
935  // Choose the value of k^2
937  "--k_sq",&GlobalParameters::K_squared);
938 
939  // Choose the value of the shift in the CSLP
941  "--alpha",&GlobalParameters::Alpha_shift);
942 
943  // Choose the value of the damping factor in the damped Jacobi solver
946 
947  // Choose the pre-smoother
949  "--presmoother",&GlobalParameters::Pre_smoother_flag);
950 
951  // Choose the post-smoother
953  "--postsmoother",&GlobalParameters::Post_smoother_flag);
954 
955  // Choose the linear solver
957  "--linear_solver",&GlobalParameters::Linear_solver_flag);
958 
959  // Decide whether or not to suppress all or some of the MG solver output
962 
963  // Decide whether or not to display convergence information
966 
967  // Decide whether or not to display convergence information
969  "--test_pml_mapping",&GlobalParameters::Enable_test_pml_mapping_flag);
970 
971  // Parse command line
973 
974  // Document what has been specified on the command line
976 
977  // Update any parameters that need to be updated
979 
980  //--------------------------------
981  // Set the documentation directory
982  //--------------------------------
983  // Set output directory
985 
986  //-------------------
987  // Set up the problem
988  //-------------------
989  // Initialise a null pointer to the class Problem
990  Problem* problem_pt=0;
991 
992  // Set the problem pointer depending on the input (defaulted to nnode_1d=2)
994  {
995  // Set up the problem with refineable 3D eight-noded elements from the
996  // QPMLHelmholtzElement family
997  typedef RefineableQPMLHelmholtzElement<3,2> ELEMENT;
998 
999  // Set the problem pointer
1000  problem_pt=new PMLStructuredCubicHelmholtz<ELEMENT>;
1001  }
1002  else if (GlobalParameters::Nnode_1d==3)
1003  {
1004  // Set up the problem with refineable 3D twenty-seven-noded elements from
1005  // the QPMLHelmholtzElement family
1006  typedef RefineableQPMLHelmholtzElement<3,3> ELEMENT;
1007 
1008  // Set the problem pointer
1009  problem_pt=new PMLStructuredCubicHelmholtz<ELEMENT>;
1010  }
1011  else if (GlobalParameters::Nnode_1d==4)
1012  {
1013  // Set up the problem with refineable 3D sixty-four-noded elements from
1014  // the QPMLHelmholtzElement family
1015  typedef RefineableQPMLHelmholtzElement<3,4> ELEMENT;
1016 
1017  // Set the problem pointer
1018  problem_pt=new PMLStructuredCubicHelmholtz<ELEMENT>;
1019  }
1020  else
1021  {
1022  // Throw an error otherwise
1023  throw OomphLibError("nnode_1d can only be 2,3 or 4.",
1026  }
1027 
1028  //------------------
1029  // Solve the problem
1030  //------------------
1031  // If the user wishes to use adaptive refinement then we use the Newton
1032  // solver with a given argument to indicate how many adaptations to use
1034  {
1035  // If the user wishes to silence everything
1037  {
1038  // Store the output stream pointer
1040 
1041  // Now set the oomph_info stream pointer to the null stream to
1042  // disable all possible output
1044  }
1045 
1046  // Keep refining until the minimum refinement level is reached
1047  for (unsigned i=0;i<GlobalParameters::Min_refinement_level;i++)
1048  {
1049  oomph_info << "\n===================="
1050  << "Initial Refinement"
1051  << "====================\n"
1052  << std::endl;
1053 
1054  // Add additional refinement
1055  problem_pt->refine_uniformly();
1056  }
1057 
1058  // If we silenced the adaptation, allow output again
1060  {
1061  // Now set the oomph_info stream pointer to the null stream to
1062  // disable all possible output
1064  }
1065 
1066  // Solve the problem
1067  problem_pt->newton_solve();
1068 
1069  // Keep refining until the minimum refinement level is reached
1070  for (unsigned i=0;i<GlobalParameters::N_adaptations;i++)
1071  {
1072  // If the user wishes to silence everything
1074  {
1075  // Store the output stream pointer
1077 
1078  // Now set the oomph_info stream pointer to the null stream to
1079  // disable all possible output
1081  }
1082 
1083  // Adapt the problem
1084  problem_pt->adapt();
1085 
1086  // If we silenced the adaptation, allow output again
1088  {
1089  // Now set the oomph_info stream pointer to the null stream to
1090  // disable all possible output
1092  }
1093 
1094  // Solve the problem
1095  problem_pt->newton_solve();
1096  }
1097  }
1098  // If the user instead wishes to use uniform refinement
1099  else
1100  {
1101  // If the user wishes to silence everything
1103  {
1104  // Store the output stream pointer
1106 
1107  // Now set the oomph_info stream pointer to the null stream to
1108  // disable all possible output
1110  }
1111 
1112  // Keep refining until the minimum refinement level is reached
1113  for (unsigned i=0;i<GlobalParameters::Min_refinement_level;i++)
1114  {
1115  oomph_info << "\n===================="
1116  << "Initial Refinement"
1117  << "====================\n"
1118  << std::endl;
1119 
1120  // Add additional refinement
1121  problem_pt->refine_uniformly();
1122  }
1123 
1124  // If we silenced the adaptation, allow output again
1126  {
1127  // Now set the oomph_info stream pointer to the null stream to
1128  // disable all possible output
1130  }
1131 
1132  // Solve the problem
1133  problem_pt->newton_solve();
1134 
1135  // Refine and solve until the additional refinements have been completed
1136  for (unsigned i=0;i<GlobalParameters::Add_refinement_level;i++)
1137  {
1138  // If the user wishes to silence everything
1140  {
1141  // Store the output stream pointer
1143 
1144  // Now set the oomph_info stream pointer to the null stream to
1145  // disable all possible output
1147  }
1148 
1149  oomph_info << "==================="
1150  << "Additional Refinement"
1151  << "==================\n"
1152  << std::endl;
1153 
1154  // Add additional refinement
1155  problem_pt->refine_uniformly();
1156 
1157  // If we silenced the adaptation, allow output again
1159  {
1160  // Now set the oomph_info stream pointer to the null stream to
1161  // disable all possible output
1163  }
1164 
1165  // Solve the problem
1166  problem_pt->newton_solve();
1167  }
1168  } // if (GlobalParameters::Use_adaptation_flag)
1169 
1170  // Delete the problem pointer
1171  delete problem_pt;
1172 
1173  // Make it a null pointer
1174  problem_pt=0;
1175 } // End of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Problem class.
Definition: structured_cubic_point_source.cc:272
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
std::ostream *& stream_pt()
Access function for the stream pointer.
Definition: oomph_definitions.h:464
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
void adapt(unsigned &n_refined, unsigned &n_unrefined)
Definition: problem.cc:13666
Definition: refineable_pml_helmholtz_elements.h:199
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
unsigned Enable_test_pml_mapping_flag
Definition: structured_cubic_point_source.cc:215
unsigned Use_adaptation_flag
Definition: structured_cubic_point_source.cc:70
unsigned Min_refinement_level
The minimum level of uniform refinement.
Definition: structured_cubic_point_source.cc:59
double Alpha_shift
Definition: structured_cubic_point_source.cc:129
unsigned Output_management_flag
Definition: structured_cubic_point_source.cc:92
unsigned Doc_convergence_flag
Definition: structured_cubic_point_source.cc:98
unsigned Pre_smoother_flag
Definition: structured_cubic_point_source.cc:75
unsigned Add_refinement_level
The additional levels of uniform refinement.
Definition: structured_cubic_point_source.cc:62
unsigned Nnode_1d
The number of nodes in one direction (default=2)
Definition: structured_cubic_point_source.cc:56
double K_squared
Square of the wavenumber.
Definition: helmholtz_point_source.cc:60
unsigned Post_smoother_flag
Definition: structured_cubic_point_source.cc:80
unsigned Linear_solver_flag
Definition: structured_cubic_point_source.cc:85
void update_parameters()
Update parameters.
Definition: extrude_with_macro_element_representation.cc:274
unsigned N_adaptations
The number of adaptations allowed by the Newton solver.
Definition: structured_cubic_point_source.cc:65
DocInfo Doc_info
Helper for documenting.
Definition: extrude_triangle_generated_mesh.cc:57
std::ostream * Stream_pt
Definition: structured_cubic_point_source.cc:104
double Omega
The value of the damping factor for the damped Jacobi smoother.
Definition: structured_cubic_point_source.cc:244
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
Nullstream oomph_nullstream
Single (global) instantiation of the Nullstream.
Definition: oomph_definitions.cc:313
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 oomph::Problem::adapt(), GlobalParameters::Add_refinement_level, GlobalParameters::Alpha_shift, GlobalParameters::Doc_convergence_flag, GlobalParameters::Doc_info, oomph::CommandLineArgs::doc_specified_flags(), GlobalParameters::Enable_test_pml_mapping_flag, i, GlobalParameters::K_squared, GlobalParameters::Linear_solver_flag, GlobalParameters::Min_refinement_level, GlobalParameters::N_adaptations, oomph::Problem::newton_solve(), GlobalParameters::Nnode_1d, Smoother_Factory_Function_Helper::Omega, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, oomph::oomph_nullstream, GlobalParameters::Output_management_flag, oomph::CommandLineArgs::parse_and_assign(), GlobalParameters::Post_smoother_flag, GlobalParameters::Pre_smoother_flag, oomph::Problem::refine_uniformly(), oomph::DocInfo::set_directory(), Flag_definition::setup(), oomph::CommandLineArgs::specify_command_line_flag(), GlobalParameters::Stream_pt, oomph::OomphInfo::stream_pt(), GlobalParameters::update_parameters(), and GlobalParameters::Use_adaptation_flag.