oomph::BiharmonicProblem< DIM > Class Template Reference

#include <biharmonic_problem.h>

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

Public Types

typedef void(* DirichletBCFctPt) (const double &s, double &u)
 
typedef void(* BiharmonicSourceFctPt) (const Vector< double > &x, double &f)
 Definition of the Source Function. More...
 
- Public Types inherited from oomph::Problem
typedef void(* SpatialErrorEstimatorFctPt) (Mesh *&mesh_pt, Vector< double > &elemental_error)
 Function pointer for spatial error estimator. More...
 
typedef void(* SpatialErrorEstimatorWithDocFctPt) (Mesh *&mesh_pt, Vector< double > &elemental_error, DocInfo &doc_info)
 Function pointer for spatial error estimator with doc. More...
 

Public Member Functions

 BiharmonicProblem ()
 Constructor. More...
 
virtual ~BiharmonicProblem ()
 Destructor. Delete the meshes. More...
 
void actions_before_newton_solve ()
 actions before solve, performs self test More...
 
void actions_after_newton_solve ()
 action after solve More...
 
void doc_solution (DocInfo &doc_info, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt=0)
 
Meshbulk_element_mesh_pt ()
 Access function to the bulk element mesh pt. More...
 
- Public Member Functions inherited from oomph::Problem
virtual void debug_hook_fct (const unsigned &i)
 
void set_analytic_dparameter (double *const &parameter_pt)
 
void unset_analytic_dparameter (double *const &parameter_pt)
 
bool is_dparameter_calculated_analytically (double *const &parameter_pt)
 
void set_analytic_hessian_products ()
 
void unset_analytic_hessian_products ()
 
bool are_hessian_products_calculated_analytically ()
 
void set_pinned_values_to_zero ()
 
bool distributed () const
 
virtual void actions_before_adapt ()
 
virtual void actions_after_adapt ()
 Actions that are to be performed after a mesh adaptation. More...
 
OomphCommunicatorcommunicator_pt ()
 access function to the oomph-lib communicator More...
 
const OomphCommunicatorcommunicator_pt () const
 access function to the oomph-lib communicator, const version More...
 
 Problem ()
 
 Problem (const Problem &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Problem &)=delete
 Broken assignment operator. More...
 
virtual ~Problem ()
 Virtual destructor to clean up memory. More...
 
Mesh *& mesh_pt ()
 Return a pointer to the global mesh. More...
 
Mesh *const & mesh_pt () const
 Return a pointer to the global mesh (const version) More...
 
Mesh *& mesh_pt (const unsigned &imesh)
 
Mesh *const & mesh_pt (const unsigned &imesh) const
 Return a pointer to the i-th submesh (const version) More...
 
unsigned nsub_mesh () const
 Return number of submeshes. More...
 
unsigned add_sub_mesh (Mesh *const &mesh_pt)
 
void flush_sub_meshes ()
 
void build_global_mesh ()
 
void rebuild_global_mesh ()
 
LinearSolver *& linear_solver_pt ()
 Return a pointer to the linear solver object. More...
 
LinearSolver *const & linear_solver_pt () const
 Return a pointer to the linear solver object (const version) More...
 
LinearSolver *& mass_matrix_solver_for_explicit_timestepper_pt ()
 
LinearSolvermass_matrix_solver_for_explicit_timestepper_pt () const
 
EigenSolver *& eigen_solver_pt ()
 Return a pointer to the eigen solver object. More...
 
EigenSolver *const & eigen_solver_pt () const
 Return a pointer to the eigen solver object (const version) More...
 
Time *& time_pt ()
 Return a pointer to the global time object. More...
 
Timetime_pt () const
 Return a pointer to the global time object (const version). More...
 
doubletime ()
 Return the current value of continuous time. More...
 
double time () const
 Return the current value of continuous time (const version) More...
 
TimeStepper *& time_stepper_pt ()
 
const TimeSteppertime_stepper_pt () const
 
TimeStepper *& time_stepper_pt (const unsigned &i)
 Return a pointer to the i-th timestepper. More...
 
ExplicitTimeStepper *& explicit_time_stepper_pt ()
 Return a pointer to the explicit timestepper. More...
 
unsigned long set_timestepper_for_all_data (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data=false)
 
virtual void shift_time_values ()
 Shift all values along to prepare for next timestep. More...
 
AssemblyHandler *& assembly_handler_pt ()
 Return a pointer to the assembly handler object. More...
 
AssemblyHandler *const & assembly_handler_pt () const
 Return a pointer to the assembly handler object (const version) More...
 
doubleminimum_dt ()
 Access function to min timestep in adaptive timestepping. More...
 
doublemaximum_dt ()
 Access function to max timestep in adaptive timestepping. More...
 
unsignedmax_newton_iterations ()
 Access function to max Newton iterations before giving up. More...
 
void problem_is_nonlinear (const bool &prob_lin)
 Access function to Problem_is_nonlinear. More...
 
doublemax_residuals ()
 
booltime_adaptive_newton_crash_on_solve_fail ()
 Access function for Time_adaptive_newton_crash_on_solve_fail. More...
 
doublenewton_solver_tolerance ()
 
