FSIDrivenCavityProblem< ELEMENT > Class Template Reference

Problem class. More...

#include <fsi_driven_cavity_problem.h>

+ Inheritance diagram for FSIDrivenCavityProblem< ELEMENT >:

Public Member Functions

 FSIDrivenCavityProblem (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, const double &gap_fraction, const double &period)
 Constructor for the collapsible channel problem. More...
 
 ~FSIDrivenCavityProblem ()
 Destructor. More...
 
AlgebraicFSIDrivenCavityMesh< ELEMENT > * bulk_mesh_pt ()
 Access function for the specific bulk (fluid) mesh. More...
 
OneDLagrangianMesh< FSIHermiteBeamElement > * wall_mesh_pt ()
 Access function for the wall mesh. More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem after solve (empty) More...
 
void actions_before_implicit_timestep ()
 Update the velocity boundary condition on the moving lid. More...
 
void actions_before_newton_convergence_check ()
 
void doc_solution (DocInfo &doc_info, ofstream &trace_file)
 Doc the solution. More...
 
void set_initial_condition ()
 Apply initial conditions. 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 ()
 

Protected Attributes

unsigned Nx
 Number of elements in the x direction. More...
 
unsigned Ny
 Number of elements in the y direction. More...
 
double Lx
 Width of domain. More...
 
double Ly
 Height of domain. More...
 
double T
 Period of oscillation. More...
 
AlgebraicFSIDrivenCavityMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
OneDLagrangianMesh< FSIHermiteBeamElement > * Wall_mesh_pt
 Pointer to the "wall" mesh. More...
 
NodeLeft_node_pt
 Pointer to the left control node. More...
 
NodeRight_node_pt
 Pointer to right control node. More...
 
NodeWall_node_pt
 Pointer to control node on the wall. More...
 
MeshAsGeomObjectWall_geom_object_pt
 
- 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
 

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_step ()
 
virtual void actions_after_newton_step ()
 
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 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 ()
 
- 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 FSIDrivenCavityProblem< ELEMENT >

Problem class.

Constructor & Destructor Documentation

◆ FSIDrivenCavityProblem()

template<class ELEMENT >
FSIDrivenCavityProblem< ELEMENT >::FSIDrivenCavityProblem ( const unsigned nx,
const unsigned ny,
const double lx,
const double ly,
const double gap_fraction,
const double period 
)

Constructor for the collapsible channel problem.

Constructor: The arguments are the number of elements, the lengths of the domain, the fractional height of the gap next to the moving lid and the period of the lid's oscillation

