oomph::PMLHelmholtzPowerElement< ELEMENT > Class Template Reference

#include <pml_helmholtz_flux_elements.h>

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

Public Member Functions

 PMLHelmholtzPowerElement (FiniteElement *const &bulk_el_pt, const int &face_index)
 
 PMLHelmholtzPowerElement ()
 Broken empty constructor. More...
 
 PMLHelmholtzPowerElement (const PMLHelmholtzPowerElement &dummy)=delete
 Broken copy constructor. More...
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 Broken assignment operator. More...
 
virtual std::complex< unsignedu_index_helmholtz () const
 
double global_power_contribution ()
 
double global_power_contribution (std::ofstream &outfile)
 
std::complex< doubleglobal_flux_contribution ()
 
std::complex< doubleglobal_flux_contribution (std::ofstream &outfile)
 
- Public Member Functions inherited from oomph::FaceElement
 FaceElement ()
 Constructor: Initialise all appropriate member data. More...
 
virtual ~FaceElement ()
 Empty virtual destructor. More...
 
 FaceElement (const FaceElement &)=delete
 Broken copy constructor. More...
 
const unsignedboundary_number_in_bulk_mesh () const
 Broken assignment operator. More...
 
void set_boundary_number_in_bulk_mesh (const unsigned &b)
 Set function for the boundary number in bulk mesh. More...
 
double J_eulerian (const Vector< double > &s) const
 
double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
double interpolated_x (const Vector< double > &s, const unsigned &i) const
 
double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 
void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
intnormal_sign ()
 
int normal_sign () const
 
intface_index ()
 
int face_index () const
 
const Vector< double > * tangent_direction_pt () const
 Public access function for the tangent direction pointer. More...
 
void set_tangent_direction (Vector< double > *tangent_direction_pt)
 Set the tangent direction vector. More...
 
void turn_on_warning_for_discontinuous_tangent ()
 
void turn_off_warning_for_discontinuous_tangent ()
 
void continuous_tangent_and_outer_unit_normal (const Vector< double > &s, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void continuous_tangent_and_outer_unit_normal (const unsigned &ipt, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void outer_unit_normal (const Vector< double > &s, Vector< double > &unit_normal) const
 Compute outer unit normal at the specified local coordinate. More...
 
void outer_unit_normal (const unsigned &ipt, Vector< double > &unit_normal) const
 Compute outer unit normal at ipt-th integration point. More...
 
FiniteElement *& bulk_element_pt ()
 Pointer to higher-dimensional "bulk" element. More...
 
FiniteElementbulk_element_pt () const
 Pointer to higher-dimensional "bulk" element (const version) More...
 
CoordinateMappingFctPtface_to_bulk_coordinate_fct_pt ()
 
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt () const
 
BulkCoordinateDerivativesFctPtbulk_coordinate_derivatives_fct_pt ()
 
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt () const
 
Vector< doublelocal_coordinate_in_bulk (const Vector< double > &s) const
 
void get_local_coordinate_in_bulk (const Vector< double > &s, Vector< double > &s_bulk) const
 
void get_ds_bulk_ds_face (const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction) const
 
unsignedbulk_position_type (const unsigned &i)
 
const unsignedbulk_position_type (const unsigned &i) const
 
void bulk_node_number_resize (const unsigned &i)
 Resize the storage for the bulk node numbers. More...
 
unsignedbulk_node_number (const unsigned &n)
 
const unsignedbulk_node_number (const unsigned &n) const
 
void bulk_position_type_resize (const unsigned &i)
 Resize the storage for bulk_position_type to i entries. More...
 
unsignednbulk_value (const unsigned &n)
 
unsigned nbulk_value (const unsigned &n) const
 
void nbulk_value_resize (const unsigned &i)
 
void resize_nodes (Vector< unsigned > &nadditional_data_values)
 
void output_zeta (std::ostream &outfile, const unsigned &nplot)
 Output boundary coordinate zeta. 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
 
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...
 
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)
 
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 (std::ostream &outfile)
 
virtual void output (std::ostream &outfile, const unsigned &n_plot)
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output (FILE *file_pt)
 
virtual void output (FILE *file_pt, const unsigned &n_plot)
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output an exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 Output a time-dependent exact solution over the element. More...
 
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, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, 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 self_test ()
 
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 Attributes

std::complex< unsignedU_index_helmholtz
 
unsigned Dim
 The spatial dimension of the problem. More...
 
- Protected Attributes inherited from oomph::FaceElement
unsigned Boundary_number_in_bulk_mesh
 The boundary number in the bulk mesh to which this element is attached. More...
 
FiniteElementBulk_element_pt
 Pointer to the associated higher-dimensional "bulk" element. More...
 
Vector< unsignedBulk_node_number
 
Vector< unsignedNbulk_value
 
Vector< double > * Tangent_direction_pt
 
- 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

- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 
- 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
 
- Protected Member Functions inherited from oomph::FaceElement
void add_additional_values (const Vector< unsigned > &nadditional_values, const unsigned &id)
 
- 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)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
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)
 
