oomph::PolarNavierStokesTractionElement< ELEMENT > Class Template Reference

#include <polar_fluid_traction_elements.h>

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

Public Member Functions

const doublealpha () const
 Alpha. More...
 
double *& alpha_pt ()
 Pointer to Alpha. More...
 
void set_external_pressure_data (Data *pext_data_pt)
 Function for setting up external pressure. More...
 
const int boundary () const
 Boundary. More...
 
void set_boundary (int bound)
 Function to set boundary. More...
 
const double get_eta () const
 Eta. More...
 
void set_eta (double eta)
 Function to set Eta. More...
 
 PolarNavierStokesTractionElement (FiniteElement *const &element_pt, const int &face_index)
 
 ~PolarNavierStokesTractionElement ()
 Destructor should not delete anything. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 This function returns just the residuals. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 This function returns the residuals and the jacobian. More...
 
void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
void output (std::ostream &outfile)
 Overload the output function. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 Output function: x,y,[z],u,v,[w],p in tecplot format. More...
 
double u (const unsigned &l, const unsigned &i)
 local velocities More...
 
double x (const unsigned &l, const unsigned &i)
 local position More...
 
- Public Member Functions inherited from oomph::FaceElement
 FaceElement ()
 Constructor: Initialise all appropriate member data. More...
 
virtual ~FaceElement ()
 Empty virtual destructor. More...
 
 FaceElement (const FaceElement &)=delete
 Broken copy constructor. More...
 
const unsignedboundary_number_in_bulk_mesh () const
 Broken assignment operator. More...
 
void set_boundary_number_in_bulk_mesh (const unsigned &b)
 Set function for the boundary number in bulk mesh. More...
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double J_eulerian (const Vector< double > &s) const
 
double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
double interpolated_x (const Vector< double > &s, const unsigned &i) const
 
double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 
void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
intnormal_sign ()
 
int normal_sign () const
 
intface_index ()
 
int face_index () const
 
const Vector< double > * tangent_direction_pt () const
 Public access function for the tangent direction pointer. More...
 
void set_tangent_direction (Vector< double > *tangent_direction_pt)
 Set the tangent direction vector. More...
 
void turn_on_warning_for_discontinuous_tangent ()
 
void turn_off_warning_for_discontinuous_tangent ()
 