292 {
293 
294  // Store problem parameters
295  Nx=nx;
296  Ny=ny;
297  Lx=lx;
298  Ly=ly;
299 
300  // Period of lid oscillation
301  T=period;
302 
303 
304  // Allow for crappy initial guess
306  Max_residuals = 1.0e8;
307 
308  // Allocate the timestepper -- this constructs the Problem's
309  // time object with a sufficient amount of storage to store the
310  // previous timsteps.
311  BDF<2>* fluid_time_stepper_pt=new BDF<2>;
312  add_time_stepper_pt(fluid_time_stepper_pt);
313 
314  // Create the timestepper for the solid
315  Newmark<2>* solid_time_stepper_pt=new Newmark<2>;
316  //Steady<2>* solid_time_stepper_pt=new Steady<2>;
317  add_time_stepper_pt(solid_time_stepper_pt);
318 
319  // Geometric object that represents the undeformed wall:
320  // A straight line at height y=ly; starting at x=0.
321  UndeformedWall* undeformed_wall_pt=new UndeformedWall(0.0,ly);
322 
323  //Create the "wall" mesh with FSI Hermite elements
325  (nx,lx,undeformed_wall_pt,solid_time_stepper_pt);
326 
327 
328  // Build a geometric object (one Lagrangian, two Eulerian coordinates)
329  // from the wall mesh
332 
333  //Build bulk (fluid) mesh
335  (nx, ny, lx, ly, gap_fraction, Wall_geom_object_pt,fluid_time_stepper_pt);
336 
337 
338  // Add the sub meshes to the problem
341 
342  // Combine all submeshes into a single Mesh
344 
345 
346  // Complete build of fluid mesh
347  //-----------------------------
348 
349  // Loop over the elements to set up element-specific
350  // things that cannot be handled by constructor
351  unsigned n_element=Bulk_mesh_pt->nelement();
352  for(unsigned e=0;e<n_element;e++)
353  {
354  // Upcast from GeneralisedElement to the present element
355  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
356 
357  //Set the Reynolds number
358  el_pt->re_pt() = &Global_Physical_Variables::Re;
359 
360  // Set the Womersley number
361  el_pt->re_st_pt() = &Global_Physical_Variables::ReSt;
362 
363  } // end loop over elements
364 
365 
366 
367  // Apply boundary conditions for fluid
368  //------------------------------------
369 
370  // Pin the velocity on all boundaries apart from 1 and 5
371  // (the gaps above the driven lid)
372  for (unsigned ibound=0;ibound<6;ibound++)
373  {
374  if ((ibound!=1)&&(ibound!=5))
375  {
376  unsigned num_nod= bulk_mesh_pt()->nboundary_node(ibound);
377  for (unsigned inod=0;inod<num_nod;inod++)
378  {
379  for(unsigned i=0;i<2;i++)
380  {
381  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->pin(i);
382  }
383  }
384  }
385  }//end of pin_velocity
386 
387 
388  // Complete build of wall elements
389  //--------------------------------
390 
391  //Loop over the elements to set physical parameters etc.
392  n_element = wall_mesh_pt()->nelement();
393  for(unsigned e=0;e<n_element;e++)
394  {
395  // Upcast to the specific element type
396  FSIHermiteBeamElement *elem_pt =
397  dynamic_cast<FSIHermiteBeamElement*>(wall_mesh_pt()->element_pt(e));
398 
399  // Set physical parameters for each element:
400  elem_pt->h_pt() = &Global_Physical_Variables::H;
401 
402  // Function that specifies the load ratios
403  elem_pt->q_pt() = &Global_Physical_Variables::Q;
404 
405  // Set the undeformed shape for each element
406  elem_pt->undeformed_beam_pt() = undeformed_wall_pt;
407 
408  // The normal on the wall elements as computed by the FSIHermiteElements
409  // points away from the fluid rather than into the fluid (as assumed
410  // by default)
412 
413  // Timescale ratio for solid
415 
416  } // end of loop over elements
417 
418 
419  // Boundary conditions for wall mesh
420  //----------------------------------
421 
422  // Set the boundary conditions: Each end of the beam is fixed in space
423  // Loop over the boundaries (ends of the beam)
424  for(unsigned b=0;b<2;b++)
425  {
426  // Pin displacements in both x and y directions
429  }
430 
431 
432 
433 
434  //Choose control nodes
435  //---------------------
436 
437  // Left boundary
438  unsigned ibound=5;
439  unsigned num_nod= bulk_mesh_pt()->nboundary_node(ibound);
440  unsigned control_nod=num_nod/2;
441  Left_node_pt= bulk_mesh_pt()->boundary_node_pt(ibound, control_nod);
442 
443  // Right boundary
444  ibound=1;
445  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
446  control_nod=num_nod/2;
447  Right_node_pt= bulk_mesh_pt()->boundary_node_pt(ibound, control_nod);
448 
449 
450  // Set the pointer to the control node on the wall
451  num_nod= wall_mesh_pt()->nnode();
453 
454 
455 
456 
457  // Setup FSI
458  //----------
459 
460  // The velocity of the fluid nodes on the wall (fluid mesh boundary 3)
461  // is set by the wall motion -- hence the no-slip condition must be
462  // re-applied whenever a node update is performed for these nodes.
463  // Such tasks may be performed automatically by the auxiliary node update
464  // function specified by a function pointer:
465  ibound=3;
466  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
467  for (unsigned inod=0;inod<num_nod;inod++)
468  {
469  static_cast<AlgebraicNode*>(
470  bulk_mesh_pt()->boundary_node_pt(ibound, inod))->
471  set_auxiliary_node_update_fct_pt(
473  }
474 
475 
476  // Work out which fluid dofs affect the residuals of the wall elements:
477  // We pass the boundary between the fluid and solid meshes and
478  // pointers to the meshes. The interaction boundary is boundary 3 of the
479  // 2D fluid mesh.
480  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
481  (this,3,Bulk_mesh_pt,Wall_mesh_pt);
482 
483  // Setup equation numbering scheme
484  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
485 
486 
487 }//end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Scalar * b
Definition: benchVecAdd.cpp:17
MeshAsGeomObject * Wall_geom_object_pt
Definition: fsi_driven_cavity_problem.h:274
OneDLagrangianMesh< FSIHermiteBeamElement > * Wall_mesh_pt
Pointer to the "wall" mesh.
Definition: fsi_driven_cavity_problem.h:261
unsigned Ny
Number of elements in the y direction.
Definition: fsi_driven_cavity_problem.h:246
Node * Wall_node_pt
Pointer to control node on the wall.
Definition: fsi_driven_cavity_problem.h:270
OneDLagrangianMesh< FSIHermiteBeamElement > * wall_mesh_pt()
Access function for the wall mesh.
Definition: fsi_driven_cavity_problem.h:193
double Lx
Width of domain.
Definition: fsi_driven_cavity_problem.h:249
double Ly
Height of domain.
Definition: fsi_driven_cavity_problem.h:252
Node * Left_node_pt
Pointer to the left control node.
Definition: fsi_driven_cavity_problem.h:264
AlgebraicFSIDrivenCavityMesh< ELEMENT > * bulk_mesh_pt()
Access function for the specific bulk (fluid) mesh.
Definition: fsi_driven_cavity_problem.h:182
Node * Right_node_pt
Pointer to right control node.
Definition: fsi_driven_cavity_problem.h:267
AlgebraicFSIDrivenCavityMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: fsi_driven_cavity_problem.h:258
unsigned Nx
Number of elements in the x direction.
Definition: fsi_driven_cavity_problem.h:243
Definition: fsi_chan_problem.h:156
Definition: fsi_driven_cavity_mesh.template.h:161
Definition: algebraic_elements.h:55
Definition: beam_elements.h:352
void set_normal_pointing_out_of_fluid()
Definition: beam_elements.h:384
double *& q_pt()
Definition: fsi.h:248
double *& h_pt()
Return a pointer to non-dim. wall thickness.
Definition: beam_elements.h:240
double *& lambda_sq_pt()
Return a pointer to timescale ratio (nondim density)
Definition: beam_elements.h:246
GeomObject *& undeformed_beam_pt()
Definition: beam_elements.h:253
Definition: mesh_as_geometric_object.h:93
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
unsigned Max_newton_iterations
Maximum number of Newton iterations.
Definition: problem.h:599
void build_global_mesh()
Definition: problem.cc:1493
double Max_residuals
Definition: problem.h:610
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
SolidNode * node_pt(const unsigned long &n)
Return a pointer to the n-th global SolidNode.
Definition: mesh.h:2594
SolidNode * boundary_node_pt(const unsigned &b, const unsigned &n)
Return n-th SolidNodes on b-th boundary.
Definition: mesh.h:2612
void pin_position(const unsigned &i)
Pin the nodal position.
Definition: nodes.h:1816
double ReSt
Womersley number.
Definition: rayleigh_instability.cc:56
double Q
Ratio of scales.
Definition: elastic_bretherton.cc:131
double Lambda_sq
Pseudo-solid mass density.
Definition: pressure_driven_torus.cc:62
double Re
Reynolds number.
Definition: fibre.cc:55
double H
Nondim thickness.
Definition: steady_ring.cc:50
const double ly
Definition: ConstraintElementsUnitTest.cpp:34
const double lx
Definition: ConstraintElementsUnitTest.cpp:33
const unsigned nx
Definition: ConstraintElementsUnitTest.cpp:30
const unsigned ny
Definition: ConstraintElementsUnitTest.cpp:31
void apply_no_slip_on_moving_wall(Node *node_pt)
Definition: fsi.cc:48

References oomph::FSI_functions::apply_no_slip_on_moving_wall(), b, e(), Global_Physical_Variables::H, oomph::KirchhoffLoveBeamEquations::h_pt(), i, Global_Physical_Variables::Lambda_sq, oomph::KirchhoffLoveBeamEquations::lambda_sq_pt(), Mesh_Parameters::lx, Global_Parameters::Lx, Mesh_Parameters::ly, Global_Parameters::Ly, oomph::Locate_zeta_helpers::Max_newton_iterations, Mesh_Parameters::nx, GlobalParameters::Nx, Mesh_Parameters::ny, GlobalParameters::Ny, Global_Physical_Variables::Q, oomph::FSIWallElement::q_pt(), Global_Physical_Variables::Re, Global_Physical_Variables::ReSt, oomph::FSIHermiteBeamElement::set_normal_pointing_out_of_fluid(), and oomph::KirchhoffLoveBeamEquations::undeformed_beam_pt().

◆ ~FSIDrivenCavityProblem()

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

Destructor.

176  {
177  // Mesh gets killed in general problem destructor
178  }

Member Function Documentation

◆ actions_after_newton_solve()

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

Update the problem after solve (empty)

Reimplemented from oomph::Problem.

Reimplemented in SegregatedFSIDrivenCavityProblem< ELEMENT >.

203 {}

◆ actions_before_implicit_timestep()

template<class ELEMENT >
void FSIDrivenCavityProblem< ELEMENT >::actions_before_implicit_timestep ( )
inlinevirtual

Update the velocity boundary condition on the moving lid.

Reimplemented from oomph::Problem.

207  {
208  // Oscillating lid
209  unsigned ibound=0;
210  unsigned num_nod=bulk_mesh_pt()->nboundary_node(ibound);
211  for (unsigned inod=0;inod<num_nod;inod++)
212  {
213  // Which node are we dealing with?
214  Node* node_pt=bulk_mesh_pt()->boundary_node_pt(ibound,inod);
215 
216  // Set velocity
217  double veloc=1.0-cos(2.0*MathematicalConstants::Pi*time_pt()->time()/T);
218 
219  // Apply no slip
220  node_pt->set_value(0,veloc);
221  }
222  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Definition: nodes.h:906
Time *& time_pt()
Return a pointer to the global time object.
Definition: problem.h:1504
double & time()
Return the current value of continuous time.
Definition: problem.cc:11531
double Pi
Definition: two_d_biharmonic.cc:235

References cos(), BiharmonicTestFunctions2::Pi, and oomph::Data::set_value().

◆ actions_before_newton_convergence_check()

template<class ELEMENT >
void FSIDrivenCavityProblem< ELEMENT >::actions_before_newton_convergence_check ( )
inlinevirtual

Update before checking Newton convergence: Update the nodal positions in the fluid mesh in response to possible changes in the wall shape

Reimplemented from oomph::Problem.

Reimplemented in SegregatedFSIDrivenCavityProblem< ELEMENT >.

229  {
230  Bulk_mesh_pt->node_update();
231  }

◆ actions_before_newton_solve()

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

Reimplemented in SegregatedFSIDrivenCavityProblem< ELEMENT >.

200 {}

◆ bulk_mesh_pt()

template<class ELEMENT >
AlgebraicFSIDrivenCavityMesh<ELEMENT>* FSIDrivenCavityProblem< ELEMENT >::bulk_mesh_pt ( )
inline

Access function for the specific bulk (fluid) mesh.

183  {
184  // Upcast from pointer to the Mesh base class to the specific
185  // element type that we're using here.
186  return dynamic_cast<
188  (Bulk_mesh_pt);
189  }

◆ doc_solution()

template<class ELEMENT >
void FSIDrivenCavityProblem< ELEMENT >::doc_solution ( DocInfo doc_info,
ofstream &  trace_file 
)

Doc the solution.

498 {
499 
500  // Doc fsi
501  if (CommandLineArgs::Argc>1)
502  {
503  FSI_functions::doc_fsi<AlgebraicNode>(Bulk_mesh_pt,Wall_mesh_pt,doc_info);
504  }
505 
506  ofstream some_file;
507  char filename[100];
508 
509  // Number of plot points
510  unsigned npts;
511  npts=5;
512 
513  // Output fluid solution
514  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
515  doc_info.number());
516  some_file.open(filename);
517  bulk_mesh_pt()->output(some_file,npts);
518  some_file.close();
519 
520  // Document the wall shape
521  sprintf(filename,"%s/beam%i.dat",doc_info.directory().c_str(),
522  doc_info.number());
523  some_file.open(filename);
524  wall_mesh_pt()->output(some_file,npts);
525  some_file.close();
526 
527 
528  // Write trace file
529  trace_file << time_pt()->time() << " "
530  << Wall_node_pt->x(1) << " "
531  << Left_node_pt->value(0) << " "
532  << Right_node_pt->value(0) << " "
533  << std::endl;
534 
535 } // end_of_doc_solution
std::string directory() const
Output directory.
Definition: oomph_utilities.h:524
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
double value(const unsigned &i) const
Definition: nodes.cc:2408
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
string filename
Definition: MergeRestartFiles.py:39
int Argc
Number of arguments + 1.
Definition: oomph_utilities.cc:407

References oomph::CommandLineArgs::Argc, oomph::DocInfo::directory(), MergeRestartFiles::filename, and oomph::DocInfo::number().

◆ set_initial_condition()

template<class ELEMENT >
void FSIDrivenCavityProblem< ELEMENT >::set_initial_condition
virtual

Apply initial conditions.

Reimplemented from oomph::Problem.

545 {
546  // Impulsive start for wall
548 
549  // Update the mesh
550  bulk_mesh_pt()->node_update();
551 
552  // Loop over the nodes to set initial guess everywhere
553  unsigned num_nod = bulk_mesh_pt()->nnode();
554  for (unsigned n=0;n<num_nod;n++)
555  {
556  // Get nodal coordinates
557  Vector<double> x(2);
558  x[0]=bulk_mesh_pt()->node_pt(n)->x(0);
559  x[1]=bulk_mesh_pt()->node_pt(n)->x(1);
560 
561  // Assign initial condition: Zero flow
562  bulk_mesh_pt()->node_pt(n)->set_value(0,0.0);
563  bulk_mesh_pt()->node_pt(n)->set_value(1,0.0);
564  }
565 
566  // Assign initial values for an impulsive start
567  bulk_mesh_pt()->assign_initial_values_impulsive();
568 
569 
570 } // end of set_initial_condition
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
void assign_initial_values_impulsive()
Assign initial values for an impulsive start.
Definition: mesh.cc:2288
list x
Definition: plotDoE.py:28

References n, and plotDoE::x.

◆ wall_mesh_pt()

template<class ELEMENT >
OneDLagrangianMesh<FSIHermiteBeamElement>* FSIDrivenCavityProblem< ELEMENT >::wall_mesh_pt ( )
inline

Access function for the wall mesh.

194  {
195  return Wall_mesh_pt;
196  } // end of access to wall mesh

Member Data Documentation

◆ Bulk_mesh_pt

template<class ELEMENT >
AlgebraicFSIDrivenCavityMesh<ELEMENT>* FSIDrivenCavityProblem< ELEMENT >::Bulk_mesh_pt
protected

◆ Left_node_pt

template<class ELEMENT >
Node* FSIDrivenCavityProblem< ELEMENT >::Left_node_pt
protected

Pointer to the left control node.

◆ Lx

template<class ELEMENT >
double FSIDrivenCavityProblem< ELEMENT >::Lx
protected

Width of domain.

◆ Ly

template<class ELEMENT >
double FSIDrivenCavityProblem< ELEMENT >::Ly
protected

Height of domain.

◆ Nx

template<class ELEMENT >
unsigned FSIDrivenCavityProblem< ELEMENT >::Nx
protected

Number of elements in the x direction.

◆ Ny

template<class ELEMENT >
unsigned FSIDrivenCavityProblem< ELEMENT >::Ny
protected

Number of elements in the y direction.

◆ Right_node_pt

template<class ELEMENT >
Node* FSIDrivenCavityProblem< ELEMENT >::Right_node_pt
protected

Pointer to right control node.

◆ T

template<class ELEMENT >
double FSIDrivenCavityProblem< ELEMENT >::T
protected

Period of oscillation.

◆ Wall_geom_object_pt

template<class ELEMENT >
MeshAsGeomObject* FSIDrivenCavityProblem< ELEMENT >::Wall_geom_object_pt
protected

Pointer to geometric object (one Lagrangian, two Eulerian coordinates) that will be built from the wall mesh

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

◆ Wall_mesh_pt

template<class ELEMENT >
OneDLagrangianMesh<FSIHermiteBeamElement>* FSIDrivenCavityProblem< ELEMENT >::Wall_mesh_pt
protected

Pointer to the "wall" mesh.

◆ Wall_node_pt

template<class ELEMENT >
Node* FSIDrivenCavityProblem< ELEMENT >::Wall_node_pt
protected

Pointer to control node on the wall.


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