void add_time_stepper_pt (TimeStepper *const &time_stepper_pt)
 
void set_explicit_time_stepper_pt (ExplicitTimeStepper *const &explicit_time_stepper_pt)
 
void initialise_dt (const double &dt)
 
void initialise_dt (const Vector< double > &dt)
 
Data *& global_data_pt (const unsigned &i)
 Return a pointer to the the i-th global data object. More...
 
void add_global_data (Data *const &global_data_pt)
 
void flush_global_data ()
 
LinearAlgebraDistribution *const & dof_distribution_pt () const
 Return the pointer to the dof distribution (read-only) More...
 
unsigned long ndof () const
 Return the number of dofs. More...
 
unsigned ntime_stepper () const
 Return the number of time steppers. More...
 
unsigned nglobal_data () const
 Return the number of global data values. More...
 
unsigned self_test ()
 Self-test: Check meshes and global data. Return 0 for OK. More...
 
void enable_store_local_dof_pt_in_elements ()
 
void disable_store_local_dof_pt_in_elements ()
 
unsigned long assign_eqn_numbers (const bool &assign_local_eqn_numbers=true)
 
void describe_dofs (std::ostream &out= *(oomph_info.stream_pt())) const
 
void enable_discontinuous_formulation ()
 
void disable_discontinuous_formulation ()
 
void get_dofs (DoubleVector &dofs) const
 
void get_dofs (const unsigned &t, DoubleVector &dofs) const
 Return vector of the t'th history value of all dofs. More...
 
void set_dofs (const DoubleVector &dofs)
 Set the values of the dofs. More...
 
void set_dofs (const unsigned &t, DoubleVector &dofs)
 Set the history values of the dofs. More...
 
void set_dofs (const unsigned &t, Vector< double * > &dof_pt)
 
void add_to_dofs (const double &lambda, const DoubleVector &increment_dofs)
 Add lambda x incremenet_dofs[l] to the l-th dof. More...
 
doubleglobal_dof_pt (const unsigned &i)
 
doubledof (const unsigned &i)
 i-th dof in the problem More...
 
double dof (const unsigned &i) const
 i-th dof in the problem (const version) More...
 
double *& dof_pt (const unsigned &i)
 Pointer to i-th dof in the problem. More...
 
doubledof_pt (const unsigned &i) const
 Pointer to i-th dof in the problem (const version) More...
 
virtual void get_inverse_mass_matrix_times_residuals (DoubleVector &Mres)
 
virtual void get_dvaluesdt (DoubleVector &f)
 
virtual void get_residuals (DoubleVector &residuals)
 Get the total residuals Vector for the problem. More...
 
virtual void get_jacobian (DoubleVector &residuals, DenseDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CRDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CCDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, SumOfMatrices &jacobian)
 
void get_fd_jacobian (DoubleVector &residuals, DenseMatrix< double > &jacobian)
 Get the full Jacobian by finite differencing. More...
 
void get_derivative_wrt_global_parameter (double *const &parameter_pt, DoubleVector &result)
 
void get_hessian_vector_products (DoubleVectorWithHaloEntries const &Y, Vector< DoubleVectorWithHaloEntries > const &C, Vector< DoubleVectorWithHaloEntries > &product)
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &steady=true)
 Solve the eigenproblem. More...
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true)
 
virtual void get_eigenproblem_matrices (CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0)
 
void assign_eigenvector_to_dofs (DoubleVector &eigenvector)
 Assign the eigenvector passed to the function to the dofs. More...
 
void add_eigenvector_to_dofs (const double &epsilon, const DoubleVector &eigenvector)
 
void store_current_dof_values ()
 Store the current values of the degrees of freedom. More...
 
void restore_dof_values ()
 Restore the stored values of the degrees of freedom. More...
 
void enable_jacobian_reuse ()
 
void disable_jacobian_reuse ()
 Disable recycling of Jacobian in Newton iteration. More...
 
bool jacobian_reuse_is_enabled ()
 Is recycling of Jacobian in Newton iteration enabled? More...
 
booluse_predictor_values_as_initial_guess ()
 
void newton_solve ()
 Use Newton method to solve the problem. More...
 
void enable_globally_convergent_newton_method ()
 enable globally convergent Newton method More...
 
void disable_globally_convergent_newton_method ()
 disable globally convergent Newton method More...
 
void newton_solve (unsigned const &max_adapt)
 
void steady_newton_solve (unsigned const &max_adapt=0)
 
void copy (Problem *orig_problem_pt)
 
virtual Problemmake_copy ()
 
virtual void read (std::ifstream &restart_file, bool &unsteady_restart)
 
virtual void read (std::ifstream &restart_file)
 
virtual void dump (std::ofstream &dump_file) const
 
void dump (const std::string &dump_file_name) const
 
void delete_all_external_storage ()
 
virtual void symmetrise_eigenfunction_for_adaptive_pitchfork_tracking ()
 
doublebifurcation_parameter_pt () const
 
void get_bifurcation_eigenfunction (Vector< DoubleVector > &eigenfunction)
 
void activate_fold_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation, const bool &block_solve=true)
 
void activate_pitchfork_tracking (double *const &parameter_pt, const DoubleVector &symmetry_vector, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const double &omega, const DoubleVector &null_real, const DoubleVector &null_imag, const bool &block_solve=true)
 
void deactivate_bifurcation_tracking ()
 
void reset_assembly_handler_to_default ()
 Reset the system to the standard non-augemented state. More...
 
double arc_length_step_solve (double *const &parameter_pt, const double &ds, const unsigned &max_adapt=0)
 
double arc_length_step_solve (Data *const &data_pt, const unsigned &data_index, const double &ds, const unsigned &max_adapt=0)
 
void reset_arc_length_parameters ()
 
intsign_of_jacobian ()
 
void explicit_timestep (const double &dt, const bool &shift_values=true)
 Take an explicit timestep of size dt. More...
 
void unsteady_newton_solve (const double &dt)
 
void unsteady_newton_solve (const double &dt, const bool &shift_values)
 
void unsteady_newton_solve (const double &dt, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const unsigned &suppress_resolve_after_spatial_adapt_flag, const bool &first, const bool &shift=true)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon, const bool &shift_values)
 
void assign_initial_values_impulsive ()
 
void assign_initial_values_impulsive (const double &dt)
 
void calculate_predictions ()
 Calculate predictions. More...
 
void enable_mass_matrix_reuse ()
 
void disable_mass_matrix_reuse ()
 
bool mass_matrix_reuse_is_enabled ()
 Return whether the mass matrix is being reused. More...
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly (DocInfo &doc_info)
 
void refine_uniformly_and_prune (DocInfo &doc_info)
 
void refine_uniformly ()
 
void refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform refinement for submesh i_mesh with documentation. More...
 
