![]() |
|
#include <poisson_elements.h>
Public Types | |
typedef void(* | PoissonSourceFctPt) (const Vector< double > &x, double &f) |
typedef void(* | PoissonSourceFctGradientPt) (const Vector< double > &x, Vector< double > &gradient) |
![]() | |
typedef void(* | SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &) |
typedef void(* | UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &) |
Public Member Functions | |
PoissonEquations () | |
Constructor (must initialise the Source_fct_pt to null) More... | |
PoissonEquations (const PoissonEquations &dummy)=delete | |
Broken copy constructor. More... | |
void | operator= (const PoissonEquations &)=delete |
Broken assignment operator. More... | |
virtual unsigned | u_index_poisson () const |
void | point_output_data (const Vector< double > &s, Vector< double > &data) |
unsigned | nscalar_paraview () const |
void | scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const |
std::string | scalar_name_paraview (const unsigned &i) const |
void | scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const |
void | output (std::ostream &outfile) |
Output with default number of plot points. More... | |
void | output (std::ostream &outfile, const unsigned &n_plot) |
void | output (FILE *file_pt) |
C_style output with default number of plot points. More... | |
void | output (FILE *file_pt, const unsigned &n_plot) |
void | output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) |
virtual void | output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) |
void | compute_norm (double &norm) |
Compute norm of solution: square of the L2 norm. More... | |
void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
Get error against and norm of exact solution. More... | |
void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm) |
Dummy, time dependent error checker. More... | |
PoissonSourceFctPt & | source_fct_pt () |
Access function: Pointer to source function. More... | |
PoissonSourceFctPt | source_fct_pt () const |
Access function: Pointer to source function. Const version. More... | |
PoissonSourceFctGradientPt & | source_fct_gradient_pt () |
Access function: Pointer to gradient of source function. More... | |
PoissonSourceFctGradientPt | source_fct_gradient_pt () const |
Access function: Pointer to gradient source function. Const version. More... | |
virtual void | get_source_poisson (const unsigned &ipt, const Vector< double > &x, double &source) const |
virtual void | get_source_gradient_poisson (const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient) const |
void | get_flux (const Vector< double > &s, Vector< double > &flux) const |
Get flux: flux[i] = du/dx_i. More... | |
void | get_dflux_dnodal_u (const Vector< double > &s, Vector< Vector< double >> &dflux_dnodal_u) const |
void | fill_in_contribution_to_residuals (Vector< double > &residuals) |
Add the element's contribution to its residual vector (wrapper) More... | |
void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
virtual double | interpolated_u_poisson (const Vector< double > &s) const |
virtual void | get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates) |
unsigned | self_test () |
Self-test: Return 0 for OK. More... | |
![]() | |
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) |
MacroElement * | macro_elem_pt () |
Access function to pointer to macro element. More... | |
void | get_x (const Vector< double > &s, Vector< double > &x) const |
void | get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) |
virtual void | get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const |
virtual void | get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x) |
virtual void | set_integration_scheme (Integral *const &integral_pt) |
Set the spatial integration scheme. More... | |
Integral *const & | integral_pt () const |
Return the pointer to the integration scheme (const version) More... | |
virtual void | shape (const Vector< double > &s, Shape &psi) const =0 |
virtual void | shape_at_knot (const unsigned &ipt, Shape &psi) const |
virtual void | dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const |
virtual void | dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const |
virtual void | d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const |
virtual void | d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const |
virtual double | J_eulerian (const Vector< double > &s) const |
virtual double | J_eulerian_at_knot (const unsigned &ipt) const |
void | check_J_eulerian_at_knots (bool &passed) const |
void | check_jacobian (const double &jacobian) const |
double | dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const |
virtual double | dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const |
virtual double | dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const |
double | d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const |
virtual double | d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const |
virtual void | assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt) |
virtual void | describe_local_dofs (std::ostream &out, const std::string ¤t_string) const |
virtual void | describe_nodal_local_dofs (std::ostream &out, const std::string ¤t_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 | 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 Node * | vertex_node_pt (const unsigned &j) const |
virtual Node * | construct_node (const unsigned &n) |
virtual Node * | construct_node (const unsigned &n, TimeStepper *const &time_stepper_pt) |
virtual Node * | construct_boundary_node (const unsigned &n) |
virtual Node * | construct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt) |
int | get_node_number (Node *const &node_pt) const |
virtual Node * | get_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 |
Data * | geom_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 |
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 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 void | output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const |
virtual void | output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const |
Output a time-dependent exact solution over the element. More... | |
virtual void | get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const |
virtual std::string | tecplot_zone_string (const unsigned &nplot) const |
virtual void | write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const |
virtual void | write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const |
virtual unsigned | nplot_points (const unsigned &nplot) const |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
Calculate the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm) |
Calculate the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error) |
void | integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral) |
Integrate Vector-valued function over element. More... | |
void | integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral) |
Integrate Vector-valued time-dep function over element. More... | |
virtual void | build_face_element (const int &face_index, FaceElement *face_element_pt) |
virtual unsigned | get_bulk_node_number (const int &face_index, const unsigned &i) const |
virtual int | face_outer_unit_normal_sign (const int &face_index) const |
Get the sign of the outer unit normal on the face given by face_index. More... | |
virtual unsigned | nnode_on_face () const |
void | face_node_number_error_check (const unsigned &i) const |
Range check for face node numbers. More... | |
virtual CoordinateMappingFctPt | face_to_bulk_coordinate_fct_pt (const int &face_index) const |
virtual BulkCoordinateDerivativesFctPt | bulk_coordinate_derivatives_fct_pt (const int &face_index) const |
![]() | |
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 ¤t_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 ¶meter_pt, Vector< double > &dres_dparam) |
virtual void | get_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
virtual void | get_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam) |
virtual void | get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
virtual void | get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product) |
virtual void | get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector) |
virtual void | compute_norm (Vector< double > &norm) |
virtual unsigned | ndof_types () const |
virtual void | get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const |
![]() | |
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 () |
TimeStepper * | time_stepper_pt () const |
virtual void | position (const double &t, const Vector< double > &zeta, Vector< double > &r) const |
virtual void | dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const |
virtual void | d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const |
virtual void | d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const |
Protected Member Functions | |
virtual double | dshape_and_dtest_eulerian_poisson (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0 |
virtual double | dshape_and_dtest_eulerian_at_knot_poisson (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0 |
virtual double | dshape_and_dtest_eulerian_at_knot_poisson (const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const =0 |
virtual void | fill_in_generic_residual_contribution_poisson (Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag) |
![]() | |
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 |
![]() | |
unsigned | add_internal_data (Data *const &data_pt, const bool &fd=true) |
bool | internal_data_fd (const unsigned &i) const |
void | exclude_internal_data_fd (const unsigned &i) |
void | include_internal_data_fd (const unsigned &i) |
void | clear_global_eqn_numbers () |
void | add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt) |
virtual void | assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt) |
virtual void | assign_additional_local_eqn_numbers () |
int | internal_local_eqn (const unsigned &i, const unsigned &j) const |
int | external_local_eqn (const unsigned &i, const unsigned &j) |
void | fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
virtual void | update_before_internal_fd () |
virtual void | reset_after_internal_fd () |
virtual void | update_in_internal_fd (const unsigned &i) |
virtual void | reset_in_internal_fd (const unsigned &i) |
virtual void | update_before_external_fd () |
virtual void | reset_after_external_fd () |
virtual void | update_in_external_fd (const unsigned &i) |
virtual void | reset_in_external_fd (const unsigned &i) |
virtual void | fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix) |
virtual void | fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
virtual void | fill_in_contribution_to_dresiduals_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam) |
virtual void | fill_in_contribution_to_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
virtual void | fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_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 | |
PoissonSourceFctPt | Source_fct_pt |
Pointer to source function: More... | |
PoissonSourceFctGradientPt | Source_fct_gradient_pt |
Pointer to gradient of source function. More... | |
![]() | |
MacroElement * | Macro_elem_pt |
Pointer to the element's macro element (NULL by default) More... | |
![]() | |
unsigned | NLagrangian |
Number of Lagrangian (intrinsic) coordinates. More... | |
unsigned | Ndim |
Number of Eulerian coordinates. More... | |
TimeStepper * | Geom_object_time_stepper_pt |
Additional Inherited Members | |
![]() | |
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 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 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 DenseMatrix< double > | Dummy_matrix |
static std::deque< double * > | Dof_pt_deque |
A class for all isoparametric elements that solve the Poisson equations.
\[ \frac{\partial^2 u}{\partial x_i^2} = f(x_j) \]
This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class.
typedef void(* oomph::PoissonEquations< DIM >::PoissonSourceFctGradientPt) (const Vector< double > &x, Vector< double > &gradient) |
Function pointer to gradient of source function fct(x,g(x)) – x is a Vector!
typedef void(* oomph::PoissonEquations< DIM >::PoissonSourceFctPt) (const Vector< double > &x, double &f) |
Function pointer to source function fct(x,f(x)) – x is a Vector!
|
inline |
Constructor (must initialise the Source_fct_pt to null)
|
delete |
Broken copy constructor.
|
virtual |
Get error against and norm of exact solution.
Validate against exact solution
Solution is provided via function pointer. Plot error at a given number of plot points.
Reimplemented from oomph::FiniteElement.
References DIM, calibrate::error, s, and plotDoE::x.
|
inlinevirtual |
Dummy, time dependent error checker.
Reimplemented from oomph::FiniteElement.
References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
|
virtual |
Compute norm of solution: square of the L2 norm.
Compute norm of the solution.
Reimplemented from oomph::GeneralisedElement.
|
protectedpure virtual |
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of mapping (J). Also compute derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
Implemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
|
protectedpure virtual |
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of mapping
Implemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
|
protectedpure virtual |
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping
Implemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
|
inlinevirtual |
Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)
Reimplemented from oomph::FiniteElement.
References oomph::PoissonEquations< DIM >::fill_in_generic_residual_contribution_poisson().
|
inlinevirtual |
Add the element's contribution to its residual vector (wrapper)
Reimplemented from oomph::GeneralisedElement.
References oomph::GeneralisedElement::Dummy_matrix, and oomph::PoissonEquations< DIM >::fill_in_generic_residual_contribution_poisson().
|
protectedvirtual |
Compute element residual Vector only (if flag=and/or element Jacobian matrix
Compute element residual Vector and/or element Jacobian matrix
flag=1: compute both flag=0: compute only residual Vector
Pure version without hanging nodes
Reimplemented in oomph::RefineablePoissonEquations< DIM >.
References DIM, i, J, j, k, TestProblem::source(), Eigen::test, w, and oomph::QuadTreeNames::W.
Referenced by oomph::PoissonEquations< DIM >::fill_in_contribution_to_jacobian(), and oomph::PoissonEquations< DIM >::fill_in_contribution_to_residuals().
|
inline |
Get derivative of flux w.r.t. to nodal values: dflux_dnodal_u[i][j] = d ( du/dx_i ) / dU_j
References DIM, oomph::FiniteElement::dshape_eulerian(), i, j, oomph::FiniteElement::nnode(), and s.
|
virtual |
Compute derivatives of elemental residual vector with respect to nodal coordinates. Overwrites default implementation in FiniteElement base class. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij}
Compute derivatives of elemental residual vector with respect to nodal coordinates (fully analytical). dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij} Overloads the FD-based version in the FE base class.
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::RefineablePoissonEquations< DIM >, and ModalPRefineableQPoissonElement< DIM >.
References DIM, J, TestProblem::source(), Eigen::test, and w.
|
inline |
Get flux: flux[i] = du/dx_i.
References DIM, oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, and oomph::PoissonEquations< DIM >::u_index_poisson().
Referenced by oomph::PRefineableQPoissonElement< DIM >::compute_energy_error(), oomph::RefineablePoissonEquations< DIM >::get_Z2_flux(), and oomph::TPoissonElement< DIM, NNODE_1D >::get_Z2_flux().
|
inlinevirtual |
Get gradient of source term at (Eulerian) position x. This function is virtual to allow overloading in multi-physics problems where the strength of the source function might be determined by another system of equations. Computed via function pointer (if set) or by finite differencing (default)
References oomph::GeneralisedElement::Default_fd_jacobian_step, DIM, oomph::PoissonEquations< DIM >::get_source_poisson(), i, TestProblem::source(), oomph::PoissonEquations< DIM >::Source_fct_gradient_pt, and plotDoE::x.
|
inlinevirtual |
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-physics problems where the strength of the source function might be determined by another system of equations.
References TestProblem::source(), oomph::PoissonEquations< DIM >::Source_fct_pt, and plotDoE::x.
Referenced by oomph::PoissonEquations< DIM >::get_source_gradient_poisson().
|
inlinevirtual |
Return FE representation of function value u_poisson(s) at local coordinate s
References oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::FiniteElement::shape(), and oomph::PoissonEquations< DIM >::u_index_poisson().
Referenced by oomph::PoissonEquations< DIM >::point_output_data(), and oomph::PoissonEquations< DIM >::scalar_value_paraview().
|
inlinevirtual |
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
Reimplemented from oomph::FiniteElement.
|
delete |
Broken assignment operator.
|
inlinevirtual |
C_style output with default number of plot points.
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
References oomph::PoissonEquations< DIM >::output().
|
virtual |
C-style output FE representation of soln: x,y,u or x,y,z,u at n_plot^DIM plot points
C-style output function:
x,y,u or x,y,z,u
nplot points in each coordinate direction
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
|
inlinevirtual |
Output with default number of plot points.
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
Referenced by oomph::PoissonEquations< DIM >::output(), oomph::QPoissonElement< DIM, NNODE_1D >::output(), oomph::QSpectralPoissonElement< DIM, NNODE_1D >::output(), and oomph::TPoissonElement< DIM, NNODE_1D >::output().
|
virtual |
Output FE representation of soln: x,y,u or x,y,z,u at n_plot^DIM plot points
Output function:
x,y,u or x,y,z,u
nplot points in each coordinate direction
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
|
inlinevirtual |
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points (dummy time-dependent version to keep intel compiler happy)
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
|
virtual |
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points
Output exact solution
Solution is provided via function pointer. Plot at a given number of plot points.
x,y,u_exact or x,y,z,u_exact
Reimplemented from oomph::FiniteElement.
Reimplemented in oomph::TPoissonElement< DIM, NNODE_1D >, oomph::QSpectralPoissonElement< DIM, NNODE_1D >, oomph::QPoissonElement< DIM, NNODE_1D >, and oomph::QPoissonElement< DIM, 2 >.
References DIM, ProblemParameters::exact_soln(), i, s, and plotDoE::x.
Referenced by oomph::QPoissonElement< DIM, NNODE_1D >::output_fct(), oomph::QSpectralPoissonElement< DIM, NNODE_1D >::output_fct(), and oomph::TPoissonElement< DIM, NNODE_1D >::output_fct().
|
inlinevirtual |
Output solution in data vector at local cordinates s: x,y [,z], u
Reimplemented from oomph::FiniteElement.
References data, oomph::FiniteElement::dim(), i, oomph::PoissonEquations< DIM >::interpolated_u_poisson(), oomph::FiniteElement::interpolated_x(), and s.
|
inlinevirtual |
Name of the i-th scalar field. Default implementation returns V1 for the first one, V2 for the second etc. Can (should!) be overloaded with more meaningful names in specific elements.
Reimplemented from oomph::FiniteElement.
References i, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
|
inlinevirtual |
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specific element type.
Reimplemented from oomph::FiniteElement.
References DIM, ProblemParameters::exact_soln(), oomph::FiniteElement::get_s_plot(), i, oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, s, and plotDoE::x.
|
inlinevirtual |
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specific element type.
Reimplemented from oomph::FiniteElement.
References DIM, oomph::FiniteElement::get_s_plot(), i, oomph::PoissonEquations< DIM >::interpolated_u_poisson(), j, oomph::FiniteElement::nplot_points_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and s.
|
virtual |
Self-test: Return 0 for OK.
Reimplemented from oomph::FiniteElement.
References oomph::FiniteElement::self_test().
|
inline |
Access function: Pointer to gradient of source function.
References oomph::PoissonEquations< DIM >::Source_fct_gradient_pt.
|
inline |
Access function: Pointer to gradient source function. Const version.
References oomph::PoissonEquations< DIM >::Source_fct_gradient_pt.
|
inline |
Access function: Pointer to source function.
References oomph::PoissonEquations< DIM >::Source_fct_pt.
Referenced by oomph::RefineablePoissonEquations< DIM >::further_build().
|
inline |
Access function: Pointer to source function. Const version.
References oomph::PoissonEquations< DIM >::Source_fct_pt.
|
inlinevirtual |
Return the index at which the unknown value is stored. The default value, 0, is appropriate for single-physics problems, when there is only one variable, the value that satisfies the poisson equation. In derived multi-physics elements, this function should be overloaded to reflect the chosen storage scheme. Note that these equations require that the unknown is always stored at the same index at each node.
Referenced by oomph::PoissonEquations< DIM >::get_flux(), oomph::RefineablePoissonEquations< DIM >::get_interpolated_values(), oomph::PoissonEquations< DIM >::interpolated_u_poisson(), and oomph::PoissonFluxElement< ELEMENT >::PoissonFluxElement().
|
protected |
Pointer to gradient of source function.
Referenced by oomph::PoissonEquations< DIM >::get_source_gradient_poisson(), and oomph::PoissonEquations< DIM >::source_fct_gradient_pt().
|
protected |
Pointer to source function:
Referenced by oomph::RefineablePoissonEquations< DIM >::further_build(), oomph::PoissonEquations< DIM >::get_source_poisson(), and oomph::PoissonEquations< DIM >::source_fct_pt().