oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT > Class Template Reference

#include <fourier_decomposed_helmholtz_bc_elements.h>

+ Inheritance diagram for oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >:

Public Member Functions

 FourierDecomposedHelmholtzDtNMesh (const double &outer_radius, const unsigned &n_terms)
 
void setup_gamma ()
 of the mesh's constituent elements More...
 
Vector< std::complex< double > > & gamma_at_gauss_point (FiniteElement *el_pt)
 
Vector< std::map< unsigned, std::complex< double > > > & d_gamma_at_gauss_point (FiniteElement *el_pt)
 
doubleouter_radius ()
 The outer radius. More...
 
unsignedn_terms ()
 
- Public Member Functions inherited from oomph::Mesh
 Mesh ()
 Default constructor. More...
 
 Mesh (const Vector< Mesh * > &sub_mesh_pt)
 
void merge_meshes (const Vector< Mesh * > &sub_mesh_pt)
 
virtual void setup_boundary_element_info ()
 
virtual void setup_boundary_element_info (std::ostream &outfile)
 
virtual void reset_boundary_element_info (Vector< unsigned > &ntmp_boundary_elements, Vector< Vector< unsigned >> &ntmp_boundary_elements_in_region, Vector< FiniteElement * > &deleted_elements)
 Virtual function to perform the reset boundary elements info rutines. More...
 
template<class BULK_ELEMENT >
void doc_boundary_coordinates (const unsigned &b, std::ofstream &the_file)
 
virtual void scale_mesh (const double &factor)
 
 Mesh (const Mesh &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Mesh &)=delete
 Broken assignment operator. More...
 
virtual ~Mesh ()
 Virtual Destructor to clean up all memory. More...
 
void flush_element_and_node_storage ()
 
void flush_element_storage ()
 
void flush_node_storage ()
 
Node *& node_pt (const unsigned long &n)
 Return pointer to global node n. More...
 
Nodenode_pt (const unsigned long &n) const
 Return pointer to global node n (const version) More...
 
GeneralisedElement *& element_pt (const unsigned long &e)
 Return pointer to element e. More...
 
GeneralisedElementelement_pt (const unsigned long &e) const
 Return pointer to element e (const version) More...
 
const Vector< GeneralisedElement * > & element_pt () const
 Return reference to the Vector of elements. More...
 
Vector< GeneralisedElement * > & element_pt ()
 Return reference to the Vector of elements. More...
 
FiniteElementfinite_element_pt (const unsigned &e) const
 
Node *& boundary_node_pt (const unsigned &b, const unsigned &n)
 Return pointer to node n on boundary b. More...
 
Nodeboundary_node_pt (const unsigned &b, const unsigned &n) const
 Return pointer to node n on boundary b. More...
 
void set_nboundary (const unsigned &nbound)
 Set the number of boundaries in the mesh. More...
 
void remove_boundary_nodes ()
 Clear all pointers to boundary nodes. More...
 
void remove_boundary_nodes (const unsigned &b)
 
void remove_boundary_node (const unsigned &b, Node *const &node_pt)
 Remove a node from the boundary b. More...
 
void add_boundary_node (const unsigned &b, Node *const &node_pt)
 Add a (pointer to) a node to the b-th boundary. More...
 
void copy_boundary_node_data_from_nodes ()
 
bool boundary_coordinate_exists (const unsigned &i) const
 Indicate whether the i-th boundary has an intrinsic coordinate. More...
 
unsigned long nelement () const
 Return number of elements in the mesh. More...
 
unsigned long nnode () const
 Return number of nodes in the mesh. More...
 
unsigned ndof_types () const
 Return number of dof types in mesh. More...
 
unsigned elemental_dimension () const
 Return number of elemental dimension in mesh. More...
 
unsigned nodal_dimension () const
 Return number of nodal dimension in mesh. More...
 
void add_node_pt (Node *const &node_pt)
 Add a (pointer to a) node to the mesh. More...
 
void add_element_pt (GeneralisedElement *const &element_pt)
 Add a (pointer to) an element to the mesh. More...
 
virtual void node_update (const bool &update_all_solid_nodes=false)
 
virtual void reorder_nodes (const bool &use_old_ordering=true)
 
virtual void get_node_reordering (Vector< Node * > &reordering, const bool &use_old_ordering=true) const
 
template<class BULK_ELEMENT , template< class > class FACE_ELEMENT>
void build_face_mesh (const unsigned &b, Mesh *const &face_mesh_pt)
 