void refine_uniformly (const unsigned &i_mesh)
 Do uniform refinement for submesh i_mesh without documentation. More...
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly (DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (DocInfo &doc_info)
 
void p_refine_uniformly ()
 
void p_refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-refinement for submesh i_mesh with documentation. More...
 
void p_refine_uniformly (const unsigned &i_mesh)
 Do uniform p-refinement for submesh i_mesh without documentation. More...
 
void refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void refine_selected_elements (const Vector< Vector< RefineableElement * >> &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< Vector< PRefineableElement * >> &elements_to_be_refined_pt)
 
unsigned unrefine_uniformly ()
 
unsigned unrefine_uniformly (const unsigned &i_mesh)
 
void p_unrefine_uniformly (DocInfo &doc_info)
 
void p_unrefine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-unrefinement for submesh i_mesh without documentation. More...
 
void adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void adapt ()
 
void p_adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void p_adapt ()
 
void adapt_based_on_error_estimates (unsigned &n_refined, unsigned &n_unrefined, Vector< Vector< double >> &elemental_error)
 
void adapt_based_on_error_estimates (Vector< Vector< double >> &elemental_error)
 
void get_all_error_estimates (Vector< Vector< double >> &elemental_error)
 
void doc_errors (DocInfo &doc_info)
 Get max and min error for all elements in submeshes. More...
 
void doc_errors ()
 Get max and min error for all elements in submeshes. More...
 
void enable_info_in_newton_solve ()
 
void disable_info_in_newton_solve ()
 Disable the output of information when in the newton solver. More...
 
- Public Member Functions inherited from oomph::ExplicitTimeSteppableObject
 ExplicitTimeSteppableObject ()
 Empty constructor. More...
 
 ExplicitTimeSteppableObject (const ExplicitTimeSteppableObject &)=delete
 Broken copy constructor. More...
 
void operator= (const ExplicitTimeSteppableObject &)=delete
 Broken assignment operator. More...
 
virtual ~ExplicitTimeSteppableObject ()
 Empty destructor. More...
 
virtual void actions_before_explicit_stage ()
 
virtual void actions_after_explicit_stage ()
 

Public Attributes

MeshBulk_element_mesh_pt
 
MeshFace_element_mesh_pt
 mesh for face elements More...
 
- Public Attributes inherited from oomph::Problem
bool Shut_up_in_newton_solve
 

Protected Member Functions

void build_bulk_mesh (const unsigned n_x, const unsigned n_y, TopologicallyRectangularDomain *domain_pt, HermiteQuadMesh< BiharmonicElement< 2 >>::MeshSpacingFnPtr spacing_fn=0)
 
void build_global_mesh_and_assign_eqn_numbers ()
 Build global mesh and assign equation numbers. More...
 
void set_source_function (const BiharmonicSourceFctPt source_fct_pt)
 
void set_dirichlet_boundary_condition (const unsigned &b, DirichletBCFctPt u_fn=0, DirichletBCFctPt dudn_fn=0)
 
void set_neumann_boundary_condition (const unsigned &b, BiharmonicFluxElement< 2 >::FluxFctPt flux0_fct_pt, BiharmonicFluxElement< 2 >::FluxFctPt flux1_fct_pt=0)
 
- Protected Member Functions inherited from oomph::Problem
unsigned setup_element_count_per_dof ()
 
virtual void sparse_assemble_row_or_column_compressed (Vector< int * > &column_or_row_index, Vector< int * > &row_or_column_start, Vector< double * > &value, Vector< unsigned > &nnz, Vector< double * > &residual, bool compressed_row_flag)
 
virtual void actions_before_newton_convergence_check ()
 
virtual void actions_before_newton_step ()
 
virtual void actions_after_newton_step ()
 
virtual void actions_before_implicit_timestep ()
 
virtual void actions_after_implicit_timestep ()
 
virtual void actions_after_implicit_timestep_and_error_estimation ()
 
virtual void actions_before_explicit_timestep ()
 Actions that should be performed before each explicit time step. More...
 
virtual void actions_after_explicit_timestep ()
 Actions that should be performed after each explicit time step. More...
 
virtual void actions_before_read_unstructured_meshes ()
 
virtual void actions_after_read_unstructured_meshes ()
 
virtual void actions_after_change_in_global_parameter (double *const &parameter_pt)
 
virtual void actions_after_change_in_bifurcation_parameter ()
 
virtual void actions_after_parameter_increase (double *const &parameter_pt)
 
doubledof_derivative (const unsigned &i)
 
doubledof_current (const unsigned &i)
 
virtual void set_initial_condition ()
 
virtual double global_temporal_error_norm ()
 
unsigned newton_solve_continuation (double *const &parameter_pt)
 
unsigned newton_solve_continuation (double *const &parameter_pt, DoubleVector &z)
 
void calculate_continuation_derivatives (double *const &parameter_pt)
 
void calculate_continuation_derivatives (const DoubleVector &z)
 
void calculate_continuation_derivatives_fd (double *const &parameter_pt)
 
bool does_pointer_correspond_to_problem_data (double *const &parameter_pt)
 
void set_consistent_pinned_values_for_continuation ()
 

Additional Inherited Members

- Static Public Attributes inherited from oomph::Problem
static bool Suppress_warning_about_actions_before_read_unstructured_meshes
 
- Protected Types inherited from oomph::Problem
enum  Assembly_method {
  Perform_assembly_using_vectors_of_pairs , Perform_assembly_using_two_vectors , Perform_assembly_using_maps , Perform_assembly_using_lists ,
  Perform_assembly_using_two_arrays
}
 Enumerated flags to determine which sparse assembly method is used. More...
 
- Protected Attributes inherited from oomph::Problem
Vector< Problem * > Copy_of_problem_pt
 
std::map< double *, boolCalculate_dparameter_analytic
 
bool Calculate_hessian_products_analytic
 
LinearAlgebraDistributionDof_distribution_pt
 
Vector< double * > Dof_pt
 Vector of pointers to dofs. More...
 
DoubleVectorWithHaloEntries Element_count_per_dof
 
double Relaxation_factor
 
double Newton_solver_tolerance
 
unsigned Max_newton_iterations
 Maximum number of Newton iterations. More...
 
unsigned Nnewton_iter_taken
 
Vector< doubleMax_res
 Maximum residuals at start and after each newton iteration. More...
 
double Max_residuals
 
bool Time_adaptive_newton_crash_on_solve_fail
 
bool Jacobian_reuse_is_enabled
 Is re-use of Jacobian in Newton iteration enabled? Default: false. More...
 
bool Jacobian_has_been_computed
 
bool Problem_is_nonlinear
 
bool Pause_at_end_of_sparse_assembly
 
bool Doc_time_in_distribute
 
unsigned Sparse_assembly_method
 
unsigned Sparse_assemble_with_arrays_initial_allocation
 
unsigned Sparse_assemble_with_arrays_allocation_increment
 
Vector< Vector< unsigned > > Sparse_assemble_with_arrays_previous_allocation
 
double Numerical_zero_for_sparse_assembly
 
double FD_step_used_in_get_hessian_vector_products
 
bool Mass_matrix_reuse_is_enabled
 
bool Mass_matrix_has_been_computed
 
bool Discontinuous_element_formulation
 
double Minimum_dt
 Minimum desired dt: if dt falls below this value, exit. More...
 
double Maximum_dt
 Maximum desired dt. More...
 
double DTSF_max_increase
 
double DTSF_min_decrease
 
double Minimum_dt_but_still_proceed
 
bool Scale_arc_length
 Boolean to control whether arc-length should be scaled. More...
 
double Desired_proportion_of_arc_length
 Proportion of the arc-length to taken by the parameter. More...
 
double Theta_squared
 
int Sign_of_jacobian
 Storage for the sign of the global Jacobian. More...
 
double Continuation_direction
 
double Parameter_derivative
 Storage for the derivative of the global parameter wrt arc-length. More...
 
double Parameter_current
 Storage for the present value of the global parameter. More...
 
bool Use_continuation_timestepper
 Boolean to control original or new storage of dof stuff. More...
 
unsigned Dof_derivative_offset
 
unsigned Dof_current_offset
 
Vector< doubleDof_derivative
 Storage for the derivative of the problem variables wrt arc-length. More...
 
Vector< doubleDof_current
 Storage for the present values of the variables. More...
 
double Ds_current
 Storage for the current step value. More...
 
unsigned Desired_newton_iterations_ds
 
double Minimum_ds
 Minimum desired value of arc-length. More...
 
bool Bifurcation_detection
 Boolean to control bifurcation detection via determinant of Jacobian. More...
 
bool Bisect_to_find_bifurcation
 Boolean to control wheter bisection is used to located bifurcation. More...
 
bool First_jacobian_sign_change
 Boolean to indicate whether a sign change has occured in the Jacobian. More...
 
bool Arc_length_step_taken
 Boolean to indicate whether an arc-length step has been taken. More...
 
bool Use_finite_differences_for_continuation_derivatives
 
OomphCommunicatorCommunicator_pt
 The communicator for this problem. More...
 
bool Always_take_one_newton_step
 
double Timestep_reduction_factor_after_nonconvergence
 
bool Keep_temporal_error_below_tolerance
 
- Static Protected Attributes inherited from oomph::Problem
static ContinuationStorageScheme Continuation_time_stepper
 Storage for the single static continuation timestorage object. More...
 

Detailed Description

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

Biharmonic Plate Problem Class - for problems where the load can be assumed to be acting normal to the surface of the plate and the deflections are small relative to the thickness of the plate. Developed for the topologically rectangular Hermite Element Mesh. Contains functions allowing the following boundary conditions to be applied (on a given edge):

  • clamped : u and du/dn imposed
  • simply supported : u and laplacian(u) imposed
  • free : laplacian(u) and dlaplacian(u)/dn imposed

Member Typedef Documentation

◆ BiharmonicSourceFctPt

template<unsigned DIM>
typedef void(* oomph::BiharmonicProblem< DIM >::BiharmonicSourceFctPt) (const Vector< double > &x, double &f)

Definition of the Source Function.

◆ DirichletBCFctPt

template<unsigned DIM>
typedef void(* oomph::BiharmonicProblem< DIM >::DirichletBCFctPt) (const double &s, double &u)

Definition of a dirichlet boundary condition function pointer. Takes the position along a boundary (s) in the macro element coordinate scheme and returns the value of the boundary condition at that point (u).

Constructor & Destructor Documentation

◆ BiharmonicProblem()

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

Constructor.

81  {
84  }
Mesh * Bulk_element_mesh_pt
Definition: biharmonic_problem.h:201
Mesh * Face_element_mesh_pt
mesh for face elements
Definition: biharmonic_problem.h:204

References oomph::BiharmonicProblem< DIM >::Bulk_element_mesh_pt, and oomph::BiharmonicProblem< DIM >::Face_element_mesh_pt.

◆ ~BiharmonicProblem()

template<unsigned DIM>
virtual oomph::BiharmonicProblem< DIM >::~BiharmonicProblem ( )
inlinevirtual

Member Function Documentation

◆ actions_after_newton_solve()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::actions_after_newton_solve ( )
inlinevirtual

action after solve

Reimplemented from oomph::Problem.

109 {}

◆ actions_before_newton_solve()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::actions_before_newton_solve ( )
inlinevirtual

actions before solve, performs self test

Reimplemented from oomph::Problem.

95  {
96 #ifdef PARANOID
97  if (0 == self_test())
98  {
99  oomph_info << "self test passed" << std::endl;
100  }
101  else
102  {
103  oomph_info << "self test failed" << std::endl;
104  }
105 #endif
106  }
unsigned self_test()
Self-test: Check meshes and global data. Return 0 for OK.
Definition: problem.cc:13276
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::oomph_info, and oomph::Problem::self_test().

◆ build_bulk_mesh()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::build_bulk_mesh ( const unsigned  n_x,
const unsigned  n_y,
TopologicallyRectangularDomain domain_pt,
HermiteQuadMesh< BiharmonicElement< 2 >>::MeshSpacingFnPtr  spacing_fn = 0 
)
inlineprotected

builds the bulk mesh on a prescribed domain with a node spacing defined by spacing fn with n_x by n_y elements

131  {
132  if (spacing_fn == 0)
133  {
135  new HermiteQuadMesh<Hijacked<BiharmonicElement<2>>>(
136  n_x, n_y, domain_pt);
137  }
138  else
139  {
141  new HermiteQuadMesh<Hijacked<BiharmonicElement<2>>>(
142  n_x, n_y, domain_pt, spacing_fn);
143  }
144  // add_sub_mesh(Bulk_element_mesh_pt);
145  // build_global_mesh();
146  }

References oomph::BiharmonicProblem< DIM >::Bulk_element_mesh_pt.

◆ build_global_mesh_and_assign_eqn_numbers()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::build_global_mesh_and_assign_eqn_numbers ( )
inlineprotected

Build global mesh and assign equation numbers.

151  {
153  if (Face_element_mesh_pt != 0)
154  {
156  }
159  }
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989

References oomph::Problem::add_sub_mesh(), oomph::Problem::assign_eqn_numbers(), oomph::Problem::build_global_mesh(), oomph::BiharmonicProblem< DIM >::Bulk_element_mesh_pt, and oomph::BiharmonicProblem< DIM >::Face_element_mesh_pt.

◆ bulk_element_mesh_pt()

template<unsigned DIM>
Mesh* oomph::BiharmonicProblem< DIM >::bulk_element_mesh_pt ( )
inline

Access function to the bulk element mesh pt.

119  {
120  return Bulk_element_mesh_pt;
121  }

References oomph::BiharmonicProblem< DIM >::Bulk_element_mesh_pt.

◆ doc_solution()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::doc_solution ( DocInfo doc_info,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt = 0 
)

documents the solution, and if an exact solution is provided, then the error between the numerical and exact solution is presented

343  {
344  std::ofstream some_file;
345  std::ostringstream filename;
346 
347  // Number of plot points: npts x npts
348  unsigned npts = 5;
349 
350  // Output solution
351  filename << doc_info.directory() << "/soln_" << doc_info.number() << ".dat";
352  some_file.open(filename.str().c_str());
353  Bulk_element_mesh_pt->output(some_file, npts);
354  some_file.close();
355 
356  // if the exact solution is not provided
357  if (exact_soln_pt != 0)
358  {
359  // Output exact solution
360  filename.str("");
361  filename << doc_info.directory() << "/exact_soln_" << doc_info.number()
362  << ".dat";
363  some_file.open(filename.str().c_str());
364  Bulk_element_mesh_pt->output_fct(some_file, npts, exact_soln_pt);
365  some_file.close();
366 
367  // Doc error and return of the square of the L2 error
368  double error, norm;
369  filename.str("");
370  filename << doc_info.directory() << "/error_" << doc_info.number()
371  << ".dat";
372  some_file.open(filename.str().c_str());
374  some_file, exact_soln_pt, error, norm);
375  some_file.close();
376 
377  // Doc L2 error and norm of solution
378  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
379  oomph_info << "Norm of solution: " << sqrt(norm) << std::endl
380  << std::endl;
381  }
382  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
Output a given Vector function at f(n_plot) points in each element.
Definition: mesh.cc:2199
virtual void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: mesh.h:1140
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
string filename
Definition: MergeRestartFiles.py:39
int error
Definition: calibrate.py:297

References oomph::DocInfo::directory(), calibrate::error, MergeRestartFiles::filename, oomph::DocInfo::number(), oomph::oomph_info, and sqrt().

◆ set_dirichlet_boundary_condition()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::set_dirichlet_boundary_condition ( const unsigned b,
DirichletBCFctPt  u_fn = 0,
DirichletBCFctPt  dudn_fn = 0 
)
protected

Imposes the prescribed dirichlet BCs u (u_fn) and du/dn (dudn_fn) dirichlet BCs by 'pinning'

Impose a Clamped Edge. Imposes the prescribed dirichlet BCs u and du/dn dirichlet BCs by pinning

44  {
45  // number of nodes on boundary b
46  unsigned n_node = Bulk_element_mesh_pt->nboundary_node(b);
47 
48  // fixed faced index for boundary
49  int face_index = Bulk_element_mesh_pt->face_index_at_boundary(b, 0);
50 
51  // Need to get the s_fixed_index
52  unsigned s_fixed_index = 0;
53  switch (face_index)
54  {
55  case -1:
56  case 1:
57  s_fixed_index = 0;
58  break;
59 
60  case -2:
61  case 2:
62  s_fixed_index = 1;
63  break;
64 
65  default:
66  throw OomphLibError("Face Index not +/-1 or +/-2: Need 2D QElements",
69  }
70 
71  // node position along edge b in macro element boundary representation
72  // [-1,1]
73  Vector<double> s(1);
74 
75  // Set the edge sign
76  int edge_sign = 0;
77  switch (face_index)
78  {
79  case -1:
80  case 2:
81  edge_sign = -1;
82  break;
83 
84  case 1:
85  case -2:
86  edge_sign = 1;
87  break;
88  }
89 
90  // finite difference step
91  const double h = 10e-8;
92 
93  // node position along edge b in macro element boundary representation
94  // [-1,1]
95  Vector<double> m(2);
96 
97  // if u is prescribed then impose
98  if (u_fn != 0)
99  {
100  // loop over nodes on boundary b
101  for (unsigned n = 0; n < n_node; n++)
102  {
103  // find node position along edge [-1,1] in macro element representation
106 
107  // get u at node
108  double u;
109  (*u_fn)(m[0], u);
110 
111  // finite difference is used to compute dudm_t
112  // u0 and u1 store values of u left/below and right/above node n
113  double u_L, u_R;
114 
115  // if left/lower corner node
116  if (n == 0)
117  {
118  (*u_fn)(m[0], u_L);
119  (*u_fn)(m[0] + h, u_R);
120  }
121  // if right/upper corner node
122  else if (n == n_node - 1)
123  {
124  (*u_fn)(m[0] - h, u_L);
125  (*u_fn)(m[0], u_R);
126  }
127  // if other node
128  else
129  {
130  (*u_fn)(m[0] - 0.5 * h, u_L);
131  (*u_fn)(m[0] + 0.5 * h, u_R);
132  }
133 
134  // compute dudm_t
135  double dudm_t = (u_R - u_L) / h;
136 
137  // compute duds_t
138  double duds_t = m[1] * dudm_t;
139 
140  // pin and set u type dof
143 
144  // pin and set duds_t type degree of freedom
145  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(2 - s_fixed_index);
147  2 - s_fixed_index, duds_t);
148  }
149  }
150 
151  // if dudn is prescribed then impose
152  if (dudn_fn != 0)
153  {
154  // loop over nodes on boundary b
155  for (unsigned n = 0; n < n_node; n++)
156  {
157  // vectors for dx_i/ds_n and dx_i/ds_t
158  Vector<double> dxds_n(2);
159  Vector<double> dxds_t(2);
161  1 + s_fixed_index, 0);
163  1 + s_fixed_index, 1);
165  2 - s_fixed_index, 0);
167  2 - s_fixed_index, 1);
168 
169  // vector for d2xi/ds_n ds_t
170  Vector<double> d2xds_nds_t(2);
171  d2xds_nds_t[0] =
173  d2xds_nds_t[1] =
175 
176  // compute dn/ds_n
177  double dnds_n =
178  ((dxds_n[0] * dxds_t[1]) - (dxds_n[1] * dxds_t[0])) /
179  (sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]) * edge_sign);
180 
181  // compute dt/ds_n
182  double dtds_n = ((dxds_n[0] * dxds_t[0]) + (dxds_n[1] * dxds_t[1])) /
183  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
184 
185  // compute dt/ds_t
186  double dtds_t = sqrt(pow(dxds_t[0], 2) + pow(dxds_t[1], 2));
187 
188  // compute ds_n/dn
189  double ds_ndn =
190  -1.0 * (edge_sign * sqrt(pow(dxds_t[0], 2) + pow(dxds_t[1], 2))) /
191  (dxds_t[0] * dxds_n[1] - dxds_n[0] * dxds_t[1]);
192 
193  // compute ds_t/dt
194  double ds_tdt = 1 / sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
195 
196  // compute d2t/ds_nds_t
197  double d2tds_nds_t =
198  (dxds_t[0] * d2xds_nds_t[0] + dxds_t[1] * d2xds_nds_t[1]) /
199  sqrt(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1]);
200 
201  // compute d2s_t/ds_ndt
202  double d2s_tds_ndt =
203  (dxds_t[0] * d2xds_nds_t[0] + dxds_t[1] * d2xds_nds_t[1]) /
204  pow(dxds_t[0] * dxds_t[0] + dxds_t[1] * dxds_t[1], 3.0 / 2.0);
205 
206  // get m_t and dm_t/ds_t for this node
207  Vector<double> m_N(2);
210 
211  // compute d2u/dm_t2 and d/dm_t(dudn) by finite difference
212  double d2udm_t2 = 0;
213  double ddm_tdudn = 0;
214  double u_2L, u_N, u_2R, dudn_L, dudn_R;
215  // evaluate u_fn and dudn_fn for current node
216  if (n == 0)
217  {
218  (*u_fn)(m_N[0], u_2L);
219  (*u_fn)(m_N[0] + h, u_N);
220  (*u_fn)(m_N[0] - 2 * h, u_2R);
221  (*dudn_fn)(m_N[0], dudn_L);
222  (*dudn_fn)(m_N[0] + h, dudn_R);
223  }
224  else if (n == n_node - 1)
225  {
226  (*u_fn)(m_N[0] - 2 * h, u_2L);
227  (*u_fn)(m_N[0] - h, u_N);
228  (*u_fn)(m_N[0], u_2R);
229  (*dudn_fn)(m_N[0] - h, dudn_L);
230  (*dudn_fn)(m_N[0], dudn_R);
231  }
232  else
233  {
234  (*u_fn)(m_N[0] - h, u_2L);
236  (*u_fn)(m_N[0] + h, u_2R);
237  (*dudn_fn)(m_N[0] - 0.5 * h, dudn_L);
238  (*dudn_fn)(m_N[0] + 0.5 * h, dudn_R);
239  }
240  // compute
241  d2udm_t2 = (u_2L + u_2R - 2 * u_N) / (h * h);
242  ddm_tdudn = (dudn_R - dudn_L) / h;
243 
244  // get dudn at the node
245  double dudn;
246  (*dudn_fn)(m_N[0], dudn);
247 
248  // compute d2u/ds_t2
249  double d2uds_t2 = m_N[1] * m_N[1] * d2udm_t2;
250 
251  // get duds_t
252  double duds_t = Bulk_element_mesh_pt->boundary_node_pt(b, n)->value(
253  2 - s_fixed_index);
254 
255  // get du/dt
256  double dudt = ds_tdt * duds_t;
257 
258  // compute d2u/dndt
259  double d2udndt = ds_tdt = dtds_t * m_N[1] * ddm_tdudn;
260 
261  // compute d2udt2
262  double dtds_nd2udt2 =
263  edge_sign * (dxds_t[0] * dxds_n[1] - dxds_n[0] * dxds_t[1]) *
264  (ds_tdt *
265  (d2udndt - ds_ndn * (d2s_tds_ndt * dudt + ds_tdt * d2uds_t2)));
266 
267  // compute dds_n(dudt)
268  double dds_ndudt = dtds_nd2udt2 + dnds_n * d2udndt;
269 
270  // compute du/ds_n
271  double duds_n = dnds_n * dudn + dtds_n * ds_tdt * duds_t;
272 
273  // compute d2u/ds_nds_t
274  double d2uds_nds_t = d2tds_nds_t * dudt + dtds_t * dds_ndudt;
275 
276  // pin du/ds_n dof and set value
277  Bulk_element_mesh_pt->boundary_node_pt(b, n)->pin(1 + s_fixed_index);
279  1 + s_fixed_index, duds_n);
280 
281  // pin d2u/ds_nds_t dof and set value
283  Bulk_element_mesh_pt->boundary_node_pt(b, n)->set_value(3, d2uds_nds_t);
284  }
285  }
286  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Scalar * b
Definition: benchVecAdd.cpp:17
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
unsigned long nboundary_node(const unsigned &ibound) const
Return number of nodes on a particular boundary.
Definition: mesh.h:833
int face_index_at_boundary(const unsigned &b, const unsigned &e) const
Definition: mesh.h:896
Node *& boundary_node_pt(const unsigned &b, const unsigned &n)
Return pointer to node n on boundary b.
Definition: mesh.h:493
double & x_gen(const unsigned &k, const unsigned &i)
Definition: nodes.h:1126
virtual void get_coordinates_on_boundary(const unsigned &b, const unsigned &k, Vector< double > &boundary_zeta)
Definition: nodes.cc:2379
double value(const unsigned &i) const
Definition: nodes.cc:2408
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
void u_N(const double &s, double &u)
Definition: two_d_biharmonic.cc:293
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References b, e(), m, n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and s.

