UnstructuredFvKProblem< ELEMENT > Class Template Reference

Class definition. More...

+ Inheritance diagram for UnstructuredFvKProblem< ELEMENT >:

Public Member Functions

 UnstructuredFvKProblem (double element_area=0.2)
 Constructor. More...
 
 ~UnstructuredFvKProblem ()
 Destructor. More...
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void doc_solution ()
 Doc the solution. More...
 
 UnstructuredFvKProblem (double element_area=0.2)
 Constructor. More...
 
 ~UnstructuredFvKProblem ()
 Destructor. More...
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void doc_solution ()
 Doc the solution. More...
 
 UnstructuredFvKProblem (double element_area=0.2)
 Constructor. More...
 
 ~UnstructuredFvKProblem ()
 Destructor. More...
 
void actions_before_adapt ()
 
void actions_after_adapt ()
 
void actions_after_newton_solve ()
 Update after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve: Re-apply boundary conditons. More...
 
void doc_solution (const std::string &comment="")
 Doc the solution. 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
 
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 ()
 

Private Types

enum  {
  Outer_boundary0 = 0 , Outer_boundary1 = 1 , Inner_boundary0 = 2 , Inner_boundary1 = 3 ,
  Inner_boundary2 = 4 , Inner_boundary3 = 5 , Inner_boundary4 = 6 , Inner_boundary5 = 7
}
 
enum  {
  Outer_boundary0 = 0 , Outer_boundary1 = 1 , Inner_boundary0 = 2 , Inner_boundary1 = 3 ,
  Inner_boundary2 = 4 , Inner_boundary3 = 5 , Inner_boundary4 = 6 , Inner_boundary5 = 7
}
 
enum  {
  Outer_boundary0 = 0 , Outer_boundary1 = 1 , Inner_boundary0 = 2 , Inner_boundary1 = 3 ,
  Inner_boundary2 = 4 , Inner_boundary3 = 5 , Inner_boundary4 = 6 , Inner_boundary5 = 7
}
 

Private Member Functions

void complete_problem_setup ()
 
void complete_problem_setup ()
 
void apply_boundary_conditions ()
 Helper function to apply boundary conditions. More...
 
void complete_problem_setup ()
 

Private Attributes

DocInfo Doc_info
 Doc info object for labeling output. More...
 
TriangleMesh< ELEMENT > * My_mesh_pt
 Pointers to specific mesh. More...
 
ofstream Trace_file
 Trace file to document norm of solution. More...
 
double Element_area
 Element area. More...
 
MeshAsGeomObjectMesh_as_geom_object_pt
 Geom object made of mesh. More...
 
Vector< std::pair< GeomObject *, Vector< double > > > Radial_sample_point_pt
 Plot points along radial line are located in this element/local coord. More...
 
RefineableTriangleMesh< ELEMENT > * My_mesh_pt
 Pointers to specific mesh. More...
 
MeshVolume_constraint_mesh_pt
 Mesh for the volume constraint element. More...
 
FoepplvonKarmanVolumeConstraintElement< ELEMENT, RefineableTriangleMesh > * Volume_constraint_element_pt
 Single volume constraint element instance. More...
 
doubleTemp_pressure_pt
 
Vector< unsignedBubble_regions
 

Additional Inherited Members

- 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 Attributes inherited from oomph::Problem
bool Shut_up_in_newton_solve
 
- 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 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 ()
 
- 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<class ELEMENT>
class UnstructuredFvKProblem< ELEMENT >

Class definition.

//////////////////////////////////////////////////////// //////////////////////////////////////////////////////// ////////////////////////////////////////////////////////

Member Enumeration Documentation

◆ anonymous enum

template<class ELEMENT >
anonymous enum
private
Enumerator
Outer_boundary0 
Outer_boundary1 
Inner_boundary0 
Inner_boundary1 
Inner_boundary2 
Inner_boundary3 
Inner_boundary4 
Inner_boundary5 
127  {
128  Outer_boundary0 = 0,
129  Outer_boundary1 = 1,
130  Inner_boundary0 = 2,
131  Inner_boundary1 = 3,
132  Inner_boundary2 = 4,
133  Inner_boundary3 = 5,
134  Inner_boundary4 = 6,
135  Inner_boundary5 = 7
136  };
@ Outer_boundary1
Definition: all_foeppl_von_karman/circular_disk.cc:129
@ Inner_boundary0
Definition: all_foeppl_von_karman/circular_disk.cc:130
@ Outer_boundary0
Definition: all_foeppl_von_karman/circular_disk.cc:128
@ Inner_boundary5
Definition: all_foeppl_von_karman/circular_disk.cc:135
@ Inner_boundary3
Definition: all_foeppl_von_karman/circular_disk.cc:133
@ Inner_boundary4
Definition: all_foeppl_von_karman/circular_disk.cc:134
@ Inner_boundary2
Definition: all_foeppl_von_karman/circular_disk.cc:132
@ Inner_boundary1
Definition: all_foeppl_von_karman/circular_disk.cc:131

