OrrSommerfeldEquations< DIM > Class Template Referenceabstract
+ Inheritance diagram for OrrSommerfeldEquations< DIM >:

Public Member Functions

 OrrSommerfeldEquations ()
 Constructor. More...
 
virtual double u (const unsigned &n, const unsigned &i) const =0
 
virtual double p (const unsigned &n, const unsigned &i) const =0
 
virtual unsigned npres () const =0
 
virtual void pshape (const Vector< double > &s, Shape &psi) const =0
 
void get_base_flow (const Vector< double > &s, double &U, double &Uz)
 
double *& re_pt ()
 
double *& a_real_pt ()
 
double *& a_imag_pt ()
 
void output (ostream &outfile)
 Output with default number of plot points. More...
 
void output (ostream &outfile, const unsigned &nplot)
 
void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 Assemble the contributions to the jacobian and mass matrices. More...
 
double interpolated_p (const Vector< double > &s, const unsigned &i) const
 Return FE representation of the pressure. More...
 
double interpolated_u (const Vector< double > &s, const unsigned &i) const
 Return FE representation of function value u(s) at local coordinate s. 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 &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 (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 Member Functions

virtual double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const =0
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const =0
 
- 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_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

DenseMatrix< intU_local_eqn
 Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC) More...
 
DenseMatrix< intP_local_eqn
 
- 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
 

Private Attributes

doubleRe_pt
 
doubleA_Real_pt
 
doubleA_Imag_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
 
- 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 OrrSommerfeldEquations< DIM >

A class for all elements that solve the simple orr-sommerfeld equations for a given base flow. The formulation is in primitive variables and so contains quite a large pressure NULL space. The equations are derived from the Navier–Stokes equations by posing a perturbation to an axially uniform base-flow of the form:

\[ u = U(z) + \hat{u}, \quad v = \hat{v}, \quad w = \hat{w} \]

where the perturbation is of the form

\[ \hat{u} = \mbox{e}^{\lambda t + \mbox{i}\alpha x} \]

The equations are split into real and imaginary parts and is currently only implemented for the two-dimensional equations. This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class.

Constructor & Destructor Documentation

◆ OrrSommerfeldEquations()

template<unsigned DIM>
OrrSommerfeldEquations< DIM >::OrrSommerfeldEquations ( )
inline

Constructor.

72 {}

Member Function Documentation

◆ a_imag_pt()

template<unsigned DIM>
double* & OrrSommerfeldEquations< DIM >::a_imag_pt ( )
inline
98 {return A_Imag_pt;}
double * A_Imag_pt
Definition: orr_sommerfeld.cc:556

◆ a_real_pt()

template<unsigned DIM>
double* & OrrSommerfeldEquations< DIM >::a_real_pt ( )
inline
96 {return A_Real_pt;}
double * A_Real_pt
Definition: orr_sommerfeld.cc:555

◆ dshape_eulerian()

template<unsigned DIM>
virtual double OrrSommerfeldEquations< DIM >::dshape_eulerian ( const Vector< double > &  s,
Shape psi,
DShape dpsidx 
) const
protectedpure virtual

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

◆ dshape_eulerian_at_knot()

template<unsigned DIM>
virtual double OrrSommerfeldEquations< DIM >::dshape_eulerian_at_knot ( const unsigned ipt,
Shape psi,
DShape dpsidx 
) const
protectedpure virtual

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

Reimplemented from oomph::FiniteElement.

◆ fill_in_contribution_to_jacobian_and_mass_matrix()

template<unsigned DIM>
void OrrSommerfeldEquations< DIM >::fill_in_contribution_to_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
)
inlinevirtual

Assemble the contributions to the jacobian and mass matrices.

Reimplemented from oomph::GeneralisedElement.