unsigned self_test ()
 Self-test: Check elements and nodes. Return 0 for OK. More...
 
void max_and_min_element_size (double &max_size, double &min_size)
 
double total_size ()
 
void check_inverted_elements (bool &mesh_has_inverted_elements, std::ofstream &inverted_element_file)
 
void check_inverted_elements (bool &mesh_has_inverted_elements)
 
unsigned check_for_repeated_nodes (const double &epsilon=1.0e-12)
 
Vector< Node * > prune_dead_nodes ()
 
unsigned nboundary () const
 Return number of boundaries. More...
 
unsigned long nboundary_node (const unsigned &ibound) const
 Return number of nodes on a particular boundary. More...
 
FiniteElementboundary_element_pt (const unsigned &b, const unsigned &e) const
 Return pointer to e-th finite element on boundary b. More...
 
Nodeget_some_non_boundary_node () const
 
unsigned nboundary_element (const unsigned &b) const
 Return number of finite elements that are adjacent to boundary b. More...
 
int face_index_at_boundary (const unsigned &b, const unsigned &e) const
 
virtual void dump (std::ofstream &dump_file, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
void dump (const std::string &dump_file_name, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
virtual void read (std::ifstream &restart_file)
 Read solution from restart file. More...
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
void output (std::ostream &outfile)
 Output for all elements. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output (FILE *file_pt)
 Output for all elements (C-style output) More...
 
void output (FILE *file_pt, const unsigned &nplot)
 Output at f(n_plot) points in each element (C-style output) More...
 
void output (const std::string &output_filename)
 Output for all elements. More...
 
void output (const std::string &output_filename, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
 Output a given Vector function at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt)
 
void output_boundaries (std::ostream &outfile)
 Output the nodes on the boundaries (into separate tecplot zones) More...
 
void output_boundaries (const std::string &output_filename)
 
void assign_initial_values_impulsive ()
 Assign initial values for an impulsive start. More...
 
void shift_time_values ()
 
void calculate_predictions ()
 
void set_nodal_and_elemental_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_consistent_pinned_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 Set consistent values for pinned data in continuation. More...
 
bool does_pointer_correspond_to_mesh_data (double *const &parameter_pt)
 Does the double pointer correspond to any mesh data. More...
 
void set_nodal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set the timestepper associated with the nodal data in the mesh. More...
 
void set_elemental_internal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void compute_norm (double &norm)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Returns the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
bool is_mesh_distributed () const
 Boolean to indicate if Mesh has been distributed. More...
 
OomphCommunicatorcommunicator_pt () const
 
void delete_all_external_storage ()
 Wipe the storage for all externally-based elements. More...
 

Private Attributes

double Outer_radius
 Outer radius. More...
 
unsigned N_terms
 Number of terms used in the Gamma computation. More...
 
std::map< FiniteElement *, Vector< std::complex< double > > > Gamma_at_gauss_point
 
std::map< FiniteElement *, Vector< std::map< unsigned, std::complex< double > > > > D_Gamma_at_gauss_point
 

Additional Inherited Members

- Public Types inherited from oomph::Mesh
typedef void(FiniteElement::* SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln)
 
typedef void(FiniteElement::* UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln)
 
- Static Public Attributes inherited from oomph::Mesh
static Steady< 0 > Default_TimeStepper
 The Steady Timestepper. More...
 
static bool Suppress_warning_about_empty_mesh_level_time_stepper_function
 Static boolean flag to control warning about mesh level timesteppers. More...
 
- Protected Member Functions inherited from oomph::Mesh
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign (global) equation numbers to the nodes. More...
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign local equation numbers in all elements. More...
 
void convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
 
void convert_to_boundary_node (Node *&node_pt)
 
- Protected Attributes inherited from oomph::Mesh
Vector< Vector< Node * > > Boundary_node_pt
 
bool Lookup_for_elements_next_boundary_is_setup
 
Vector< Vector< FiniteElement * > > Boundary_element_pt
 
Vector< Vector< int > > Face_index_at_boundary
 
Vector< Node * > Node_pt
 Vector of pointers to nodes. More...
 
Vector< GeneralisedElement * > Element_pt
 Vector of pointers to generalised elements. More...
 
std::vector< boolBoundary_coordinate_exists
 

Detailed Description

template<class ELEMENT>
class oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >

/////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// ================================================================= Mesh for DtN boundary condition elements – provides functionality to apply Sommerfeld radiation condtion

via DtN BC.

Constructor & Destructor Documentation

◆ FourierDecomposedHelmholtzDtNMesh()

template<class ELEMENT >
oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::FourierDecomposedHelmholtzDtNMesh ( const double outer_radius,
const unsigned n_terms 
)
inline

Constructor: Specify radius of outer boundary and number of terms used in the computation of the gamma integral

358  {
359  }
double & outer_radius()
The outer radius.
Definition: fourier_decomposed_helmholtz_bc_elements.h:381
unsigned N_terms
Number of terms used in the Gamma computation.
Definition: fourier_decomposed_helmholtz_bc_elements.h:398
double Outer_radius
Outer radius.
Definition: fourier_decomposed_helmholtz_bc_elements.h:395
unsigned & n_terms()
Definition: fourier_decomposed_helmholtz_bc_elements.h:388

Member Function Documentation

◆ d_gamma_at_gauss_point()

template<class ELEMENT >
Vector<std::map<unsigned, std::complex<double> > >& oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::d_gamma_at_gauss_point ( FiniteElement el_pt)
inline

Derivative of Gamma integral w.r.t global unknown, evaluated at Gauss points for specified element

376  {
377  return D_Gamma_at_gauss_point[el_pt];
378  }
std::map< FiniteElement *, Vector< std::map< unsigned, std::complex< double > > > > D_Gamma_at_gauss_point
Definition: fourier_decomposed_helmholtz_bc_elements.h:409

References oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::D_Gamma_at_gauss_point.

◆ gamma_at_gauss_point()

template<class ELEMENT >
Vector<std::complex<double> >& oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::gamma_at_gauss_point ( FiniteElement el_pt)
inline

Gamma integral evaluated at Gauss points for specified element

368  {
369  return Gamma_at_gauss_point[el_pt];
370  }
std::map< FiniteElement *, Vector< std::complex< double > > > Gamma_at_gauss_point
Definition: fourier_decomposed_helmholtz_bc_elements.h:403

References oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::Gamma_at_gauss_point.

◆ n_terms()

template<class ELEMENT >
unsigned& oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::n_terms ( )
inline

Number of terms used in the computation of the gamma integral

389  {
390  return N_terms;
391  }

References oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::N_terms.

◆ outer_radius()

template<class ELEMENT >
double& oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::outer_radius ( )
inline

The outer radius.

382  {
383  return Outer_radius;
384  }

References oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::Outer_radius.

◆ setup_gamma()

template<class ELEMENT >
void oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::setup_gamma

of the mesh's constituent elements

Compute and store the gamma integral at all integration points of the constituent elements.

////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// ================================================================ Compute and store the gamma integral and derivates

Imaginary unit

937  {
938 #ifdef PARANOID
939  {
940  // Loop over elements e
941  unsigned nel = this->nelement();
942  for (unsigned e = 0; e < nel; e++)
943  {
944  FiniteElement* fe_pt = finite_element_pt(e);
945  unsigned nnod = fe_pt->nnode();
946  for (unsigned j = 0; j < nnod; j++)
947  {
948  Node* nod_pt = fe_pt->node_pt(j);
949 
950  // Extract nodal coordinates from node:
951  Vector<double> x(2);
952  x[0] = nod_pt->x(0);
953  x[1] = nod_pt->x(1);
954 
955  // Evaluate the radial distance
956  double r = sqrt(x[0] * x[0] + x[1] * x[1]);
957 
958  // Check
959  if (Outer_radius == 0.0)
960  {
961  throw OomphLibError("Outer radius for DtN BC must not be zero!",
964  }
965 
966  if (std::fabs((r - this->Outer_radius) / Outer_radius) >
968  {
969  std::ostringstream error_stream;
970  error_stream
971  << "Node at " << x[0] << " " << x[1] << " has radius " << r
972  << " which does not "
973  << " agree with \nspecified outer radius " << this->Outer_radius
974  << " within relative tolerance "
976  << ".\nYou can adjust the tolerance via\n"
977  << "ToleranceForFourierDecomposedHelmholtzOuterBoundary::Tol\n"
978  << "or recompile without PARANOID.\n";
979  throw OomphLibError(error_stream.str(),
982  }
983  }
984  }
985  }
986 #endif
987 
988 
989  // Get parameters from first element
990  FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>* el_pt =
991  dynamic_cast<FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>*>(
992  this->element_pt(0));
993  double k =
994  sqrt(dynamic_cast<ELEMENT*>(el_pt->bulk_element_pt())->k_squared());
995  int n_fourier_decomposed =
996  dynamic_cast<ELEMENT*>(el_pt->bulk_element_pt())->fourier_wavenumber();
997  double n_hankel_order_max = double(N_terms) + 0.5;
998  double n_hankel_order_tmp = 0.0;
999 
1001  std::complex<double> I(0.0, 1.0);
1002 
1003  // Precompute factors in sum
1004  Vector<std::complex<double>> h_a(N_terms + 1), hp_a(N_terms + 1),
1005  q(N_terms + 1, std::complex<double>(0.0, 0.0));
1006  Vector<double> jv(N_terms + 1), yv(N_terms + 1), djv(N_terms + 1),
1007  dyv(N_terms + 1);
1008 
1009  // Get Bessel functions
1010  CRBond_Bessel::bessjyv(n_hankel_order_max,
1011  Outer_radius * k,
1012  n_hankel_order_tmp,
1013  &jv[0],
1014  &yv[0],
1015  &djv[0],
1016  &dyv[0]);
1017 
1018  // Assemble Hankel functions and their derivatives
1019  for (unsigned j = 0; j < N_terms + 1; j++)
1020  {
1021  h_a[j] = jv[j] + I * yv[j];
1022  hp_a[j] = djv[j] + I * dyv[j];
1023  }
1024 
1025  // Precompute relevant terms in sum
1026  for (unsigned i = n_fourier_decomposed; i < N_terms; i++)
1027  {
1028  double n_1 =
1029  Legendre_functions_helper::factorial(i - n_fourier_decomposed);
1030  double n_2 =
1031  Legendre_functions_helper::factorial(i + n_fourier_decomposed);
1032 
1033  q[i] = k * sqrt(MathematicalConstants::Pi / (2.0 * k * Outer_radius)) *
1034  (hp_a[i] - h_a[i] / (2.0 * k * Outer_radius)) *
1035  (2.0 * double(i) + 1.0) /
1036  (2.0 * sqrt(MathematicalConstants::Pi / (2.0 * k * Outer_radius)) *
1037  h_a[i]) *
1038  (n_1 / n_2);
1039  }
1040 
1041  // first loop over elements e
1042  unsigned nel = this->nelement();
1043  for (unsigned e = 0; e < nel; e++)
1044  {
1045  // Get a pointer to element
1046  FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>* el_pt =
1047  dynamic_cast<FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>*>(
1048  this->element_pt(e));
1049 
1050  // Set the value of n_intpt
1051  const unsigned n_intpt = el_pt->integral_pt()->nweight();
1052 
1053  // initialise gamma integral and its derivatives
1054  Vector<std::complex<double>> gamma_vector(n_intpt,
1055  std::complex<double>(0.0, 0.0));
1056  Vector<std::map<unsigned, std::complex<double>>> d_gamma_vector(n_intpt);
1057 
1058  // Loop over the integration points
1059  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1060  {
1061  // Allocate and initialise coordinate
1062  Vector<double> x(el_pt->dim() + 1, 0.0);
1063 
1064  // Set the Vector to hold local coordinates
1065  Vector<double> s(el_pt->dim(), 0.0);
1066  for (unsigned i = 0; i < el_pt->dim(); i++)
1067  {
1068  s[i] = el_pt->integral_pt()->knot(ipt, i);
1069  }
1070 
1071  // Get the coordinates of the integration point
1072  el_pt->interpolated_x(s, x);
1073 
1074  // Polar angle
1075  double theta = atan2(x[0], x[1]);
1076 
1077  // Elemental contribution to gamma integral and its derivative
1078  std::complex<double> gamma_con(0.0, 0.0);
1079  std::map<unsigned, std::complex<double>> d_gamma_con;
1080 
1081  // loop over the Fourier terms
1082  for (unsigned nn = n_fourier_decomposed; nn < N_terms; nn++)
1083  {
1084  // Second loop over the elements
1085  // to evaluate the complete integral
1086  for (unsigned ee = 0; ee < nel; ee++)
1087  {
1088  FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>* eel_pt =
1089  dynamic_cast<
1090  FourierDecomposedHelmholtzDtNBoundaryElement<ELEMENT>*>(
1091  this->element_pt(ee));
1092 
1093  // contribution of the positive term in the sum
1094  eel_pt->compute_gamma_contribution(
1095  theta, nn, gamma_con, d_gamma_con);
1096 
1097  unsigned n_node = eel_pt->nnode();
1098 
1099  gamma_vector[ipt] += q[nn] * gamma_con;
1100  for (unsigned l = 0; l < n_node; l++)
1101  {
1102  // Add the contribution of the real local data
1103  int local_unknown_real = eel_pt->nodal_local_eqn(
1104  l, eel_pt->u_index_fourier_decomposed_helmholtz().real());
1105 
1106  if (local_unknown_real >= 0)
1107  {
1108  int global_eqn_real = eel_pt->eqn_number(local_unknown_real);
1109  d_gamma_vector[ipt][global_eqn_real] +=
1110  q[nn] * d_gamma_con[global_eqn_real];
1111  }
1112 
1113  // Add the contribution of the imag local data
1114  int local_unknown_imag = eel_pt->nodal_local_eqn(
1115  l, eel_pt->u_index_fourier_decomposed_helmholtz().imag());
1116 
1117  if (local_unknown_imag >= 0)
1118  {
1119  int global_eqn_imag = eel_pt->eqn_number(local_unknown_imag);
1120  d_gamma_vector[ipt][global_eqn_imag] +=
1121  q[nn] * d_gamma_con[global_eqn_imag];
1122  }
1123  } // end of loop over the node
1124  } // End of second loop over the elements
1125  } // End of loop over Fourier terms
1126  } // end of loop over integration point
1127 
1128  // Store it in map
1129  Gamma_at_gauss_point[el_pt] = gamma_vector;
1130  D_Gamma_at_gauss_point[el_pt] = d_gamma_vector;
1131 
1132  } // end of first loop over element
1133  }
AnnoyingScalar atan2(const AnnoyingScalar &y, const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:139
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
Definition: mesh.h:460
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
#define I
Definition: main.h:127
double theta
Definition: two_d_biharmonic.cc:236
int bessjyv(double v, double x, double &vm, double *jv, double *yv, double *djv, double *dyv)
Definition: crbond_bessel.cc:1050
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
r
Definition: UniformPSDSelfTest.py:20
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
double factorial(const unsigned &l)
Factorial.
Definition: fourier_decomposed_helmholtz_elements.cc:40
const double Pi
50 digits from maple
Definition: oomph_utilities.h:157
double Tol
Definition: fourier_decomposed_helmholtz_bc_elements.h:921
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References atan2(), CRBond_Bessel::bessjyv(), oomph::FaceElement::bulk_element_pt(), oomph::FourierDecomposedHelmholtzDtNBoundaryElement< ELEMENT >::compute_gamma_contribution(), oomph::FiniteElement::dim(), e(), oomph::GeneralisedElement::eqn_number(), boost::multiprecision::fabs(), oomph::Legendre_functions_helper::factorial(), i, I, oomph::FiniteElement::integral_pt(), oomph::FaceElement::interpolated_x(), j, k, oomph::Integral::knot(), PlanarWave::N_terms, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::FiniteElement::node_pt(), oomph::Integral::nweight(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, GlobalParameters::Outer_radius, oomph::MathematicalConstants::Pi, Eigen::numext::q, UniformPSDSelfTest::r, s, sqrt(), BiharmonicTestFunctions2::theta, oomph::ToleranceForFourierDecomposedHelmholtzOuterBoundary::Tol, oomph::FourierDecomposedHelmholtzBCElementBase< ELEMENT >::u_index_fourier_decomposed_helmholtz(), plotDoE::x, and oomph::Node::x().

Member Data Documentation

◆ D_Gamma_at_gauss_point

template<class ELEMENT >
std::map<FiniteElement*, Vector<std::map<unsigned, std::complex<double> > > > oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::D_Gamma_at_gauss_point
private

Container to store the derivate of Gamma integral w.r.t global unknown evaluated at Gauss points for specified element

Referenced by oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::d_gamma_at_gauss_point().

◆ Gamma_at_gauss_point

template<class ELEMENT >
std::map<FiniteElement*, Vector<std::complex<double> > > oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::Gamma_at_gauss_point
private

Container to store the gamma integral for given Gauss point and element

Referenced by oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::gamma_at_gauss_point().

◆ N_terms

template<class ELEMENT >
unsigned oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::N_terms
private

Number of terms used in the Gamma computation.

Referenced by oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::n_terms().

◆ Outer_radius

template<class ELEMENT >
double oomph::FourierDecomposedHelmholtzDtNMesh< ELEMENT >::Outer_radius
private

The documentation for this class was generated from the following file: