oomph::LinearWaveEquations< DIM > Class Template Referenceabstract

#include <linear_wave_elements.h>

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

Public Types

typedef void(* LinearWaveSourceFctPt) (const double &time, const Vector< double > &x, double &u)
 
- 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

 LinearWaveEquations ()
 Constructor (must initialise the Source_fct_pt to null) More...
 
 LinearWaveEquations (const LinearWaveEquations &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const LinearWaveEquations &)=delete
 Broken assignment operator. More...
 
virtual unsigned u_index_lin_wave () const
 
double du_dt_lin_wave (const unsigned &n) const
 
double d2u_dt2_lin_wave (const unsigned &n) const
 
void output (std::ostream &outfile)
 Output with default number of plot points. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void output (FILE *file_pt)
 Output with default number of plot points. More...
 
void output (FILE *file_pt, const unsigned &nplot)
 
void output_fct (std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt 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)
 Get error against and norm of exact solution. More...
 
LinearWaveSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
LinearWaveSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
void get_source_lin_wave (const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
 
void get_flux (const Vector< double > &s, Vector< double > &flux) const
 Get flux: flux[i] = du/dx_i. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Compute element residual Vector (wrapper) More...
 
virtual void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Compute element residual Vector and element Jacobian matrix (wrapper) More...
 
double interpolated_u_lin_wave (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
double interpolated_du_dt_lin_wave (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
double interpolated_d2u_dt2_lin_wave (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual 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 void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

virtual double dshape_and_dtest_eulerian_lin_wave (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_lin_wave (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_lin_wave (Vector< double > &residuals, DenseMatrix< double > &jacobian, 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

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

Additional Inherited Members

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

Detailed Description

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

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

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

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

Member Typedef Documentation

◆ LinearWaveSourceFctPt

template<unsigned DIM>
typedef void(* oomph::LinearWaveEquations< DIM >::LinearWaveSourceFctPt) (const double &time, const Vector< double > &x, double &u)

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

Constructor & Destructor Documentation

◆ LinearWaveEquations() [1/2]

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

Constructor (must initialise the Source_fct_pt to null)

62 : Source_fct_pt(0) {}
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
Definition: linear_wave_elements.h:380

◆ LinearWaveEquations() [2/2]

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

Broken copy constructor.

Member Function Documentation

◆ compute_error() [1/2]

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

Get error against and norm of exact solution.

Validate against exact solution

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

Reimplemented from oomph::FiniteElement.

419  {
420  // Initialise
421  error = 0.0;
422  norm = 0.0;
423 
424  // Vector of local coordinates
425  Vector<double> s(DIM);
426 
427  // Vector for coordintes
428  Vector<double> x(DIM);
429 
430 
431  // Find out how many nodes there are in the element
432  unsigned n_node = nnode();
433 
434  Shape psi(n_node);
435 
436  // Set the value of n_intpt
437  unsigned n_intpt = integral_pt()->nweight();
438 
439  // Tecplot header info
440  outfile << "ZONE" << std::endl;
441 
442  // Exact solution Vector (value, and first and second time derivs)
443  Vector<double> exact_soln(3);
444 
445  // Loop over the integration points
446  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
447  {
448  // Assign values of s
449  for (unsigned i = 0; i < DIM; i++)
450  {
451  s[i] = integral_pt()->knot(ipt, i);
452  }
453 
454  // Get the integral weight
455  double w = integral_pt()->weight(ipt);
456 
457  // Get jacobian of mapping
458  double J = J_eulerian(s);
459 
460  // Premultiply the weights and the Jacobian
461  double W = w * J;
462 
463  // Get x position as Vector
464  interpolated_x(s, x);
465 
466  // Get FE function value
467  double u_fe = interpolated_u_lin_wave(s);
468 
469  // Get exact solution at this point
470  (*exact_soln_pt)(x, exact_soln);
471 
472  // Output x,y,...,error
473  for (unsigned i = 0; i < DIM; i++)
474  {
475  outfile << x[i] << " ";
476  }
477  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
478 
479  // Add to error and norm
480  norm += exact_soln[0] * exact_soln[0] * W;
481  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
482  }
483  }
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.
double interpolated_u_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:274
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

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

◆ compute_error() [2/2]

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

Get error against and norm of exact solution.

Validate against exact solution at time t.

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

Reimplemented from oomph::FiniteElement.

501  {
502  // Initialise
503  error = 0.0;
504  norm = 0.0;
505 
506  // Vector of local coordinates
507  Vector<double> s(DIM);
508 
509  // Vector for coordintes
510  Vector<double> x(DIM);
511 
512  // Find out how many nodes there are in the element
513  unsigned n_node = nnode();
514 
515  Shape psi(n_node);
516 
517  // Set the value of n_intpt
518  unsigned n_intpt = integral_pt()->nweight();
519 
520  // Tecplot header info
521  outfile << "ZONE" << std::endl;
522 
523  // Exact solution Vector (value, and first and second time derivs)
524  Vector<double> exact_soln(3);
525 
526  // Loop over the integration points
527  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
528  {
529  // Assign values of s
530  for (unsigned i = 0; i < DIM; i++)
531  {
532  s[i] = integral_pt()->knot(ipt, i);
533  }
534 
535  // Get the integral weight
536  double w = integral_pt()->weight(ipt);
537 
538  // Get jacobian of mapping
539  double J = J_eulerian(s);
540 
541  // Premultiply the weights and the Jacobian
542  double W = w * J;
543 
544  // Get x position as Vector
545  interpolated_x(s, x);
546 
547  // Get FE function value
548  double u_fe = interpolated_u_lin_wave(s);
549 
550  // Get exact solution at this point
551  (*exact_soln_pt)(time, x, exact_soln);
552 
553  // Output x,y,...,error
554  for (unsigned i = 0; i < DIM; i++)
555  {
556  outfile << x[i] << " ";
557  }
558  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
559 
560  // Add to error and norm
561  norm += exact_soln[0] * exact_soln[0] * W;
562  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
563  }
564  }

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

◆ d2u_dt2_lin_wave()

template<unsigned DIM>
double oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave ( const unsigned n) const
inline

d^2u/dt^2 at local node n. Uses suitably interpolated value for hanging nodes.

113  {
114  // Get the data's timestepper
115  TimeStepper* time_stepper_pt = node_pt(n)->time_stepper_pt();
116 
117  // Initialise d^2u/dt^2
118  double ddudt = 0.0;
119  // Loop over the timesteps, if there is a non steady timestepper
120  if (!time_stepper_pt->is_steady())
121  {
122  // Find the index at which the linear wave unknown is stored
123  const unsigned u_nodal_index = u_index_lin_wave();
124 
125  // Number of timsteps (past & present)
126  const unsigned n_time = time_stepper_pt->ntstorage();
127 
128  // Add the contributions to the derivative
129  for (unsigned t = 0; t < n_time; t++)
130  {
131  ddudt +=
132  time_stepper_pt->weight(2, t) * nodal_value(t, n, u_nodal_index);
133  }
134  }
135  return ddudt;
136  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
virtual unsigned u_index_lin_wave() const
Definition: linear_wave_elements.h:77
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
t
Definition: plotPSD.py:36

References oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), oomph::LinearWaveEquations< DIM >::u_index_lin_wave(), and oomph::TimeStepper::weight().

Referenced by oomph::LinearWaveEquations< DIM >::interpolated_d2u_dt2_lin_wave().

◆ dshape_and_dtest_eulerian_at_knot_lin_wave()

template<unsigned DIM>
virtual double oomph::LinearWaveEquations< DIM >::dshape_and_dtest_eulerian_at_knot_lin_wave ( 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::QLinearWaveElement< DIM, NNODE_1D >.

◆ dshape_and_dtest_eulerian_lin_wave()

template<unsigned DIM>
virtual double oomph::LinearWaveEquations< DIM >::dshape_and_dtest_eulerian_lin_wave ( 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::QLinearWaveElement< DIM, NNODE_1D >.

◆ du_dt_lin_wave()

template<unsigned DIM>
double oomph::LinearWaveEquations< DIM >::du_dt_lin_wave ( const unsigned n) const
inline

du/dt at local node n. Uses suitably interpolated value for hanging nodes.

85  {
86  // Get the data's timestepper
87  TimeStepper* time_stepper_pt = node_pt(n)->time_stepper_pt();
88 
89  // Initialise d^2u/dt^2
90  double dudt = 0.0;
91  // Loop over the timesteps, if there is a non steady timestepper
92  if (!time_stepper_pt->is_steady())
93  {
94  // Find the index at which the linear wave unknown is stored
95  const unsigned u_nodal_index = u_index_lin_wave();
96 
97  // Number of timsteps (past & present)
98  const unsigned n_time = time_stepper_pt->ntstorage();
99 
100  // Add the contributions to the derivative
101  for (unsigned t = 0; t < n_time; t++)
102  {
103  dudt +=
104  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
105  }
106  }
107  return dudt;
108  }

References oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), oomph::LinearWaveEquations< DIM >::u_index_lin_wave(), and oomph::TimeStepper::weight().

Referenced by oomph::LinearWaveEquations< DIM >::interpolated_du_dt_lin_wave().

◆ fill_in_contribution_to_jacobian()

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

Compute element residual Vector and element Jacobian matrix (wrapper)

Reimplemented from oomph::FiniteElement.

267  {
268  // Call the generic routine with the flag set to 1
269  fill_in_generic_residual_contribution_lin_wave(residuals, jacobian, 1);
270  }
virtual void fill_in_generic_residual_contribution_lin_wave(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: linear_wave_elements.cc:79

References oomph::LinearWaveEquations< DIM >::fill_in_generic_residual_contribution_lin_wave().

◆ fill_in_contribution_to_residuals()

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

Compute element residual Vector (wrapper)

Reimplemented from oomph::GeneralisedElement.

256  {
257  // Call the generic residuals function with flag set to 0
258  // using the dummy matrix argument
260  residuals, GeneralisedElement::Dummy_matrix, 0);
261  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::LinearWaveEquations< DIM >::fill_in_generic_residual_contribution_lin_wave().

◆ fill_in_generic_residual_contribution_lin_wave()

template<unsigned DIM>
void oomph::LinearWaveEquations< DIM >::fill_in_generic_residual_contribution_lin_wave ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
unsigned  flag 
)
protectedvirtual

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

Compute element residual Vector and/or element Jacobian matrix

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

Pure version without hanging nodes

Reimplemented in oomph::RefineableLinearWaveEquations< DIM >.

81  {
82  // Find out how many nodes there are
83  unsigned n_node = nnode();
84 
85  // Get continuous time from timestepper of first node
86  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
87 
88  // Find the index at which the linear wave unknown is stored
89  unsigned u_nodal_index = u_index_lin_wave();
90 
91  // Set up memory for the shape and test functions
92  Shape psi(n_node), test(n_node);
93  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
94 
95  // Set the value of n_intpt
96  unsigned n_intpt = integral_pt()->nweight();
97 
98  // Set the Vector to hold local coordinates
99  Vector<double> s(DIM);
100 
101  // Integers to hold the local equation and unknowns
102  int local_eqn = 0, local_unknown = 0;
103 
104  // Loop over the integration points
105  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
106  {
107  // Assign values of s
108  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
109 
110  // Get the integral weight
111  double w = integral_pt()->weight(ipt);
112 
113  // Call the derivatives of the shape and test functions
115  ipt, psi, dpsidx, test, dtestdx);
116 
117  // Premultiply the weights and the Jacobian
118  double W = w * J;
119 
120  // Allocate memory for local quantities and initialise to zero
121  double interpolated_u = 0.0;
122  double ddudt = 0.0;
123  Vector<double> interpolated_x(DIM, 0.0);
124  Vector<double> interpolated_dudx(DIM, 0.0);
125 
126  // Calculate function value and derivatives:
127  // Loop over nodes
128  for (unsigned l = 0; l < n_node; l++)
129  {
130  // Get the nodal value
131  double u_value = raw_nodal_value(l, u_nodal_index);
132  interpolated_u += u_value * psi(l);
133  ddudt += d2u_dt2_lin_wave(l) * psi(l);
134  // Loop over directions
135  for (unsigned j = 0; j < DIM; j++)
136  {
137  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
138  interpolated_dudx[j] += u_value * dpsidx(l, j);
139  }
140  }
141 
142 
143  // Get source function
144  //-------------------
145  double source;
147 
148  // Assemble residuals and Jacobian
149  //--------------------------------
150 
151  // Loop over the test functions
152  for (unsigned l = 0; l < n_node; l++)
153  {
154  local_eqn = nodal_local_eqn(l, u_nodal_index);
155  /*IF it's not a boundary condition*/
156  if (local_eqn >= 0)
157  {
158  // Add body force/source term and time derivative
159  residuals[local_eqn] += (ddudt + source) * test(l) * W;
160 
161  // Laplace operator
162  for (unsigned k = 0; k < DIM; k++)
163  {
164  residuals[local_eqn] += interpolated_dudx[k] * dtestdx(l, k) * W;
165  }
166 
167 
168  // Calculate the jacobian
169  //-----------------------
170  if (flag)
171  {
172  // Loop over the velocity shape functions again
173  for (unsigned l2 = 0; l2 < n_node; l2++)
174  {
175  local_unknown = nodal_local_eqn(l2, u_nodal_index);
176  // If at a non-zero degree of freedom add in the entry
177  if (local_unknown >= 0)
178  {
179  // Mass matrix
180  jacobian(local_eqn, local_unknown) +=
181  test(l) * psi(l2) *
182  node_pt(l2)->time_stepper_pt()->weight(2, 0) * W;
183 
184  // Laplace operator
185  for (unsigned i = 0; i < DIM; i++)
186  {
187  jacobian(local_eqn, local_unknown) +=
188  dpsidx(l2, i) * dtestdx(l, i) * W;
189  }
190  }
191  }
192  }
193  }
194  }
195 
196 
197  } // End of loop over integration points
198  }
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
void get_source_lin_wave(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: linear_wave_elements.h:202
virtual double dshape_and_dtest_eulerian_at_knot_lin_wave(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
double d2u_dt2_lin_wave(const unsigned &n) const
Definition: linear_wave_elements.h:112
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
char char char int int * k
Definition: level2_impl.h:374
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

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

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

◆ get_flux()

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

Get flux: flux[i] = du/dx_i.

222  {
223  // Find out how many nodes there are in the element
224  unsigned n_node = nnode();
225 
226  // Find the index at which the linear wave unknown is stored
227  unsigned u_nodal_index = u_index_lin_wave();
228 
229  // Set up memory for the shape and test functions
230  Shape psi(n_node);
231  DShape dpsidx(n_node, DIM);
232 
233  // Call the derivatives of the shape and test functions
234  dshape_eulerian(s, psi, dpsidx);
235 
236  // Initialise to zero
237  for (unsigned j = 0; j < DIM; j++)
238  {
239  flux[j] = 0.0;
240  }
241 
242  // Loop over nodes
243  for (unsigned l = 0; l < n_node; l++)
244  {
245  // Loop over derivative directions
246  for (unsigned j = 0; j < DIM; j++)
247  {
248  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
249  }
250  }
251  }
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

References DIM, oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, and oomph::LinearWaveEquations< DIM >::u_index_lin_wave().

Referenced by oomph::RefineableLinearWaveEquations< DIM >::get_Z2_flux().

◆ get_source_lin_wave()

template<unsigned DIM>
void oomph::LinearWaveEquations< DIM >::get_source_lin_wave ( const double t,
const unsigned ipt,
const Vector< double > &  x,
double source 
) const
inline

Get source term at continous time t and (Eulerian) position x. Virtual so it can be overloaded in derived multiphysics elements.

206  {
207  // If no source function has been set, return zero
208  if (Source_fct_pt == 0)
209  {
210  source = 0.0;
211  }
212  else
213  {
214  // Get source strength
215  (*Source_fct_pt)(t, x, source);
216  }
217  }

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

◆ interpolated_d2u_dt2_lin_wave()

template<unsigned DIM>
double oomph::LinearWaveEquations< DIM >::interpolated_d2u_dt2_lin_wave ( const Vector< double > &  s) const
inline

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

327  {
328  // Find number of nodes
329  unsigned n_node = nnode();
330 
331  // Local shape function
332  Shape psi(n_node);
333 
334  // Find values of shape function
335  shape(s, psi);
336 
337  // Initialise value of u
338  double interpolated_d2u_dt2 = 0.0;
339 
340  // Loop over the local nodes and sum
341  for (unsigned l = 0; l < n_node; l++)
342  {
343  interpolated_d2u_dt2 += d2u_dt2_lin_wave(l) * psi[l];
344  }
345 
346  return (interpolated_d2u_dt2);
347  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::FiniteElement::nnode(), s, and oomph::FiniteElement::shape().

◆ interpolated_du_dt_lin_wave()

template<unsigned DIM>
double oomph::LinearWaveEquations< DIM >::interpolated_du_dt_lin_wave ( const Vector< double > &  s) const
inline

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

303  {
304  // Find number of nodes
305  unsigned n_node = nnode();
306 
307  // Local shape function
308  Shape psi(n_node);
309 
310  // Find values of shape function
311  shape(s, psi);
312 
313  // Initialise value of u
314  double interpolated_du_dt = 0.0;
315 
316  // Loop over the local nodes and sum
317  for (unsigned l = 0; l < n_node; l++)
318  {
319  interpolated_du_dt += du_dt_lin_wave(l) * psi[l];
320  }
321 
322  return (interpolated_du_dt);
323  }
double du_dt_lin_wave(const unsigned &n) const
Definition: linear_wave_elements.h:84

References oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::FiniteElement::nnode(), s, and oomph::FiniteElement::shape().

◆ interpolated_u_lin_wave()

template<unsigned DIM>
double oomph::LinearWaveEquations< DIM >::interpolated_u_lin_wave ( const Vector< double > &  s) const
inline

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

275  {
276  // Find number of nodes
277  unsigned n_node = nnode();
278 
279  // Find the index at which the linear wave unknown is stored
280  unsigned u_nodal_index = u_index_lin_wave();
281 
282  // Local shape function
283  Shape psi(n_node);
284 
285  // Find values of shape function
286  shape(s, psi);
287 
288  // Initialise value of u
289  double interpolated_u = 0.0;
290 
291  // Loop over the local nodes and sum
292  for (unsigned l = 0; l < n_node; l++)
293  {
294  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
295  }
296 
297  return (interpolated_u);
298  }

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

◆ operator=()

template<unsigned DIM>
void oomph::LinearWaveEquations< DIM >::operator= ( const LinearWaveEquations< DIM > &  )
delete

Broken assignment operator.

◆ output() [1/4]

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

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

152  {
153  unsigned nplot = 5;
154  output(file_pt, nplot);
155  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: linear_wave_elements.h:139

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

◆ output() [2/4]

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

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

C-style output function:

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

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

272  {
273  // Vector of local coordinates
274  Vector<double> s(DIM);
275 
276  // Tecplot header info
277  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
278 
279  // Loop over plot points
280  unsigned num_plot_points = nplot_points(nplot);
281  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
282  {
283  // Get local coordinates of plot point
284  get_s_plot(iplot, nplot, s);
285 
286  for (unsigned i = 0; i < DIM; i++)
287  {
288  // outfile << interpolated_x(s,i) << " ";
289  fprintf(file_pt, "%g ", interpolated_x(s, i));
290  }
291  // outfile << interpolated_u(s) << std::endl;
292  fprintf(file_pt, "%g \n", interpolated_u_lin_wave(s));
293  }
294 
295  // Write tecplot footer (e.g. FE connectivity lists)
296  write_tecplot_zone_footer(file_pt, nplot);
297  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References DIM, i, and s.

◆ output() [3/4]

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

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

140  {
141  unsigned nplot = 5;
142  output(outfile, nplot);
143  }

Referenced by oomph::LinearWaveEquations< DIM >::output(), and oomph::QLinearWaveElement< DIM, NNODE_1D >::output().

◆ output() [4/4]

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

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

Output function:

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

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

237  {
238  // Vector of local coordinates
239  Vector<double> s(DIM);
240 
241  // Tecplot header info
242  outfile << tecplot_zone_string(nplot);
243 
244  // Loop over plot points
245  unsigned num_plot_points = nplot_points(nplot);
246  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
247  {
248  // Get local coordinates of plot point
249  get_s_plot(iplot, nplot, s);
250  for (unsigned i = 0; i < DIM; i++)
251  {
252  outfile << interpolated_x(s, i) << " ";
253  }
254  outfile << interpolated_u_lin_wave(s) << " ";
255  outfile << interpolated_du_dt_lin_wave(s) << " ";
256  outfile << interpolated_d2u_dt2_lin_wave(s) << std::endl;
257  }
258 
259  // Write tecplot footer (e.g. FE connectivity lists)
260  write_tecplot_zone_footer(outfile, nplot);
261  }
double interpolated_du_dt_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:302
double interpolated_d2u_dt2_lin_wave(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: linear_wave_elements.h:326

References DIM, i, and s.

◆ output_fct() [1/2]

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

Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points (time-dependent version)

Output exact solution at time t

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

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

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

366  {
367  // Vector of local coordinates
368  Vector<double> s(DIM);
369 
370  // Vector for coordintes
371  Vector<double> x(DIM);
372 
373  // Tecplot header info
374  outfile << tecplot_zone_string(nplot);
375 
376  // Exact solution Vector u, dudt, d2udt2
377  Vector<double> exact_soln(3);
378 
379  // Loop over plot points
380  unsigned num_plot_points = nplot_points(nplot);
381  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
382  {
383  // Get local coordinates of plot point
384  get_s_plot(iplot, nplot, s);
385 
386  // Get x position as Vector
387  interpolated_x(s, x);
388 
389  // Get exact solution at this point
390  (*exact_soln_pt)(time, x, exact_soln);
391 
392  // Output x,y,...,u_exact
393  for (unsigned i = 0; i < DIM; i++)
394  {
395  outfile << x[i] << " ";
396  }
397  outfile << exact_soln[0] << " " << exact_soln[1] << " " << exact_soln[2]
398  << " " << std::endl;
399  }
400 
401  // Write tecplot footer (e.g. FE connectivity lists)
402  write_tecplot_zone_footer(outfile, nplot);
403  }

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

◆ output_fct() [2/2]

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

Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points.

Output exact solution

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

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

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QLinearWaveElement< DIM, NNODE_1D >.

313  {
314  // Vector of local coordinates
315  Vector<double> s(DIM);
316 
317  // Vector for coordintes
318  Vector<double> x(DIM);
319 
320  // Tecplot header info
321  outfile << tecplot_zone_string(nplot);
322 
323  // Exact solution Vector (value and first and second time derivs)
324  Vector<double> exact_soln(3);
325 
326  // Loop over plot points
327  unsigned num_plot_points = nplot_points(nplot);
328  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
329  {
330  // Get local coordinates of plot point
331  get_s_plot(iplot, nplot, s);
332 
333  // Get x position as Vector
334  interpolated_x(s, x);
335 
336  // Get exact solution at this point
337  (*exact_soln_pt)(x, exact_soln);
338 
339  // Output x,y,...,u_exact
340  for (unsigned i = 0; i < DIM; i++)
341  {
342  outfile << x[i] << " ";
343  }
344  outfile << exact_soln[0] << std::endl;
345  }
346 
347  // Write tecplot footer (e.g. FE connectivity lists)
348  write_tecplot_zone_footer(outfile, nplot);
349  }

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

Referenced by oomph::QLinearWaveElement< DIM, NNODE_1D >::output_fct().

◆ self_test()

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

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

206  {
207  bool passed = true;
208 
209  // Check lower-level stuff
210  if (FiniteElement::self_test() != 0)
211  {
212  passed = false;
213  }
214 
215  // Return verdict
216  if (passed)
217  {
218  return 0;
219  }
220  else
221  {
222  return 1;
223  }
224  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

template<unsigned DIM>
LinearWaveSourceFctPt& oomph::LinearWaveEquations< DIM >::source_fct_pt ( )
inline

Access function: Pointer to source function.

189  {
190  return Source_fct_pt;
191  }

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

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

◆ source_fct_pt() [2/2]

template<unsigned DIM>
LinearWaveSourceFctPt oomph::LinearWaveEquations< DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

195  {
196  return Source_fct_pt;
197  }

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

◆ u_index_lin_wave()

template<unsigned DIM>
virtual unsigned oomph::LinearWaveEquations< DIM >::u_index_lin_wave ( ) const
inlinevirtual

Return the index at which the unknown value is stored. The default value, 0, is appropriate for single-physics problems, when there is only one variable, the value that satisfies the linear wave equation. In derived multi-physics elements, this function should be overloaded to reflect the chosen storage scheme. Note that these equations require that the unknown is always stored at the same index at each node.

78  {
79  return 0;
80  }

Referenced by oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::LinearWaveEquations< DIM >::get_flux(), oomph::RefineableLinearWaveEquations< DIM >::get_interpolated_values(), oomph::LinearWaveEquations< DIM >::interpolated_u_lin_wave(), and oomph::LinearWaveFluxElement< ELEMENT >::LinearWaveFluxElement().

Member Data Documentation

◆ Source_fct_pt


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