void continuous_tangent_and_outer_unit_normal (const Vector< double > &s, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void continuous_tangent_and_outer_unit_normal (const unsigned &ipt, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void outer_unit_normal (const Vector< double > &s, Vector< double > &unit_normal) const
 Compute outer unit normal at the specified local coordinate. More...
 
void outer_unit_normal (const unsigned &ipt, Vector< double > &unit_normal) const
 Compute outer unit normal at ipt-th integration point. More...
 
FiniteElement *& bulk_element_pt ()
 Pointer to higher-dimensional "bulk" element. More...
 
FiniteElementbulk_element_pt () const
 Pointer to higher-dimensional "bulk" element (const version) More...
 
CoordinateMappingFctPtface_to_bulk_coordinate_fct_pt ()
 
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt () const
 
BulkCoordinateDerivativesFctPtbulk_coordinate_derivatives_fct_pt ()
 
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt () const
 
Vector< doublelocal_coordinate_in_bulk (const Vector< double > &s) const
 
void get_local_coordinate_in_bulk (const Vector< double > &s, Vector< double > &s_bulk) const
 
void get_ds_bulk_ds_face (const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction) const
 
unsignedbulk_position_type (const unsigned &i)
 
const unsignedbulk_position_type (const unsigned &i) const
 
void bulk_node_number_resize (const unsigned &i)
 Resize the storage for the bulk node numbers. More...
 
unsignedbulk_node_number (const unsigned &n)
 
const unsignedbulk_node_number (const unsigned &n) const
 
void bulk_position_type_resize (const unsigned &i)
 Resize the storage for bulk_position_type to i entries. More...
 
unsignednbulk_value (const unsigned &n)
 
unsigned nbulk_value (const unsigned &n) const
 
void nbulk_value_resize (const unsigned &i)
 
void resize_nodes (Vector< unsigned > &nadditional_data_values)
 
void output_zeta (std::ostream &outfile, const unsigned &nplot)
 Output boundary coordinate zeta. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (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
 

Public Attributes

void(*&)(const double &t, const Vector< double > &x, Vector< double > &result) traction_fct_pt ()
 

Protected Member Functions

virtual int u_local_eqn (const unsigned &n, const unsigned &i)
 
double shape_and_test_at_knot (const unsigned &ipt, Shape &psi, Shape &test) const
 
void get_traction (double time, const Vector< double > &x, Vector< double > &result)
 Function to calculate the traction applied to the fluid. More...
 
void fill_in_generic_residual_contribution (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
 
- Protected Member Functions inherited from oomph::FaceElement
void add_additional_values (const Vector< unsigned > &nadditional_values, const unsigned &id)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
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_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

doubleAlpha_pt
 Pointer to the angle alpha. More...
 
DataPext_data_pt
 Pointer to the Data item that stores the external pressure. More...
 
unsigned External_data_number_of_Pext
 
int Boundary
 
double Eta
 
- Protected Attributes inherited from oomph::FaceElement
unsigned Boundary_number_in_bulk_mesh
 The boundary number in the bulk mesh to which this element is attached. More...
 
FiniteElementBulk_element_pt
 Pointer to the associated higher-dimensional "bulk" element. More...
 
Vector< unsignedBulk_node_number
 
Vector< unsignedNbulk_value
 
Vector< double > * Tangent_direction_pt
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Private Attributes

void(* Traction_fct_pt )(const double &time, const Vector< double > &x, Vector< double > &result)
 Pointer to an imposed traction function. More...
 
unsigned Dim
 The highest dimension of the problem. More...
 

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<class ELEMENT>
class oomph::PolarNavierStokesTractionElement< ELEMENT >

A class for elements that allow the imposition of an applied traction to the Navier–Stokes equations The geometrical information can be read from the FaceGeometery<ELEMENT> class and and thus, we can be generic enough without the need to have a separate equations class

Constructor & Destructor Documentation

◆ PolarNavierStokesTractionElement()

template<class ELEMENT >
oomph::PolarNavierStokesTractionElement< ELEMENT >::PolarNavierStokesTractionElement ( FiniteElement *const &  element_pt,
const int face_index 
)
inline

Constructor, which takes a "bulk" element and the value of the index and its limit

194  : FaceGeometry<ELEMENT>(), FaceElement()
195  {
196  // Attach the geometrical information to the element. N.B. This function
197  // also assigns nbulk_value from the required_nvalue of the bulk element
198  element_pt->build_face_element(face_index, this);
199 
200 #ifdef PARANOID
201  {
202  // Check that the element is not a refineable 3d element
203  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(element_pt);
204  // If it's three-d
205  if (elem_pt->dim() == 3)
206  {
207  // Is it refineable
208  RefineableElement* ref_el_pt =
209  dynamic_cast<RefineableElement*>(elem_pt);
210  if (ref_el_pt != 0)
211  {
212  if (this->has_hanging_nodes())
213  {
214  throw OomphLibError("This flux element will not work correctly "
215  "if nodes are hanging\n",
218  }
219  }
220  }
221  }
222 #endif
223 
224  // Set the body force function pointer to zero
225  Traction_fct_pt = 0;
226 
227  // Set the external pressure pointer to be zero
228  Pext_data_pt = 0;
229 
230  // Set the dimension from the dimension of the first node
231  Dim = this->node_pt(0)->ndim();
232 
233  // Set Eta to one by default
234  Eta = 1.0;
235  }
int & face_index()
Definition: elements.h:4626
FaceElement()
Constructor: Initialise all appropriate member data.
Definition: elements.h:4446
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to an imposed traction function.
Definition: polar_fluid_traction_elements.h:57
unsigned Dim
The highest dimension of the problem.
Definition: polar_fluid_traction_elements.h:62
Data * Pext_data_pt
Pointer to the Data item that stores the external pressure.
Definition: polar_fluid_traction_elements.h:127
double Eta
Definition: polar_fluid_traction_elements.h:140
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FiniteElement::build_face_element(), oomph::PolarNavierStokesTractionElement< ELEMENT >::Dim, oomph::FaceElement::face_index(), oomph::FiniteElement::has_hanging_nodes(), oomph::Node::ndim(), oomph::FiniteElement::node_pt(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::PolarNavierStokesTractionElement< ELEMENT >::Pext_data_pt, and oomph::PolarNavierStokesTractionElement< ELEMENT >::Traction_fct_pt.

◆ ~PolarNavierStokesTractionElement()

template<class ELEMENT >
oomph::PolarNavierStokesTractionElement< ELEMENT >::~PolarNavierStokesTractionElement ( )
inline

Destructor should not delete anything.

238 {}

Member Function Documentation

◆ alpha()

template<class ELEMENT >
const double& oomph::PolarNavierStokesTractionElement< ELEMENT >::alpha ( ) const
inline

Alpha.

145  {
146  return *Alpha_pt;
147  }
double * Alpha_pt
Pointer to the angle alpha.
Definition: polar_fluid_traction_elements.h:124

References oomph::PolarNavierStokesTractionElement< ELEMENT >::Alpha_pt.

◆ alpha_pt()

template<class ELEMENT >
double*& oomph::PolarNavierStokesTractionElement< ELEMENT >::alpha_pt ( )
inline

◆ boundary()

template<class ELEMENT >
const int oomph::PolarNavierStokesTractionElement< ELEMENT >::boundary ( ) const
inline

Boundary.

168  {
169  return Boundary;
170  }
int Boundary
Definition: polar_fluid_traction_elements.h:137

References oomph::PolarNavierStokesTractionElement< ELEMENT >::Boundary.

◆ fill_in_contribution_to_jacobian()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

This function returns the residuals and the jacobian.

Reimplemented from oomph::FiniteElement.

262  {
263  // Call the generic routine with the flag set to 1
265  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
266  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: polar_fluid_traction_elements.h:317

References oomph::GeneralisedElement::Dummy_matrix, and oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution().

◆ fill_in_contribution_to_jacobian_and_mass_matrix()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
)
inlinevirtual

Compute the element's residual Vector and the jacobian matrix Plus the mass matrix especially for eigenvalue problems

Reimplemented from oomph::GeneralisedElement.

274  {
275  // Call the generic routine with the flag set to 2
277  residuals, jacobian, GeneralisedElement::Dummy_matrix, 2);
278  }

References oomph::GeneralisedElement::Dummy_matrix, and oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution().

◆ fill_in_contribution_to_residuals()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

This function returns just the residuals.

Reimplemented from oomph::GeneralisedElement.

250  {
251  // Call the generic residuals function with flag set to 0
252  // using a dummy matrix argument
256  0);
257  }

References oomph::GeneralisedElement::Dummy_matrix, and oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution().

◆ fill_in_generic_residual_contribution()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix,
unsigned  flag 
)
protected

This function returns the residuals for the traction function. flag=1(or 0): do (or don't) compute the Jacobian as well.

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Function that returns the residuals for the imposed traction Navier_Stokes equations

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

/////////////////////////////////////NEW//////////////////////////////////////////

The additional residual for the mass flux (ie. the extra equation for pext) This is an integral equation along the whole boundary It lies outside the loop over shape functions above

/////////////////////////////////////NEW//////////////////////////////////////////

321  {
322  // Find out how many nodes there are
323  unsigned n_node = nnode();
324 
325  // Get continuous time from timestepper of first node
326  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
327 
328  // Set up memory for the shape and test functions
329  Shape psif(n_node), testf(n_node);
330 
331  // Set the value of n_intpt
332  unsigned n_intpt = integral_pt()->nweight();
333 
334  // Get Alpha
335  const double Alpha = alpha();
336 
337  // Storage for external pressure
338  double pext = 0.0;
339 
340  // Get boundary multiplier
341  // This is necessary because the sign of the traction is
342  // dependent on the boundary
343  const int multiplier = boundary();
344 
345  // Get the homotopy parameter (if necessary)
346  const double eta = get_eta();
347 
348  // Integers to store local equation numbers
349  int local_eqn = 0, local_unknown = 0, pext_local_eqn = 0,
350  pext_local_unknown = 0;
351 
353 
354  // Get local equation number of external pressure
355  // Note that if we have not passed an external pressure pointer to this
356  // element (and at the same time added it to the element's external data)
357  // than this will be -1 to indicate that it is not a degree of freedom here
358  if (Pext_data_pt == 0)
359  {
360  pext_local_eqn = -1;
361  }
362  else
363  {
364  // If at a non-zero degree of freedom add in the entry
366 
367  // Get external pressure
368  pext = Pext_data_pt->value(0);
369  }
370 
371  // The local unkown number of pext will be the same
372  pext_local_unknown = pext_local_eqn;
373 
375 
376  // Loop over the integration points
377  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
378  {
379  // Get the integral weight
380  double w = integral_pt()->weight(ipt);
381 
382  // Find the shape and test functions and return the Jacobian
383  // of the mapping
384  double J = shape_and_test_at_knot(ipt, psif, testf);
385 
386  // Premultiply the weights and the Jacobian
387  double W = w * J;
388 
389  // Need to find position to feed into Traction function
390  Vector<double> interpolated_x(Dim);
391  Vector<double> interpolated_u(Dim);
392 
393  // Initialise to zero
394  for (unsigned i = 0; i < Dim; i++)
395  {
396  interpolated_x[i] = 0.0;
397  interpolated_u[i] = 0.0;
398  }
399 
400  // Calculate velocities and derivatives:
401  // Loop over nodes
402  for (unsigned l = 0; l < n_node; l++)
403  {
404  // Loop over directions
405  for (unsigned i = 0; i < Dim; i++)
406  {
407  // Get the nodal value
408  interpolated_u[i] += this->nodal_value(l, i) * psif[l];
409  interpolated_x[i] += this->nodal_position(l, i) * psif[l];
410  }
411  }
412 
413  // Get the user-defined traction terms
414  Vector<double> traction(Dim);
415  get_traction(time, interpolated_x, traction);
416 
417  // Now add to the appropriate equations
418 
419  // Loop over the test functions
420  for (unsigned l = 0; l < n_node; l++)
421  {
422  // Only alter u velocity component
423  {
424  unsigned i = 0;
425  local_eqn = u_local_eqn(l, i);
426  /*IF it's not a boundary condition*/
427  if (local_eqn >= 0)
428  {
429  // Add the user-defined traction terms
430  residuals[local_eqn] -= multiplier * eta * 3.0 *
431  (interpolated_u[i] / interpolated_x[0]) *
432  testf[l] * interpolated_x[0] * Alpha * W;
433 
435 
436  // Plus additional external pressure contribution at inlet
437  // This is zero if we haven't passed a Pext_data_pt to the element
438  residuals[local_eqn] +=
439  pext * testf[l] * interpolated_x[0] * Alpha * W;
440 
442 
443  // CALCULATE THE JACOBIAN
444  if (flag)
445  {
446  // Loop over the velocity shape functions again
447  for (unsigned l2 = 0; l2 < n_node; l2++)
448  {
449  // We only have an i2=0 contribution
450  unsigned i2 = 0;
451  {
452  // If at a non-zero degree of freedom add in the entry
453  local_unknown = u_local_eqn(l2, i2);
454  if (local_unknown >= 0)
455  {
456  // Add contribution to Elemental Matrix
457  jacobian(local_eqn, local_unknown) -=
458  multiplier * eta * 3.0 * (psif[l2] / interpolated_x[0]) *
459  testf[l] * interpolated_x[0] * Alpha * W;
460 
461  } // End of (Jacobian's) if not boundary condition statement
462  } // End of i2=0 section
463  } // End of l2 loop
464 
466  // Add pext's contribution to these residuals
467  // This only needs to be done once hence why it is outside the l2
468  // loop
469  if (pext_local_unknown >= 0)
470  {
471  // Add contribution to Elemental Matrix
472  jacobian(local_eqn, pext_local_unknown) +=
473  testf[l] * interpolated_x[0] * Alpha * W;
474  }
476 
477  } /*End of Jacobian calculation*/
478 
479  } // end of if not boundary condition statement
480  } // End of i=0 section
481 
482  } // End of loop over shape functions
483 
485 
490  {
491  /*IF it's not a boundary condition*/
492  if (pext_local_eqn >= 0)
493  {
494  // Add the user-defined traction terms
495  residuals[pext_local_eqn] +=
496  interpolated_u[0] * interpolated_x[0] * Alpha * W;
497 
498  // No longer necessary due to my FluxCosntraint element
499  // Now take off a fraction of the desired mass flux
500  // Divided by number of elements and number of int points in each
501  // HACK
502  // residuals[pext_local_eqn] -= (1.0/(30.*3.));
503  // HACK
504 
505  // CALCULATE THE JACOBIAN
506  if (flag)
507  {
508  // Loop over the velocity shape functions again
509  for (unsigned l2 = 0; l2 < n_node; l2++)
510  {
511  // We only have an i2=0 contribution
512  unsigned i2 = 0;
513  {
514  // If at a non-zero degree of freedom add in the entry
515  local_unknown = u_local_eqn(l2, i2);
516  if (local_unknown >= 0)
517  {
518  // Add contribution to Elemental Matrix
519  jacobian(pext_local_eqn, local_unknown) +=
520  psif[l2] * interpolated_x[0] * Alpha * W;
521 
522  } // End of (Jacobian's) if not boundary condition statement
523  } // End of i2=0 section
524 
525  } // End of l2 loop
526  } /*End of Jacobian calculation*/
527 
528  } // end of if not boundary condition statement
529  } // End of additional residual for the mass flux
530 
532 
533  } // End of loop over integration points
534 
535  } // End of fill_in_generic_residual_contribution
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
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
double value(const unsigned &i) const
Definition: nodes.h:293
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
int external_local_eqn(const unsigned &i, const unsigned &j)
Definition: elements.h:311
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.
const int boundary() const
Boundary.
Definition: polar_fluid_traction_elements.h:167
void get_traction(double time, const Vector< double > &x, Vector< double > &result)
Function to calculate the traction applied to the fluid.
Definition: polar_fluid_traction_elements.h:96
unsigned External_data_number_of_Pext
Definition: polar_fluid_traction_elements.h:131
const double get_eta() const
Eta.
Definition: polar_fluid_traction_elements.h:179
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: polar_fluid_traction_elements.h:77
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Definition: polar_fluid_traction_elements.h:70
const double & alpha() const
Alpha.
Definition: polar_fluid_traction_elements.h:144
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
double multiplier(const Vector< double > &xi)
Definition: disk_oscillation.cc:63
double Alpha
Parameter for steepness of step.
Definition: two_d_adv_diff_adapt.cc:53
double eta
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:45
@ W
Definition: quadtree.h:63

References alpha, TanhSolnForAdvectionDiffusion::Alpha, Global_Variables::Dim, TestSoln::eta, i, J, Global_Physical_Variables::multiplier(), Global_Physical_Variables::Pext_data_pt, oomph::Data::value(), w, and oomph::QuadTreeNames::W.

Referenced by oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_jacobian_and_mass_matrix(), and oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_residuals().

◆ get_eta()

template<class ELEMENT >
const double oomph::PolarNavierStokesTractionElement< ELEMENT >::get_eta ( ) const
inline

Eta.

180  {
181  return Eta;
182  }

References oomph::PolarNavierStokesTractionElement< ELEMENT >::Eta.

◆ get_traction()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::get_traction ( double  time,
const Vector< double > &  x,
Vector< double > &  result 
)
inlineprotected

Function to calculate the traction applied to the fluid.

99  {
100  // If the function pointer is zero return zero
101  if (Traction_fct_pt == 0)
102  {
103  // Loop over dimensions and set body forces to zero
104  for (unsigned i = 0; i < Dim; i++)
105  {
106  result[i] = 0.0;
107  }
108  }
109  // Otherwise call the function
110  else
111  {
112  (*Traction_fct_pt)(time, x, result);
113  }
114  }
double x(const unsigned &l, const unsigned &i)
local position
Definition: polar_fluid_traction_elements.h:299

References oomph::PolarNavierStokesTractionElement< ELEMENT >::Dim, i, oomph::PolarNavierStokesTractionElement< ELEMENT >::Traction_fct_pt, and oomph::PolarNavierStokesTractionElement< ELEMENT >::x().

◆ output() [1/2]

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::output ( std::ostream &  outfile)
inlinevirtual

Overload the output function.

Reimplemented from oomph::FiniteElement.

282  {
283  FiniteElement::output(outfile);
284  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/2]

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned nplot 
)
inlinevirtual

Output function: x,y,[z],u,v,[w],p in tecplot format.

Reimplemented from oomph::FiniteElement.

288  {
289  FiniteElement::output(outfile, nplot);
290  }

References oomph::FiniteElement::output().

◆ set_boundary()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::set_boundary ( int  bound)
inline

Function to set boundary.

174  {
175  Boundary = bound;
176  }

References oomph::PolarNavierStokesTractionElement< ELEMENT >::Boundary.

Referenced by oomph::jh_mesh< ELEMENT >::make_traction_elements().

◆ set_eta()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::set_eta ( double  eta)
inline

Function to set Eta.

186  {
187  Eta = eta;
188  }

References TestSoln::eta, and oomph::PolarNavierStokesTractionElement< ELEMENT >::Eta.

Referenced by oomph::jh_mesh< ELEMENT >::make_traction_elements().

◆ set_external_pressure_data()

template<class ELEMENT >
void oomph::PolarNavierStokesTractionElement< ELEMENT >::set_external_pressure_data ( Data pext_data_pt)
inline

Function for setting up external pressure.

157  {
158  // Set external pressure pointer
159  Pext_data_pt = pext_data_pt;
160 
161  // Add to the element's external data so it gets included
162  // in the black-box local equation numbering scheme
163  External_data_number_of_Pext = this->add_external_data(pext_data_pt);
164  }
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307

References oomph::GeneralisedElement::add_external_data(), oomph::PolarNavierStokesTractionElement< ELEMENT >::External_data_number_of_Pext, and oomph::PolarNavierStokesTractionElement< ELEMENT >::Pext_data_pt.

◆ shape_and_test_at_knot()

template<class ELEMENT >
double oomph::PolarNavierStokesTractionElement< ELEMENT >::shape_and_test_at_knot ( const unsigned ipt,
Shape psi,
Shape test 
) const
inlineprotected

Function to compute the shape and test functions and to return the Jacobian of mapping

80  {
81  // Find number of nodes
82  unsigned n_node = nnode();
83  // Calculate the shape functions
84  shape_at_knot(ipt, psi);
85  // Set the test functions to be the same as the shape functions
86  for (unsigned i = 0; i < n_node; i++)
87  {
88  test[i] = psi[i];
89  }
90  // Return the value of the jacobian
91  return J_eulerian_at_knot(ipt);
92  }
double J_eulerian_at_knot(const unsigned &ipt) const
Definition: elements.cc:5328
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
Definition: indexed_view.cpp:20

References i, oomph::FaceElement::J_eulerian_at_knot(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::shape_at_knot().

◆ u()

template<class ELEMENT >
double oomph::PolarNavierStokesTractionElement< ELEMENT >::u ( const unsigned l,
const unsigned i 
)
inline

local velocities

294  {
295  return nodal_value(l, i);
296  }

References i, and oomph::FiniteElement::nodal_value().

◆ u_local_eqn()

template<class ELEMENT >
virtual int oomph::PolarNavierStokesTractionElement< ELEMENT >::u_local_eqn ( const unsigned n,
const unsigned i 
)
inlineprotectedvirtual

Access function that returns the local equation numbers for velocity components. u_local_eqn(n,i) = local equation number or < 0 if pinned. The default is to asssume that n is the local node number and the i-th velocity component is the i-th unknown stored at the node.

71  {
72  return nodal_local_eqn(n, i);
73  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432

References i, n, and oomph::FiniteElement::nodal_local_eqn().

◆ x()

template<class ELEMENT >
double oomph::PolarNavierStokesTractionElement< ELEMENT >::x ( const unsigned l,
const unsigned i 
)
inline

local position

300  {
301  return nodal_position(l, i);
302  }

References i, and oomph::FiniteElement::nodal_position().

Referenced by oomph::PolarNavierStokesTractionElement< ELEMENT >::get_traction().

Member Data Documentation

◆ Alpha_pt

◆ Boundary

◆ Dim

◆ Eta

◆ External_data_number_of_Pext

template<class ELEMENT >
unsigned oomph::PolarNavierStokesTractionElement< ELEMENT >::External_data_number_of_Pext
protected

The Data that contains the traded pressure is stored as external Data for the element. Which external Data item is it?

Referenced by oomph::PolarNavierStokesTractionElement< ELEMENT >::set_external_pressure_data().

◆ Pext_data_pt

template<class ELEMENT >
Data* oomph::PolarNavierStokesTractionElement< ELEMENT >::Pext_data_pt
protected

◆ traction_fct_pt

template<class ELEMENT >
void(*&)(const double& t, const Vector<double>& x, Vector<double>& result) oomph::PolarNavierStokesTractionElement< ELEMENT >::traction_fct_pt()
inline
244  {
245  return Traction_fct_pt;
246  }

Referenced by oomph::jh_mesh< ELEMENT >::make_traction_elements().

◆ Traction_fct_pt

template<class ELEMENT >
void(* oomph::PolarNavierStokesTractionElement< ELEMENT >::Traction_fct_pt) (const double &time, const Vector< double > &x, Vector< double > &result)
private

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