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

Public Types

typedef void(* PoissonSourceFctPt) (const Vector< double > &x, double &f)
 
typedef void(* PoissonDiffFctPt) (const Vector< double > &x, double &f)
 Function pointer to a diffusion function fxt(x,f(x)) More...
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 

Public Member Functions

 ModalPoissonEquations ()
 Constructor (must initialise the Source_fct_pt to null) More...
 
virtual double u (const unsigned &n) const =0
 
virtual unsigned nbasis () const =0
 Number of basis functions. More...
 
virtual void basis (const Vector< double > &s, Shape &basis) const =0
 
void output (ostream &outfile)
 Output with default number of plot points. More...
 
void output (ostream &outfile, const unsigned &nplot)
 
void output_fct (ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points. More...
 
virtual void output_fct (ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
void compute_error (ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Dummy, time dependent error checker. More...
 
PoissonSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
PoissonSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
virtual void get_source_poisson (const Vector< double > &x, double &source) const
 
PoissonDiffFctPtdiff_fct_pt ()
 Access function: Pointer to diffusivity function. More...
 
PoissonDiffFctPt diff_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
virtual void get_diff (const Vector< double > &x, double &diff) const
 
void get_flux (const Vector< double > &s, Vector< double > &flux, unsigned p_order) const
 Get flux: flux[i] = du/dx_i. More...
 
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)
 
double interpolated_u (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output (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 get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

virtual double dshape_dbasis_and_dtest_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &basis, DShape &dbasisdx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_dbasis_and_dtest_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &basis, DShape &dbasisdx, Shape &test, DShape &dtestdx) const =0
 
virtual void add_generic_residual_contribution (Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 

Protected Attributes

PoissonSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
PoissonDiffFctPt Diff_fct_pt
 Pointer to the diffusivity function. More...
 
Vector< intLocal_eqn
 Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC) More...
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Additional Inherited Members

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

Detailed Description

template<unsigned DIM>
class ModalPoissonEquations< DIM >

A class for all 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.

Member Typedef Documentation

◆ PoissonDiffFctPt

template<unsigned DIM>
typedef void(* ModalPoissonEquations< DIM >::PoissonDiffFctPt) (const Vector< double > &x, double &f)

Function pointer to a diffusion function fxt(x,f(x))

◆ PoissonSourceFctPt

template<unsigned DIM>
typedef void(* ModalPoissonEquations< DIM >::PoissonSourceFctPt) (const Vector< double > &x, double &f)

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

Constructor & Destructor Documentation

◆ ModalPoissonEquations()

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

Constructor (must initialise the Source_fct_pt to null)

101 : Source_fct_pt(0), Diff_fct_pt(0) {}
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
Definition: one_d_poisson_hp_adapt.cc:680
PoissonDiffFctPt Diff_fct_pt
Pointer to the diffusivity function.
Definition: one_d_poisson_hp_adapt.cc:683

Member Function Documentation

◆ add_generic_residual_contribution()

template<unsigned DIM>
virtual void ModalPoissonEquations< DIM >::add_generic_residual_contribution ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
unsigned  flag 
)
inlineprotectedvirtual

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

558  {
559  //Find out how many nodes there are
560  unsigned n_node = nnode();
561  //Find out the total number of basis and test functions there are
562  unsigned n_basis = nbasis();
563 
564  //Set up memory for the shape and test functions
565  Shape psi(n_node), basis(n_basis), test(n_basis);
566  DShape dpsidx(n_node,DIM), dbasisdx(n_basis,DIM), dtestdx(n_basis,DIM);
567 
568  //Set the value of n_intpt
569  unsigned n_intpt = integral_pt()->nweight();
570 
571  //Set the Vector to hold local coordinates
573 
574  //Integers to store the local equation and unknown numbers
575  int local_eqn=0, local_unknown=0;
576 
577  //Loop over the integration points
578  for(unsigned ipt=0;ipt<n_intpt;ipt++)
579  {
580  //Assign values of s
581  for(unsigned i=0;i<DIM;i++) s[i] = integral_pt()->knot(ipt,i);
582 
583  //Get the integral weight
584  double w = integral_pt()->weight(ipt);
585 
586  //Call the derivatives of the shape and test functions
588  psi,dpsidx,
589  basis,dbasisdx,
590  test,dtestdx);
591 
592  //Premultiply the weights and the Jacobian
593  double W = w*J;
594 
595  //Calculate local values of the pressure and velocity components
596  //Allocate and initialise to zero
597  double interpolated_u=0.0;
599  Vector<double> interpolated_dudx(DIM,0.0);
600 
601  //Calculate function value and derivatives:
602  //-----------------------------------------
603  // Loop over nodes
604  for(unsigned l=0;l<n_node;l++)
605  {
606  // Loop over directions
607  for(unsigned j=0;j<DIM;j++)
608  {
609  interpolated_x[j] += nodal_position(l,j)*psi[l];
610  }
611  }
612 
613 
614  for(unsigned l=0;l<n_basis;l++)
615  {
616  interpolated_u += u(l)*basis[l];
617  for(unsigned j=0;j<DIM;j++)
618  {
619  interpolated_dudx[j] += u(l)*dbasisdx(l,j);
620  }
621  }
622 
623  //Get source function
624  //-------------------
625  double source;
627 
628  //Get the diffusivity
629  double diff;
630  get_diff(interpolated_x,diff);
631 
632  // Assemble residuals and Jacobian
633  //--------------------------------
634 
635  // Loop over the test functions
636  for(unsigned l=0;l<n_basis;l++)
637  {
638  //Get the local equation
639  local_eqn = Local_eqn[l];
640  //IF it's not a boundary condition
641  if(local_eqn >= 0)
642  {
643  // Add body force/source term here
644  residuals[local_eqn] += source*test[l]*W;
645 
646  // The Poisson bit itself
647  for(unsigned k=0;k<DIM;k++)
648  {
649  residuals[local_eqn] += diff*interpolated_dudx[k]*dtestdx(l,k)*W;
650  }
651 
652  // Calculate the jacobian
653  //-----------------------
654  if(flag)
655  {
656  //Loop over the velocity shape functions again
657  for(unsigned l2=0;l2<n_basis;l2++)
658  {
659  local_unknown = Local_eqn[l2];
660  //If at a non-zero degree of freedom add in the entry
661  if(local_unknown >= 0)
662  {
663  //Add contribution to Elemental Matrix
664  for(unsigned i=0;i<DIM;i++)
665  {
666  jacobian(local_eqn,local_unknown)
667  += diff*dbasisdx(l2,i)*dtestdx(l,i)*W;
668  }
669  }
670  }
671  }
672  }
673  }
674 
675  } // End of loop over integration points
676 
677  }
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
double interpolated_u(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: one_d_poisson_hp_adapt.cc:504
Vector< int > Local_eqn
Array to hold local eqn numbers: Local_eqn[n] (=-1 for BC)
Definition: one_d_poisson_hp_adapt.cc:686
virtual void basis(const Vector< double > &s, Shape &basis) const =0
virtual unsigned nbasis() const =0
Number of basis functions.
virtual double dshape_dbasis_and_dtest_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &basis, DShape &dbasisdx, Shape &test, DShape &dtestdx) const =0
virtual double u(const unsigned &n) const =0
virtual void get_source_poisson(const Vector< double > &x, double &source) const
Definition: one_d_poisson_hp_adapt.cc:405
virtual void get_diff(const Vector< double > &x, double &diff) const
Definition: one_d_poisson_hp_adapt.cc:434
Definition: shape.h:278
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
Definition: shape.h:76
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
#define DIM
Definition: linearised_navier_stokes_elements.h:44
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
@ W
Definition: quadtree.h:63
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

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

◆ basis()

template<unsigned DIM>
virtual void ModalPoissonEquations< DIM >::basis ( const Vector< double > &  s,
Shape basis 
) const
pure virtual

◆ compute_error() [1/2]

template<unsigned DIM>
void ModalPoissonEquations< DIM >::compute_error ( ostream &  outfile,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt,
double error,
double norm 
)
inline

Get error against and norm of exact solution.

222  {
223 
224  // Initialise
225  error=0.0;
226  norm=0.0;
227 
228  //Vector of local coordinates
230 
231  // Vector for coordintes
233 
234  //Set the value of n_intpt
235  unsigned n_intpt = integral_pt()->nweight();
236 
237  // Setup output structure: Conversion is fishy but it's only output...
238  unsigned nplot;
239  if (DIM==1)
240  {
241  nplot=n_intpt;
242  }
243  else
244  {
245  nplot=unsigned(pow(n_intpt,1.0/double(DIM)));
246  }
247 
248  // Tecplot header info
249  outfile << tecplot_zone_string(nplot);
250 
251  // Exact solution Vector (here a scalar)
253 
254  //Loop over the integration points
255  for(unsigned ipt=0;ipt<n_intpt;ipt++)
256  {
257 
258  //Assign values of s
259  for(unsigned i=0;i<DIM;i++)
260  {
261  s[i] = integral_pt()->knot(ipt,i);
262  }
263 
264  //Get the integral weight
265  double w = integral_pt()->weight(ipt);
266 
267  // Get jacobian of mapping
268  double J=J_eulerian(s);
269 
270  //Premultiply the weights and the Jacobian
271  double W = w*J;
272 
273  // Get x position as Vector
274  interpolated_x(s,x);
275 
276  // Get FE function value
277  double u_fe=interpolated_u(s);
278 
279  // Get exact solution at this point
280  (*exact_soln_pt)(x,exact_soln);
281 
282  //Output x,y,...,error
283  for(unsigned i=0;i<DIM;i++)
284  {
285  outfile << x[i] << " ";
286  }
287  outfile << exact_soln[0] << " " << exact_soln[0]-u_fe << std::endl;
288 
289  // Add to error and norm
290  norm+=exact_soln[0]*exact_soln[0]*W;
291  error+=(exact_soln[0]-u_fe)*(exact_soln[0]-u_fe)*W;
292 
293  }
294  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
list x
Definition: plotDoE.py:28

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

Referenced by ModalPRefineableQPoissonElement< DIM >::compute_error().

◆ compute_error() [2/2]

template<unsigned DIM>
void ModalPoissonEquations< DIM >::compute_error ( ostream &  outfile,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt,
const double time,
double error,
double norm 
)
inline

Dummy, time dependent error checker.

380  {
381  oomph_info << "No time-dep. compute_error() for Poisson elements " << std::endl;
382  exit(1);
383  }
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::oomph_info.

◆ diff_fct_pt() [1/2]

template<unsigned DIM>
PoissonDiffFctPt& ModalPoissonEquations< DIM >::diff_fct_pt ( )
inline

Access function: Pointer to diffusivity function.

418  {
419  return Diff_fct_pt;
420  }

References oomph::Diff_fct_pt.

◆ diff_fct_pt() [2/2]

template<unsigned DIM>
PoissonDiffFctPt ModalPoissonEquations< DIM >::diff_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

425  {
426  return Diff_fct_pt;
427  }

References oomph::Diff_fct_pt.

◆ dshape_dbasis_and_dtest_eulerian()

template<unsigned DIM>
virtual double ModalPoissonEquations< DIM >::dshape_dbasis_and_dtest_eulerian ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape basis,
DShape dbasisdx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in ModalPRefineableQPoissonElement< DIM >.

◆ dshape_dbasis_and_dtest_eulerian_at_knot()

template<unsigned DIM>
virtual double ModalPoissonEquations< DIM >::dshape_dbasis_and_dtest_eulerian_at_knot ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape basis,
DShape dbasisdx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in ModalPRefineableQPoissonElement< DIM >.

◆ fill_in_contribution_to_jacobian()

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

Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)

Reimplemented from oomph::GeneralisedElement.

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

497  {
498  //Call the generic routine with the flag set to 1
499  add_generic_residual_contribution(residuals,jacobian,1);
500  }
virtual void add_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: one_d_poisson_hp_adapt.cc:555

Referenced by ModalPRefineableQPoissonElement< DIM >::fill_in_contribution_to_jacobian().

◆ fill_in_contribution_to_residuals()

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

Add the element's contribution to its residual vector (wrapper)

Reimplemented from oomph::GeneralisedElement.

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

484  {
485  //Create a dummy matrix
486  DenseMatrix<double> dummy(1);
487 
488  //Call the generic residuals function with flag set to 0
489  add_generic_residual_contribution(residuals,dummy,0);
490  }

Referenced by ModalPRefineableQPoissonElement< DIM >::fill_in_contribution_to_residuals().

◆ get_diff()

template<unsigned DIM>
virtual void ModalPoissonEquations< DIM >::get_diff ( const Vector< double > &  x,
double diff 
) const
inlinevirtual

Get diffusivity term at (Eulerian) position x. This function is virtual to allow overloading in multi-physics problems where the strength of the diffusivity function might be determined by another system of equations

435  {
436  //If no source function has been set, return one
437  if(Diff_fct_pt==0) {diff = 1.0;}
438  else
439  {
440  // Get source strength
441  (*Diff_fct_pt)(x,diff);
442  }
443  }

References oomph::Diff_fct_pt, and plotDoE::x.

◆ get_flux()

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

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

449  {
450  //Find out how many nodes there are
451  unsigned n_node = nnode();
452  //Find out the total number of basis and test functions there are
453  unsigned n_basis = nbasis();
454 
455  //Set up memory for the shape and test functions
456  Shape psi(n_node), basis(n_basis), test(n_basis);
457  DShape dpsidx(n_node,DIM), dbasisdx(n_basis,DIM), dtestdx(n_basis,DIM);
458 
459  //Call the derivatives of the shape and test functions
460  //dshape_eulerian(s,psi,dpsidx);
461  dshape_dbasis_and_dtest_eulerian(s, psi, dpsidx, basis, dbasisdx, test, dtestdx);
462 
463  //Initialise to zero
464  for(unsigned j=0;j<DIM;j++)
465  {
466  flux[j] = 0.0;
467  }
468 
469  // Loop over modes
470  for(unsigned l=0;l<p_order;l++)
471  {
472  //Loop over derivative directions
473  for(unsigned j=0;j<DIM;j++)
474  {
475  //flux[j] += u(l)*dpsidx(l,j);
476  flux[j] += u(l)*dbasisdx(l,j);
477  }
478  }
479  }
virtual double dshape_dbasis_and_dtest_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &basis, DShape &dbasisdx, Shape &test, DShape &dtestdx) const =0
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59