virtual void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
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)
 
- 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<class ELEMENT>
class oomph::PMLHelmholtzPowerElement< ELEMENT >

A class for elements that allow the post-processing of radiated power and flux on the boundaries of PMLHelmholtz elements. The element geometry is obtained from the FaceGeometry<ELEMENT> policy class.

Constructor & Destructor Documentation

◆ PMLHelmholtzPowerElement() [1/3]

template<class ELEMENT >
oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement ( FiniteElement *const &  bulk_el_pt,
const int face_index 
)

Constructor, takes the pointer to the "bulk" element and the index of the face to which the element is attached.

/////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// Constructor, takes the pointer to the "bulk" element, the index of the fixed local coordinate and its value represented by an integer (+/- 1), indicating that the face is located at the max. or min. value of the "fixed" local coordinate in the bulk element.

396  : FaceGeometry<ELEMENT>(), FaceElement()
397  {
398 #ifdef PARANOID
399  {
400  // Check that the element is not a refineable 3d element
401  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(bulk_el_pt);
402  // If it's three-d
403  if (elem_pt->dim() == 3)
404  {
405  // Is it refineable
406  RefineableElement* ref_el_pt =
407  dynamic_cast<RefineableElement*>(elem_pt);
408  if (ref_el_pt != 0)
409  {
410  if (this->has_hanging_nodes())
411  {
412  throw OomphLibError("This flux element will not work correctly if "
413  "nodes are hanging\n",
416  }
417  }
418  }
419  }
420 #endif
421 
422  // Let the bulk element build the FaceElement, i.e. setup the pointers
423  // to its nodes (by referring to the appropriate nodes in the bulk
424  // element), etc.
425  bulk_el_pt->build_face_element(face_index, this);
426 
427  // Extract the dimension of the problem from the dimension of
428  // the first node
429  Dim = this->node_pt(0)->ndim();
430 
431  // Set up U_index_helmholtz. Initialise to zero, which probably won't change
432  // in most cases, oh well, the price we pay for generality
433  U_index_helmholtz = std::complex<unsigned>(0, 1);
434 
435  // Cast to the appropriate PMLHelmholtzEquation so that we can
436  // find the index at which the variable is stored
437  // We assume that the dimension of the full problem is the same
438  // as the dimension of the node, if this is not the case you will have
439  // to write custom elements, sorry
440  switch (Dim)
441  {
442  // One dimensional problem
443  case 1:
444  {
445  PMLHelmholtzEquations<1>* eqn_pt =
446  dynamic_cast<PMLHelmholtzEquations<1>*>(bulk_el_pt);
447  // If the cast has failed die
448  if (eqn_pt == 0)
449  {
450  std::string error_string =
451  "Bulk element must inherit from PMLHelmholtzEquations.";
452  error_string +=
453  "Nodes are one dimensional, but cannot cast the bulk element to\n";
454  error_string += "PMLHelmholtzEquations<1>\n.";
455  error_string += "If you desire this functionality, you must "
456  "implement it yourself\n";
457 
458  throw OomphLibError(
460  }
461  // Otherwise read out the value
462  else
463  {
464  // Read the index from the (cast) bulk element
465  U_index_helmholtz = eqn_pt->u_index_helmholtz();
466  }
467  }
468  break;
469 
470  // Two dimensional problem
471  case 2:
472  {
473  PMLHelmholtzEquations<2>* eqn_pt =
474  dynamic_cast<PMLHelmholtzEquations<2>*>(bulk_el_pt);
475  // If the cast has failed die
476  if (eqn_pt == 0)
477  {
478  std::string error_string =
479  "Bulk element must inherit from PMLHelmholtzEquations.";
480  error_string +=
481  "Nodes are two dimensional, but cannot cast the bulk element to\n";
482  error_string += "PMLHelmholtzEquations<2>\n.";
483  error_string += "If you desire this functionality, you must "
484  "implement it yourself\n";
485 
486  throw OomphLibError(
488  }
489  else
490  {
491  // Read the index from the (cast) bulk element
492  U_index_helmholtz = eqn_pt->u_index_helmholtz();
493  }
494  }
495 
496  break;
497 
498  // Three dimensional problem
499  case 3:
500  {
501  PMLHelmholtzEquations<3>* eqn_pt =
502  dynamic_cast<PMLHelmholtzEquations<3>*>(bulk_el_pt);
503  // If the cast has failed die
504  if (eqn_pt == 0)
505  {
506  std::string error_string =
507  "Bulk element must inherit from PMLHelmholtzEquations.";
508  error_string += "Nodes are three dimensional, but cannot cast the "
509  "bulk element to\n";
510  error_string += "PMLHelmholtzEquations<3>\n.";
511  error_string += "If you desire this functionality, you must "
512  "implement it yourself\n";
513 
514  throw OomphLibError(
516  }
517  else
518  {
519  // Read the index from the (cast) bulk element
520  U_index_helmholtz = eqn_pt->u_index_helmholtz();
521  }
522  }
523  break;
524 
525  // Any other case is an error
526  default:
527  std::ostringstream error_stream;
528  error_stream << "Dimension of node is " << Dim
529  << ". It should be 1,2, or 3!" << std::endl;
530 
531  throw OomphLibError(
532  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
533  break;
534  }
535  }
int & face_index()
Definition: elements.h:4626
FaceElement()
Constructor: Initialise all appropriate member data.
Definition: elements.h:4446
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
std::complex< unsigned > U_index_helmholtz
Definition: pml_helmholtz_flux_elements.h:375
unsigned Dim
The spatial dimension of the problem.
Definition: pml_helmholtz_flux_elements.h:378
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FiniteElement::build_face_element(), oomph::PMLHelmholtzPowerElement< ELEMENT >::Dim, oomph::FaceElement::face_index(), oomph::FiniteElement::has_hanging_nodes(), oomph::Node::ndim(), oomph::FiniteElement::node_pt(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Global_string_for_annotation::string(), oomph::PMLHelmholtzEquations< DIM >::u_index_helmholtz(), and oomph::PMLHelmholtzPowerElement< ELEMENT >::U_index_helmholtz.

◆ PMLHelmholtzPowerElement() [2/3]

template<class ELEMENT >
oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement ( )
inline

Broken empty constructor.

62  {
63  throw OomphLibError(
64  "Don't call empty constructor for PMLHelmholtzPowerElement",
67  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ PMLHelmholtzPowerElement() [3/3]

template<class ELEMENT >
oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement ( const PMLHelmholtzPowerElement< ELEMENT > &  dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ global_flux_contribution() [1/2]

template<class ELEMENT >
std::complex<double> oomph::PMLHelmholtzPowerElement< ELEMENT >::global_flux_contribution ( )
inline

Compute the element's contribution to the time-averaged flux

248  {
249  // Dummy output file
250  std::ofstream outfile;
251  return global_flux_contribution(outfile);
252  }
std::complex< double > global_flux_contribution()
Definition: pml_helmholtz_flux_elements.h:247

◆ global_flux_contribution() [2/2]

template<class ELEMENT >
std::complex<double> oomph::PMLHelmholtzPowerElement< ELEMENT >::global_flux_contribution ( std::ofstream &  outfile)
inline

Compute the element's contribution to the integral of the flux over the artificial boundary. Also output the flux in the specified output file if it's open.

259  {
260  // pointer to the corresponding bulk element
261  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(this->bulk_element_pt());
262 
263  // Number of nodes in bulk element
264  unsigned nnode_bulk = bulk_elem_pt->nnode();
265  const unsigned n_node_local = nnode();
266 
267  // get the dim of the bulk and local nodes
268  const unsigned bulk_dim = bulk_elem_pt->dim();
269  const unsigned local_dim = this->dim();
270 
271  // Set up memory for the shape and test functions
272  Shape psi(n_node_local);
273 
274  // Set up memory for the shape functions
275  Shape psi_bulk(nnode_bulk);
276  DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
277 
278  // Set up memory for the outer unit normal
279  Vector<double> unit_normal(bulk_dim);
280 
281  // Set the value of n_intpt
282  const unsigned n_intpt = integral_pt()->nweight();
283 
284  // Set the Vector to hold local coordinates
285  Vector<double> s(local_dim);
286  std::complex<double> flux(0.0, 0.0);
287 
288  // Output?
289  if (outfile.is_open())
290  {
291  outfile << "ZONE\n";
292  }
293 
294  // Loop over the integration points
295  //--------------------------------
296  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
297  {
298  // Assign values of s
299  for (unsigned i = 0; i < local_dim; i++)
300  {
301  s[i] = integral_pt()->knot(ipt, i);
302  }
303  // get the outer_unit_ext vector
304  this->outer_unit_normal(s, unit_normal);
305 
306  // Get the integral weight
307  double w = integral_pt()->weight(ipt);
308 
309  // Get jacobian of mapping
310  double J = J_eulerian(s);
311 
312  // Premultiply the weights and the Jacobian
313  double W = w * J;
314 
315  // Get local coordinates in bulk element by copy construction
316  Vector<double> s_bulk(local_coordinate_in_bulk(s));
317 
318  // Call the derivatives of the shape functions
319  // in the bulk -- must do this via s because this point
320  // is not an integration point the bulk element!
321  (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
322  this->shape(s, psi);
323 
324  // Derivs of Eulerian coordinates w.r.t. local coordinates
325  std::complex<double> dphi_dn(0.0, 0.0);
326  Vector<std::complex<double>> interpolated_dphidx(
327  bulk_dim, std::complex<double>(0.0, 0.0));
328  Vector<double> x(bulk_dim);
329 
330  // Calculate function value and derivatives:
331  //-----------------------------------------
332  // Loop over nodes
333  for (unsigned l = 0; l < nnode_bulk; l++)
334  {
335  // Get the nodal value of the helmholtz unknown
336  const std::complex<double> phi_value(
337  bulk_elem_pt->nodal_value(l,
338  bulk_elem_pt->u_index_helmholtz().real()),
339  bulk_elem_pt->nodal_value(
340  l, bulk_elem_pt->u_index_helmholtz().imag()));
341 
342  // Loop over directions
343  for (unsigned i = 0; i < bulk_dim; i++)
344  {
345  interpolated_dphidx[i] += phi_value * dpsi_bulk_dx(l, i);
346  }
347  } // End of loop over the bulk_nodes
348 
349 
350  // define dphi_dn
351  for (unsigned i = 0; i < bulk_dim; i++)
352  {
353  dphi_dn += interpolated_dphidx[i] * unit_normal[i];
354  }
355 
356  // Output?
357  if (outfile.is_open())
358  {
359  interpolated_x(s, x);
360  outfile << x[0] << " " << x[1] << " " << dphi_dn.real() << " "
361  << dphi_dn.imag() << "\n";
362  }
363 
364  // ...add to integral
365  flux += dphi_dn * W;
366  }
367 
368  return flux;
369  }
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
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
Vector< double > local_coordinate_in_bulk(const Vector< double > &s) const
Definition: elements.cc:6353
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:5242
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned dim() const
Definition: elements.h:2611
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 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
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

References oomph::FaceElement::bulk_element_pt(), oomph::FiniteElement::dim(), ProblemParameters::flux(), i, oomph::FiniteElement::integral_pt(), oomph::FaceElement::interpolated_x(), J, oomph::FaceElement::J_eulerian(), oomph::Integral::knot(), oomph::FaceElement::local_coordinate_in_bulk(), oomph::FiniteElement::nnode(), oomph::Integral::nweight(), oomph::FaceElement::outer_unit_normal(), s, oomph::FiniteElement::shape(), w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and plotDoE::x.

◆ global_power_contribution() [1/2]

template<class ELEMENT >
double oomph::PMLHelmholtzPowerElement< ELEMENT >::global_power_contribution ( )
inline

Compute the element's contribution to the time-averaged radiated power over the artificial boundary

106  {
107  // Dummy output file
108  std::ofstream outfile;
109  return global_power_contribution(outfile);
110  }
double global_power_contribution()
Definition: pml_helmholtz_flux_elements.h:105

◆ global_power_contribution() [2/2]

template<class ELEMENT >
double oomph::PMLHelmholtzPowerElement< ELEMENT >::global_power_contribution ( std::ofstream &  outfile)
inline

Compute the element's contribution to the time-averaged radiated power over the artificial boundary. Also output the power density in the specified output file if it's open.

117  {
118  // pointer to the corresponding bulk element
119  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(this->bulk_element_pt());
120 
121  // Number of nodes in bulk element
122  unsigned nnode_bulk = bulk_elem_pt->nnode();
123  const unsigned n_node_local = nnode();
124 
125  // get the dim of the bulk and local nodes
126  const unsigned bulk_dim = bulk_elem_pt->dim();
127  const unsigned local_dim = this->dim();
128 
129  // Set up memory for the shape and test functions
130  Shape psi(n_node_local);
131 
132  // Set up memory for the shape functions
133  Shape psi_bulk(nnode_bulk);
134  DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
135 
136  // Set up memory for the outer unit normal
137  Vector<double> unit_normal(bulk_dim);
138 
139  // Set the value of n_intpt
140  const unsigned n_intpt = integral_pt()->nweight();
141 
142  // Set the Vector to hold local coordinates
143  Vector<double> s(local_dim);
144  double power = 0.0;
145 
146  // Output?
147  if (outfile.is_open())
148  {
149  outfile << "ZONE\n";
150  }
151 
152  // Loop over the integration points
153  //--------------------------------
154  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
155  {
156  // Assign values of s
157  for (unsigned i = 0; i < local_dim; i++)
158  {
159  s[i] = integral_pt()->knot(ipt, i);
160  }
161  // get the outer_unit_ext vector
162  this->outer_unit_normal(s, unit_normal);
163 
164  // Get the integral weight
165  double w = integral_pt()->weight(ipt);
166 
167  // Get jacobian of mapping
168  double J = J_eulerian(s);
169 
170  // Premultiply the weights and the Jacobian
171  double W = w * J;
172 
173  // Get local coordinates in bulk element by copy construction
174  Vector<double> s_bulk(local_coordinate_in_bulk(s));
175 
176  // Call the derivatives of the shape functions
177  // in the bulk -- must do this via s because this point
178  // is not an integration point the bulk element!
179  (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
180  this->shape(s, psi);
181 
182  // Derivs of Eulerian coordinates w.r.t. local coordinates
183  std::complex<double> dphi_dn(0.0, 0.0);
184  Vector<std::complex<double>> interpolated_dphidx(
185  bulk_dim, std::complex<double>(0.0, 0.0));
186  std::complex<double> interpolated_phi(0.0, 0.0);
187  Vector<double> x(bulk_dim);
188 
189  // Calculate function value and derivatives:
190  //-----------------------------------------
191  // Loop over nodes
192  for (unsigned l = 0; l < nnode_bulk; l++)
193  {
194  // Get the nodal value of the helmholtz unknown
195  const std::complex<double> phi_value(
196  bulk_elem_pt->nodal_value(l,
197  bulk_elem_pt->u_index_helmholtz().real()),
198  bulk_elem_pt->nodal_value(
199  l, bulk_elem_pt->u_index_helmholtz().imag()));
200 
201  // Loop over directions
202  for (unsigned i = 0; i < bulk_dim; i++)
203  {
204  interpolated_dphidx[i] += phi_value * dpsi_bulk_dx(l, i);
205  }
206  } // End of loop over the bulk_nodes
207 
208  for (unsigned l = 0; l < n_node_local; l++)
209  {
210  // Get the nodal value of the helmholtz unknown
211  const std::complex<double> phi_value(
214 
215  interpolated_phi += phi_value * psi(l);
216  }
217 
218  // define dphi_dn
219  for (unsigned i = 0; i < bulk_dim; i++)
220  {
221  dphi_dn += interpolated_dphidx[i] * unit_normal[i];
222  }
223 
224  // Power density
225  double integrand = 0.5 * (interpolated_phi.real() * dphi_dn.imag() -
226  interpolated_phi.imag() * dphi_dn.real());
227 
228  // Output?
229  if (outfile.is_open())
230  {
231  interpolated_x(s, x);
232  double phi = atan2(x[1], x[0]);
233  outfile << x[0] << " " << x[1] << " " << phi << " " << integrand
234  << "\n";
235  }
236 
237  // ...add to integral
238  power += integrand * W;
239  }
240 
241  return power;
242  }
AnnoyingScalar atan2(const AnnoyingScalar &y, const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:139
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual std::complex< unsigned > u_index_helmholtz() const
Definition: pml_helmholtz_flux_elements.h:96
float real
Definition: datatypes.h:10

References atan2(), oomph::FaceElement::bulk_element_pt(), oomph::FiniteElement::dim(), i, imag(), oomph::FiniteElement::integral_pt(), oomph::FaceElement::interpolated_x(), J, oomph::FaceElement::J_eulerian(), oomph::Integral::knot(), oomph::FaceElement::local_coordinate_in_bulk(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), oomph::Integral::nweight(), oomph::FaceElement::outer_unit_normal(), s, oomph::FiniteElement::shape(), oomph::PMLHelmholtzPowerElement< ELEMENT >::u_index_helmholtz(), w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and plotDoE::x.

◆ u_index_helmholtz()

template<class ELEMENT >
virtual std::complex<unsigned> oomph::PMLHelmholtzPowerElement< ELEMENT >::u_index_helmholtz ( ) const
inlinevirtual

Return the index at which the unknown value is stored.

97  {
98  return std::complex<unsigned>(U_index_helmholtz.real(),
99  U_index_helmholtz.imag());
100  }

References oomph::PMLHelmholtzPowerElement< ELEMENT >::U_index_helmholtz.

Referenced by oomph::PMLHelmholtzPowerElement< ELEMENT >::global_power_contribution().

◆ zeta_nodal()

template<class ELEMENT >
double oomph::PMLHelmholtzPowerElement< ELEMENT >::zeta_nodal ( const unsigned n,
const unsigned k,
const unsigned i 
) const
inlinevirtual

Broken assignment operator.

Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default (needed to break indeterminacy if bulk element is SolidElement)

Reimplemented from oomph::FaceElement.

89  {
90  return FaceElement::zeta_nodal(n, k, i);
91  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
char char char int int * k
Definition: level2_impl.h:374

References i, k, n, and oomph::FaceElement::zeta_nodal().

Member Data Documentation

◆ Dim

template<class ELEMENT >
unsigned oomph::PMLHelmholtzPowerElement< ELEMENT >::Dim
protected

The spatial dimension of the problem.

Referenced by oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement().

◆ U_index_helmholtz

template<class ELEMENT >
std::complex<unsigned> oomph::PMLHelmholtzPowerElement< ELEMENT >::U_index_helmholtz
protected

The index at which the real and imag part of the unknown is stored at the nodes

Referenced by oomph::PMLHelmholtzPowerElement< ELEMENT >::PMLHelmholtzPowerElement(), and oomph::PMLHelmholtzPowerElement< ELEMENT >::u_index_helmholtz().


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