oomph::UnsteadyHeatEquations< DIM > Class Template Referenceabstract

#include <unsteady_heat_elements.h>

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

Public Types

typedef void(* UnsteadyHeatSourceFctPt) (const double &time, const Vector< double > &x, double &u)
 
- Public Types inherited from oomph::UnsteadyHeatEquationsBase
typedef void(* UnsteadyHeatSourceFctPt) (const double &time, const Vector< double > &x, double &u)
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 

Public Member Functions

 UnsteadyHeatEquations ()
 
 UnsteadyHeatEquations (const UnsteadyHeatEquations &dummy)=delete
 Broken copy constructor. More...
 
virtual unsigned u_index_ust_heat () const
 Broken assignment operator. More...
 
double du_dt_ust_heat (const unsigned &n) const
 
void disable_ALE ()
 
void enable_ALE ()
 
void compute_norm (double &norm)
 Compute norm of fe solution. More...
 
void output (std::ostream &outfile)
 Output with default number of plot points. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void output (FILE *file_pt)
 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 &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
UnsteadyHeatSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
UnsteadyHeatSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
virtual void get_source_ust_heat (const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
 
const doublealpha () const
 Alpha parameter (thermal inertia) More...
 
double *& alpha_pt ()
 Pointer to Alpha parameter (thermal inertia) More...
 
const doublebeta () const
 Beta parameter (thermal conductivity) More...
 
double *& beta_pt ()
 Pointer to Beta parameter (thermal conductivity) More...
 
void get_flux (const Vector< double > &s, Vector< double > &flux) const
 Get flux: flux[i] = du/dx_i. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Compute element residual Vector (wrapper) More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Compute element residual Vector and element Jacobian matrix (wrapper) More...
 
double interpolated_u_ust_heat (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
double interpolated_u_ust_heat (const unsigned &t, const Vector< double > &s) const
 
double interpolated_du_dt_ust_heat (const Vector< double > &s) const
 
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 unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

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

Protected Attributes

UnsteadyHeatSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
bool ALE_is_disabled
 
doubleAlpha_pt
 Pointer to Alpha parameter (thermal inertia) More...
 
doubleBeta_pt
 Pointer to Beta parameter (thermal conductivity) 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
 

Static Private Attributes

static double Default_alpha_parameter = 1.0
 2D UnsteadyHeat elements More...
 
static double Default_beta_parameter = 1.0
 Default value for Beta parameter (thermal conductivity) More...
 

Additional Inherited Members

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

Detailed Description

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

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

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

This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class. Note that this class assumes an isoparametric formulation, i.e. that the scalar unknown is interpolated using the same shape funcitons as the position.

Member Typedef Documentation

◆ UnsteadyHeatSourceFctPt

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

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

Constructor & Destructor Documentation

◆ UnsteadyHeatEquations() [1/2]

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

Constructor: Initialises the Source_fct_pt to null and sets flag to use ALE formulation of the equations. Also set Alpha (thermal inertia) and Beta (thermal conductivity) parameters to defaults (both one for natural scaling)

85  : Source_fct_pt(0), ALE_is_disabled(false)
86  {
87  // Set Alpha and Beta parameter to default (one for natural scaling of
88  // time)
91  }
UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
Definition: unsteady_heat_elements.h:446
static double Default_beta_parameter
Default value for Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:466
bool ALE_is_disabled
Definition: unsteady_heat_elements.h:451
double * Beta_pt
Pointer to Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:457
static double Default_alpha_parameter
2D UnsteadyHeat elements
Definition: unsteady_heat_elements.h:462
double * Alpha_pt
Pointer to Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:454

References oomph::UnsteadyHeatEquations< DIM >::Alpha_pt, oomph::UnsteadyHeatEquations< DIM >::Beta_pt, oomph::UnsteadyHeatEquations< DIM >::Default_alpha_parameter, and oomph::UnsteadyHeatEquations< DIM >::Default_beta_parameter.

◆ UnsteadyHeatEquations() [2/2]

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

Broken copy constructor.

Member Function Documentation

◆ alpha()

template<unsigned DIM>
const double& oomph::UnsteadyHeatEquations< DIM >::alpha ( ) const
inline

Alpha parameter (thermal inertia)

256  {
257  return *Alpha_pt;
258  }

References oomph::UnsteadyHeatEquations< DIM >::Alpha_pt.

◆ alpha_pt()

template<unsigned DIM>
double*& oomph::UnsteadyHeatEquations< DIM >::alpha_pt ( )
inline

Pointer to Alpha parameter (thermal inertia)

262  {
263  return Alpha_pt;
264  }

References oomph::UnsteadyHeatEquations< DIM >::Alpha_pt.

◆ beta()

template<unsigned DIM>
const double& oomph::UnsteadyHeatEquations< DIM >::beta ( ) const
inline

Beta parameter (thermal conductivity)

269  {
270  return *Beta_pt;
271  }

References oomph::UnsteadyHeatEquations< DIM >::Beta_pt.

◆ beta_pt()

template<unsigned DIM>
double*& oomph::UnsteadyHeatEquations< DIM >::beta_pt ( )
inline

Pointer to Beta parameter (thermal conductivity)

275  {
276  return Beta_pt;
277  }

References oomph::UnsteadyHeatEquations< DIM >::Beta_pt.

◆ compute_error() [1/2]

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

Get error against and norm of exact solution.

Validate against exact solution

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

Reimplemented from oomph::FiniteElement.

477  {
478  // Initialise
479  error = 0.0;
480  norm = 0.0;
481 
482  // Vector of local coordinates
483  Vector<double> s(DIM);
484 
485  // Vector for coordintes
486  Vector<double> x(DIM);
487 
488  // Find out how many nodes there are in the element
489  unsigned n_node = nnode();
490 
491  Shape psi(n_node);
492 
493  // Set the value of n_intpt
494  unsigned n_intpt = integral_pt()->nweight();
495 
496  // Tecplot header info
497  outfile << "ZONE" << std::endl;
498 
499  // Exact solution Vector (here a scalar)
500  Vector<double> exact_soln(1);
501 
502  // Loop over the integration points
503  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
504  {
505  // Assign values of s
506  for (unsigned i = 0; i < DIM; i++)
507  {
508  s[i] = integral_pt()->knot(ipt, i);
509  }
510 
511  // Get the integral weight
512  double w = integral_pt()->weight(ipt);
513 
514  // Get jacobian of mapping
515  double J = J_eulerian(s);
516 
517  // Premultiply the weights and the Jacobian
518  double W = w * J;
519 
520  // Get x position as Vector
521  interpolated_x(s, x);
522 
523  // Get FE function value
524  double u_fe = interpolated_u_ust_heat(s);
525 
526  // Get exact solution at this point
527  (*exact_soln_pt)(x, exact_soln);
528 
529  // Output x,y,...,error
530  for (unsigned i = 0; i < DIM; i++)
531  {
532  outfile << x[i] << " ";
533  }
534  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
535 
536  // Add to error and norm
537  norm += exact_soln[0] * exact_soln[0] * W;
538  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
539  }
540  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double interpolated_u_ust_heat(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: unsteady_heat_elements.h:333
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

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

◆ compute_error() [2/2]

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

Get error against and norm of exact solution.

Validate against exact solution at time t.

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

Reimplemented from oomph::FiniteElement.

558  {
559  // Initialise
560  error = 0.0;
561  norm = 0.0;
562 
563  // Vector of local coordinates
564  Vector<double> s(DIM);
565 
566  // Vector for coordintes
567  Vector<double> x(DIM);
568 
569 
570  // Find out how many nodes there are in the element
571  unsigned n_node = nnode();
572 
573  Shape psi(n_node);
574 
575  // Set the value of n_intpt
576  unsigned n_intpt = integral_pt()->nweight();
577 
578  // Tecplot
579  outfile << "ZONE" << std::endl;
580 
581  // Exact solution Vector (here a scalar)
582  Vector<double> exact_soln(1);
583 
584  // Loop over the integration points
585  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
586  {
587  // Assign values of s
588  for (unsigned i = 0; i < DIM; i++)
589  {
590  s[i] = integral_pt()->knot(ipt, i);
591  }
592 
593  // Get the integral weight
594  double w = integral_pt()->weight(ipt);
595 
596  // Get jacobian of mapping
597  double J = J_eulerian(s);
598 
599  // Premultiply the weights and the Jacobian
600  double W = w * J;
601 
602  // Get x position as Vector
603  interpolated_x(s, x);
604 
605  // Get FE function value
606  double u_fe = interpolated_u_ust_heat(s);
607 
608  // Get exact solution at this point
609  (*exact_soln_pt)(time, x, exact_soln);
610 
611  // Output x,y,...,error
612  for (unsigned i = 0; i < DIM; i++)
613  {
614  outfile << x[i] << " ";
615  }
616  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
617 
618  // Add to error and norm
619  norm += exact_soln[0] * exact_soln[0] * W;
620  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
621  }
622  }

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

◆ compute_norm()

template<unsigned DIM>
void oomph::UnsteadyHeatEquations< DIM >::compute_norm ( double norm)
virtual

Compute norm of fe solution.

Reimplemented from oomph::GeneralisedElement.

216  {
217  // Initialise
218  norm = 0.0;
219 
220  // Vector of local coordinates
221  Vector<double> s(DIM);
222 
223  // Vector for coordintes
224  Vector<double> x(DIM);
225 
226  // Find out how many nodes there are in the element
227  unsigned n_node = nnode();
228 
229  Shape psi(n_node);
230 
231  // Set the value of n_intpt
232  unsigned n_intpt = integral_pt()->nweight();
233 
234  // Loop over the integration points
235  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
236  {
237  // Assign values of s
238  for (unsigned i = 0; i < DIM; i++)
239  {
240  s[i] = integral_pt()->knot(ipt, i);
241  }
242 
243  // Get the integral weight
244  double w = integral_pt()->weight(ipt);
245 
246  // Get jacobian of mapping
247  double J = J_eulerian(s);
248 
249  // Premultiply the weights and the Jacobian
250  double W = w * J;
251 
252  // Get FE function value
253  double u = interpolated_u_ust_heat(s);
254 
255  // Add to norm
256  norm += u * u * W;
257  }
258  }

References DIM, i, J, s, w, oomph::QuadTreeNames::W, and plotDoE::x.

Referenced by QThermalPVDElement< DIM >::compute_norm().

◆ disable_ALE()

template<unsigned DIM>
void oomph::UnsteadyHeatEquations< DIM >::disable_ALE ( )
inlinevirtual

Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!

Reimplemented from oomph::FiniteElement.

149  {
150  ALE_is_disabled = true;
151  }

References oomph::UnsteadyHeatEquations< DIM >::ALE_is_disabled.

◆ dshape_and_dtest_eulerian_at_knot_ust_heat()

template<unsigned DIM>
virtual double oomph::UnsteadyHeatEquations< DIM >::dshape_and_dtest_eulerian_at_knot_ust_heat ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, oomph::QUnsteadyHeatElement< DIM, 3 >, and oomph::TUnsteadyHeatElement< DIM, NNODE_1D >.

◆ dshape_and_dtest_eulerian_ust_heat()

template<unsigned DIM>
virtual double oomph::UnsteadyHeatEquations< DIM >::dshape_and_dtest_eulerian_ust_heat ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, oomph::QUnsteadyHeatElement< DIM, 3 >, and oomph::TUnsteadyHeatElement< DIM, NNODE_1D >.

◆ du_dt_ust_heat()

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat ( const unsigned n) const
inline

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

120  {
121  // Get the data's timestepper
122  TimeStepper* time_stepper_pt = this->node_pt(n)->time_stepper_pt();
123 
124  // Initialise dudt
125  double dudt = 0.0;
126 
127  // Loop over the timesteps, if there is a non Steady timestepper
128  if (!time_stepper_pt->is_steady())
129  {
130  // Find the index at which the variable is stored
131  const unsigned u_nodal_index = u_index_ust_heat();
132 
133  // Number of timsteps (past & present)
134  const unsigned n_time = time_stepper_pt->ntstorage();
135 
136  // Add the contributions to the time derivative
137  for (unsigned t = 0; t < n_time; t++)
138  {
139  dudt +=
140  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
141  }
142  }
143  return dudt;
144  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
Definition: unsteady_heat_elements.h:112
t
Definition: plotPSD.py:36

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

Referenced by oomph::UnsteadyHeatEquations< DIM >::interpolated_du_dt_ust_heat().

◆ enable_ALE()

template<unsigned DIM>
void oomph::UnsteadyHeatEquations< DIM >::enable_ALE ( )
inlinevirtual

(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative. Note: By default, ALE is enabled, at the expense of possibly creating unnecessary work in problems where the mesh is, in fact, stationary.

Reimplemented from oomph::FiniteElement.

159  {
160  ALE_is_disabled = false;
161  }

References oomph::UnsteadyHeatEquations< DIM >::ALE_is_disabled.

◆ fill_in_contribution_to_jacobian()

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

Compute element residual Vector and element Jacobian matrix (wrapper)

Reimplemented from oomph::FiniteElement.

326  {
327  // Call the generic routine with the flag set to 1
328  fill_in_generic_residual_contribution_ust_heat(residuals, jacobian, 1);
329  }
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: unsteady_heat_elements.cc:60

References oomph::UnsteadyHeatEquations< DIM >::fill_in_generic_residual_contribution_ust_heat().

◆ fill_in_contribution_to_residuals()

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

Compute element residual Vector (wrapper)

Reimplemented from oomph::GeneralisedElement.

315  {
316  // Call the generic residuals function with flag set to 0
317  // using a dummy matrix argument
319  residuals, GeneralisedElement::Dummy_matrix, 0);
320  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::UnsteadyHeatEquations< DIM >::fill_in_generic_residual_contribution_ust_heat().

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

◆ fill_in_generic_residual_contribution_ust_heat()

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

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

Compute element residual Vector and/or element Jacobian matrix

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

Pure version without hanging nodes

Reimplemented in oomph::RefineableUnsteadyHeatEquations< DIM >.

62  {
63  // Find out how many nodes there are
64  unsigned n_node = nnode();
65 
66  // Get continuous time from timestepper of first node
67  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
68 
69  // Find the index at which the variable is stored
70  unsigned u_nodal_index = u_index_ust_heat();
71 
72  // Set up memory for the shape and test functions
73  Shape psi(n_node), test(n_node);
74  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
75 
76  // Set the value of n_intpt
77  unsigned n_intpt = integral_pt()->nweight();
78 
79  // Set the Vector to hold local coordinates
80  Vector<double> s(DIM);
81 
82  // Get Alpha and beta parameters number
83  double alpha_local = alpha();
84  double beta_local = beta();
85 
86  // Integers to hold the local equation and unknowns
87  int local_eqn = 0, local_unknown = 0;
88 
89  // Loop over the integration points
90  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
91  {
92  // Assign values of s
93  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
94 
95  // Get the integral weight
96  double w = integral_pt()->weight(ipt);
97 
98  // Call the derivatives of the shape and test functions
100  ipt, psi, dpsidx, test, dtestdx);
101 
102  // Premultiply the weights and the Jacobian
103  double W = w * J;
104 
105  // Allocate memory for local quantities and initialise to zero
106  double interpolated_u = 0.0;
107  double dudt = 0.0;
108  Vector<double> interpolated_x(DIM, 0.0);
109  Vector<double> interpolated_dudx(DIM, 0.0);
110  Vector<double> mesh_velocity(DIM, 0.0);
111 
112  // Calculate function value and derivatives:
113  // Loop over nodes
114  for (unsigned l = 0; l < n_node; l++)
115  {
116  // Calculate the value at the nodes
117  double u_value = raw_nodal_value(l, u_nodal_index);
118  interpolated_u += u_value * psi(l);
119  dudt += du_dt_ust_heat(l) * psi(l);
120  // Loop over directions
121  for (unsigned j = 0; j < DIM; j++)
122  {
123  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
124  interpolated_dudx[j] += u_value * dpsidx(l, j);
125  }
126  }
127 
128  // Mesh velocity?
129  if (!ALE_is_disabled)
130  {
131  for (unsigned l = 0; l < n_node; l++)
132  {
133  for (unsigned j = 0; j < DIM; j++)
134  {
135  mesh_velocity[j] += raw_dnodal_position_dt(l, j) * psi(l);
136  }
137  }
138  }
139 
140  // Get source function
141  //-------------------
142  double source = 0.0;
144 
145  // Assemble residuals and Jacobian
146  //--------------------------------
147 
148  // Loop over the test functions
149  for (unsigned l = 0; l < n_node; l++)
150  {
151  local_eqn = nodal_local_eqn(l, u_nodal_index);
152  /*IF it's not a boundary condition*/
153  if (local_eqn >= 0)
154  {
155  // Add body force/source term and time derivative
156  residuals[local_eqn] += (alpha_local * dudt + source) * test(l) * W;
157 
158  // The mesh velocity bit
159  if (!ALE_is_disabled)
160  {
161  for (unsigned k = 0; k < DIM; k++)
162  {
163  residuals[local_eqn] -= alpha_local * mesh_velocity[k] *
164  interpolated_dudx[k] * test(l) * W;
165  }
166  }
167 
168  // Laplace operator
169  for (unsigned k = 0; k < DIM; k++)
170  {
171  residuals[local_eqn] +=
172  beta_local * interpolated_dudx[k] * dtestdx(l, k) * W;
173  }
174 
175 
176  // Calculate the jacobian
177  //-----------------------
178  if (flag)
179  {
180  // Loop over the velocity shape functions again
181  for (unsigned l2 = 0; l2 < n_node; l2++)
182  {
183  local_unknown = nodal_local_eqn(l2, u_nodal_index);
184  // If at a non-zero degree of freedom add in the entry
185  if (local_unknown >= 0)
186  {
187  // Mass matrix
188  jacobian(local_eqn, local_unknown) +=
189  alpha_local * test(l) * psi(l2) *
190  node_pt(l2)->time_stepper_pt()->weight(1, 0) * W;
191 
192  // Laplace operator & mesh velocity bit
193  for (unsigned i = 0; i < DIM; i++)
194  {
195  double tmp = beta_local * dtestdx(l, i);
196  if (!ALE_is_disabled)
197  tmp -= alpha_local * mesh_velocity[i] * test(l);
198  jacobian(local_eqn, local_unknown) += dpsidx(l2, i) * tmp * W;
199  }
200  }
201  }
202  }
203  }
204  }
205 
206 
207  } // End of loop over integration points
208  }
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
double raw_dnodal_position_dt(const unsigned &n, const unsigned &i) const
Definition: elements.h:2256
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.cc:1686
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
const double & alpha() const
Alpha parameter (thermal inertia)
Definition: unsteady_heat_elements.h:255
virtual void get_source_ust_heat(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: unsteady_heat_elements.h:237
virtual double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
const double & beta() const
Beta parameter (thermal conductivity)
Definition: unsteady_heat_elements.h:268
double du_dt_ust_heat(const unsigned &n) const
Definition: unsteady_heat_elements.h:119
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

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

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

◆ get_flux()

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

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

281  {
282  // Find out how many nodes there are in the element
283  unsigned n_node = nnode();
284 
285  // Find the index at which the variable is stored
286  unsigned u_nodal_index = u_index_ust_heat();
287 
288  // Set up memory for the shape and test functions
289  Shape psi(n_node);
290  DShape dpsidx(n_node, DIM);
291 
292  // Call the derivatives of the shape and test functions
293  dshape_eulerian(s, psi, dpsidx);
294 
295  // Initialise to zero
296  for (unsigned j = 0; j < DIM; j++)
297  {
298  flux[j] = 0.0;
299  }
300 
301  // Loop over nodes
302  for (unsigned l = 0; l < n_node; l++)
303  {
304  // Loop over derivative directions
305  for (unsigned j = 0; j < DIM; j++)
306  {
307  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
308  }
309  }
310  }
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59

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

Referenced by oomph::RefineableUnsteadyHeatEquations< DIM >::get_Z2_flux(), and oomph::TUnsteadyHeatElement< DIM, NNODE_1D >::get_Z2_flux().

◆ get_source_ust_heat()

template<unsigned DIM>
virtual void oomph::UnsteadyHeatEquations< DIM >::get_source_ust_heat ( const double t,
const unsigned ipt,
const Vector< double > &  x,
double source 
) const
inlinevirtual

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

241  {
242  // If no source function has been set, return zero
243  if (Source_fct_pt == 0)
244  {
245  source = 0.0;
246  }
247  else
248  {
249  // Get source strength
250  (*Source_fct_pt)(t, x, source);
251  }
252  }

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

◆ interpolated_du_dt_ust_heat()

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::interpolated_du_dt_ust_heat ( const Vector< double > &  s) const
inline

Return FE representation of function value du/dt(s) at local coordinate s

393  {
394  // Find number of nodes
395  unsigned n_node = nnode();
396 
397  // Local shape function
398  Shape psi(n_node);
399 
400  // Find values of shape function
401  shape(s, psi);
402 
403  // Initialise value of du/dt
404  double interpolated_dudt = 0.0;
405 
406  // Loop over the local nodes and sum
407  for (unsigned l = 0; l < n_node; l++)
408  {
409  interpolated_dudt += du_dt_ust_heat(l) * psi[l];
410  }
411 
412  return (interpolated_dudt);
413  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::FiniteElement::nnode(), s, and oomph::FiniteElement::shape().

◆ interpolated_u_ust_heat() [1/2]

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::interpolated_u_ust_heat ( const unsigned t,
const Vector< double > &  s 
) const
inline

Return FE representation of function value u(s) at local coordinate s at previous time t (t=0: present)

364  {
365  // Find number of nodes
366  unsigned n_node = nnode();
367 
368  // Find the index at which the variable is stored
369  unsigned u_nodal_index = u_index_ust_heat();
370 
371  // Local shape function
372  Shape psi(n_node);
373 
374  // Find values of shape function
375  shape(s, psi);
376 
377  // Initialise value of u
378  double interpolated_u = 0.0;
379 
380  // Loop over the local nodes and sum
381  for (unsigned l = 0; l < n_node; l++)
382  {
383  interpolated_u += nodal_value(t, l, u_nodal_index) * psi[l];
384  }
385 
386  return (interpolated_u);
387  }

References oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::FiniteElement::shape(), plotPSD::t, and oomph::UnsteadyHeatEquations< DIM >::u_index_ust_heat().

◆ interpolated_u_ust_heat() [2/2]

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::interpolated_u_ust_heat ( const Vector< double > &  s) const
inline

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

334  {
335  // Find number of nodes
336  unsigned n_node = nnode();
337 
338  // Find the index at which the variable is stored
339  unsigned u_nodal_index = u_index_ust_heat();
340 
341  // Local shape function
342  Shape psi(n_node);
343 
344  // Find values of shape function
345  shape(s, psi);
346 
347  // Initialise value of u
348  double interpolated_u = 0.0;
349 
350  // Loop over the local nodes and sum
351  for (unsigned l = 0; l < n_node; l++)
352  {
353  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
354  }
355 
356  return (interpolated_u);
357  }

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

◆ output() [1/4]

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

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, and oomph::QUnsteadyHeatElement< DIM, 3 >.

180  {
181  unsigned n_plot = 5;
182  output(file_pt, n_plot);
183  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: unsteady_heat_elements.h:167

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

◆ output() [2/4]

template<unsigned DIM>
void oomph::UnsteadyHeatEquations< DIM >::output ( FILE *  file_pt,
const unsigned nplot 
)
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::QUnsteadyHeatElement< DIM, NNODE_1D >, and oomph::QUnsteadyHeatElement< DIM, 3 >.

331  {
332  // Vector of local coordinates
333  Vector<double> s(DIM);
334 
335  // Tecplot header info
336  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
337 
338  // Loop over plot points
339  unsigned num_plot_points = nplot_points(nplot);
340  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
341  {
342  // Get local coordinates of plot point
343  get_s_plot(iplot, nplot, s);
344 
345  for (unsigned i = 0; i < DIM; i++)
346  {
347  fprintf(file_pt, "%g ", interpolated_x(s, i));
348  }
349  fprintf(file_pt, "%g \n", interpolated_u_ust_heat(s));
350  }
351 
352  // Write tecplot footer (e.g. FE connectivity lists)
353  write_tecplot_zone_footer(file_pt, nplot);
354  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References DIM, i, and s.

◆ output() [3/4]

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

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, and oomph::QUnsteadyHeatElement< DIM, 3 >.

168  {
169  unsigned nplot = 5;
170  output(outfile, nplot);
171  }

Referenced by oomph::UnsteadyHeatEquations< DIM >::output(), oomph::TUnsteadyHeatElement< DIM, NNODE_1D >::output(), and oomph::QUnsteadyHeatElement< DIM, NNODE_1D >::output().

◆ output() [4/4]

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

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

Output function:

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

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, and oomph::QUnsteadyHeatElement< DIM, 3 >.

296  {
297  // Vector of local coordinates
298  Vector<double> s(DIM);
299 
300  // Tecplot header info
301  outfile << tecplot_zone_string(nplot);
302 
303  // Loop over plot points
304  unsigned num_plot_points = nplot_points(nplot);
305  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
306  {
307  // Get local coordinates of plot point
308  get_s_plot(iplot, nplot, s);
309 
310  for (unsigned i = 0; i < DIM; i++)
311  {
312  outfile << interpolated_x(s, i) << " ";
313  }
314  outfile << interpolated_u_ust_heat(s) << std::endl;
315  }
316 
317  // Write tecplot footer (e.g. FE connectivity lists)
318  write_tecplot_zone_footer(outfile, nplot);
319  }

References DIM, i, and s.

◆ output_fct() [1/2]

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

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

Output exact solution at time t

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

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

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TLinearHeatAndElasticityElement< DIM, NNODE_1D >, oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, oomph::QUnsteadyHeatElement< DIM, 3 >, and oomph::TUnsteadyHeatElement< DIM, NNODE_1D >.

424  {
425  // Vector of local coordinates
426  Vector<double> s(DIM);
427 
428  // Vector for coordintes
429  Vector<double> x(DIM);
430 
431 
432  // Tecplot header info
433  outfile << tecplot_zone_string(nplot);
434 
435  // Exact solution Vector (here a scalar)
436  Vector<double> exact_soln(1);
437 
438  // Loop over plot points
439  unsigned num_plot_points = nplot_points(nplot);
440  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
441  {
442  // Get local coordinates of plot point
443  get_s_plot(iplot, nplot, s);
444 
445  // Get x position as Vector
446  interpolated_x(s, x);
447 
448  // Get exact solution at this point
449  (*exact_soln_pt)(time, x, exact_soln);
450 
451  // Output x,y,...,u_exact
452  for (unsigned i = 0; i < DIM; i++)
453  {
454  outfile << x[i] << " ";
455  }
456  outfile << exact_soln[0] << std::endl;
457  }
458 
459  // Write tecplot footer (e.g. FE connectivity lists)
460  write_tecplot_zone_footer(outfile, nplot);
461  }

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

◆ output_fct() [2/2]

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

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

Output exact solution

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

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

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QUnsteadyHeatElement< DIM, NNODE_1D >, and oomph::QUnsteadyHeatElement< DIM, 3 >.

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

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

Referenced by oomph::TUnsteadyHeatElement< DIM, NNODE_1D >::output_fct(), and oomph::QUnsteadyHeatElement< DIM, NNODE_1D >::output_fct().

◆ self_test()

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

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

265  {
266  bool passed = true;
267 
268  // Check lower-level stuff
269  if (FiniteElement::self_test() != 0)
270  {
271  passed = false;
272  }
273 
274  // Return verdict
275  if (passed)
276  {
277  return 0;
278  }
279  else
280  {
281  return 1;
282  }
283  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

template<unsigned DIM>
UnsteadyHeatSourceFctPt& oomph::UnsteadyHeatEquations< DIM >::source_fct_pt ( )
inlinevirtual

Access function: Pointer to source function.

Implements oomph::UnsteadyHeatEquationsBase.

223  {
224  return Source_fct_pt;
225  }

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

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

◆ source_fct_pt() [2/2]

template<unsigned DIM>
UnsteadyHeatSourceFctPt oomph::UnsteadyHeatEquations< DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

230  {
231  return Source_fct_pt;
232  }

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

◆ u_index_ust_heat()

template<unsigned DIM>
virtual unsigned oomph::UnsteadyHeatEquations< DIM >::u_index_ust_heat ( ) const
inlinevirtual

Broken assignment operator.

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 unsteady heat 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.

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

113  {
114  return 0;
115  }

Referenced by oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::UnsteadyHeatEquations< DIM >::get_flux(), oomph::RefineableUnsteadyHeatEquations< DIM >::get_interpolated_values(), oomph::HeatedPenetratorFluxElement< ELEMENT >::HeatedPenetratorFluxElement(), oomph::UnsteadyHeatEquations< DIM >::interpolated_u_ust_heat(), oomph::UnsteadyHeatBaseFaceElement< ELEMENT >::UnsteadyHeatBaseFaceElement(), oomph::UnsteadyHeatFluxElement< ELEMENT >::UnsteadyHeatFluxElement(), and oomph::UnsteadyHeatFluxPseudoMeltElement< ELEMENT >::UnsteadyHeatFluxPseudoMeltElement().

Member Data Documentation

◆ ALE_is_disabled

template<unsigned DIM>
bool oomph::UnsteadyHeatEquations< DIM >::ALE_is_disabled
protected

Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed. Only set to true if you're sure that the mesh is stationary.

Referenced by oomph::UnsteadyHeatEquations< DIM >::disable_ALE(), oomph::UnsteadyHeatEquations< DIM >::enable_ALE(), and oomph::RefineableUnsteadyHeatEquations< DIM >::further_build().

◆ Alpha_pt

◆ Beta_pt

◆ Default_alpha_parameter

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::Default_alpha_parameter = 1.0
staticprivate

2D UnsteadyHeat elements

Static default value for the Alpha parameter: (thermal inertia): One for natural scaling

Default value for Alpha parameter (thermal inertia)

Referenced by oomph::UnsteadyHeatEquations< DIM >::UnsteadyHeatEquations().

◆ Default_beta_parameter

template<unsigned DIM>
double oomph::UnsteadyHeatEquations< DIM >::Default_beta_parameter = 1.0
staticprivate

Default value for Beta parameter (thermal conductivity)

Static default value for the Beta parameter (thermal conductivity): One for natural scaling

Referenced by oomph::UnsteadyHeatEquations< DIM >::UnsteadyHeatEquations().

◆ Source_fct_pt


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