References DIM, ProblemParameters::flux(), j, s, and Eigen::test.

Referenced by ModalPRefineableQPoissonElement< DIM >::compute_energy_error().

◆ get_source_poisson()

template<unsigned DIM>
virtual void ModalPoissonEquations< DIM >::get_source_poisson ( const Vector< double > &  x,
double source 
) const
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

406  {
407  //If no source function has been set, return zero
408  if(Source_fct_pt==0) {source = 0.0;}
409  else
410  {
411  // Get source strength
412  (*Source_fct_pt)(x,source);
413  }
414  }

References TestProblem::source(), oomph::Source_fct_pt, and plotDoE::x.

◆ interpolated_u()

template<unsigned DIM>
double ModalPoissonEquations< DIM >::interpolated_u ( const Vector< double > &  s) const
inline

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

505  {
506  //Find number of basis functions
507  unsigned n_basis = nbasis();
508 
509  //Local basis functions
510  Shape psi(n_basis);
511 
512  //Find values of basis function
513  basis(s,psi);
514 
515  //Initialise value of u
516  double interpolated_u = 0.0;
517 
518  //Loop over the local nodes and sum
519  for(unsigned l=0;l<n_basis;l++)
520  {
521  interpolated_u+=u(l)*psi(l);
522  }
523 
524  return(interpolated_u);
525  }