148  {
149  //Find out how many nodes there are
150  unsigned n_node = nnode();
151  //Number of pressure degrees of freedom
152  unsigned n_pres = npres();
153 
154  //Set up memory for the shape functions
155  Shape psi(n_node);
156  DShape dpsidx(n_node,DIM);
157  Shape psip(n_pres);
158 
159  //Set the value of n_intpt
160  unsigned n_intpt = integral_pt()->nweight();
161 
162  //Set the Vector to hold local coordinates
164 
165  //Integers to store the local equation and unknown numbers
166  int local_eqn=0, local_unknown=0;
167 
168  //Storage for the base flow
169  double Ubar=0.0, Ubardz=0.0;
170 
171  //Get the Reynolds number
172  double Re = *Re_pt;
173  double A_real = *A_Real_pt;
174  double A_imag = *A_Imag_pt;
175 
176  //Loop over the integration points
177  for(unsigned ipt=0;ipt<n_intpt;ipt++)
178  {
179  //Assign values of s
180  for(unsigned i=0;i<DIM;i++) s[i] = integral_pt()->knot(ipt,i);
181 
182  //Get the integral weight
183  double w = integral_pt()->weight(ipt);
184 
185  //Call the derivatives of the shape and test functions
186  double J = dshape_eulerian_at_knot(ipt,
187  psi,dpsidx);
188 
189  //Get the pressure shape functions
190  pshape(s,psip);
191 
192  //Get the base velocity field
193  get_base_flow(s,Ubar,Ubardz);
194 
195  //Calculate the real and imaginary parts of our variable B
196  double B_real = (A_real*A_real - A_imag*A_imag)/Re - A_imag*Ubar;
197  double B_imag = 2.0*A_real*A_imag/Re + A_real*Ubar;
198 
199  //Premultiply the weights and the Jacobian
200  double W = w*J;
201 
202  //Assemble the contributions to the matrices
203 
204  // Loop over the test functions
205  for(unsigned l=0;l<n_node;l++)
206  {
207  //REAL PART OF THE U EQUATION
208  //Get the local equation
209  local_eqn = U_local_eqn(l,0);
210  /*IF it's not a boundary condition*/
211  if(local_eqn >= 0)
212  {
213  //Loop over the shape functions
214  for(unsigned l2=0;l2<n_node;l2++)
215  {
216  //Get the contribution from the real part of u
217  local_unknown = U_local_eqn(l2,0);
218  //If at a non-zero degree of freedom add in the entry
219  if(local_unknown >= 0)
220  {
221  //Mass matrix
222  mass_matrix(local_eqn,local_unknown) += psi(l)*psi(l2)*W;
223  //Jacobian terms
224  jacobian(local_eqn,local_unknown) -=
225  B_real*psi(l2)*psi(l)*W;
226  jacobian(local_eqn,local_unknown) -=
227  (1.0/Re)*dpsidx(l2,0)*dpsidx(l,0)*W;
228  }
229 
230  //Get the contribution from the imaginary part of u
231  local_unknown = U_local_eqn(l2,1);
232  if(local_unknown >= 0)
233  {
234  jacobian(local_eqn,local_unknown) += B_imag*psi(l2)*psi(l)*W;
235  }
236 
237  //Get the contribution from the real part of w
238  local_unknown = U_local_eqn(l2,2);
239  if(local_unknown >= 0)
240  {
241  jacobian(local_eqn,local_unknown) -= Ubardz*psi(l2)*psi(l)*W;
242  }
243 
244  }
245 
246  //Loop over the pressure shape functions
247  for(unsigned l2=0;l2<n_pres;l2++)
248  {
249  //Get the real pressure contribution to the equation
250  local_unknown = P_local_eqn(l2,0);
251  if(local_unknown >= 0)
252  {
253  jacobian(local_eqn,local_unknown) +=
254  A_imag*psip(l2)*psi(l)*W;
255  }
256  //Get the imaginary pressure contribution to the equation
257  local_unknown = P_local_eqn(l2,1);
258  if(local_unknown >= 0)
259  {
260  jacobian(local_eqn,local_unknown) +=
261  A_real*psip(l2)*psi(l)*W;
262  }
263  }
264  }
265 
266  //IMAGINARY PART OF THE U EQUATION
267  //Get the local equation
268  local_eqn = U_local_eqn(l,1);
269  /*IF it's not a boundary condition*/
270  if(local_eqn >= 0)
271  {
272  //Loop over the shape functions
273  for(unsigned l2=0;l2<n_node;l2++)
274  {
275  //Get the contribution from the real part of u
276  local_unknown = U_local_eqn(l2,0);
277  //If at a non-zero degree of freedom add in the entry
278  if(local_unknown >= 0)
279  {
280  jacobian(local_eqn,local_unknown) -= B_imag*psi(l2)*psi(l)*W;
281  }
282 
283  //Get the contribution from the imaginary part of u
284  local_unknown = U_local_eqn(l2,1);
285  if(local_unknown >= 0)
286  {
287  //Mass matrix
288  mass_matrix(local_eqn,local_unknown) += psi(l)*psi(l2)*W;
289  //Jacobian terms
290  jacobian(local_eqn,local_unknown) -= B_real*psi(l2)*psi(l)*W;
291  jacobian(local_eqn,local_unknown) -=
292  (1.0/Re)*dpsidx(l2,0)*dpsidx(l,0)*W;
293  }
294 
295  //Get the contibution from the imaginary part of w
296  local_unknown = U_local_eqn(l2,3);
297  if(local_unknown >= 0)
298  {
299  jacobian(local_eqn,local_unknown) -= Ubardz*psi(l2)*psi(l)*W;
300  }
301  }
302 
303  //Loop over the pressure shape functions
304  for(unsigned l2=0;l2<n_pres;l2++)
305  {
306  //Get the real pressure contribution to the equation
307  local_unknown = P_local_eqn(l2,0);
308  if(local_unknown >= 0)
309  {
310  jacobian(local_eqn,local_unknown) -=
311  A_real*psip(l2)*psi(l)*W;
312  }
313  //Get the imaginary pressure contribution to the equation
314  local_unknown = P_local_eqn(l2,1);
315  if(local_unknown >= 0)
316  {
317  jacobian(local_eqn,local_unknown) +=
318  A_imag*psip(l2)*psi(l)*W;
319  }
320  }
321  }
322 
323 
324  //REAL PART OF THE W EQUATION
325  //Get the local equation
326  local_eqn = U_local_eqn(l,2);
327  /*IF it's not a boundary condition*/
328  if(local_eqn >= 0)
329  {
330  //Loop over the shape functions
331  for(unsigned l2=0;l2<n_node;l2++)
332  {
333  //Get the contribution from the real part of w
334  local_unknown = U_local_eqn(l2,2);
335  //If at a non-zero degree of freedom add in the entry
336  if(local_unknown >= 0)
337  {
338  //Mass matrix
339  mass_matrix(local_eqn,local_unknown) += psi(l)*psi(l2)*W;
340  //Jacobian terms
341  jacobian(local_eqn,local_unknown) -= B_real*psi(l2)*psi(l)*W;
342  jacobian(local_eqn,local_unknown) -=
343  (1.0/Re)*dpsidx(l2,0)*dpsidx(l,0)*W;
344  }
345 
346  //Get the contribution from the imaginary part of w
347  local_unknown = U_local_eqn(l2,3);
348  if(local_unknown >= 0)
349  {
350  jacobian(local_eqn,local_unknown) += B_imag*psi(l2)*psi(l)*W;
351  }
352 
353  } //End of loop over the shape functions
354 
355 
356  //Loop over the pressure shape functions
357  for(unsigned l2=0;l2<n_pres;l2++)
358  {
359  //Get the real pressure contribution to the equation
360  local_unknown = P_local_eqn(l2,0);
361  if(local_unknown >= 0)
362  {
363  jacobian(local_eqn,local_unknown) +=
364  psip(l2)*dpsidx(l,0)*W;
365  }
366  }
367  }
368 
369  //IMAGINARY PART OF THE W EQUATION
370  //Get the local equation
371  local_eqn = U_local_eqn(l,3);
372  /*IF it's not a boundary condition*/
373  if(local_eqn >= 0)
374  {
375  //Loop over the shape functions
376  for(unsigned l2=0;l2<n_node;l2++)
377  {
378  //Get the contribution from the real part of w
379  local_unknown = U_local_eqn(l2,2);
380  //If at a non-zero degree of freedom add in the entry
381  if(local_unknown >= 0)
382  {
383  jacobian(local_eqn,local_unknown) -= B_imag*psi(l2)*psi(l)*W;
384  }
385 
386  //Get the contribution from the imaginary part of w
387  local_unknown = U_local_eqn(l2,3);
388  if(local_unknown >= 0)
389  {
390  //Mass matrix
391  mass_matrix(local_eqn,local_unknown) += psi(l)*psi(l2)*W;
392  //Jacobian terms
393  jacobian(local_eqn,local_unknown) -= B_real*psi(l2)*psi(l)*W;
394  jacobian(local_eqn,local_unknown) -=
395  (1.0/Re)*dpsidx(l2,0)*dpsidx(l,0)*W;
396  }
397  }
398 
399  //Loop over the pressure shape functions
400  for(unsigned l2=0;l2<n_pres;l2++)
401  {
402  //Get the imaginary pressure contribution to the equation
403  local_unknown = P_local_eqn(l2,1);
404  if(local_unknown >= 0)
405  {
406  jacobian(local_eqn,local_unknown) +=
407  psip(l2)*dpsidx(l,0)*W;
408  }
409  }
410  }
411  }
412 
413  //Loop over the pressure shape functions
414  for(unsigned l=0;l<n_pres;l++)
415  {
416  //Real part of the continuity equation
417  local_eqn = P_local_eqn(l,0);
418  if(local_eqn >= 0)
419  {
420  //Loop over the shape functions
421  for(unsigned l2=0;l2<n_node;l2++)
422  {
423  //REAL U conributions
424  local_unknown = U_local_eqn(l2,0);
425  if(local_unknown >= 0)
426  {
427  jacobian(local_eqn,local_unknown) -=
428  A_imag*psi(l2)*psip(l)*W;
429  }
430 
431  //IMAGINARY U contributions
432  local_unknown = U_local_eqn(l2,1);
433  if(local_unknown >= 0)
434  {
435  jacobian(local_eqn,local_unknown) -=
436  A_real*psi(l2)*psip(l)*W;
437  }
438 
439  //REAL W contribution
440  local_unknown = U_local_eqn(l2,2);
441  if(local_unknown >= 0)
442  {
443  jacobian(local_eqn,local_unknown) += dpsidx(l2,0)*psip(l)*W;
444  }
445  }
446  }
447 
448 
449  //Imaginary part of the continuity equation
450  local_eqn = P_local_eqn(l,1);
451  if(local_eqn >= 0)
452  {
453  //Loop over the shape functions
454  for(unsigned l2=0;l2<n_node;l2++)
455  {
456  //REAL U conributions
457  local_unknown = U_local_eqn(l2,0);
458  if(local_unknown >= 0)
459  {
460  jacobian(local_eqn,local_unknown) +=
461  A_real*psi(l2)*psip(l)*W;
462  }
463 
464  //IMAGINARY U contributions
465  local_unknown = U_local_eqn(l2,1);
466  if(local_unknown >= 0)
467  {
468  jacobian(local_eqn,local_unknown) -=
469  A_imag*psi(l2)*psip(l)*W;
470  }
471 
472  //IMAGINARY W contribution
473  local_unknown = U_local_eqn(l2,3);
474  if(local_unknown >= 0)
475  {
476  jacobian(local_eqn,local_unknown) += dpsidx(l2,0)*psip(l)*W;
477  }
478  }
479  }
480  } //End of loop over pressure shape functions
481  }
482  }
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 void pshape(const Vector< double > &s, Shape &psi) const =0
void get_base_flow(const Vector< double > &s, double &U, double &Uz)
Definition: orr_sommerfeld.cc:84
DenseMatrix< int > U_local_eqn
Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC)
Definition: orr_sommerfeld.cc:547
double * Re_pt
Definition: orr_sommerfeld.cc:554
DenseMatrix< int > P_local_eqn
Definition: orr_sommerfeld.cc:550
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const =0
virtual unsigned npres() const =0
Definition: shape.h:278
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.
Definition: shape.h:76
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
double A_imag
Definition: orr_sommerfeld.cc:44
double Re
Definition: orr_sommerfeld.cc:43
@ W
Definition: quadtree.h:63