◆ anonymous enum

template<class ELEMENT >
anonymous enum
private
Enumerator
Outer_boundary0 
Outer_boundary1 
Inner_boundary0 
Inner_boundary1 
Inner_boundary2 
Inner_boundary3 
Inner_boundary4 
Inner_boundary5 
127  {
128  Outer_boundary0 = 0,
129  Outer_boundary1 = 1,
130  Inner_boundary0 = 2,
131  Inner_boundary1 = 3,
132  Inner_boundary2 = 4,
133  Inner_boundary3 = 5,
134  Inner_boundary4 = 6,
135  Inner_boundary5 = 7
136  };

◆ anonymous enum

template<class ELEMENT >
anonymous enum
private
Enumerator
Outer_boundary0 
Outer_boundary1 
Inner_boundary0 
Inner_boundary1 
Inner_boundary2 
Inner_boundary3 
Inner_boundary4 
Inner_boundary5 
186  {
187  Outer_boundary0 = 0,
188  Outer_boundary1 = 1,
189  Inner_boundary0 = 2,
190  Inner_boundary1 = 3,
191  Inner_boundary2 = 4,
192  Inner_boundary3 = 5,
193  Inner_boundary4 = 6,
194  Inner_boundary5 = 7
195  };

Constructor & Destructor Documentation