References s.

◆ nbasis()

template<unsigned DIM>
virtual unsigned ModalPoissonEquations< DIM >::nbasis ( ) const
pure virtual

Number of basis functions.

Implemented in ModalPRefineableQPoissonElement< DIM >.

◆ output() [1/2]

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

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

115  {
116  unsigned nplot=5;
117  output(outfile,nplot);
118  }
void output(ostream &outfile)
Output with default number of plot points.
Definition: one_d_poisson_hp_adapt.cc:114

References output().

Referenced by ModalPRefineableQPoissonElement< DIM >::output().

◆ output() [2/2]

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

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

Reimplemented from oomph::FiniteElement.

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

123  {
124  //Vector of local coordinates
126 
127  //Find out how many nodes there are
128  unsigned n_node = nnode();
129  //Find out the total number of basis and test functions there are
130  unsigned n_basis = nbasis();
131 
132  //Set up memory for the shape and test functions
133  Shape psi(n_node), basis(n_basis), test(n_basis);
134  DShape dpsidx(n_node,DIM), dbasisdx(n_basis,DIM), dtestdx(n_basis,DIM);
135 
136  // Tecplot header info
137  outfile << tecplot_zone_string(nplot);
138 
139  // Loop over plot points
140  unsigned num_plot_points=nplot_points(nplot);
141  for (unsigned iplot=0;iplot<num_plot_points;iplot++)
142  {
143 
144  // Get local coordinates of plot point
145  get_s_plot(iplot,nplot,s);
146 
147  //get_basis(s,basis);
149  dbasisdx,test,dtestdx);
150 
151  for(unsigned i=0;i<DIM;i++)
152  {outfile << interpolated_x(s,i) << " ";}
153 
154  outfile << interpolated_u(s) << std::endl;
155  }
156 
157  }
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186

