oomph::HelmholtzEquations< DIM > Class Template Referenceabstract

#include <helmholtz_elements.h>

+ Inheritance diagram for oomph::HelmholtzEquations< DIM >:

Public Types

typedef void(* HelmholtzSourceFctPt) (const Vector< double > &x, std::complex< double > &f)
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 

Public Member Functions

 HelmholtzEquations ()
 Constructor (must initialise the Source_fct_pt to null) More...
 
 HelmholtzEquations (const HelmholtzEquations &dummy)=delete
 Broken copy constructor. More...
 
virtual std::complex< unsignedu_index_helmholtz () const
 Broken assignment operator. More...
 
double *& k_squared_pt ()
 Get pointer to square of wavenumber. More...
 
double k_squared ()
 Get the square of wavenumber. More...
 
unsigned nscalar_paraview () const
 
void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
std::string scalar_name_paraview (const unsigned &i) const
 
void output (std::ostream &outfile)
 Output with default number of plot points. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 
void output_real (std::ostream &outfile, const double &phi, const unsigned &n_plot)
 
void output (FILE *file_pt)
 C_style output with default number of plot points. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
void output_real_fct (std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Dummy, time dependent error checker. More...
 
HelmholtzSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
HelmholtzSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
virtual void get_source_helmholtz (const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
 
void get_flux (const Vector< double > &s, Vector< std::complex< double >> &flux) const
 Get flux: flux[i] = du/dx_i for real and imag part. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the element's contribution to its residual vector (wrapper) More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
std::complex< doubleinterpolated_u_helmholtz (const Vector< double > &s) const
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (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 (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

virtual double dshape_and_dtest_eulerian_helmholtz (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_helmholtz (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_helmholtz (Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 

Protected Attributes

HelmholtzSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
doubleK_squared_pt
 Pointer to square of wavenumber. More...
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Additional Inherited Members

- Static Public Attributes inherited from oomph::FiniteElement
static double Tolerance_for_singular_jacobian = 1.0e-16
 Tolerance below which the jacobian is considered singular. More...
 
static bool Accept_negative_jacobian = false
 
static bool Suppress_output_while_checking_for_inverted_elements
 
- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Static Protected Attributes inherited from oomph::FiniteElement
static const unsigned Default_Initial_Nvalue = 0
 Default value for the number of values at a node. More...
 
static const double Node_location_tolerance = 1.0e-14
 
static const unsigned N2deriv [] = {0, 1, 3, 6}
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

template<unsigned DIM>
class oomph::HelmholtzEquations< DIM >

A class for all isoparametric elements that solve the Helmholtz equations.

\[ \frac{\partial^2 u}{\partial x_i^2} + k^2 u = f(x_j) \]

This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class.

Member Typedef Documentation

◆ HelmholtzSourceFctPt

template<unsigned DIM>
typedef void(* oomph::HelmholtzEquations< DIM >::HelmholtzSourceFctPt) (const Vector< double > &x, std::complex< double > &f)

Function pointer to source function fct(x,f(x)) – x is a Vector!

Constructor & Destructor Documentation

◆ HelmholtzEquations() [1/2]

template<unsigned DIM>
oomph::HelmholtzEquations< DIM >::HelmholtzEquations ( )
inline

Constructor (must initialise the Source_fct_pt to null)

65 : Source_fct_pt(0), K_squared_pt(0) {}
double * K_squared_pt
Pointer to square of wavenumber.
Definition: helmholtz_elements.h:437
HelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: helmholtz_elements.h:434

◆ HelmholtzEquations() [2/2]

template<unsigned DIM>
oomph::HelmholtzEquations< DIM >::HelmholtzEquations ( const HelmholtzEquations< DIM > &  dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ compute_error() [1/2]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::compute_error ( std::ostream &  outfile,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt,
double error,
double norm 
)
virtual

Get error against and norm of exact solution.

Validate against exact solution

Solution is provided via function pointer. Plot error at a given number of plot points.

Reimplemented from oomph::FiniteElement.

487  {
488  // Initialise
489  error = 0.0;
490  norm = 0.0;
491 
492  // Vector of local coordinates
493  Vector<double> s(DIM);
494 
495  // Vector for coordintes
496  Vector<double> x(DIM);
497 
498  // Find out how many nodes there are in the element
499  unsigned n_node = nnode();
500 
501  Shape psi(n_node);
502 
503  // Set the value of n_intpt
504  unsigned n_intpt = integral_pt()->nweight();
505 
506  // Tecplot
507  outfile << "ZONE" << std::endl;
508 
509  // Exact solution Vector
510  Vector<double> exact_soln(2);
511 
512  // Loop over the integration points
513  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
514  {
515  // Assign values of s
516  for (unsigned i = 0; i < DIM; i++)
517  {
518  s[i] = integral_pt()->knot(ipt, i);
519  }
520 
521  // Get the integral weight
522  double w = integral_pt()->weight(ipt);
523 
524  // Get jacobian of mapping
525  double J = J_eulerian(s);
526 
527  // Premultiply the weights and the Jacobian
528  double W = w * J;
529 
530  // Get x position as Vector
531  interpolated_x(s, x);
532 
533  // Get FE function value
534  std::complex<double> u_fe = interpolated_u_helmholtz(s);
535 
536  // Get exact solution at this point
537  (*exact_soln_pt)(x, exact_soln);
538 
539  // Output x,y,...,error
540  for (unsigned i = 0; i < DIM; i++)
541  {
542  outfile << x[i] << " ";
543  }
544  outfile << exact_soln[0] << " " << exact_soln[1] << " "
545  << exact_soln[0] - u_fe.real() << " "
546  << exact_soln[1] - u_fe.imag() << std::endl;
547 
548  // Add to error and norm
549  norm +=
550  (exact_soln[0] * exact_soln[0] + exact_soln[1] * exact_soln[1]) * W;
551  error += ((exact_soln[0] - u_fe.real()) * (exact_soln[0] - u_fe.real()) +
552  (exact_soln[1] - u_fe.imag()) * (exact_soln[1] - u_fe.imag())) *
553  W;
554  }
555  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
std::complex< double > interpolated_u_helmholtz(const Vector< double > &s) const
Definition: helmholtz_elements.h:369
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

References DIM, calibrate::error, ProblemParameters::exact_soln(), i, J, s, w, oomph::QuadTreeNames::W, and plotDoE::x.

◆ compute_error() [2/2]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::compute_error ( std::ostream &  outfile,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt,
const double time,
double error,
double norm 
)
inlinevirtual

Dummy, time dependent error checker.

Reimplemented from oomph::FiniteElement.

268  {
269  throw OomphLibError(
270  "There is no time-dependent compute_error() for Helmholtz elements",
273  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dshape_and_dtest_eulerian_at_knot_helmholtz()

template<unsigned DIM>
virtual double oomph::HelmholtzEquations< DIM >::dshape_and_dtest_eulerian_at_knot_helmholtz ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of mapping

Implemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

◆ dshape_and_dtest_eulerian_helmholtz()

template<unsigned DIM>
virtual double oomph::HelmholtzEquations< DIM >::dshape_and_dtest_eulerian_helmholtz ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping

Implemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

◆ fill_in_contribution_to_jacobian()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)

Reimplemented from oomph::FiniteElement.

361  {
362  // Call the generic routine with the flag set to 1
363  fill_in_generic_residual_contribution_helmholtz(residuals, jacobian, 1);
364  }
virtual void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: helmholtz_elements.cc:49

References oomph::HelmholtzEquations< DIM >::fill_in_generic_residual_contribution_helmholtz().

◆ fill_in_contribution_to_residuals()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Add the element's contribution to its residual vector (wrapper)

Reimplemented from oomph::GeneralisedElement.

349  {
350  // Call the generic residuals function with flag set to 0
351  // using a dummy matrix argument
353  residuals, GeneralisedElement::Dummy_matrix, 0);
354  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::HelmholtzEquations< DIM >::fill_in_generic_residual_contribution_helmholtz().

◆ fill_in_generic_residual_contribution_helmholtz()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::fill_in_generic_residual_contribution_helmholtz ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const unsigned flag 
)
protectedvirtual

Compute element residual Vector only (if flag=and/or element Jacobian matrix

Compute element residual Vector and/or element Jacobian matrix

flag=1: compute both flag=0: compute only residual Vector

Pure version without hanging nodes

Reimplemented in oomph::RefineableHelmholtzEquations< DIM >.

53  {
54  // Find out how many nodes there are
55  const unsigned n_node = nnode();
56 
57  // Set up memory for the shape and test functions
58  Shape psi(n_node), test(n_node);
59  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
60 
61  // Set the value of n_intpt
62  const unsigned n_intpt = integral_pt()->nweight();
63 
64  // Integers to store the local equation and unknown numbers
65  int local_eqn_real = 0, local_unknown_real = 0;
66  int local_eqn_imag = 0, local_unknown_imag = 0;
67 
68  // Loop over the integration points
69  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
70  {
71  // Get the integral weight
72  double w = integral_pt()->weight(ipt);
73 
74  // Call the derivatives of the shape and test functions
76  ipt, psi, dpsidx, test, dtestdx);
77 
78  // Premultiply the weights and the Jacobian
79  double W = w * J;
80 
81  // Calculate local values of unknown
82  // Allocate and initialise to zero
83  std::complex<double> interpolated_u(0.0, 0.0);
84  Vector<double> interpolated_x(DIM, 0.0);
85  Vector<std::complex<double>> interpolated_dudx(DIM);
86 
87  // Calculate function value and derivatives:
88  //-----------------------------------------
89  // Loop over nodes
90  for (unsigned l = 0; l < n_node; l++)
91  {
92  // Loop over directions
93  for (unsigned j = 0; j < DIM; j++)
94  {
95  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
96  }
97 
98  // Get the nodal value of the helmholtz unknown
99  const std::complex<double> u_value(
102 
103  // Add to the interpolated value
104  interpolated_u += u_value * psi(l);
105 
106  // Loop over directions
107  for (unsigned j = 0; j < DIM; j++)
108  {
109  interpolated_dudx[j] += u_value * dpsidx(l, j);
110  }
111  }
112 
113  // Get source function
114  //-------------------
115  std::complex<double> source(0.0, 0.0);
117 
118  // Assemble residuals and Jacobian
119  //--------------------------------
120 
121  // Loop over the test functions
122  for (unsigned l = 0; l < n_node; l++)
123  {
124  // first, compute the real part contribution
125  //-------------------------------------------
126 
127  // Get the local equation
128  local_eqn_real = nodal_local_eqn(l, u_index_helmholtz().real());
129 
130  /*IF it's not a boundary condition*/
131  if (local_eqn_real >= 0)
132  {
133  // Add body force/source term and Helmholtz bit
134  residuals[local_eqn_real] +=
135  (source.real() - k_squared() * interpolated_u.real()) * test(l) * W;
136 
137  // The Helmholtz bit itself
138  for (unsigned k = 0; k < DIM; k++)
139  {
140  residuals[local_eqn_real] +=
141  interpolated_dudx[k].real() * dtestdx(l, k) * W;
142  }
143 
144  // Calculate the jacobian
145  //-----------------------
146  if (flag)
147  {
148  // Loop over the velocity shape functions again
149  for (unsigned l2 = 0; l2 < n_node; l2++)
150  {
151  local_unknown_real =
153  // If at a non-zero degree of freedom add in the entry
154  if (local_unknown_real >= 0)
155  {
156  // Add contribution to Elemental Matrix
157  for (unsigned i = 0; i < DIM; i++)
158  {
159  jacobian(local_eqn_real, local_unknown_real) +=
160  dpsidx(l2, i) * dtestdx(l, i) * W;
161  }
162  // Add the helmholtz contribution
163  jacobian(local_eqn_real, local_unknown_real) -=
164  k_squared() * psi(l2) * test(l) * W;
165  } // end of local_unknown
166  }
167  }
168  }
169 
170  // Second, compute the imaginary part contribution
171  //------------------------------------------------
172 
173  // Get the local equation
174  local_eqn_imag = nodal_local_eqn(l, u_index_helmholtz().imag());
175 
176  /*IF it's not a boundary condition*/
177  if (local_eqn_imag >= 0)
178  {
179  // Add body force/source term and Helmholtz bit
180  residuals[local_eqn_imag] +=
181  (source.imag() - k_squared() * interpolated_u.imag()) * test(l) * W;
182 
183  // The Helmholtz bit itself
184  for (unsigned k = 0; k < DIM; k++)
185  {
186  residuals[local_eqn_imag] +=
187  interpolated_dudx[k].imag() * dtestdx(l, k) * W;
188  }
189 
190  // Calculate the jacobian
191  //-----------------------
192  if (flag)
193  {
194  // Loop over the velocity shape functions again
195  for (unsigned l2 = 0; l2 < n_node; l2++)
196  {
197  local_unknown_imag =
199  // If at a non-zero degree of freedom add in the entry
200  if (local_unknown_imag >= 0)
201  {
202  // Add contribution to Elemental Matrix
203  for (unsigned i = 0; i < DIM; i++)
204  {
205  jacobian(local_eqn_imag, local_unknown_imag) +=
206  dpsidx(l2, i) * dtestdx(l, i) * W;
207  }
208  // Add the helmholtz contribution
209  jacobian(local_eqn_imag, local_unknown_imag) -=
210  k_squared() * psi(l2) * test(l) * W;
211  }
212  }
213  }
214  }
215  }
216  } // End of loop over integration points
217  }
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.cc:1686
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
Definition: helmholtz_elements.h:80
virtual double dshape_and_dtest_eulerian_at_knot_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void get_source_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Definition: helmholtz_elements.h:292
double k_squared()
Get the square of wavenumber.
Definition: helmholtz_elements.h:94
float real
Definition: datatypes.h:10
char char char int int * k
Definition: level2_impl.h:374
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References DIM, i, imag(), J, j, k, TestProblem::source(), Eigen::test, w, and oomph::QuadTreeNames::W.

Referenced by oomph::HelmholtzEquations< DIM >::fill_in_contribution_to_jacobian(), and oomph::HelmholtzEquations< DIM >::fill_in_contribution_to_residuals().

◆ get_flux()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::get_flux ( const Vector< double > &  s,
Vector< std::complex< double >> &  flux 
) const
inline

Get flux: flux[i] = du/dx_i for real and imag part.

312  {
313  // Find out how many nodes there are in the element
314  const unsigned n_node = nnode();
315 
316 
317  // Set up memory for the shape and test functions
318  Shape psi(n_node);
319  DShape dpsidx(n_node, DIM);
320 
321  // Call the derivatives of the shape and test functions
322  dshape_eulerian(s, psi, dpsidx);
323 
324  // Initialise to zero
325  const std::complex<double> zero(0.0, 0.0);
326  for (unsigned j = 0; j < DIM; j++)
327  {
328  flux[j] = zero;
329  }
330 
331  // Loop over nodes
332  for (unsigned l = 0; l < n_node; l++)
333  {
334  // Cache the complex value of the unknown
335  const std::complex<double> u_value(
336  this->nodal_value(l, u_index_helmholtz().real()),
337  this->nodal_value(l, u_index_helmholtz().imag()));
338  // Loop over derivative directions
339  for (unsigned j = 0; j < DIM; j++)
340  {
341  flux[j] += u_value * dpsidx(l, j);
342  }
343  }
344  }
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
EIGEN_DONT_INLINE Scalar zero()
Definition: svd_common.h:232

References DIM, oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), imag(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::HelmholtzEquations< DIM >::u_index_helmholtz(), and zero().

Referenced by oomph::RefineableHelmholtzEquations< DIM >::get_Z2_flux(), and oomph::THelmholtzElement< DIM, NNODE_1D >::get_Z2_flux().

◆ get_source_helmholtz()

template<unsigned DIM>
virtual void oomph::HelmholtzEquations< DIM >::get_source_helmholtz ( const unsigned ipt,
const Vector< double > &  x,
std::complex< double > &  source 
) const
inlinevirtual

Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-physics problems where the strength of the source function might be determined by another system of equations.

295  {
296  // If no source function has been set, return zero
297  if (Source_fct_pt == 0)
298  {
299  source = std::complex<double>(0.0, 0.0);
300  }
301  else
302  {
303  // Get source strength
304  (*Source_fct_pt)(x, source);
305  }
306  }

References TestProblem::source(), oomph::HelmholtzEquations< DIM >::Source_fct_pt, and plotDoE::x.

◆ interpolated_u_helmholtz()

template<unsigned DIM>
std::complex<double> oomph::HelmholtzEquations< DIM >::interpolated_u_helmholtz ( const Vector< double > &  s) const
inline

Return FE representation of function value u_helmholtz(s) at local coordinate s

371  {
372  // Find number of nodes
373  const unsigned n_node = nnode();
374 
375  // Local shape function
376  Shape psi(n_node);
377 
378  // Find values of shape function
379  shape(s, psi);
380 
381  // Initialise value of u
382  std::complex<double> interpolated_u(0.0, 0.0);
383 
384  // Get the index at which the helmholtz unknown is stored
385  const unsigned u_nodal_index_real = u_index_helmholtz().real();
386  const unsigned u_nodal_index_imag = u_index_helmholtz().imag();
387 
388  // Loop over the local nodes and sum
389  for (unsigned l = 0; l < n_node; l++)
390  {
391  // Make a temporary complex number from the stored data
392  const std::complex<double> u_value(
393  this->nodal_value(l, u_nodal_index_real),
394  this->nodal_value(l, u_nodal_index_imag));
395  // Add to the interpolated value
396  interpolated_u += u_value * psi[l];
397  }
398  return interpolated_u;
399  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::FiniteElement::shape(), and oomph::HelmholtzEquations< DIM >::u_index_helmholtz().

Referenced by oomph::HelmholtzEquations< DIM >::scalar_value_paraview().

◆ k_squared()

template<unsigned DIM>
double oomph::HelmholtzEquations< DIM >::k_squared ( )
inline

Get the square of wavenumber.

95  {
96 #ifdef PARANOID
97  if (K_squared_pt == 0)
98  {
99  throw OomphLibError(
100  "Please set pointer to k_squared using access fct to pointer!",
103  }
104 #endif
105  return *K_squared_pt;
106  }

References oomph::HelmholtzEquations< DIM >::K_squared_pt, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ k_squared_pt()

template<unsigned DIM>
double*& oomph::HelmholtzEquations< DIM >::k_squared_pt ( )
inline

Get pointer to square of wavenumber.

88  {
89  return K_squared_pt;
90  }

References oomph::HelmholtzEquations< DIM >::K_squared_pt.

Referenced by oomph::RefineableHelmholtzEquations< DIM >::further_build().

◆ nscalar_paraview()

template<unsigned DIM>
unsigned oomph::HelmholtzEquations< DIM >::nscalar_paraview ( ) const
inlinevirtual

Number of scalars/fields output by this element. Reimplements broken virtual function in base class.

Reimplemented from oomph::FiniteElement.

112  {
113  return 2;
114  }

◆ output() [1/4]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output ( FILE *  file_pt)
inlinevirtual

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

214  {
215  const unsigned n_plot = 5;
216  output(file_pt, n_plot);
217  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: helmholtz_elements.h:193

References oomph::HelmholtzEquations< DIM >::output().

◆ output() [2/4]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output ( FILE *  file_pt,
const unsigned nplot 
)
virtual

C-style output FE representation of soln: x,y,u_re,u_im or x,y,z,u_re,u_im at n_plot^DIM plot points

C-style output function:

x,y,u or x,y,z,u

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

332  {
333  // Vector of local coordinates
334  Vector<double> s(DIM);
335 
336  // Tecplot header info
337  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
338 
339  // Loop over plot points
340  unsigned num_plot_points = nplot_points(nplot);
341  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
342  {
343  // Get local coordinates of plot point
344  get_s_plot(iplot, nplot, s);
345  std::complex<double> u(interpolated_u_helmholtz(s));
346 
347  for (unsigned i = 0; i < DIM; i++)
348  {
349  fprintf(file_pt, "%g ", interpolated_x(s, i));
350  }
351 
352  for (unsigned i = 0; i < DIM; i++)
353  {
354  fprintf(file_pt, "%g ", interpolated_x(s, i));
355  }
356  fprintf(file_pt, "%g ", u.real());
357  fprintf(file_pt, "%g \n", u.imag());
358  }
359 
360  // Write tecplot footer (e.g. FE connectivity lists)
361  write_tecplot_zone_footer(file_pt, nplot);
362  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References DIM, i, and s.

◆ output() [3/4]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output ( std::ostream &  outfile)
inlinevirtual

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

194  {
195  const unsigned n_plot = 5;
196  output(outfile, n_plot);
197  }

Referenced by oomph::HelmholtzEquations< DIM >::output(), oomph::QHelmholtzElement< DIM, NNODE_1D >::output(), and oomph::THelmholtzElement< DIM, NNODE_1D >::output().

◆ output() [4/4]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output ( std::ostream &  outfile,
const unsigned nplot 
)
virtual

Output FE representation of soln: x,y,u_re,u_im or x,y,z,u_re,u_im at n_plot^DIM plot points

Output function:

x,y,u_re,u_imag or x,y,z,u_re,u_imag

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

256  {
257  // Vector of local coordinates
258  Vector<double> s(DIM);
259 
260  // Tecplot header info
261  outfile << tecplot_zone_string(nplot);
262 
263  // Loop over plot points
264  unsigned num_plot_points = nplot_points(nplot);
265  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
266  {
267  // Get local coordinates of plot point
268  get_s_plot(iplot, nplot, s);
269  std::complex<double> u(interpolated_u_helmholtz(s));
270  for (unsigned i = 0; i < DIM; i++)
271  {
272  outfile << interpolated_x(s, i) << " ";
273  }
274  outfile << u.real() << " " << u.imag() << std::endl;
275  }
276 
277  // Write tecplot footer (e.g. FE connectivity lists)
278  write_tecplot_zone_footer(outfile, nplot);
279  }

References DIM, i, and s.

◆ output_fct() [1/2]

template<unsigned DIM>
virtual void oomph::HelmholtzEquations< DIM >::output_fct ( std::ostream &  outfile,
const unsigned n_plot,
const double time,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt 
)
inlinevirtual

Output exact soln: (dummy time-dependent version to keep intel compiler happy)

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

236  {
237  throw OomphLibError(
238  "There is no time-dependent output_fct() for Helmholtz elements ",
241  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ output_fct() [2/2]

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output_fct ( std::ostream &  outfile,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
virtual

Output exact soln: x,y,u_re_exact,u_im_exact or x,y,z,u_re_exact,u_im_exact at n_plot^DIM plot points

Output exact solution

Solution is provided via function pointer. Plot at a given number of plot points.

x,y,u_exact or x,y,z,u_exact

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::THelmholtzElement< DIM, NNODE_1D >, and oomph::QHelmholtzElement< DIM, NNODE_1D >.

378  {
379  // Vector of local coordinates
380  Vector<double> s(DIM);
381 
382  // Vector for coordintes
383  Vector<double> x(DIM);
384 
385  // Tecplot header info
386  outfile << tecplot_zone_string(nplot);
387 
388  // Exact solution Vector
389  Vector<double> exact_soln(2);
390 
391  // Loop over plot points
392  unsigned num_plot_points = nplot_points(nplot);
393  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
394  {
395  // Get local coordinates of plot point
396  get_s_plot(iplot, nplot, s);
397 
398  // Get x position as Vector
399  interpolated_x(s, x);
400 
401  // Get exact solution at this point
402  (*exact_soln_pt)(x, exact_soln);
403 
404  // Output x,y,...,u_exact
405  for (unsigned i = 0; i < DIM; i++)
406  {
407  outfile << x[i] << " ";
408  }
409  outfile << exact_soln[0] << " " << exact_soln[1] << std::endl;
410  }
411 
412  // Write tecplot footer (e.g. FE connectivity lists)
413  write_tecplot_zone_footer(outfile, nplot);
414  }

References DIM, ProblemParameters::exact_soln(), i, s, and plotDoE::x.

Referenced by oomph::QHelmholtzElement< DIM, NNODE_1D >::output_fct(), and oomph::THelmholtzElement< DIM, NNODE_1D >::output_fct().

◆ output_real()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output_real ( std::ostream &  outfile,
const double phi,
const unsigned nplot 
)

Output function for real part of full time-dependent solution u = Re( (u_r +i u_i) exp(-i omega t) at phase angle omega t = phi. x,y,u or x,y,z,u at n_plot plot points in each coordinate direction

Output function for real part of full time-dependent solution

u = Re( (u_r +i u_i) exp(-i omega t)

at phase angle omega t = phi.

x,y,u or x,y,z,u

Output at nplot points in each coordinate direction

297  {
298  // Vector of local coordinates
299  Vector<double> s(DIM);
300 
301  // Tecplot header info
302  outfile << tecplot_zone_string(nplot);
303 
304  // Loop over plot points
305  unsigned num_plot_points = nplot_points(nplot);
306  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
307  {
308  // Get local coordinates of plot point
309  get_s_plot(iplot, nplot, s);
310  std::complex<double> u(interpolated_u_helmholtz(s));
311  for (unsigned i = 0; i < DIM; i++)
312  {
313  outfile << interpolated_x(s, i) << " ";
314  }
315  outfile << u.real() * cos(phi) + u.imag() * sin(phi) << std::endl;
316  }
317 
318  // Write tecplot footer (e.g. FE connectivity lists)
319  write_tecplot_zone_footer(outfile, nplot);
320  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137

References cos(), DIM, i, s, and sin().

Referenced by oomph::QHelmholtzElement< DIM, NNODE_1D >::output_real().

◆ output_real_fct()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::output_real_fct ( std::ostream &  outfile,
const double phi,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)

Output function for real part of full time-dependent fct u = Re( (u_r +i u_i) exp(-i omega t) at phase angle omega t = phi. x,y,u or x,y,z,u at n_plot plot points in each coordinate direction

Output function for real part of full time-dependent fct

u = Re( (u_r +i u_i) exp(-i omega t)

at phase angle omega t = phi.

x,y,u or x,y,z,u

Output at nplot points in each coordinate direction

434  {
435  // Vector of local coordinates
436  Vector<double> s(DIM);
437 
438  // Vector for coordintes
439  Vector<double> x(DIM);
440 
441  // Tecplot header info
442  outfile << tecplot_zone_string(nplot);
443 
444  // Exact solution Vector
445  Vector<double> exact_soln(2);
446 
447  // Loop over plot points
448  unsigned num_plot_points = nplot_points(nplot);
449  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
450  {
451  // Get local coordinates of plot point
452  get_s_plot(iplot, nplot, s);
453 
454  // Get x position as Vector
455  interpolated_x(s, x);
456 
457  // Get exact solution at this point
458  (*exact_soln_pt)(x, exact_soln);
459 
460  // Output x,y,...,u_exact
461  for (unsigned i = 0; i < DIM; i++)
462  {
463  outfile << x[i] << " ";
464  }
465  outfile << exact_soln[0] * cos(phi) + exact_soln[1] * sin(phi)
466  << std::endl;
467  }
468 
469  // Write tecplot footer (e.g. FE connectivity lists)
470  write_tecplot_zone_footer(outfile, nplot);
471  }

References cos(), DIM, ProblemParameters::exact_soln(), i, s, sin(), and plotDoE::x.

Referenced by oomph::QHelmholtzElement< DIM, NNODE_1D >::output_real_fct().

◆ scalar_name_paraview()

template<unsigned DIM>
std::string oomph::HelmholtzEquations< DIM >::scalar_name_paraview ( const unsigned i) const
inlinevirtual

Name of the i-th scalar field. Default implementation returns V1 for the first one, V2 for the second etc. Can (should!) be overloaded with more meaningful names in specific elements.

Reimplemented from oomph::FiniteElement.

165  {
166  switch (i)
167  {
168  case 0:
169  return "Real part";
170  break;
171 
172  case 1:
173  return "Imaginary part";
174  break;
175 
176  // Never get here
177  default:
178  std::stringstream error_stream;
179  error_stream
180  << "Helmholtz elements only store 2 fields so i must be 0 or 1"
181  << std::endl;
182  throw OomphLibError(error_stream.str(),
185 
186  // Dummy return for the default statement
187  return " ";
188  break;
189  }
190  }

References i, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ scalar_value_paraview()

template<unsigned DIM>
void oomph::HelmholtzEquations< DIM >::scalar_value_paraview ( std::ofstream &  file_out,
const unsigned i,
const unsigned nplot 
) const
inlinevirtual

Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specific element type.

Reimplemented from oomph::FiniteElement.

121  {
122  // Vector of local coordinates
123  Vector<double> s(DIM);
124 
125  // Loop over plot points
126  unsigned num_plot_points = nplot_points_paraview(nplot);
127  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
128  {
129  // Get local coordinates of plot point
130  get_s_plot(iplot, nplot, s);
131  std::complex<double> u(interpolated_u_helmholtz(s));
132 
133  // Paraview need to ouput the fileds seperatly so it loops through all
134  // the elements twice
135  switch (i)
136  {
137  // Real part first
138  case 0:
139  file_out << u.real() << std::endl;
140  break;
141 
142  // Imaginary part second
143  case 1:
144  file_out << u.imag() << std::endl;
145  break;
146 
147  // Never get here
148  default:
149  std::stringstream error_stream;
150  error_stream
151  << "Helmholtz elements only store 2 fields so i must be 0 or 1"
152  << std::endl;
153  throw OomphLibError(error_stream.str(),
156  break;
157  }
158  } // end of plotpoint loop
159  } // end scalar_value_paraview
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862

References DIM, oomph::FiniteElement::get_s_plot(), i, oomph::HelmholtzEquations< DIM >::interpolated_u_helmholtz(), oomph::FiniteElement::nplot_points_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and s.

◆ self_test()

template<unsigned DIM>
unsigned oomph::HelmholtzEquations< DIM >::self_test
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

225  {
226  bool passed = true;
227 
228  // Check lower-level stuff
229  if (FiniteElement::self_test() != 0)
230  {
231  passed = false;
232  }
233 
234  // Return verdict
235  if (passed)
236  {
237  return 0;
238  }
239  else
240  {
241  return 1;
242  }
243  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

template<unsigned DIM>
HelmholtzSourceFctPt& oomph::HelmholtzEquations< DIM >::source_fct_pt ( )
inline

Access function: Pointer to source function.

277  {
278  return Source_fct_pt;
279  }

References oomph::HelmholtzEquations< DIM >::Source_fct_pt.

Referenced by oomph::RefineableHelmholtzEquations< DIM >::further_build().

◆ source_fct_pt() [2/2]

template<unsigned DIM>
HelmholtzSourceFctPt oomph::HelmholtzEquations< DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

283  {
284  return Source_fct_pt;
285  }

References oomph::HelmholtzEquations< DIM >::Source_fct_pt.

◆ u_index_helmholtz()

Member Data Documentation

◆ K_squared_pt

◆ Source_fct_pt


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