oomph::PMLFourierDecomposedHelmholtzEquations Class Referenceabstract

#include <pml_fourier_decomposed_helmholtz_elements.h>

+ Inheritance diagram for oomph::PMLFourierDecomposedHelmholtzEquations:

Public Types

typedef void(* PMLFourierDecomposedHelmholtzSourceFctPt) (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

 PMLFourierDecomposedHelmholtzEquations ()
 Constructor. More...
 
 PMLFourierDecomposedHelmholtzEquations (const PMLFourierDecomposedHelmholtzEquations &dummy)=delete
 Broken copy constructor. More...
 
virtual std::complex< unsignedu_index_pml_fourier_decomposed_helmholtz () const
 Broken assignment operator. More...
 
double *& k_squared_pt ()
 Get pointer to frequency. More...
 
double k_squared ()
 Get k squared. More...
 
double *& alpha_pt ()
 Get pointer to complex shift. More...
 
double alpha ()
 Get complex shift. More...
 
int *& pml_fourier_wavenumber_pt ()
 Get pointer to Fourier wavenumber. More...
 
int pml_fourier_wavenumber ()
 Get the Fourier wavenumber. More...
 
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...
 
void compute_norm (double &norm)
 Compute norm of fe solution. More...
 
PMLFourierDecomposedHelmholtzSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
PMLFourierDecomposedHelmholtzSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
virtual void get_source_pml_fourier_decomposed_helmholtz (const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
 
void values_to_be_pinned_on_outer_pml_boundary (Vector< unsigned > &values_to_pin)
 
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_pml_fourier_decomposed_helmholtz (const Vector< double > &s) const
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
- Public Member Functions inherited from oomph::PMLElementBase< 2 >
 PMLElementBase ()
 Constructor. More...
 
virtual ~PMLElementBase ()
 Virtual destructor. More...
 
void disable_pml ()
 
void enable_pml (const int &direction, const double &interface_border_value, const double &outer_domain_border_value)
 
- 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 unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) 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 std::string scalar_name_paraview (const unsigned &i) 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 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

void compute_pml_coefficients (const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_laplace_factor, std::complex< double > &pml_k_squared_factor)
 
void compute_complex_r (const unsigned &ipt, const Vector< double > &x, std::complex< double > &complex_r)
 
PMLMappingAndTransformedCoordinate *& pml_mapping_and_transformed_coordinate_pt ()
 Return a pointer to the PML Mapping object. More...
 
PMLMappingAndTransformedCoordinate *const & pml_mapping_and_transformed_coordinate_pt () const
 Return a pointer to the PML Mapping object (const version) More...
 
virtual double dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_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

PMLFourierDecomposedHelmholtzSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
doubleK_squared_pt
 Pointer to k^2 (wavenumber squared) More...
 
PMLMappingAndTransformedCoordinatePml_mapping_and_transformed_coordinate_pt
 
doubleAlpha_pt
 Pointer to wavenumber complex shift. More...
 
intN_pml_fourier_pt
 Pointer to Fourier wave number. More...
 
- Protected Attributes inherited from oomph::PMLElementBase< 2 >
bool Pml_is_enabled
 Boolean indicating if element is used in pml mode. More...
 
std::vector< boolPml_direction_active
 
Vector< doublePml_inner_boundary
 
Vector< doublePml_outer_boundary
 
- 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
 

Static Protected Attributes

static BermudezPMLMappingAndTransformedCoordinate Default_pml_mapping_and_transformed_coordinate
 
static double Default_Physical_Constant_Value
 Static default value for the physical constants (initialised to zero) More...
 
- 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
 

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
 

Detailed Description

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// A class for all isoparametric elements that solve the Helmholtz equations with pml capabilities. in Fourier decomposed form (cylindrical polars):

\[ U(r,\varphi,z) = \Re( u^{(n)}(r,z) \exp(-i n \varphi)) \]

We are solving for \( u^{(n)}(r,z)\) for given parameters \( k^2 \) and \( n \) . This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class.

Member Typedef Documentation

◆ PMLFourierDecomposedHelmholtzSourceFctPt

typedef void(* oomph::PMLFourierDecomposedHelmholtzEquations::PMLFourierDecomposedHelmholtzSourceFctPt) (const Vector< double > &x, std::complex< double > &f)

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

Constructor & Destructor Documentation

◆ PMLFourierDecomposedHelmholtzEquations() [1/2]

oomph::PMLFourierDecomposedHelmholtzEquations::PMLFourierDecomposedHelmholtzEquations ( )
inline

Constructor.

172  {
173  // Provide pointer to static method (Save memory)
177 
179  }
double * K_squared_pt
Pointer to k^2 (wavenumber squared)
Definition: pml_fourier_decomposed_helmholtz_elements.h:653
static double Default_Physical_Constant_Value
Static default value for the physical constants (initialised to zero)
Definition: pml_fourier_decomposed_helmholtz_elements.h:664
double * Alpha_pt
Pointer to wavenumber complex shift.
Definition: pml_fourier_decomposed_helmholtz_elements.h:661
int * N_pml_fourier_pt
Pointer to Fourier wave number.
Definition: pml_fourier_decomposed_helmholtz_elements.h:667
PMLFourierDecomposedHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: pml_fourier_decomposed_helmholtz_elements.h:650
static BermudezPMLMappingAndTransformedCoordinate Default_pml_mapping_and_transformed_coordinate
Definition: pml_fourier_decomposed_helmholtz_elements.h:620
PMLMappingAndTransformedCoordinate * Pml_mapping_and_transformed_coordinate_pt
Definition: pml_fourier_decomposed_helmholtz_elements.h:658

References Alpha_pt, Default_Physical_Constant_Value, Default_pml_mapping_and_transformed_coordinate, and Pml_mapping_and_transformed_coordinate_pt.

◆ PMLFourierDecomposedHelmholtzEquations() [2/2]

oomph::PMLFourierDecomposedHelmholtzEquations::PMLFourierDecomposedHelmholtzEquations ( const PMLFourierDecomposedHelmholtzEquations dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ alpha()

double oomph::PMLFourierDecomposedHelmholtzEquations::alpha ( )
inline

Get complex shift.

235  {
236  return *Alpha_pt;
237  }

References Alpha_pt.

◆ alpha_pt()

double*& oomph::PMLFourierDecomposedHelmholtzEquations::alpha_pt ( )
inline

Get pointer to complex shift.

228  {
229  return Alpha_pt;
230  }

References Alpha_pt.

◆ compute_complex_r()

void oomph::PMLFourierDecomposedHelmholtzEquations::compute_complex_r ( const unsigned ipt,
const Vector< double > &  x,
std::complex< double > &  complex_r 
)
inlineprotected

Compute complex variable r at position x[0] and integration point ipt

The complex r variable is only imaginary on two conditions: First, the decaying nature of the pml layers is active. Secondly, the integration point is contained in the right pml layer or the two corner pml layers.

573  {
574  // Cache current position r
575  double r = x[0];
576 
582 
583  // If the complex r variable is imaginary
584  if (this->Pml_is_enabled && (this->Pml_direction_active[0]))
585  {
586  double nu = x[0] - Pml_inner_boundary[0];
587  double pml_width = Pml_outer_boundary[0] - Pml_inner_boundary[0];
588  double k_squared_local = k_squared();
589 
590  // Determine the complex r variable
591  complex_r =
593  nu, pml_width, Pml_inner_boundary[0], k_squared_local);
594  }
595  else
596  {
597  // The complex r variable is infact purely real, and
598  // is equal to x[0]
599  complex_r = std::complex<double>(r, 0.0);
600  }
601 
602  } // end of compute_complex_r
bool Pml_is_enabled
Boolean indicating if element is used in pml mode.
Definition: pml_meshes.h:119
std::vector< bool > Pml_direction_active
Definition: pml_meshes.h:124
Vector< double > Pml_outer_boundary
Definition: pml_meshes.h:134
Vector< double > Pml_inner_boundary
Definition: pml_meshes.h:129
double k_squared()
Get k squared.
Definition: pml_fourier_decomposed_helmholtz_elements.h:212
virtual std::complex< double > transformed_coordinate(const double &nu_i, const double &pml_width_i, const double &pml_inner_boundary, const double &k_squared)=0
r
Definition: UniformPSDSelfTest.py:20
list x
Definition: plotDoE.py:28

References k_squared(), oomph::PMLElementBase< 2 >::Pml_direction_active, oomph::PMLElementBase< 2 >::Pml_inner_boundary, oomph::PMLElementBase< 2 >::Pml_is_enabled, Pml_mapping_and_transformed_coordinate_pt, oomph::PMLElementBase< 2 >::Pml_outer_boundary, UniformPSDSelfTest::r, oomph::PMLMappingAndTransformedCoordinate::transformed_coordinate(), and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ compute_error() [1/2]

void oomph::PMLFourierDecomposedHelmholtzEquations::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.

720  {
721  // Initialise
722  error = 0.0;
723  norm = 0.0;
724 
725  // Vector of local coordinates
726  Vector<double> s(2);
727 
728  // Vector for coordintes
729  Vector<double> x(2);
730 
731  // Find out how many nodes there are in the element
732  unsigned n_node = nnode();
733 
734  Shape psi(n_node);
735 
736  // Set the value of n_intpt
737  unsigned n_intpt = integral_pt()->nweight();
738 
739  // Tecplot
740  outfile << "ZONE" << std::endl;
741 
742  // Exact solution Vector
743  Vector<double> exact_soln(2);
744 
745  // Loop over the integration points
746  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
747  {
748  // Assign values of s
749  for (unsigned i = 0; i < 2; i++)
750  {
751  s[i] = integral_pt()->knot(ipt, i);
752  }
753 
754  // Get the integral weight
755  double w = integral_pt()->weight(ipt);
756 
757  // Get jacobian of mapping
758  double J = J_eulerian(s);
759 
760  // Premultiply the weights and the Jacobian
761  double W = w * J;
762 
763  // Get x position as Vector
764  interpolated_x(s, x);
765 
766  // Get FE function value
767  std::complex<double> u_fe =
769 
770  // Get exact solution at this point
771  (*exact_soln_pt)(x, exact_soln);
772 
773  // Output r,z,error
774  for (unsigned i = 0; i < 2; i++)
775  {
776  outfile << x[i] << " ";
777  }
778  outfile << exact_soln[0] << " " << exact_soln[1] << " "
779  << exact_soln[0] - u_fe.real() << " "
780  << exact_soln[1] - u_fe.imag() << std::endl;
781 
782  // Add to error and norm
783  norm +=
784  (exact_soln[0] * exact_soln[0] + exact_soln[1] * exact_soln[1]) * W;
785  error += ((exact_soln[0] - u_fe.real()) * (exact_soln[0] - u_fe.real()) +
786  (exact_soln[1] - u_fe.imag()) * (exact_soln[1] - u_fe.imag())) *
787  W;
788  }
789  }
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
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.
std::complex< double > interpolated_u_pml_fourier_decomposed_helmholtz(const Vector< double > &s) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:458
RealScalar s
Definition: level1_cplx_impl.h:130
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

References calibrate::error, ProblemParameters::exact_soln(), i, oomph::FiniteElement::integral_pt(), interpolated_u_pml_fourier_decomposed_helmholtz(), oomph::FiniteElement::interpolated_x(), J, oomph::FiniteElement::J_eulerian(), oomph::Integral::knot(), oomph::FiniteElement::nnode(), oomph::Integral::nweight(), s, w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and plotDoE::x.

◆ compute_error() [2/2]

void oomph::PMLFourierDecomposedHelmholtzEquations::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.

335  {
336  throw OomphLibError("There is no time-dependent compute_error() for "
337  "PMLFourierDecomposedHelmholtz elements",
340  }
#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.

◆ compute_norm()

void oomph::PMLFourierDecomposedHelmholtzEquations::compute_norm ( double norm)
virtual

Compute norm of fe solution.

Reimplemented from oomph::GeneralisedElement.

796  {
797  // Initialise
798  norm = 0.0;
799 
800  // Vector of local coordinates
801  Vector<double> s(2);
802 
803  // Vector for coordintes
804  Vector<double> x(2);
805 
806  // Find out how many nodes there are in the element
807  unsigned n_node = nnode();
808 
809  Shape psi(n_node);
810 
811  // Set the value of n_intpt
812  unsigned n_intpt = integral_pt()->nweight();
813 
814  // Loop over the integration points
815  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
816  {
817  // Assign values of s
818  for (unsigned i = 0; i < 2; i++)
819  {
820  s[i] = integral_pt()->knot(ipt, i);
821  }
822 
823  // Get the integral weight
824  double w = integral_pt()->weight(ipt);
825 
826  // Get jacobian of mapping
827  double J = J_eulerian(s);
828 
829  // Premultiply the weights and the Jacobian
830  double W = w * J;
831 
832  // Get FE function value
833  std::complex<double> u_fe =
835 
836  // Add to norm
837  norm += (u_fe.real() * u_fe.real() + u_fe.imag() * u_fe.imag()) * W;
838  }
839  }

References i, oomph::FiniteElement::integral_pt(), interpolated_u_pml_fourier_decomposed_helmholtz(), J, oomph::FiniteElement::J_eulerian(), oomph::Integral::knot(), oomph::FiniteElement::nnode(), oomph::Integral::nweight(), s, w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and plotDoE::x.

◆ compute_pml_coefficients()

void oomph::PMLFourierDecomposedHelmholtzEquations::compute_pml_coefficients ( const unsigned ipt,
const Vector< double > &  x,
Vector< std::complex< double >> &  pml_laplace_factor,
std::complex< double > &  pml_k_squared_factor 
)
inlineprotected

Compute pml coefficients at position x and integration point ipt. pml_laplace_factor is used in the residual contribution from the laplace operator, similarly pml_k_squared_factor is used in the contribution from the k^2 of the Helmholtz operator.

Vector which points from the inner boundary to x

Vector which points from the inner boundary to the edge of the boundary

for 2D, in order: g_y/g_x, g_x/g_y for Laplace bit and g_x*g_y for Helmholtz bit for 3D, in order: g_y*g_x/g_x, g*x*g_z/g_y, g_x*g_y/g_z for Laplace bit and g_x*g_y*g_z for Helmholtz factor

The weights all default to 1.0 as if the propagation medium is the physical domain

507  {
509  Vector<double> nu(2);
510  for (unsigned k = 0; k < 2; k++)
511  {
512  nu[k] = x[k] - this->Pml_inner_boundary[k];
513  }
514 
517  Vector<double> pml_width(2);
518  for (unsigned k = 0; k < 2; k++)
519  {
520  pml_width[k] =
521  this->Pml_outer_boundary[k] - this->Pml_inner_boundary[k];
522  }
523 
524  // Declare gamma_i vectors of complex numbers for PML weights
525  Vector<std::complex<double>> pml_gamma(2);
526 
527  if (this->Pml_is_enabled)
528  {
529  // Cache k_squared to pass into mapping function
530  double k_squared_local = k_squared();
531 
532  for (unsigned k = 0; k < 2; k++)
533  {
534  // If PML is enabled in the respective direction
535  if (this->Pml_direction_active[k])
536  {
538  nu[k], pml_width[k], k_squared_local);
539  }
540  else
541  {
542  pml_gamma[k] = 1.0;
543  }
544  }
545 
550  pml_laplace_factor[0] = pml_gamma[1] / pml_gamma[0];
551  pml_laplace_factor[1] = pml_gamma[0] / pml_gamma[1];
552  pml_k_squared_factor = pml_gamma[0] * pml_gamma[1];
553  }
554  else
555  {
558  for (unsigned k = 0; k < 2; k++)
559  {
560  pml_laplace_factor[k] = std::complex<double>(1.0, 0.0);
561  }
562 
563  pml_k_squared_factor = std::complex<double>(1.0, 0.0);
564  }
565  }
virtual std::complex< double > gamma(const double &nu_i, const double &pml_width_i, const double &k_squared)=0
char char char int int * k
Definition: level2_impl.h:374

References oomph::PMLMappingAndTransformedCoordinate::gamma(), k, k_squared(), oomph::PMLElementBase< 2 >::Pml_direction_active, oomph::PMLElementBase< 2 >::Pml_inner_boundary, oomph::PMLElementBase< 2 >::Pml_is_enabled, Pml_mapping_and_transformed_coordinate_pt, oomph::PMLElementBase< 2 >::Pml_outer_boundary, and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz()

virtual double oomph::PMLFourierDecomposedHelmholtzEquations::dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_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::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz()

virtual double oomph::PMLFourierDecomposedHelmholtzEquations::dshape_and_dtest_eulerian_pml_fourier_decomposed_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::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

◆ fill_in_contribution_to_jacobian()

void oomph::PMLFourierDecomposedHelmholtzEquations::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.

449  {
450  // Call the generic routine with the flag set to 1
452  residuals, jacobian, 1);
453  }
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: pml_fourier_decomposed_helmholtz_elements.cc:192

References fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ fill_in_contribution_to_residuals()

void oomph::PMLFourierDecomposedHelmholtzEquations::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

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

Reimplemented from oomph::GeneralisedElement.

437  {
438  // Call the generic residuals function with flag set to 0
439  // using a dummy matrix argument
441  residuals, GeneralisedElement::Dummy_matrix, 0);
442  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz()

void oomph::PMLFourierDecomposedHelmholtzEquations::fill_in_generic_residual_contribution_pml_fourier_decomposed_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

196  {
197  // Find out how many nodes there are
198  const unsigned n_node = nnode();
199 
200  // Set up memory for the shape and test functions
201  Shape psi(n_node), test(n_node);
202  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
203 
204  // Set the value of n_intpt
205  const unsigned n_intpt = integral_pt()->nweight();
206 
207  // Integers to store the local equation and unknown numbers
208  int local_eqn_real = 0, local_unknown_real = 0;
209  int local_eqn_imag = 0, local_unknown_imag = 0;
210 
211  // Loop over the integration points
212  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
213  {
214  // Get the integral weight
215  double w = integral_pt()->weight(ipt);
216 
217  // Call the derivatives of the shape and test functions
218  double J =
220  ipt, psi, dpsidx, test, dtestdx);
221 
222  // Premultiply the weights and the Jacobian
223  double W = w * J;
224 
225  // Calculate local values of unknown
226  // Allocate and initialise to zero
227  std::complex<double> interpolated_u(0.0, 0.0);
228  Vector<double> interpolated_x(2, 0.0);
229  Vector<std::complex<double>> interpolated_dudx(2);
230 
231  // Calculate function value and derivatives:
232  //-----------------------------------------
233  // Loop over nodes
234  for (unsigned l = 0; l < n_node; l++)
235  {
236  // Loop over directions
237  for (unsigned j = 0; j < 2; j++)
238  {
239  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
240  }
241 
242  // Get the nodal value of the helmholtz unknown
243  const std::complex<double> u_value(
245  raw_nodal_value(l,
247 
248  // Add to the interpolated value
249  interpolated_u += u_value * psi(l);
250 
251  // Loop over directions
252  for (unsigned j = 0; j < 2; j++)
253  {
254  interpolated_dudx[j] += u_value * dpsidx(l, j);
255  }
256  }
257 
258  // Get source function
259  //-------------------
260  std::complex<double> source(0.0, 0.0);
262 
263  double n = (double)pml_fourier_wavenumber();
264  double n_squared = n * n;
265 
266 
267  // Declare a vector of complex numbers for pml weights on the Laplace bit
268  Vector<std::complex<double>> pml_laplace_factor(2);
269  // Declare a complex number for pml weights on the mass matrix bit
270  std::complex<double> pml_k_squared_factor =
271  std::complex<double>(1.0, 0.0);
272 
273  // All the PML weights that participate in the assemby process
274  // are computed here. pml_laplace_factor will contain the entries
275  // for the Laplace bit, while pml_k_squared_factor contains the
276  // contributions to the Helmholtz bit. Both default to 1.0, should the PML
277  // not be enabled via enable_pml.
279  ipt, interpolated_x, pml_laplace_factor, pml_k_squared_factor);
280 
281 
282  // Determine the complex r variable. The variable is
283  // only complex once it enters the right pml domain or either
284  // of the two corner pml domains, otherwise it acts like the
285  // variable r.
286  std::complex<double> complex_r = std::complex<double>(1.0, 0.0);
287  compute_complex_r(ipt, interpolated_x, complex_r);
288 
289  // Calculate Jacobian
290  // ------------------
291  std::complex<double> pml_k_squared_jacobian =
292  -pml_k_squared_factor *
293  (k_squared() - n_squared / (complex_r * complex_r)) * complex_r * W;
294 
295  // Term from the Laplace operator
296  Vector<std::complex<double>> pml_laplace_jacobian(2);
297  for (unsigned k = 0; k < 2; k++)
298  {
299  pml_laplace_jacobian[k] = pml_laplace_factor[k] * complex_r * W;
300  }
301 
302  // Calculate residuals
303  // -------------------
304  // Note: it is a linear equation so residual = jacobian * u
305  std::complex<double> pml_k_squared_residual =
306  pml_k_squared_jacobian * interpolated_u;
307 
308  // Term from the Laplace operator
309  Vector<std::complex<double>> pml_laplace_residual(2);
310  for (unsigned k = 0; k < 2; k++)
311  {
312  pml_laplace_residual[k] =
313  pml_laplace_jacobian[k] * interpolated_dudx[k];
314  }
315 
316  // Assemble residuals and Jacobian
317  //--------------------------------
318  // Loop over the test functions
319  for (unsigned l = 0; l < n_node; l++)
320  {
321  // Get the equation numbers
322  local_eqn_real =
324  local_eqn_imag =
326 
327  // first, add the real contributions
328  //-------------------------------------------
329 
330  /*IF it's not a boundary condition*/
331  if (local_eqn_real >= 0)
332  {
333  // Add k squared term of equation
334  residuals[local_eqn_real] += pml_k_squared_residual.real() * test(l);
335 
336  // Add the term from the Laplace operator
337  for (unsigned k = 0; k < 2; k++)
338  {
339  residuals[local_eqn_real] +=
340  pml_laplace_residual[k].real() * dtestdx(l, k);
341  }
342 
343  // Add contributions to the jacobian
344  //----------------------------------
345  if (flag)
346  {
347  // Loop over the velocity shape functions again
348  for (unsigned l2 = 0; l2 < n_node; l2++)
349  {
350  // Get the unknown numbers
351  local_unknown_real = nodal_local_eqn(
353  local_unknown_imag = nodal_local_eqn(
355 
356  // If at a non-zero degree of freedom add in the entry
357  if (local_unknown_real >= 0)
358  {
359  // Add the helmholtz contribution to the jacobian
360  jacobian(local_eqn_real, local_unknown_real) +=
361  pml_k_squared_jacobian.real() * psi(l2) * test(l);
362 
363  // Add the term from the Laplace operator to the jacobian
364  for (unsigned k = 0; k < 2; k++)
365  {
366  jacobian(local_eqn_real, local_unknown_real) +=
367  pml_laplace_jacobian[k].real() * dpsidx(l2, k) *
368  dtestdx(l, k);
369  }
370  }
371  // If at a non-zero degree of freedom add in the entry
372  if (local_unknown_imag >= 0)
373  {
374  // Add k squared term to jacobian
375  jacobian(local_eqn_real, local_unknown_imag) +=
376  -pml_k_squared_jacobian.imag() * psi(l2) * test(l);
377 
378  // Add contribution to elemental Matrix
379  for (unsigned k = 0; k < 2; k++)
380  {
381  jacobian(local_eqn_real, local_unknown_imag) +=
382  -pml_laplace_jacobian[k].imag() * dpsidx(l2, k) *
383  dtestdx(l, k);
384  }
385  }
386 
387  } // end of loop over velocity shape functions
388  } // end of if(flag)
389  }
390 
391  // Second, add the imaginary contribution
392  //------------------------------------------------
393 
394  // IF it's not a boundary condition
395  if (local_eqn_imag >= 0)
396  {
397  // Add k squared term of equation
398  residuals[local_eqn_imag] += pml_k_squared_residual.imag() * test(l);
399 
400  // Add the term from the Laplace operator
401  for (unsigned k = 0; k < 2; k++)
402  {
403  residuals[local_eqn_imag] +=
404  pml_laplace_residual[k].imag() * dtestdx(l, k);
405  }
406 
407  // Add the contribution to the jacobian
408  //-----------------------
409  if (flag)
410  {
411  // Loop over the velocity shape functions again
412  for (unsigned l2 = 0; l2 < n_node; l2++)
413  {
414  local_unknown_imag = nodal_local_eqn(
416  local_unknown_real = nodal_local_eqn(
418 
419  // If at a non-zero degree of freedom add in the entry
420  if (local_unknown_imag >= 0)
421  {
422  // Add the helmholtz contribution
423  jacobian(local_eqn_imag, local_unknown_imag) +=
424  pml_k_squared_jacobian.real() * psi(l2) * test(l);
425 
426  // Add contribution to elemental Matrix
427  for (unsigned k = 0; k < 2; k++)
428  {
429  jacobian(local_eqn_imag, local_unknown_imag) +=
430  pml_laplace_jacobian[k].real() * dpsidx(l2, k) *
431  dtestdx(l, k);
432  }
433  }
434  // If at a non-zero degree of freedom add in the entry
435  if (local_unknown_real >= 0)
436  {
437  // Add the helmholtz contribution
438  jacobian(local_eqn_imag, local_unknown_real) +=
439  pml_k_squared_jacobian.imag() * psi(l2) * test(l);
440 
441  // Add contribution to elemental Matrix
442  for (unsigned k = 0; k < 2; k++)
443  {
444  jacobian(local_eqn_imag, local_unknown_real) +=
445  pml_laplace_jacobian[k].imag() * dpsidx(l2, k) *
446  dtestdx(l, k);
447  }
448  }
449  }
450  }
451  }
452  }
453  } // End of loop over integration points
454  }
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
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
int pml_fourier_wavenumber()
Get the Fourier wavenumber.
Definition: pml_fourier_decomposed_helmholtz_elements.h:246
void compute_complex_r(const unsigned &ipt, const Vector< double > &x, std::complex< double > &complex_r)
Definition: pml_fourier_decomposed_helmholtz_elements.h:570
virtual void get_source_pml_fourier_decomposed_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:364
void compute_pml_coefficients(const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_laplace_factor, std::complex< double > &pml_k_squared_factor)
Definition: pml_fourier_decomposed_helmholtz_elements.h:502
virtual double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
Definition: pml_fourier_decomposed_helmholtz_elements.h:197
float real
Definition: datatypes.h:10
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 compute_complex_r(), compute_pml_coefficients(), dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(), get_source_pml_fourier_decomposed_helmholtz(), imag(), oomph::FiniteElement::integral_pt(), oomph::FiniteElement::interpolated_x(), J, j, k, k_squared(), n, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::Integral::nweight(), pml_fourier_wavenumber(), oomph::FiniteElement::raw_nodal_position(), oomph::FiniteElement::raw_nodal_value(), TestProblem::source(), Eigen::test, u_index_pml_fourier_decomposed_helmholtz(), w, oomph::QuadTreeNames::W, and oomph::Integral::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().

◆ get_flux()

void oomph::PMLFourierDecomposedHelmholtzEquations::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.

398  {
399  // Find out how many nodes there are in the element
400  const unsigned n_node = nnode();
401 
402  // Set up memory for the shape and test functions
403  Shape psi(n_node);
404  DShape dpsidx(n_node, 2);
405 
406  // Call the derivatives of the shape and test functions
407  dshape_eulerian(s, psi, dpsidx);
408 
409  // Initialise to zero
410  const std::complex<double> zero(0.0, 0.0);
411  for (unsigned j = 0; j < 2; j++)
412  {
413  flux[j] = zero;
414  }
415 
416  // Loop over nodes
417  for (unsigned l = 0; l < n_node; l++)
418  {
419  // Cache the complex value of the unknown
420  const std::complex<double> u_value(
421  this->nodal_value(l,
423  this->nodal_value(l,
425 
426  // Loop over derivative directions
427  for (unsigned j = 0; j < 2; j++)
428  {
429  flux[j] += u_value * dpsidx(l, j);
430  }
431  }
432  }
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 oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), imag(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, u_index_pml_fourier_decomposed_helmholtz(), and zero().

Referenced by oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >::get_Z2_flux().

◆ get_source_pml_fourier_decomposed_helmholtz()

virtual void oomph::PMLFourierDecomposedHelmholtzEquations::get_source_pml_fourier_decomposed_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.

368  {
369  // If no source function has been set, return zero
370  if (Source_fct_pt == 0)
371  {
372  source = std::complex<double>(0.0, 0.0);
373  }
374  else
375  {
376  // Get source strength
377  (*Source_fct_pt)(x, source);
378  }
379  }

References TestProblem::source(), Source_fct_pt, and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ interpolated_u_pml_fourier_decomposed_helmholtz()

std::complex<double> oomph::PMLFourierDecomposedHelmholtzEquations::interpolated_u_pml_fourier_decomposed_helmholtz ( const Vector< double > &  s) const
inline

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

460  {
461  // Find number of nodes
462  const unsigned n_node = nnode();
463 
464  // Local shape function
465  Shape psi(n_node);
466 
467  // Find values of shape function
468  shape(s, psi);
469 
470  // Initialise value of u
471  std::complex<double> interpolated_u(0.0, 0.0);
472 
473  // Get the index at which the helmholtz unknown is stored
474  const unsigned u_nodal_index_real =
476  const unsigned u_nodal_index_imag =
478 
479  // Loop over the local nodes and sum
480  for (unsigned l = 0; l < n_node; l++)
481  {
482  // Make a temporary complex number from the stored data
483  const std::complex<double> u_value(
484  this->nodal_value(l, u_nodal_index_real),
485  this->nodal_value(l, u_nodal_index_imag));
486  // Add to the interpolated value
487  interpolated_u += u_value * psi[l];
488  }
489  return interpolated_u;
490  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::FiniteElement::shape(), and u_index_pml_fourier_decomposed_helmholtz().

Referenced by compute_error(), compute_norm(), output(), and output_real().

◆ k_squared()

double oomph::PMLFourierDecomposedHelmholtzEquations::k_squared ( )
inline

Get k squared.

213  {
214 #ifdef PARANOID
215  if (K_squared_pt == 0)
216  {
217  throw OomphLibError(
218  "Please set pointer to k_squared using access fct to pointer!",
221  }
222 #endif
223  return *K_squared_pt;
224  }

References K_squared_pt, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by compute_complex_r(), compute_pml_coefficients(), and fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ k_squared_pt()

double*& oomph::PMLFourierDecomposedHelmholtzEquations::k_squared_pt ( )
inline

Get pointer to frequency.

206  {
207  return K_squared_pt;
208  }

References K_squared_pt.

Referenced by PMLFourierDecomposedHelmholtzProblem< ELEMENT >::complete_problem_setup().

◆ output() [1/4]

void oomph::PMLFourierDecomposedHelmholtzEquations::output ( FILE *  file_pt)
inlinevirtual

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

281  {
282  const unsigned n_plot = 5;
283  output(file_pt, n_plot);
284  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: pml_fourier_decomposed_helmholtz_elements.h:260

References output().

◆ output() [2/4]

void oomph::PMLFourierDecomposedHelmholtzEquations::output ( FILE *  file_pt,
const unsigned nplot 
)
virtual

C-style output FE representation of soln: r,z,u_re,u_im or at n_plot^2 plot points

C-style output function:

r,z,u

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

567  {
568  // Vector of local coordinates
569  Vector<double> s(2);
570 
571  // Tecplot header info
572  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
573 
574  // Loop over plot points
575  unsigned num_plot_points = nplot_points(nplot);
576  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
577  {
578  // Get local coordinates of plot point
579  get_s_plot(iplot, nplot, s);
580  std::complex<double> u(
582 
583  for (unsigned i = 0; i < 2; i++)
584  {
585  fprintf(file_pt, "%g ", interpolated_x(s, i));
586  }
587 
588  for (unsigned i = 0; i < 2; i++)
589  {
590  fprintf(file_pt, "%g ", interpolated_x(s, i));
591  }
592  fprintf(file_pt, "%g ", u.real());
593  fprintf(file_pt, "%g \n", u.imag());
594  }
595 
596  // Write tecplot footer (e.g. FE connectivity lists)
597  write_tecplot_zone_footer(file_pt, nplot);
598  }
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 oomph::FiniteElement::get_s_plot(), i, interpolated_u_pml_fourier_decomposed_helmholtz(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), and oomph::FiniteElement::write_tecplot_zone_footer().

◆ output() [3/4]

void oomph::PMLFourierDecomposedHelmholtzEquations::output ( std::ostream &  outfile)
inlinevirtual

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

261  {
262  const unsigned n_plot = 5;
263  output(outfile, n_plot);
264  }

Referenced by output(), oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output(), and oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output().

◆ output() [4/4]

void oomph::PMLFourierDecomposedHelmholtzEquations::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^2 plot points

Output function:

r,z,u_re,u_imag

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

491  {
492  // Vector of local coordinates
493  Vector<double> s(2);
494 
495  // Tecplot header info
496  outfile << tecplot_zone_string(nplot);
497 
498  // Loop over plot points
499  unsigned num_plot_points = nplot_points(nplot);
500  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
501  {
502  // Get local coordinates of plot point
503  get_s_plot(iplot, nplot, s);
504  std::complex<double> u(
506  for (unsigned i = 0; i < 2; i++)
507  {
508  outfile << interpolated_x(s, i) << " ";
509  }
510  outfile << u.real() << " " << u.imag() << std::endl;
511  }
512 
513  // Write tecplot footer (e.g. FE connectivity lists)
514  write_tecplot_zone_footer(outfile, nplot);
515  }

References oomph::FiniteElement::get_s_plot(), i, interpolated_u_pml_fourier_decomposed_helmholtz(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), and oomph::FiniteElement::write_tecplot_zone_footer().

◆ output_fct() [1/2]

virtual void oomph::PMLFourierDecomposedHelmholtzEquations::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::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

303  {
304  throw OomphLibError("There is no time-dependent output_fct() for "
305  "PMLFourierDecomposedHelmholtz elements ",
308  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ output_fct() [2/2]

void oomph::PMLFourierDecomposedHelmholtzEquations::output_fct ( std::ostream &  outfile,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
virtual

Output exact soln: r,z,u_re_exact,u_im_exact at n_plot^2 plot points

Output exact solution

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

r,z,u_exact

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >, and oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >.

613  {
614  // Vector of local coordinates
615  Vector<double> s(2);
616 
617  // Vector for coordintes
618  Vector<double> x(2);
619 
620  // Tecplot header info
621  outfile << tecplot_zone_string(nplot);
622 
623  // Exact solution Vector
624  Vector<double> exact_soln(2);
625 
626  // Loop over plot points
627  unsigned num_plot_points = nplot_points(nplot);
628  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
629  {
630  // Get local coordinates of plot point
631  get_s_plot(iplot, nplot, s);
632 
633  // Get x position as Vector
634  interpolated_x(s, x);
635 
636  // Get exact solution at this point
637  (*exact_soln_pt)(x, exact_soln);
638 
639  // Output r,z,u_exact
640  for (unsigned i = 0; i < 2; i++)
641  {
642  outfile << x[i] << " ";
643  }
644  outfile << exact_soln[0] << " " << exact_soln[1] << std::endl;
645  }
646 
647  // Write tecplot footer (e.g. FE connectivity lists)
648  write_tecplot_zone_footer(outfile, nplot);
649  }

References ProblemParameters::exact_soln(), oomph::FiniteElement::get_s_plot(), i, oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

Referenced by oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output_fct(), and oomph::TPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output_fct().

◆ output_real()

void oomph::PMLFourierDecomposedHelmholtzEquations::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. r,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.

r,z,u

Output at nplot points in each coordinate direction

531  {
532  // Vector of local coordinates
533  Vector<double> s(2);
534 
535  // Tecplot header info
536  outfile << tecplot_zone_string(nplot);
537 
538  // Loop over plot points
539  unsigned num_plot_points = nplot_points(nplot);
540  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
541  {
542  // Get local coordinates of plot point
543  get_s_plot(iplot, nplot, s);
544  std::complex<double> u(
546  for (unsigned i = 0; i < 2; i++)
547  {
548  outfile << interpolated_x(s, i) << " ";
549  }
550  outfile << u.real() * cos(phi) + u.imag() * sin(phi) << std::endl;
551  }
552 
553  // Write tecplot footer (e.g. FE connectivity lists)
554  write_tecplot_zone_footer(outfile, nplot);
555  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137

References cos(), oomph::FiniteElement::get_s_plot(), i, interpolated_u_pml_fourier_decomposed_helmholtz(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, sin(), oomph::FiniteElement::tecplot_zone_string(), and oomph::FiniteElement::write_tecplot_zone_footer().

Referenced by oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output_real().

◆ output_real_fct()

void oomph::PMLFourierDecomposedHelmholtzEquations::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. r,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.

r,z,u

Output at nplot points in each coordinate direction

668  {
669  // Vector of local coordinates
670  Vector<double> s(2);
671 
672  // Vector for coordintes
673  Vector<double> x(2);
674 
675  // Tecplot header info
676  outfile << tecplot_zone_string(nplot);
677 
678  // Exact solution Vector
679  Vector<double> exact_soln(2);
680 
681  // Loop over plot points
682  unsigned num_plot_points = nplot_points(nplot);
683  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
684  {
685  // Get local coordinates of plot point
686  get_s_plot(iplot, nplot, s);
687 
688  // Get x position as Vector
689  interpolated_x(s, x);
690 
691  // Get exact solution at this point
692  (*exact_soln_pt)(x, exact_soln);
693 
694  // Output x,y,...,u_exact
695  for (unsigned i = 0; i < 2; i++)
696  {
697  outfile << x[i] << " ";
698  }
699  outfile << exact_soln[0] * cos(phi) + exact_soln[1] * sin(phi)
700  << std::endl;
701  }
702 
703  // Write tecplot footer (e.g. FE connectivity lists)
704  write_tecplot_zone_footer(outfile, nplot);
705  }

References cos(), ProblemParameters::exact_soln(), oomph::FiniteElement::get_s_plot(), i, oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, sin(), oomph::FiniteElement::tecplot_zone_string(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

Referenced by oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >::output_real_fct().

◆ pml_fourier_wavenumber()

int oomph::PMLFourierDecomposedHelmholtzEquations::pml_fourier_wavenumber ( )
inline

Get the Fourier wavenumber.

247  {
248  if (N_pml_fourier_pt == 0)
249  {
250  return 0;
251  }
252  else
253  {
254  return *N_pml_fourier_pt;
255  }
256  }

References N_pml_fourier_pt.

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz().

◆ pml_fourier_wavenumber_pt()

int*& oomph::PMLFourierDecomposedHelmholtzEquations::pml_fourier_wavenumber_pt ( )
inline

Get pointer to Fourier wavenumber.

241  {
242  return N_pml_fourier_pt;
243  }

References N_pml_fourier_pt.

Referenced by PMLFourierDecomposedHelmholtzProblem< ELEMENT >::complete_problem_setup().

◆ pml_mapping_and_transformed_coordinate_pt() [1/2]

PMLMappingAndTransformedCoordinate*& oomph::PMLFourierDecomposedHelmholtzEquations::pml_mapping_and_transformed_coordinate_pt ( )
inlineprotected

Return a pointer to the PML Mapping object.

606  {
608  }

References Pml_mapping_and_transformed_coordinate_pt.

◆ pml_mapping_and_transformed_coordinate_pt() [2/2]

PMLMappingAndTransformedCoordinate* const& oomph::PMLFourierDecomposedHelmholtzEquations::pml_mapping_and_transformed_coordinate_pt ( ) const
inlineprotected

Return a pointer to the PML Mapping object (const version)

613  {
615  }

References Pml_mapping_and_transformed_coordinate_pt.

◆ self_test()

unsigned oomph::PMLFourierDecomposedHelmholtzEquations::self_test ( )
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

461  {
462  bool passed = true;
463 
464  // Check lower-level stuff
465  if (FiniteElement::self_test() != 0)
466  {
467  passed = false;
468  }
469 
470  // Return verdict
471  if (passed)
472  {
473  return 0;
474  }
475  else
476  {
477  return 1;
478  }
479  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

PMLFourierDecomposedHelmholtzSourceFctPt& oomph::PMLFourierDecomposedHelmholtzEquations::source_fct_pt ( )
inline

Access function: Pointer to source function.

348  {
349  return Source_fct_pt;
350  }

References Source_fct_pt.

◆ source_fct_pt() [2/2]

PMLFourierDecomposedHelmholtzSourceFctPt oomph::PMLFourierDecomposedHelmholtzEquations::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

355  {
356  return Source_fct_pt;
357  }

References Source_fct_pt.

◆ u_index_pml_fourier_decomposed_helmholtz()

virtual std::complex<unsigned> oomph::PMLFourierDecomposedHelmholtzEquations::u_index_pml_fourier_decomposed_helmholtz ( ) const
inlinevirtual

Broken assignment operator.

Return the index at which the unknown value is stored: Real/imag part of index contains (real) index of real/imag part.

199  {
200  return std::complex<unsigned>(0, 1);
201  }

Referenced by fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(), get_flux(), interpolated_u_pml_fourier_decomposed_helmholtz(), oomph::PMLFourierDecomposedHelmholtzFluxElement< ELEMENT >::PMLFourierDecomposedHelmholtzFluxElement(), and oomph::PMLFourierDecomposedHelmholtzPowerMonitorElement< ELEMENT >::PMLFourierDecomposedHelmholtzPowerMonitorElement().

◆ values_to_be_pinned_on_outer_pml_boundary()

void oomph::PMLFourierDecomposedHelmholtzEquations::values_to_be_pinned_on_outer_pml_boundary ( Vector< unsigned > &  values_to_pin)
inlinevirtual

Pure virtual function in which we specify the values to be pinned (and set to zero) on the outer edge of the pml layer. All of them! Vector is resized internally.

Implements oomph::PMLElementBase< 2 >.

386  {
387  values_to_pin.resize(2);
388  for (unsigned j = 0; j < 2; j++)
389  {
390  values_to_pin[j] = j;
391  }
392  }

References j.

Member Data Documentation

◆ Alpha_pt

double* oomph::PMLFourierDecomposedHelmholtzEquations::Alpha_pt
protected

Pointer to wavenumber complex shift.

Referenced by alpha(), alpha_pt(), and PMLFourierDecomposedHelmholtzEquations().

◆ Default_Physical_Constant_Value

double oomph::PMLFourierDecomposedHelmholtzEquations::Default_Physical_Constant_Value
staticprotected
Initial value:
=
0.0

Static default value for the physical constants (initialised to zero)

PML Helmholtz equations static data, so that by default we can point to a 0

Referenced by PMLFourierDecomposedHelmholtzEquations().

◆ Default_pml_mapping_and_transformed_coordinate

BermudezPMLMappingAndTransformedCoordinate oomph::PMLFourierDecomposedHelmholtzEquations::Default_pml_mapping_and_transformed_coordinate
staticprotected

Static so that the class doesn't need to instantiate a new default everytime it uses it

Referenced by PMLFourierDecomposedHelmholtzEquations().

◆ K_squared_pt

double* oomph::PMLFourierDecomposedHelmholtzEquations::K_squared_pt
protected

Pointer to k^2 (wavenumber squared)

Referenced by k_squared(), and k_squared_pt().

◆ N_pml_fourier_pt

int* oomph::PMLFourierDecomposedHelmholtzEquations::N_pml_fourier_pt
protected

Pointer to Fourier wave number.

Referenced by pml_fourier_wavenumber(), and pml_fourier_wavenumber_pt().

◆ Pml_mapping_and_transformed_coordinate_pt

PMLMappingAndTransformedCoordinate* oomph::PMLFourierDecomposedHelmholtzEquations::Pml_mapping_and_transformed_coordinate_pt
protected

Pointer to class which holds the pml mapping function (also known as gamma) and the associated transformed coordinate

Referenced by compute_complex_r(), compute_pml_coefficients(), pml_mapping_and_transformed_coordinate_pt(), and PMLFourierDecomposedHelmholtzEquations().

◆ Source_fct_pt

PMLFourierDecomposedHelmholtzSourceFctPt oomph::PMLFourierDecomposedHelmholtzEquations::Source_fct_pt
protected

Pointer to source function:

Referenced by get_source_pml_fourier_decomposed_helmholtz(), and source_fct_pt().


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