References Param::A_imag, DIM, i, J, Param::Re, s, w, and oomph::QuadTreeNames::W.

◆ get_base_flow()

template<unsigned DIM>
void OrrSommerfeldEquations< DIM >::get_base_flow ( const Vector< double > &  s,
double U,
double Uz 
)
inline
85  {
86  //Get the coordinate
87  double z = interpolated_x(s,0);
88  //Return the value of the base flow
89  U = (1.0-z*z);
90  Uz = -2.0*z;
91  }
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
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53

References s, and RachelsAdvectionDiffusion::U.

◆ interpolated_p()

template<unsigned DIM>
double OrrSommerfeldEquations< DIM >::interpolated_p ( const Vector< double > &  s,
const unsigned i 
) const
inline

Return FE representation of the pressure.

486  {
487  unsigned n_pres = npres();
488 
489  //Pressure shape functions
490  Shape psip(n_pres);
491 
492  //Find values of shape function
493  this->pshape(s,psip);
494 
495  //Initialise value of p
496  double interpolated_p = 0.0;
497 
498  //Loop over the local nodes and sum
499  for(unsigned l=0;l<n_pres;l++)
500  {
501  interpolated_p+=p(l,i)*psip[l];
502  }
503 
504  return(interpolated_p);
505  }
double interpolated_p(const Vector< double > &s, const unsigned &i) const
Return FE representation of the pressure.
Definition: orr_sommerfeld.cc:485
virtual double p(const unsigned &n, const unsigned &i) const =0