References DIM, i, s, and Eigen::test.

◆ output_fct() [1/2]

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

Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points (dummy time-dependent version to keep intel compiler happy)

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

211  {
212 
213  oomph_info << "No time-dep. output_fct() for Poisson elements " << std::endl;
214  exit(1);
215  }

References oomph::oomph_info.

◆ output_fct() [2/2]

template<unsigned DIM>
void ModalPoissonEquations< DIM >::output_fct ( ostream &  outfile,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
inline

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

164  {
165 
166  //Vector of local coordinates
168 
169 
170  // Vector for coordintes
172 
173  // Tecplot header info
174  outfile << tecplot_zone_string(nplot);
175 
176  // Exact solution Vector (here a scalar)
178 
179  // Loop over plot points
180  unsigned num_plot_points=nplot_points(nplot);
181  for (unsigned iplot=0;iplot<num_plot_points;iplot++)
182  {
183 
184  // Get local coordinates of plot point
185  get_s_plot(iplot,nplot,s);
186 
187  // Get x position as Vector
188  interpolated_x(s,x);
189 
190  // Get exact solution at this point
191  (*exact_soln_pt)(x,exact_soln);
192 
193  //Output x,y,...,u_exact
194  for(unsigned i=0;i<DIM;i++)
195  {
196  outfile << x[i] << " ";
197  }
198  outfile << exact_soln[0] << std::endl;
199 
200  }
201  }

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

Referenced by ModalPRefineableQPoissonElement< DIM >::output_fct().

◆ self_test()

template<unsigned DIM>
unsigned ModalPoissonEquations< DIM >::self_test ( )
inlinevirtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

Reimplemented in ModalPRefineableQPoissonElement< DIM >.

529 {return 0;}

Referenced by ModalPRefineableQPoissonElement< DIM >::self_test().

◆ source_fct_pt() [1/2]

template<unsigned DIM>
PoissonSourceFctPt& ModalPoissonEquations< DIM >::source_fct_pt ( )
inline

Access function: Pointer to source function.

389  {
390  return Source_fct_pt;
391  }

References oomph::Source_fct_pt.

Referenced by ModalPRefineableQPoissonElement< DIM >::source_fct_pt().

◆ source_fct_pt() [2/2]

template<unsigned DIM>
PoissonSourceFctPt ModalPoissonEquations< DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

396  {
397  return Source_fct_pt;
398  }

References oomph::Source_fct_pt.

◆ u()

template<unsigned DIM>
virtual double ModalPoissonEquations< DIM >::u ( const unsigned n) const
pure virtual

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

Implemented in ModalPRefineableQPoissonElement< DIM >.

Member Data Documentation

◆ Diff_fct_pt

template<unsigned DIM>
PoissonDiffFctPt ModalPoissonEquations< DIM >::Diff_fct_pt
protected

Pointer to the diffusivity function.

◆ Local_eqn

template<unsigned DIM>
Vector<int> ModalPoissonEquations< DIM >::Local_eqn
protected

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

◆ Source_fct_pt

template<unsigned DIM>
PoissonSourceFctPt ModalPoissonEquations< DIM >::Source_fct_pt
protected

Pointer to source function:


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