◆ set_neumann_boundary_condition()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::set_neumann_boundary_condition ( const unsigned b,
BiharmonicFluxElement< 2 >::FluxFctPt  flux0_fct_pt,
BiharmonicFluxElement< 2 >::FluxFctPt  flux1_fct_pt = 0 
)
protected

Imposes the prescribed Neumann BCs laplacian(u) (flux0_fct_pt) and dlaplacian(u)/dn (flux1_fct_pt) with flux edge elements

Imposes a 'free' edge. Imposes the prescribed Neumann BCs laplacian(u) and dlaplacian(u)/dn with flux edge elements

298  {
299  // if the face element mesh pt does not exist then build it
300  if (Face_element_mesh_pt == 0)
301  {
302  Face_element_mesh_pt = new Mesh();
303  }
304 
305  // How many bulk elements are adjacent to boundary b?
306  unsigned n_element = Bulk_element_mesh_pt->nboundary_element(b);
307 
308  // Loop over the bulk elements adjacent to boundary b?
309  for (unsigned e = 0; e < n_element; e++)
310  {
311  // Get pointer to the bulk element that is adjacent to boundary b
312  FiniteElement* bulk_elem_pt = dynamic_cast<FiniteElement*>(
314 
315  // What is the face index along the boundary
316  int face_index = Bulk_element_mesh_pt->face_index_at_boundary(b, e);
317 
318  // Build the corresponding prescribed-flux element
319  BiharmonicFluxElement<2>* flux_element_pt =
320  new BiharmonicFluxElement<2>(bulk_elem_pt, face_index, b);
321 
322 
323  // pass the flux BC pointers to the flux elements
324  flux_element_pt->flux0_fct_pt() = flux0_fct_pt;
325  if (flux1_fct_pt != 0)
326  {
327  flux_element_pt->flux1_fct_pt() = flux1_fct_pt;
328  }
329 
330  // Add the prescribed-flux element to the mesh
331  Face_element_mesh_pt->add_element_pt(flux_element_pt);
332  }
333  }
unsigned nboundary_element(const unsigned &b) const
Return number of finite elements that are adjacent to boundary b.
Definition: mesh.h:878
FiniteElement * boundary_element_pt(const unsigned &b, const unsigned &e) const
Return pointer to e-th finite element on boundary b.
Definition: mesh.h:840
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617

References b, e(), oomph::BiharmonicFluxElement< DIM >::flux0_fct_pt(), and oomph::BiharmonicFluxElement< DIM >::flux1_fct_pt().

◆ set_source_function()

template<unsigned DIM>
void oomph::BiharmonicProblem< DIM >::set_source_function ( const BiharmonicSourceFctPt  source_fct_pt)
inlineprotected

Impose a load to the surface of the plate. note : MUST be called before neumann boundary conditions are imposed, i.e. a free edge or a simply supported edge with laplacian(u) imposed

165  {
166  // number of elements in mesh
167  unsigned n_bulk_element = Bulk_element_mesh_pt->nelement();
168 
169  // loop over bulk elements
170  for (unsigned i = 0; i < n_bulk_element; i++)
171  {
172  // upcast from generalised element to specific element
173  BiharmonicElement<2>* element_pt = dynamic_cast<BiharmonicElement<2>*>(
175 
176  // set the source function pointer
177  element_pt->source_fct_pt() = source_fct_pt;
178  }
179  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: gen_axisym_advection_diffusion_elements.h:229

References oomph::BiharmonicProblem< DIM >::Bulk_element_mesh_pt, oomph::Mesh::element_pt(), i, oomph::Mesh::nelement(), oomph::source_fct_pt(), and oomph::BiharmonicEquations< DIM >::source_fct_pt().

Member Data Documentation

◆ Bulk_element_mesh_pt

◆ Face_element_mesh_pt


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