References i, p, and s.

◆ interpolated_u()

template<unsigned DIM>
double OrrSommerfeldEquations< DIM >::interpolated_u ( const Vector< double > &  s,
const unsigned i 
) const
inline

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

511  {
512  unsigned n_node = nnode();
513 
514  //Local shape functions
515  Shape psi(n_node);
516 
517  //Find values of basis function
518  this->shape(s,psi);
519 
520  //Initialise value of u
521  double interpolated_u = 0.0;
522 
523  //Loop over the local nodes and sum
524  for(unsigned l=0;l<n_node;l++)
525  {
526  interpolated_u+=u(l,i)*psi[l];
527  }
528 
529  return(interpolated_u);
530  }
double interpolated_u(const Vector< double > &s, const unsigned &i) const
Return FE representation of function value u(s) at local coordinate s.
Definition: orr_sommerfeld.cc:510
virtual double u(const unsigned &n, const unsigned &i) const =0
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References i, and oomph::OneDimLagrange::shape().

◆ npres()

template<unsigned DIM>
virtual unsigned OrrSommerfeldEquations< DIM >::npres ( ) const
pure virtual

◆ output() [1/2]

template<unsigned DIM>
void OrrSommerfeldEquations< DIM >::output ( ostream &  outfile)
inlinevirtual

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