◆ UnstructuredFvKProblem() [1/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem ( double  element_area = 0.2)

Constructor.

Constructor: Pass in element area.

Create the initial volume constraint element

Set the problem pressure to be that of the volume constraint element

Set initial pressure guess

Create the volume constraint mesh and add the element to the mesh

Add the sub meshes and build the global mesh

155  : Element_area(element_area)
156 {
157  Vector<double> zeta(1);
158  Vector<double> posn(2);
159 
160  //Outer boundary
161  //--------------
162  Ellipse* outer_boundary_ellipse_pt = new Ellipse(1.0,1.0);
163 
164  TriangleMeshClosedCurve* outer_boundary_pt = 0;
165 
166  Vector<TriangleMeshCurveSection*> outer_curvilinear_boundary_pt(2);
167 
168  //First bit
169  double zeta_start = 0.0;
170  double zeta_end = MathematicalConstants::Pi;
171  unsigned nsegment = (int)(MathematicalConstants::Pi/sqrt(element_area));
172  outer_curvilinear_boundary_pt[0] =
173  new TriangleMeshCurviLine(outer_boundary_ellipse_pt, zeta_start,
174  zeta_end, nsegment, Outer_boundary0);
175 
176  //Second bit
177  zeta_start = MathematicalConstants::Pi;
178  zeta_end = 2.0*MathematicalConstants::Pi;
179  nsegment = (int)(MathematicalConstants::Pi/sqrt(element_area));
180  outer_curvilinear_boundary_pt[1] =
181  new TriangleMeshCurviLine(outer_boundary_ellipse_pt, zeta_start,
182  zeta_end, nsegment, Outer_boundary1);
183 
184  outer_boundary_pt =
185  new TriangleMeshClosedCurve(outer_curvilinear_boundary_pt);
186 
187 
188  //Inner boundary
189  //--------------
190  Ellipse* inner_boundary_ellipse_pt = new Ellipse(GlobalParameters::R_b,
192 
193  Vector<TriangleMeshClosedCurve*> inner_boundary_pt(1);
194  Vector<TriangleMeshCurveSection*> inner_curvilinear_boundary_pt(2);
195 
196  //First part
197  zeta_start = 0.0;
198  zeta_end = MathematicalConstants::Pi;
199  nsegment = (int)(GlobalParameters::R_b*
200  MathematicalConstants::Pi/sqrt(element_area));
201  inner_curvilinear_boundary_pt[0] =
202  new TriangleMeshCurviLine(inner_boundary_ellipse_pt, zeta_start,
203  zeta_end, nsegment, Inner_boundary0);
204 
205  //Second part
206  zeta_start = MathematicalConstants::Pi;
207  zeta_end = 2.0*MathematicalConstants::Pi;
208  nsegment = (int)(GlobalParameters::R_b*
209  MathematicalConstants::Pi/sqrt(element_area));
210  inner_curvilinear_boundary_pt[1] =
211  new TriangleMeshCurviLine(inner_boundary_ellipse_pt, zeta_start,
212  zeta_end, nsegment, Inner_boundary1);
213 
214  //Combine to internal curvilinear boundary
215  Vector<double> left_region_coords(2);
216  Vector<double> right_region_coords(2);
217  left_region_coords[0] = -0.5*GlobalParameters::R_b;
218  left_region_coords[1] = 0.0;
219 
220  right_region_coords[0] = 0.5*GlobalParameters::R_b;
221  right_region_coords[1] = 0.0;
222 
223  inner_boundary_pt[0] =
224  new TriangleMeshClosedCurve(inner_curvilinear_boundary_pt);
225 
226  //Diameter boundary line
227  Vector<TriangleMeshOpenCurve*> inner_open_boundary_pt(4);
228 
229  Vector<Vector<double> > vertex_coord(2, Vector<double>(2));
230 
231  //Section 1 of line
232  vertex_coord[0][0] = 0.0;
233  vertex_coord[0][1] = -1.0;
234  vertex_coord[1][0] = 0.0;
235  vertex_coord[1][1] = -GlobalParameters::R_b;
236 
237  TriangleMeshPolyLine *inner_open_polyline1_pt =
238  new TriangleMeshPolyLine(vertex_coord, Inner_boundary2);
239 
240  inner_open_polyline1_pt->connect_initial_vertex_to_curviline(
241  dynamic_cast<TriangleMeshCurviLine *>(outer_curvilinear_boundary_pt[1]),
243 
244  inner_open_polyline1_pt->connect_final_vertex_to_curviline(
245  dynamic_cast<TriangleMeshCurviLine *>(inner_curvilinear_boundary_pt[1]),
247 
248  Vector<TriangleMeshCurveSection *> inner_boundary_line1_pt(1);
249  inner_boundary_line1_pt[0] = inner_open_polyline1_pt;
250 
251  inner_open_boundary_pt[0] =
252  new TriangleMeshOpenCurve(inner_boundary_line1_pt);
253 
254  //Section 2 of line
255  vertex_coord[0][0] = 0.0;
256  vertex_coord[0][1] = -GlobalParameters::R_b;
257  vertex_coord[1][0] = 0.0;
258  vertex_coord[1][1] = 0.0;
259 
260  TriangleMeshPolyLine *inner_open_polyline2_pt =
261  new TriangleMeshPolyLine(vertex_coord, Inner_boundary3);
262 
263  inner_open_polyline2_pt->connect_initial_vertex_to_curviline(
264  dynamic_cast<TriangleMeshCurviLine *>(inner_curvilinear_boundary_pt[1]),
266 
267  Vector<TriangleMeshCurveSection *> inner_boundary_line2_pt(1);
268  inner_boundary_line2_pt[0] = inner_open_polyline2_pt;
269 
270  inner_open_boundary_pt[1] =
271  new TriangleMeshOpenCurve(inner_boundary_line2_pt);
272 
273  //Section 3 of line
274  vertex_coord[0][0] = 0.0;
275  vertex_coord[0][1] = 0.0;
276  vertex_coord[1][0] = 0.0;
277  vertex_coord[1][1] = GlobalParameters::R_b;
278 
279  TriangleMeshPolyLine *inner_open_polyline3_pt =
280  new TriangleMeshPolyLine(vertex_coord, Inner_boundary4);
281 
282  inner_open_polyline3_pt->connect_initial_vertex_to_polyline(
283  inner_open_polyline2_pt,1);
284 
285  inner_open_polyline3_pt->connect_final_vertex_to_curviline(
286  dynamic_cast<TriangleMeshCurviLine *>(inner_curvilinear_boundary_pt[0]),
288 
289  Vector<TriangleMeshCurveSection *> inner_boundary_line3_pt(1);
290  inner_boundary_line3_pt[0] = inner_open_polyline3_pt;
291 
292  inner_open_boundary_pt[2] =
293  new TriangleMeshOpenCurve(inner_boundary_line3_pt);
294 
295  //Section 4 of line
296  vertex_coord[0][0] = 0.0;
297  vertex_coord[0][1] = GlobalParameters::R_b;
298  vertex_coord[1][0] = 0.0;
299  vertex_coord[1][1] = 1.0;
300 
301  TriangleMeshPolyLine *inner_open_polyline4_pt =
302  new TriangleMeshPolyLine(vertex_coord, Inner_boundary5);
303 
304  inner_open_polyline4_pt->connect_initial_vertex_to_curviline(
305  dynamic_cast<TriangleMeshCurviLine *>(inner_curvilinear_boundary_pt[0]),
307 
308  inner_open_polyline4_pt->connect_final_vertex_to_curviline(
309  dynamic_cast<TriangleMeshCurviLine *>(outer_curvilinear_boundary_pt[0]),
311 
312  Vector<TriangleMeshCurveSection *> inner_boundary_line4_pt(1);
313  inner_boundary_line4_pt[0] = inner_open_polyline4_pt;
314 
315  inner_open_boundary_pt[3] =
316  new TriangleMeshOpenCurve(inner_boundary_line4_pt);
317 
318  //Create the mesh
319  //---------------
320 
321  //Create mesh parameters object
322  TriangleMeshParameters mesh_parameters(outer_boundary_pt);
323  mesh_parameters.internal_closed_curve_pt() = inner_boundary_pt;
324  mesh_parameters.internal_open_curves_pt() = inner_open_boundary_pt;
325  mesh_parameters.add_region_coordinates(1,left_region_coords);
326  mesh_parameters.add_region_coordinates(2,right_region_coords);
327  mesh_parameters.element_area() = element_area;
328 
329  // Build the bloody thing
330  My_mesh_pt = new TriangleMesh<ELEMENT>(mesh_parameters);
331 
332 
336 
337  char filename[100];
338  sprintf(filename, "RESLT/trace.dat");
339  Trace_file.open(filename);
340 
341  oomph_info << "Number of equations: "
342  << this->assign_eqn_numbers() << '\n';
343 }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
ofstream Trace_file
Trace file to document norm of solution.
Definition: all_foeppl_von_karman/circular_disk.cc:123
TriangleMesh< ELEMENT > * My_mesh_pt
Pointers to specific mesh.
Definition: all_foeppl_von_karman/circular_disk.cc:120
void complete_problem_setup()
Definition: all_foeppl_von_karman/circular_disk.cc:352
double Element_area
Element area.
Definition: all_foeppl_von_karman/circular_disk.cc:140
Definition: geom_objects.h:644
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
Base class defining a closed curve for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1339
void connect_initial_vertex_to_polyline(TriangleMeshPolyLine *polyline_pt, const unsigned &vertex_number, const double &tolerance_for_connection=1.0e-14)
Definition: unstructured_two_d_mesh_geometry_base.cc:1120
void connect_initial_vertex_to_curviline(TriangleMeshCurviLine *curviline_pt, const double &s_value, const double &tolerance_for_connection=1.0e-14)
Definition: unstructured_two_d_mesh_geometry_base.cc:1259
void connect_final_vertex_to_curviline(TriangleMeshCurviLine *curviline_pt, const double &s_value, const double &tolerance_for_connection=1.0e-14)
Definition: unstructured_two_d_mesh_geometry_base.cc:1338
Definition: unstructured_two_d_mesh_geometry_base.h:662
Definition: unstructured_two_d_mesh_geometry_base.h:1642
Definition: triangle_mesh.template.h:94
Class defining a polyline for use in Triangle Mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:868
Definition: triangle_mesh.template.h:424
return int(ret)+1
double Pi
Definition: two_d_biharmonic.cc:235
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
double R_b
The "bubble" radius.
Definition: all_foeppl_von_karman/circular_disk.cc:54
string filename
Definition: MergeRestartFiles.py:39
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::TriangleMeshParameters::add_region_coordinates(), oomph::Problem::add_sub_mesh(), oomph::Problem::assign_eqn_numbers(), oomph::Problem::build_global_mesh(), UnstructuredFvKProblem< ELEMENT >::complete_problem_setup(), oomph::TriangleMeshCurveSection::connect_final_vertex_to_curviline(), oomph::TriangleMeshCurveSection::connect_initial_vertex_to_curviline(), oomph::TriangleMeshCurveSection::connect_initial_vertex_to_polyline(), oomph::TriangleMeshParameters::element_area(), MergeRestartFiles::filename, UnstructuredFvKProblem< ELEMENT >::Inner_boundary0, UnstructuredFvKProblem< ELEMENT >::Inner_boundary1, UnstructuredFvKProblem< ELEMENT >::Inner_boundary2, UnstructuredFvKProblem< ELEMENT >::Inner_boundary3, UnstructuredFvKProblem< ELEMENT >::Inner_boundary4, UnstructuredFvKProblem< ELEMENT >::Inner_boundary5, int(), oomph::TriangleMeshParameters::internal_closed_curve_pt(), oomph::TriangleMeshParameters::internal_open_curves_pt(), UnstructuredFvKProblem< ELEMENT >::My_mesh_pt, oomph::oomph_info, UnstructuredFvKProblem< ELEMENT >::Outer_boundary0, UnstructuredFvKProblem< ELEMENT >::Outer_boundary1, BiharmonicTestFunctions2::Pi, GlobalParameters::R_b, sqrt(), UnstructuredFvKProblem< ELEMENT >::Trace_file, and Eigen::zeta().

◆ ~UnstructuredFvKProblem() [1/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::~UnstructuredFvKProblem ( )
inline

Destructor.

95  {
96  delete My_mesh_pt;
97  }

◆ UnstructuredFvKProblem() [2/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem ( double  element_area = 0.2)

Constructor.

◆ ~UnstructuredFvKProblem() [2/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::~UnstructuredFvKProblem ( )
inline

Destructor.

95  {
96  delete My_mesh_pt;
97  }

◆ UnstructuredFvKProblem() [3/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem ( double  element_area = 0.2)

Constructor.

◆ ~UnstructuredFvKProblem() [3/3]

template<class ELEMENT >
UnstructuredFvKProblem< ELEMENT >::~UnstructuredFvKProblem ( )
inline

Destructor.

102  {
103  delete My_mesh_pt;
104  };

Member Function Documentation

◆ actions_after_adapt()

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_after_adapt ( )
inlinevirtual

Actions after adapt: Setup the problem again – remember that the mesh has been completely rebuilt and its element's don't have any pointers to source fcts etc. yet Also create new volume constraint element with previous pressure value

Reimplemented from oomph::Problem.

123  {
135  }
FoepplvonKarmanVolumeConstraintElement< ELEMENT, RefineableTriangleMesh > * Volume_constraint_element_pt
Single volume constraint element instance.
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:173
Vector< unsigned > Bubble_regions
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:179
Mesh * Volume_constraint_mesh_pt
Mesh for the volume constraint element.
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:169
double * Temp_pressure_pt
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:177
Definition: foeppl_von_karman_volume_constraint_element.h:53
Data * pressure_data_pt() const
Definition: foeppl_von_karman_volume_constraint_element.h:240
void set_prescribed_volume(double *volume_pt)
Set pointer to target volume.
Definition: foeppl_von_karman_volume_constraint_element.h:233
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
void rebuild_global_mesh()
Definition: problem.cc:1533
Unstructured refineable Triangle Mesh.
Definition: triangle_mesh.template.h:2249
double prescribed_volume
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:54
Data * p_b_pt
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:50

References oomph::Mesh::add_element_pt(), TestSoln::p_b_pt, and TestSoln::prescribed_volume.

◆ actions_after_newton_solve() [1/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update after solve (empty)

Reimplemented from oomph::Problem.

101 {}

◆ actions_after_newton_solve() [2/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update after solve (empty)

Reimplemented from oomph::Problem.

101 {}

◆ actions_after_newton_solve() [3/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update after solve (empty)

Reimplemented from oomph::Problem.

139  {
140  }

◆ actions_before_adapt()

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_before_adapt ( )
inlinevirtual

Actions before adapt. Delete old volume constraint element and update global mesh

Reimplemented from oomph::Problem.

109  {
115  }
double value(const unsigned &i) const
Definition: nodes.h:293
void flush_element_and_node_storage()
Definition: mesh.h:407

◆ actions_before_newton_solve() [1/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

104 {}

◆ actions_before_newton_solve() [2/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

104 {}

◆ actions_before_newton_solve() [3/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve: Re-apply boundary conditons.

Reimplemented from oomph::Problem.

144  {
146  }
void apply_boundary_conditions()
Helper function to apply boundary conditions.
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:510

◆ apply_boundary_conditions()

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::apply_boundary_conditions
private

Helper function to apply boundary conditions.

511 {
512 
513  // Loop over all boundary nodes
514  //unsigned nbound=this->My_mesh_pt->nboundary();
515  //Just loop over outer boundary since inner boundary doesn't have boundary
516  //conditions
517  unsigned nbound = Outer_boundary1 + 1;
518 
519  for(unsigned ibound=0;ibound<nbound;ibound++)
520  {
521  unsigned num_nod=this->My_mesh_pt->nboundary_node(ibound);
522  for (unsigned inod=0;inod<num_nod;inod++)
523  {
524  // Get node
525  Node* nod_pt=this->My_mesh_pt->boundary_node_pt(ibound,inod);
526 
527  // Extract nodal coordinates from node:
528  Vector<double> x(2);
529  x[0]=nod_pt->x(0);
530  x[1]=nod_pt->x(1);
531  }
532  }
533 
534 } // end set bc
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
list x
Definition: plotDoE.py:28

References plotDoE::x, and oomph::Node::x().

◆ complete_problem_setup() [1/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::complete_problem_setup
private

Helper function to (re-)set boundary condition and complete the build of all elements

Set boundary condition exactly, and complete the build of all elements

Number of sample points

353 {
354 
355  // Set the boundary conditions for problem: All nodes are
356  // free by default -- just pin the ones that have Dirichlet
357  // conditions here.
358  unsigned nbound = Outer_boundary1 + 1;
359  for(unsigned ibound=0;ibound<nbound;ibound++)
360  {
361  unsigned num_nod=My_mesh_pt->nboundary_node(ibound);
362  for (unsigned inod=0;inod<num_nod;inod++)
363  {
364  // Get node
365  Node* nod_pt=My_mesh_pt->boundary_node_pt(ibound,inod);
366 
367  // Pin unknown values
368  nod_pt->pin(0); // pin vertical displacement; Laplacian of w unpinned
369  // imposes clamping
370  nod_pt->pin(2); // pin Airy stress fct
371  }
372  } // end loop over boundaries
373 
374 
375 
376  // Complete the build of all elements so they are fully functional
377  unsigned n_element = My_mesh_pt->nelement();
378  for(unsigned e=0;e<n_element;e++)
379  {
380  // Upcast from GeneralisedElement to the present element
381  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(My_mesh_pt->element_pt(e));
382 
383  //Set the pressure function pointers and the physical constants
384  el_pt->eta_pt() = &GlobalParameters::Eta;
385  el_pt->pressure_fct_pt() = &GlobalParameters::get_pressure;
386 
387  }
388 
389  // Create the mesh as Geom Object
391 
392 
393 
395  unsigned n_sample=50;
396  const double dr = 1.0/double(n_sample);
397 
398  // Get 'em
399  Vector<double> x(2,0.0);
400  Radial_sample_point_pt.resize(n_sample);
401  for (unsigned j=0;j<n_sample;j++)
402  {
403  Radial_sample_point_pt[j].second.resize(2);
404  x[0]=double(j)*dr;
405 
406  // Get the element and its local coordinates
408  Radial_sample_point_pt[j].first,
409  Radial_sample_point_pt[j].second);
410  }
411 }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Vector< std::pair< GeomObject *, Vector< double > > > Radial_sample_point_pt
Plot points along radial line are located in this element/local coord.
Definition: all_foeppl_von_karman/circular_disk.cc:145
MeshAsGeomObject * Mesh_as_geom_object_pt
Geom object made of mesh.
Definition: all_foeppl_von_karman/circular_disk.cc:142
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
Definition: mesh_as_geometric_object.h:93
void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: mesh_as_geometric_object.h:373
double Eta
FvK parameter.
Definition: all_foeppl_von_karman/circular_disk.cc:51
void get_pressure(const Vector< double > &x, double &pressure)
Pressure depending on the position (x,y)
Definition: all_foeppl_von_karman/circular_disk.cc:60
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References e(), GlobalParameters::Eta, GlobalParameters::get_pressure(), j, oomph::Data::pin(), and plotDoE::x.

Referenced by UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem().

◆ complete_problem_setup() [2/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::complete_problem_setup ( )
private

Helper function to (re-)set boundary condition and complete the build of all elements

◆ complete_problem_setup() [3/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::complete_problem_setup ( )
private

Helper function to (re-)set boundary condition and complete the build of all elements

◆ doc_solution() [1/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::doc_solution

Doc the solution.

419 {
420  ofstream some_file;
421  char filename[100];
422 
423  // Number of plot points
424  unsigned npts = 5;
425 
426  sprintf(filename,"RESLT/soln%i.dat",Doc_info.number());
427  some_file.open(filename);
428  this->My_mesh_pt->output(some_file,npts);
429  some_file.close();
430 
431  // Output boundaries
432  //------------------
433  sprintf(filename,"RESLT/boundaries%i.dat",Doc_info.number());
434  some_file.open(filename);
435  My_mesh_pt->output_boundaries(some_file);
436  some_file.close();
437 
438  // Find solution at r = 0
439  // ----------------------
440  unsigned n_boundary3_node = My_mesh_pt->nboundary_node(Inner_boundary3);
441  unsigned n_boundary4_node = My_mesh_pt->nboundary_node(Inner_boundary4);
442  double w_0 = 0;
443  for(unsigned inode=0; inode < n_boundary3_node; inode++)
444  {
445  for(unsigned jnode=0; jnode < n_boundary4_node; jnode++)
446  {
447  if(My_mesh_pt->boundary_node_pt(Inner_boundary3,inode)
448  == My_mesh_pt->boundary_node_pt(Inner_boundary4,jnode))
449  {
450  w_0 = My_mesh_pt->boundary_node_pt(Inner_boundary3,inode)->value(0);
451  }
452  }
453  }
454 
455  // Plot solution along radial line
456  sprintf(filename,"RESLT/soln_along_radial_line%i.dat",Doc_info.number());
457  some_file.open(filename);
458  Vector<double> x(2);
459  unsigned nplot=Radial_sample_point_pt.size();
460  for (unsigned j=0;j<nplot;j++)
461  {
462  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(Radial_sample_point_pt[j].first);
464  el_pt->interpolated_x(s,x);
465  some_file << x[0] << " "
466  << x[1] << " "
467  << el_pt->interpolated_w_fvk(s) << " " // w
468  << std::endl;
469  }
470  some_file.close();
471 
472 
473  // Doc error and return of the square of the L2 error
474  //---------------------------------------------------
475  double dummy_error,zero_norm;
476  sprintf(filename,"RESLT/norm%i.dat",Doc_info.number());
477  some_file.open(filename);
478  My_mesh_pt->compute_error(some_file,GlobalParameters::zero,
479  dummy_error,zero_norm);
480  some_file.close();
481 
482  // Doc L2 error and norm of solution
483  oomph_info << "Norm of computed solution: "
484  << sqrt(dummy_error) << std::endl;
485 
487  << w_0 << '\n';
488 
489  // Increment the doc_info number
490  Doc_info.number()++;
491 
492 } // end of doc
DocInfo Doc_info
Doc info object for labeling output.
Definition: all_foeppl_von_karman/circular_disk.cc:113
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
RealScalar s
Definition: level1_cplx_impl.h:130
double Pressure
The pressure.
Definition: all_foeppl_von_karman/circular_disk.cc:57
void zero(const Vector< double > &x, Vector< double > &u)
Definition: all_foeppl_von_karman/circular_disk.cc:67

References GlobalParameters::Doc_info, MergeRestartFiles::filename, j, oomph::DocInfo::number(), oomph::oomph_info, GlobalParameters::Pressure, s, sqrt(), oomph::Problem_Parameter::Trace_file, plotDoE::x, and GlobalParameters::zero().

◆ doc_solution() [2/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::doc_solution ( )

Doc the solution.

◆ doc_solution() [3/3]

template<class ELEMENT >
void UnstructuredFvKProblem< ELEMENT >::doc_solution ( const std::string &  comment = "")

Doc the solution.

543 {
544  ofstream some_file;
545  char filename[100];
546 
547  // Number of plot points
548  unsigned npts = 5;
549 
550  sprintf(filename,"RESLT/soln%i-%f.dat",Doc_info.number(),Element_area);
551  some_file.open(filename);
552  this->My_mesh_pt->output(some_file,npts);
553  some_file << "TEXT X = 22, Y = 92, CS=FRAME T = \""
554  << comment << "\"\n";
555  some_file.close();
556 
557  // Output boundaries
558  //------------------
559  sprintf(filename,"RESLT/boundaries%i-%f.dat",Doc_info.number(),Element_area);
560  some_file.open(filename);
561  My_mesh_pt->output_boundaries(some_file);
562  some_file.close();
563 
564 
565  // Output along radial boundary
566  // ----------------------------
567  sprintf(filename,"RESLT/radial_soln%i-%f.dat",Doc_info.number(),Element_area);
568  some_file.open(filename);
569  for(unsigned b = Inner_boundary2; b <= Inner_boundary5; b++)
570  {
571  unsigned n_boundary_el = My_mesh_pt->nboundary_element(b);
572  for(unsigned e = 0; e < n_boundary_el; e++)
573  {
574  My_mesh_pt->boundary_element_pt(b,e)->output(some_file,npts);
575  }
576  }
577  some_file.close();
578 
579  // Find solution at r = 0
580  // ----------------------
581  unsigned n_boundary3_node = My_mesh_pt->nboundary_node(Inner_boundary3);
582  unsigned n_boundary4_node = My_mesh_pt->nboundary_node(Inner_boundary4);
583  double w_0 = 0;
584  for(unsigned inode=0; inode < n_boundary3_node; inode++)
585  {
586  for(unsigned jnode=0; jnode < n_boundary4_node; jnode++)
587  {
588  if(My_mesh_pt->boundary_node_pt(Inner_boundary3,inode)
589  == My_mesh_pt->boundary_node_pt(Inner_boundary4,jnode))
590  {
591  w_0 = My_mesh_pt->boundary_node_pt(Inner_boundary3,inode)->value(0);
592  }
593  }
594  }
595 
596 
597  // Output regions
598  unsigned n_region = My_mesh_pt->nregion();
599  if (n_region > 1)
600  {
601  for (unsigned r = 0; r < n_region; r++)
602  {
603  //Attempt to output elements in different regions
604  sprintf(filename,"RESLT/region%i%i-%f.dat",r,Doc_info.number(),
605  Element_area);
606  some_file.open(filename);
607  unsigned nel = My_mesh_pt->nregion_element(r);
608  for (unsigned e = 0; e < nel; e++)
609  {
610  My_mesh_pt->region_element_pt(r,e)->output(some_file,npts);
611  }
612  some_file.close();
613  }
614  }
615 
616  // Calculate the integral (i.e. volume) over the inner region
617  double bubble_volume = 0;
618  if (n_region > 1)
619  {
620  sprintf(filename,"RESLT/bubble_volume%i.dat",Doc_info.number());
621  some_file.open(filename);
622  unsigned n_inner_el;
623  unsigned n_bubble_regions = Bubble_regions.size();
624 
625  for(unsigned r = 0;r < n_bubble_regions; r++)
626  {
627  n_inner_el = My_mesh_pt->nregion_element(Bubble_regions[r]);
628  for(unsigned e = 0; e < n_inner_el; e++)
629  {
630  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(
631  My_mesh_pt->region_element_pt(Bubble_regions[r],e));
632  if(el_pt != 0)
633  {
634  bubble_volume += el_pt->get_bounded_volume();
635  }
636  }
637  }
638  some_file << bubble_volume << '\n';
639  some_file.close();
640  }
641 
642 
643  // Doc error and return of the square of the L2 error
644  //---------------------------------------------------
645  //double error,norm,dummy_error,zero_norm;
646  double dummy_error,zero_norm;
647  sprintf(filename,"RESLT/error%i-%f.dat",Doc_info.number(),Element_area);
648  some_file.open(filename);
649 
650  My_mesh_pt->compute_error(some_file,TestSoln::zero,
651  dummy_error,zero_norm);
652  some_file.close();
653 
654  // Doc L2 error and norm of solution
655  oomph_info << "Norm of computed solution: " << sqrt(dummy_error) << std::endl;
656 
657  Trace_file << TestSoln::p_b_pt->value(0) << " "
658  << TestSoln::p_0_pt->value(0) << " "
659  << w_0 << " " << bubble_volume << '\n';
660 
661  // Increment the doc_info number
662  Doc_info.number()++;
663 
664 } // end of doc
Scalar * b
Definition: benchVecAdd.cpp:17
void zero(const Vector< double > &x, Vector< double > &u)
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:76
Data * p_0_pt
Definition: foeppl_von_karman/circular_disk/circular_disk.cc:51
r
Definition: UniformPSDSelfTest.py:20

References b, GlobalParameters::Doc_info, e(), ProblemParameters::Element_area, MergeRestartFiles::filename, oomph::DocInfo::number(), oomph::oomph_info, TestSoln::p_0_pt, TestSoln::p_b_pt, UniformPSDSelfTest::r, sqrt(), oomph::Problem_Parameter::Trace_file, oomph::Data::value(), and TestSoln::zero().

Member Data Documentation

◆ Bubble_regions

template<class ELEMENT >
Vector<unsigned> UnstructuredFvKProblem< ELEMENT >::Bubble_regions
private

◆ Doc_info

template<class ELEMENT >
DocInfo UnstructuredFvKProblem< ELEMENT >::Doc_info
private

Doc info object for labeling output.

◆ Element_area

template<class ELEMENT >
double UnstructuredFvKProblem< ELEMENT >::Element_area
private

Element area.

◆ Mesh_as_geom_object_pt

template<class ELEMENT >
MeshAsGeomObject * UnstructuredFvKProblem< ELEMENT >::Mesh_as_geom_object_pt
private

Geom object made of mesh.

◆ My_mesh_pt [1/2]

template<class ELEMENT >
TriangleMesh< ELEMENT > * UnstructuredFvKProblem< ELEMENT >::My_mesh_pt
private

Pointers to specific mesh.

Referenced by UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem().

◆ My_mesh_pt [2/2]

template<class ELEMENT >
RefineableTriangleMesh<ELEMENT>* UnstructuredFvKProblem< ELEMENT >::My_mesh_pt
private

Pointers to specific mesh.

◆ Radial_sample_point_pt

template<class ELEMENT >
Vector< std::pair< GeomObject *, Vector< double > > > UnstructuredFvKProblem< ELEMENT >::Radial_sample_point_pt
private

Plot points along radial line are located in this element/local coord.

◆ Temp_pressure_pt

template<class ELEMENT >
double* UnstructuredFvKProblem< ELEMENT >::Temp_pressure_pt
private

Temporary storage for the pressure used when switching between instances of volume constraint element

◆ Trace_file

template<class ELEMENT >
ofstream UnstructuredFvKProblem< ELEMENT >::Trace_file
private

Trace file to document norm of solution.

Referenced by UnstructuredFvKProblem< ELEMENT >::UnstructuredFvKProblem().

◆ Volume_constraint_element_pt

template<class ELEMENT >
FoepplvonKarmanVolumeConstraintElement<ELEMENT, RefineableTriangleMesh>* UnstructuredFvKProblem< ELEMENT >::Volume_constraint_element_pt
private

Single volume constraint element instance.

◆ Volume_constraint_mesh_pt

template<class ELEMENT >
Mesh* UnstructuredFvKProblem< ELEMENT >::Volume_constraint_mesh_pt
private

Mesh for the volume constraint element.


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