unstructured_sphere_scattering.cc File Reference
#include <complex>
#include <cmath>
#include "generic.h"
#include "fourier_decomposed_helmholtz.h"
#include "meshes/triangle_mesh.h"
#include "oomph_crbond_bessel.h"

Classes

class  FourierDecomposedHelmholtzProblem< ELEMENT >
 Problem class. More...
 

Namespaces

 PlanarWave
 
 ProblemParameters
 Namespace for exact solution and problem parameters.
 

Functions

std::complex< doublePlanarWave::I (0.0, 1.0)
 Imaginary unit. More...
 
void PlanarWave::get_exact_u (const Vector< double > &x, Vector< double > &u)
 Exact solution as a Vector of size 2, containing real and imag parts. More...
 
void PlanarWave::plot ()
 Plot. More...
 
std::complex< doubleProblemParameters::I (0.0, 1.0)
 Imaginary unit. More...
 
void ProblemParameters::get_exact_u (const Vector< double > &x, Vector< double > &u)
 Exact solution as a Vector of size 2, containing real and imag parts. More...
 
void ProblemParameters::exact_minus_dudr (const Vector< double > &x, std::complex< double > &flux)
 
int main (int argc, char **argv)
 Driver code for Fourier decomposed Helmholtz problem. More...
 

Variables

Vector< doubleProblemParameters::Coeff (N_terms, 1.0)
 Coefficients in the exact solution. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char **  argv 
)

Driver code for Fourier decomposed Helmholtz problem.

882 {
883  // Store command line arguments
884  CommandLineArgs::setup(argc,argv);
885 
886  // Define possible command line arguments and parse the ones that
887  // were actually specified
888 
889  // Square domain without DtN
891 
892  // Parse command line
894 
895  // Doc what has actually been specified on the command line
897 
898  // Check if the claimed representation of a planar wave in
899  // the tutorial is correct -- of course it is!
900  //PlanarWave::plot();
901 
902  // Test Bessel/Hankel functions
903  //-----------------------------
904  {
905  // Number of Bessel functions to be computed
906  unsigned n=3;
907 
908  // Offset of Bessel function order (less than 1!)
909  double bessel_offset=0.5;
910 
911  ofstream bessely_file("besselY.dat");
912  ofstream bessely_deriv_file("dbesselY.dat");
913 
914  ofstream besselj_file("besselJ.dat");
915  ofstream besselj_deriv_file("dbesselJ.dat");
916 
917  // Evaluate Bessel/Hankel functions
918  Vector<double> jv(n+1);
919  Vector<double> yv(n+1);
920  Vector<double> djv(n+1);
921  Vector<double> dyv(n+1);
922  double x_min=0.5;
923  double x_max=5.0;
924  unsigned nplot=100;
925  for (unsigned i=0;i<nplot;i++)
926  {
927  double x=x_min+(x_max-x_min)*double(i)/double(nplot-1);
928  double order_max_in=double(n)+bessel_offset;
929  double order_max_out=0;
930 
931  // This function returns vectors containing
932  // J_k(x), Y_k(x) and their derivatives
933  // up to k=order_max, with k increasing in
934  // integer increments starting with smallest
935  // positive value. So, e.g. for order_max=3.5
936  // jv[0] contains J_{1/2}(x),
937  // jv[1] contains J_{3/2}(x),
938  // jv[2] contains J_{5/2}(x),
939  // jv[3] contains J_{7/2}(x).
940  CRBond_Bessel::bessjyv(order_max_in,x,
941  order_max_out,
942  &jv[0],&yv[0],
943  &djv[0],&dyv[0]);
944  bessely_file << x << " ";
945  for (unsigned j=0;j<=n;j++)
946  {
947  bessely_file << yv[j] << " ";
948  }
949  bessely_file << std::endl;
950 
951  besselj_file << x << " ";
952  for (unsigned j=0;j<=n;j++)
953  {
954  besselj_file << jv[j] << " ";
955  }
956  besselj_file << std::endl;
957 
958  bessely_deriv_file << x << " ";
959  for (unsigned j=0;j<=n;j++)
960  {
961  bessely_deriv_file << dyv[j] << " ";
962  }
963  bessely_deriv_file << std::endl;
964 
965  besselj_deriv_file << x << " ";
966  for (unsigned j=0;j<=n;j++)
967  {
968  besselj_deriv_file << djv[j] << " ";
969  }
970  besselj_deriv_file << std::endl;
971 
972  }
973  bessely_file.close();
974  besselj_file.close();
975  bessely_deriv_file.close();
976  besselj_deriv_file.close();
977  }
978 
979 
980  // Test Legrendre Polynomials
981  //---------------------------
982  {
983  // Number of lower indices
984  unsigned n=3;
985 
986  ofstream some_file("legendre3.dat");
987  unsigned nplot=100;
988  for (unsigned i=0;i<nplot;i++)
989  {
990  double x=double(i)/double(nplot-1);
991 
992  some_file << x << " ";
993  for (unsigned j=0;j<=n;j++)
994  {
995  some_file << Legendre_functions_helper::plgndr2(n,j,x) << " ";
996  }
997  some_file << std::endl;
998  }
999  some_file.close();
1000  }
1001 
1002 
1003 #ifdef ADAPTIVE
1004 
1005  // Create the problem with 2D six-node elements from the
1006  // TFourierDecomposedHelmholtzElement family.
1009 
1010 #else
1011 
1012  // Create the problem with 2D six-node elements from the
1013  // TFourierDecomposedHelmholtzElement family.
1015  problem;
1016 
1017 #endif
1018 
1019  // Create label for output
1020  DocInfo doc_info;
1021 
1022  // Set output directory
1023  doc_info.set_directory("RESLT");
1024 
1025  // Solve for a few Fourier wavenumbers
1028  {
1029  // Step number
1031 
1032 
1033 
1034 #ifdef ADAPTIVE
1035 
1036  // Max. number of adaptations
1037  unsigned max_adapt=1;
1038 
1039  // Solve the problem with Newton's method, allowing
1040  // up to max_adapt mesh adaptations after every solve.
1041  problem.newton_solve(max_adapt);
1042 
1043 #else
1044 
1045  // Solve the problem
1046  problem.newton_solve();
1047 
1048 #endif
1049 
1050  //Output the solution
1051  problem.doc_solution(doc_info);
1052  }
1053 
1054 } //end of main
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Problem class.
Definition: sphere_scattering.cc:380
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
Fourier decomposed Helmholtz upgraded to become projectable.
Definition: fourier_decomposed_helmholtz_elements.h:647
Definition: Tfourier_decomposed_helmholtz_elements.h:65
int bessjyv(double v, double x, double &vm, double *jv, double *yv, double *djv, double *dyv)
Definition: crbond_bessel.cc:1050
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
int N_fourier
Fourier wave number.
Definition: sphere_scattering.cc:228
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
double plgndr2(const unsigned &l, const unsigned &m, const double &x)
Legendre polynomials depending on two parameters.
Definition: fourier_decomposed_helmholtz_elements.cc:97
list x
Definition: plotDoE.py:28
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References CRBond_Bessel::bessjyv(), oomph::CommandLineArgs::doc_specified_flags(), i, j, n, ProblemParameters::N_fourier, oomph::DocInfo::number(), oomph::CommandLineArgs::parse_and_assign(), oomph::Legendre_functions_helper::plgndr2(), problem, oomph::DocInfo::set_directory(), oomph::CommandLineArgs::setup(), oomph::CommandLineArgs::specify_command_line_flag(), and plotDoE::x.