102  {
103  unsigned nplot=5;
104  output(outfile,nplot);
105  }
void output(ostream &outfile)
Output with default number of plot points.
Definition: orr_sommerfeld.cc:101

References output().

Referenced by QOrrSommerfeldElement< 1, NNODE_1D >::output().

◆ output() [2/2]

template<unsigned DIM>
void OrrSommerfeldEquations< DIM >::output ( ostream &  outfile,
const unsigned nplot 
)
inlinevirtual

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

Reimplemented from oomph::FiniteElement.

110  {
111  //Vector of local coordinates
113 
114  // Tecplot header info
115  outfile << tecplot_zone_string(nplot);
116 
117  // Loop over plot points
118  unsigned num_plot_points=nplot_points(nplot);
119  for (unsigned iplot=0;iplot<num_plot_points;iplot++)
120  {
121 
122  // Get local coordinates of plot point
123  get_s_plot(iplot,nplot,s);
124 
125  for(unsigned i=0;i<DIM;i++)
126  {
127  outfile << interpolated_x(s,i) << " ";
128  }
129 
130  for(unsigned i=0;i<4;i++)
131  {
132  outfile << interpolated_u(s,i) << " ";
133  }
134 
135  outfile << interpolated_p(s,0) << " "
136  << interpolated_p(s,1) << " ";
137  outfile << std::endl;
138  }
139 
140  // Write tecplot footer (e.g. FE connectivity lists)
141  write_tecplot_zone_footer(outfile,nplot);
142  }
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.

◆ p()

template<unsigned DIM>
virtual double OrrSommerfeldEquations< DIM >::p ( const unsigned n,
const unsigned i 
) const
pure virtual

◆ pshape()

template<unsigned DIM>
virtual void OrrSommerfeldEquations< DIM >::pshape ( const Vector< double > &  s,
Shape psi 
) const
pure virtual

◆ re_pt()

template<unsigned DIM>
double* & OrrSommerfeldEquations< DIM >::re_pt ( )
inline
94 {return Re_pt;}

◆ u()

template<unsigned DIM>
virtual double OrrSommerfeldEquations< DIM >::u ( const unsigned n,
const unsigned i 
) const
pure virtual

Access function: Nodal function value at local node n Uses suitably interpolated value for hanging nodes.

Implemented in QOrrSommerfeldElement< 1, NNODE_1D >.

Member Data Documentation

◆ A_Imag_pt

template<unsigned DIM>
double* OrrSommerfeldEquations< DIM >::A_Imag_pt
private

◆ A_Real_pt

template<unsigned DIM>
double* OrrSommerfeldEquations< DIM >::A_Real_pt
private

◆ P_local_eqn

template<unsigned DIM>
DenseMatrix<int> OrrSommerfeldEquations< DIM >::P_local_eqn
protected

◆ Re_pt

template<unsigned DIM>
double* OrrSommerfeldEquations< DIM >::Re_pt
private

◆ U_local_eqn

template<unsigned DIM>
DenseMatrix<int> OrrSommerfeldEquations< DIM >::U_local_eqn
protected

Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC)


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