CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT > Class Template Reference

Coated disk FSI. More...

+ Inheritance diagram for CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >:

Public Member Functions

 CoatedDiskProblem ()
 Constructor: More...
 
void actions_before_newton_solve ()
 Update function (empty) More...
 
void actions_after_newton_solve ()
 Update function (empty) More...
 
void actions_before_newton_convergence_check ()
 Recompute gamma integral before checking Newton residuals. More...
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of traction elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of traction elements. More...
 
void doc_solution ()
 Doc the solution. More...
 
 CoatedDiskProblem ()
 Constructor: More...
 
void actions_before_newton_solve ()
 Update function (empty) More...
 
void actions_after_newton_solve ()
 Update function (empty) More...
 
void actions_before_adapt ()
 Actions before adapt: Wipe the face meshes. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the face meshes. More...
 
void doc_solution ()
 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 Member Functions

void create_fsi_traction_elements ()
 Create FSI traction elements. More...
 
void create_helmholtz_fsi_flux_elements ()
 Create Helmholtz FSI flux elements. More...
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete (face) elements in specified mesh. More...
 
void create_helmholtz_DtN_elements ()
 Create DtN face elements. More...
 
void setup_interaction ()
 Setup interaction. More...
 
void complete_problem_setup ()
 
void create_solid_traction_elements ()
 Create solid traction elements. More...
 
void create_fsi_traction_elements ()
 Create FSI traction elements. More...
 
void create_helmholtz_fsi_flux_elements ()
 Create Helmholtz FSI flux elements. More...
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete (face) elements in specified mesh. More...
 
void create_helmholtz_ABC_elements ()
 Create ABC face elements. More...
 
void setup_interaction ()
 Setup interaction. More...
 

Private Attributes

TreeBasedRefineableMeshBaseSolid_mesh_pt
 Pointer to solid mesh. More...
 
MeshFSI_traction_mesh_pt
 Pointer to mesh of FSI traction elements. More...
 
TreeBasedRefineableMeshBaseHelmholtz_mesh_pt
 Pointer to Helmholtz mesh. More...
 
MeshHelmholtz_fsi_flux_mesh_pt
 Pointer to mesh of Helmholtz FSI flux elements. More...
 
HelmholtzDtNMesh< HELMHOLTZ_ELEMENT > * Helmholtz_outer_boundary_mesh_pt
 Pointer to mesh containing the DtN elements. More...
 
DocInfo Doc_info
 DocInfo object for output. More...
 
ofstream Trace_file
 Trace file. More...
 
RefineableTriangleMesh< ELASTICITY_ELEMENT > * Solid_mesh_pt
 Pointer to refineable solid mesh. More...
 
MeshSolid_traction_mesh_pt
 Pointer to mesh of solid traction elements. More...
 
MeshHelmholtz_outer_boundary_mesh_pt
 Pointer to mesh containing the ABC elements. More...
 
unsigned Upper_symmetry_boundary_id
 Boundary ID of upper symmetry boundary. More...
 
unsigned Lower_symmetry_boundary_id
 Boundary ID of lower symmetry boundary. More...
 
unsigned Upper_inner_boundary_id
 Boundary ID of upper inner boundary. More...
 
unsigned Lower_inner_boundary_id
 Boundary ID of lower inner boundary. More...
 
unsigned Outer_boundary_id
 Boundary ID of outer boundary. More...
 
unsigned Rib_divider_boundary_id
 

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_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 ELASTICITY_ELEMENT, class HELMHOLTZ_ELEMENT>
class CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >

Coated disk FSI.

Constructor & Destructor Documentation

◆ CoatedDiskProblem() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedDiskProblem

Constructor:

Constructor.

310 {
311 
312  // The coating mesh is periodic
313  bool periodic=true;
314  double azimuthal_fraction_of_coating=1.0;
315 
316  // Solid mesh
317  //-----------
318  // Number of elements in azimuthal direction
319  unsigned ntheta_solid=10*Global_Parameters::El_multiplier;
320 
321  // Number of elements in radial direction
322  unsigned nr_solid=3*Global_Parameters::El_multiplier;
323 
324  // Innermost radius for solid mesh
325  double a=1.0-Global_Parameters::H_coating;
326 
327  // Build solid mesh
328  Solid_mesh_pt = new
330  (periodic,azimuthal_fraction_of_coating,
331  ntheta_solid,nr_solid,a,Global_Parameters::H_coating);
332 
333 
334  // Helmholtz mesh
335  //---------------
336 
337  // Number of elements in azimuthal direction in Helmholtz mesh
338  unsigned ntheta_helmholtz=11*Global_Parameters::El_multiplier;
339 
340  // Number of elements in radial direction in Helmholtz mesh
341  unsigned nr_helmholtz=3*Global_Parameters::El_multiplier;
342 
343  // Innermost radius of Helmholtz mesh
344  a=1.0;
345 
346  // Thickness of Helmholtz mesh
347  double h_thick_helmholtz=Global_Parameters::Outer_radius-a;
348 
349  // Build mesh
350  Helmholtz_mesh_pt = new
352  (periodic,azimuthal_fraction_of_coating,
353  ntheta_helmholtz,nr_helmholtz,a,h_thick_helmholtz);
354 
355 
356  // Set error estimators
359 
360 
361  // Mesh containing the Helmholtz DtN
362  // elements. Specify outer radius and number of Fourier terms to be
363  // used in gamma integral
364  unsigned nfourier=20;
367  nfourier);
368 
369  //Assign the physical properties to the elements before any refinement
370  //Loop over the elements in the solid mesh
371  unsigned n_element=Solid_mesh_pt->nelement();
372  for(unsigned i=0;i<n_element;i++)
373  {
374  //Cast to a solid element
375  ELASTICITY_ELEMENT *el_pt =
376  dynamic_cast<ELASTICITY_ELEMENT*>(Solid_mesh_pt->element_pt(i));
377 
378  // Set the constitutive law
379  el_pt->elasticity_tensor_pt() = &Global_Parameters::E;
380 
381  // Square of non-dim frequency
382  el_pt->omega_sq_pt()= &Global_Parameters::Omega_sq;
383  }
384 
385 
386  // Same for Helmholtz mesh
387  n_element =Helmholtz_mesh_pt->nelement();
388  for(unsigned i=0;i<n_element;i++)
389  {
390  //Cast to a solid element
391  HELMHOLTZ_ELEMENT *el_pt =
392  dynamic_cast<HELMHOLTZ_ELEMENT*>(Helmholtz_mesh_pt->element_pt(i));
393 
394  //Set the pointer to square of Helmholtz wavenumber
395  el_pt->k_squared_pt() = &Global_Parameters::K_squared;
396  }
397 
398 
399  // Output meshes and their boundaries so far so we can double
400  // check the boundary enumeration
401  Solid_mesh_pt->output("solid_mesh.dat");
402  Helmholtz_mesh_pt->output("helmholtz_mesh.dat");
403  Solid_mesh_pt->output_boundaries("solid_mesh_boundary.dat");
404  Helmholtz_mesh_pt->output_boundaries("helmholtz_mesh_boundary.dat");
405 
406 
407  // Create FaceElement meshes for boundary conditions
408  //---------------------------------------------------
409 
410  // Construct the fsi traction element mesh
413 
414  // Construct the Helmholtz fsi flux element mesh
417 
418  // Create DtN elements on outer boundary of Helmholtz mesh
420 
421 
422  // Combine sub meshes
423  //-------------------
424 
425  // Solid mesh is first sub-mesh
427 
428  // Add traction sub-mesh
430 
431  // Add Helmholtz mesh
433 
434  // Add Helmholtz FSI flux mesh
436 
437  // Add Helmholtz DtN mesh
439 
440  // Build combined "global" mesh
442 
443 
444  // Solid boundary conditions:
445  //---------------------------
446  // Pin displacements on innermost boundary (boundary 0) of solid mesh
447  unsigned n_node = Solid_mesh_pt->nboundary_node(0);
449  Vector<double> x(2);
450  for(unsigned i=0;i<n_node;i++)
451  {
452  Node* nod_pt=Solid_mesh_pt->boundary_node_pt(0,i);
453  nod_pt->pin(0);
454  nod_pt->pin(1);
455  nod_pt->pin(2);
456  nod_pt->pin(3);
457 
458  // Assign displacements
459  x[0]=nod_pt->x(0);
460  x[1]=nod_pt->x(1);
462 
463  // Real part of x-displacement
464  nod_pt->set_value(0,u[0].real());
465 
466  // Imag part of x-displacement
467  nod_pt->set_value(1,u[1].real());
468 
469  // Real part of y-displacement
470  nod_pt->set_value(2,u[0].imag());
471 
472  //Imag part of y-displacement
473  nod_pt->set_value(3,u[1].imag());
474  }
475 
476 
477  // Setup fluid-structure interaction
478  //----------------------------------
480 
481  // Assign equation numbers
482  oomph_info << "Number of unknowns: " << assign_eqn_numbers() << std::endl;
483 
484  // Set output directory
486 
487  // Open trace file
488  char filename[100];
489  sprintf(filename,"%s/trace.dat",Doc_info.directory().c_str());
490  Trace_file.open(filename);
491 
492 
493 } //end of constructor
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
int i
Definition: BiCGSTAB_step_by_step.cpp:9
HelmholtzDtNMesh< HELMHOLTZ_ELEMENT > * Helmholtz_outer_boundary_mesh_pt
Pointer to mesh containing the DtN elements.
Definition: acoustic_fsi.cc:294
void create_fsi_traction_elements()
Create FSI traction elements.
Definition: acoustic_fsi.cc:575
Mesh * Helmholtz_fsi_flux_mesh_pt
Pointer to mesh of Helmholtz FSI flux elements.
Definition: acoustic_fsi.cc:291
void create_helmholtz_DtN_elements()
Create DtN face elements.
Definition: acoustic_fsi.cc:665
void create_helmholtz_fsi_flux_elements()
Create Helmholtz FSI flux elements.
Definition: acoustic_fsi.cc:621
ofstream Trace_file
Trace file.
Definition: acoustic_fsi.cc:300
TreeBasedRefineableMeshBase * Solid_mesh_pt
Pointer to solid mesh.
Definition: acoustic_fsi.cc:282
DocInfo Doc_info
DocInfo object for output.
Definition: acoustic_fsi.cc:297
void setup_interaction()
Setup interaction.
Definition: acoustic_fsi.cc:706
TreeBasedRefineableMeshBase * Helmholtz_mesh_pt
Pointer to Helmholtz mesh.
Definition: acoustic_fsi.cc:288
Mesh * FSI_traction_mesh_pt
Pointer to mesh of FSI traction elements.
Definition: acoustic_fsi.cc:285
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
std::string directory() const
Output directory.
Definition: oomph_utilities.h:524
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
Definition: mesh.h:67
unsigned long nboundary_node(const unsigned &ibound) const
Return number of nodes on a particular boundary.
Definition: mesh.h:833
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
Node *& boundary_node_pt(const unsigned &b, const unsigned &n)
Return pointer to node n on boundary b.
Definition: mesh.h:493
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
void output_boundaries(std::ostream &outfile)
Output the nodes on the boundaries (into separate tecplot zones)
Definition: mesh.cc:1064
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
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
ErrorEstimator *& spatial_error_estimator_pt()
Access to spatial error estimator.
Definition: refineable_mesh.h:143
Definition: annular_mesh.template.h:108
Definition: oomph-lib/src/generic/Vector.h:58
Definition: error_estimator.h:266
float real
Definition: datatypes.h:10
const Scalar * a
Definition: level2_cplx_impl.h:32
string Directory
Output directory.
Definition: acoustic_fsi.cc:101
unsigned El_multiplier
Multiplier for number of elements.
Definition: acoustic_fsi.cc:104
double Outer_radius
Radius of outer boundary of Helmholtz domain.
Definition: acoustic_fsi.cc:55
double K_squared
Square of wavenumber for the Helmholtz equation.
Definition: acoustic_fsi.cc:52
void solid_boundary_displacement(const Vector< double > &x, Vector< std::complex< double > > &u)
Displacement field on inner boundary of solid.
Definition: acoustic_fsi.cc:86
double E
Define the non-dimensional Young's modulus.
Definition: axisym_linear_elasticity/cylinder/cylinder.cc:49
double H_coating
Non-dim thickness of elastic coating.
Definition: acoustic_fsi.cc:61
double Omega_sq
Square of the frequency of the time dependence.
Definition: axisym_linear_elasticity/cylinder/cylinder.cc:56
string filename
Definition: MergeRestartFiles.py:39
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
list x
Definition: plotDoE.py:28

References a, Global_Parameters::Directory, oomph::DocInfo::directory(), GlobalParameters::Doc_info, Global_Parameters::E, Global_Parameters::El_multiplier, MergeRestartFiles::filename, Global_Parameters::H_coating, i, imag(), Global_Parameters::K_squared, Global_Parameters::Omega_sq, oomph::oomph_info, Global_Parameters::Outer_radius, oomph::Data::pin(), oomph::DocInfo::set_directory(), oomph::Data::set_value(), Global_Parameters::solid_boundary_displacement(), oomph::Problem_Parameter::Trace_file, plotDoE::x, and oomph::Node::x().

◆ CoatedDiskProblem() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedDiskProblem ( )

Constructor:

Member Function Documentation

◆ actions_after_adapt() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_after_adapt
virtual

Actions after adapt: Rebuild the mesh of traction elements.

Actions after adapt: Rebuild the meshes of face elements.

Reimplemented from oomph::Problem.

525 {
526  // Create fsi traction elements from all elements that are
527  // adjacent to FSI boundaries and add them to surface meshes
529 
530  // Create Helmholtz fsi flux elements
532 
533  // Create DtN elements from all elements that are
534  // adjacent to the outer boundary of Helmholtz mesh
536 
537  // Setup interaction
539 
540  // Rebuild the Problem's global mesh from its various sub-meshes
542 
543 }// end of actions_after_adapt
void rebuild_global_mesh()
Definition: problem.cc:1533

◆ actions_after_adapt() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the face meshes.

Reimplemented from oomph::Problem.

◆ actions_after_newton_solve() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update function (empty)

Reimplemented from oomph::Problem.

247 {}

◆ actions_after_newton_solve() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update function (empty)

Reimplemented from oomph::Problem.

191 {}

◆ actions_before_adapt() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_before_adapt
virtual

Actions before adapt: Wipe the mesh of traction elements.

Actions before adapt: Wipe the meshes face elements.

Reimplemented from oomph::Problem.

502 {
503  // Kill the fsi traction elements and wipe surface mesh
505 
506  // Kill Helmholtz FSI flux elements
508 
509  // Kill Helmholtz BC elements
511 
512  // Rebuild the Problem's global mesh from its various sub-meshes
514 
515 }// end of actions_before_adapt
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete (face) elements in specified mesh.
Definition: acoustic_fsi.cc:551

◆ actions_before_adapt() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the face meshes.

Reimplemented from oomph::Problem.

◆ actions_before_newton_convergence_check()

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_before_newton_convergence_check ( )
inlinevirtual

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

251  {
253  }
void setup_gamma()
of the mesh's constituent elements
Definition: helmholtz_bc_elements.h:1674

◆ actions_before_newton_solve() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update function (empty)

Reimplemented from oomph::Problem.

244 {}

◆ actions_before_newton_solve() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update function (empty)

Reimplemented from oomph::Problem.

188 {}

◆ complete_problem_setup()

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::complete_problem_setup
private

Complete problem setup: Apply boundary conditions and set physical properties

771 {
772 
773  // Solid boundary conditions:
774  //---------------------------
775  // Pin real and imag part of horizontal displacement components
776  //-------------------------------------------------------------
777  // on vertical boundaries
778  //-----------------------
779  {
780  //Loop over the nodes to pin and assign boundary displacements on
781  //solid boundary
783  for(unsigned i=0;i<n_node;i++)
784  {
786 
787  // Real part of x-displacement
788  nod_pt->pin(0);
789  nod_pt->set_value(0,0.0);
790 
791  // Imag part of x-displacement
792  nod_pt->pin(2);
793  nod_pt->set_value(2,0.0);
794  }
795  }
796  {
797  //Loop over the nodes to pin and assign boundary displacements on
798  //solid boundary
800  for(unsigned i=0;i<n_node;i++)
801  {
803 
804  // Real part of x-displacement
805  nod_pt->pin(0);
806  nod_pt->set_value(0,0.0);
807 
808  // Imag part of x-displacement
809  nod_pt->pin(2);
810  nod_pt->set_value(2,0.0);
811  }
812  }
813 
814 
815 
816  //Assign the physical properties to the elements
817  //----------------------------------------------
818  unsigned nreg=Solid_mesh_pt->nregion();
819  for (unsigned r=0;r<nreg;r++)
820  {
821  unsigned nel=Solid_mesh_pt->nregion_element(r);
822  for (unsigned e=0;e<nel;e++)
823  {
824  //Cast to a solid element
825  ELASTICITY_ELEMENT *el_pt =
826  dynamic_cast<ELASTICITY_ELEMENT*>(Solid_mesh_pt->
827  region_element_pt(r,e));
828 
829  // Set the constitutive law
830  el_pt->elasticity_tensor_pt() = Global_Parameters::E_pt[r];
831 
832  // Square of non-dim frequency
833  el_pt->omega_sq_pt()= &Global_Parameters::Omega_sq[r];
834  }
835  }
836 }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
unsigned Lower_symmetry_boundary_id
Boundary ID of lower symmetry boundary.
Definition: unstructured_acoustic_fsi.cc:247
unsigned Upper_symmetry_boundary_id
Boundary ID of upper symmetry boundary.
Definition: unstructured_acoustic_fsi.cc:244
Vector< TimeHarmonicIsotropicElasticityTensor * > E_pt
The elasticity tensors for the two regions.
Definition: unstructured_acoustic_fsi.cc:135
r
Definition: UniformPSDSelfTest.py:20

References e(), Global_Parameters::E_pt, i, Global_Parameters::Omega_sq, oomph::Data::pin(), UniformPSDSelfTest::r, and oomph::Data::set_value().

◆ create_fsi_traction_elements() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_fsi_traction_elements
private

Create FSI traction elements.

Create fsi traction elements.

576 {
577  // We're on boundary 2 of the solid mesh
578  unsigned b=2;
579 
580  // How many bulk elements are adjacent to boundary b?
581  unsigned n_element = Solid_mesh_pt->nboundary_element(b);
582 
583  // Loop over the bulk elements adjacent to boundary b
584  for(unsigned e=0;e<n_element;e++)
585  {
586  // Get pointer to the bulk element that is adjacent to boundary b
587  ELASTICITY_ELEMENT* bulk_elem_pt = dynamic_cast<ELASTICITY_ELEMENT*>(
589 
590  //Find the index of the face of element e along boundary b
591  int face_index = Solid_mesh_pt->face_index_at_boundary(b,e);
592 
593  // Create element
595  <ELASTICITY_ELEMENT,HELMHOLTZ_ELEMENT>* el_pt=
597  <ELASTICITY_ELEMENT,HELMHOLTZ_ELEMENT>(bulk_elem_pt,
598  face_index);
599  // Add to mesh
601 
602  // Associate element with bulk boundary (to allow it to access
603  // the boundary coordinates in the bulk mesh)
604  el_pt->set_boundary_number_in_bulk_mesh(b);
605 
606  // Set FSI parameter
607  el_pt->q_pt()=&Global_Parameters::Q;
608  }
609 
610 } // end of create_traction_elements
Scalar * b
Definition: benchVecAdd.cpp:17
int face_index_at_boundary(const unsigned &b, const unsigned &e) const
Definition: mesh.h:896
unsigned nboundary_element(const unsigned &b) const
Return number of finite elements that are adjacent to boundary b.
Definition: mesh.h:878
FiniteElement * boundary_element_pt(const unsigned &b, const unsigned &e) const
Return pointer to e-th finite element on boundary b.
Definition: mesh.h:840
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
Definition: helmholtz_time_harmonic_linear_elasticity_interaction.h:49
double Q
FSI parameter.
Definition: acoustic_fsi.cc:58

References b, e(), Global_Parameters::Q, and oomph::FaceElement::set_boundary_number_in_bulk_mesh().

◆ create_fsi_traction_elements() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_fsi_traction_elements ( )
private

Create FSI traction elements.

◆ create_helmholtz_ABC_elements()

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_helmholtz_ABC_elements
private

Create ABC face elements.

Create ABC elements on the outer boundary of the Helmholtz mesh

1022 {
1023  // We're on boundary 2 of the Helmholtz mesh
1024  unsigned b=2;
1025 
1026  // How many bulk elements are adjacent to boundary b?
1027  unsigned n_element = Helmholtz_mesh_pt->nboundary_element(b);
1028 
1029  // Loop over the bulk elements adjacent to boundary b?
1030  for(unsigned e=0;e<n_element;e++)
1031  {
1032  // Get pointer to the bulk element that is adjacent to boundary b
1033  HELMHOLTZ_ELEMENT* bulk_elem_pt = dynamic_cast<HELMHOLTZ_ELEMENT*>(
1035 
1036  //Find the index of the face of element e along boundary b
1037  int face_index = Helmholtz_mesh_pt->face_index_at_boundary(b,e);
1038 
1039  // Build the corresponding ABC element
1040  HelmholtzAbsorbingBCElement<HELMHOLTZ_ELEMENT>* flux_element_pt = new
1041  HelmholtzAbsorbingBCElement<HELMHOLTZ_ELEMENT>(bulk_elem_pt,face_index);
1042 
1043  // Set pointer to outer radius of artificial boundary
1045 
1046  // Set order of absorbing boundary condition
1047  flux_element_pt->abc_order_pt()=&Global_Parameters::ABC_order;
1048 
1049  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
1051  }
1052 
1053 } // end of create_helmholtz_ABC_elements
Definition: helmholtz_bc_elements.h:639
double *& outer_radius_pt()
Get pointer to radius of outer boundary (must be a cirle)
Definition: helmholtz_bc_elements.h:682
unsigned *& abc_order_pt()
Pointer to order of absorbing boundary condition.
Definition: helmholtz_bc_elements.h:655
unsigned ABC_order
Order of absorbing/appproximate boundary condition.
Definition: unstructured_acoustic_fsi.cc:117

References Global_Parameters::ABC_order, oomph::HelmholtzAbsorbingBCElement< ELEMENT >::abc_order_pt(), b, e(), Global_Parameters::Outer_radius, and oomph::HelmholtzAbsorbingBCElement< ELEMENT >::outer_radius_pt().

◆ create_helmholtz_DtN_elements()

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_helmholtz_DtN_elements
private

Create DtN face elements.

Create DtN elements on the outer boundary of the Helmholtz mesh

666 {
667  // We're on boundary 2 of the Helmholtz mesh
668  unsigned b=2;
669 
670  // How many bulk elements are adjacent to boundary b?
671  unsigned n_element = Helmholtz_mesh_pt->nboundary_element(b);
672 
673  // Loop over the bulk elements adjacent to boundary b?
674  for(unsigned e=0;e<n_element;e++)
675  {
676  // Get pointer to the bulk element that is adjacent to boundary b
677  HELMHOLTZ_ELEMENT* bulk_elem_pt = dynamic_cast<HELMHOLTZ_ELEMENT*>(
679 
680  //Find the index of the face of element e along boundary b
681  int face_index = Helmholtz_mesh_pt->face_index_at_boundary(b,e);
682 
683  // Build the corresponding DtN element
684  HelmholtzDtNBoundaryElement<HELMHOLTZ_ELEMENT>* flux_element_pt = new
685  HelmholtzDtNBoundaryElement<HELMHOLTZ_ELEMENT>(bulk_elem_pt,face_index);
686 
687  // Set pointer to the mesh that contains all the boundary condition
688  // elements on this boundary
689  flux_element_pt->set_outer_boundary_mesh_pt(
691 
692  //Add the DtN element to the mesh
694 
695  }
696 
697 } // end of create_helmholtz_DtN_elements
Definition: helmholtz_bc_elements.h:1121
void set_outer_boundary_mesh_pt(HelmholtzDtNMesh< ELEMENT > *mesh_pt)
Set mesh of all DtN boundary condition elements.
Definition: helmholtz_bc_elements.h:1169

References b, e(), and oomph::HelmholtzDtNBoundaryElement< ELEMENT >::set_outer_boundary_mesh_pt().

◆ create_helmholtz_fsi_flux_elements() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_helmholtz_fsi_flux_elements
private

Create Helmholtz FSI flux elements.

Create Helmholtz fsi flux elements.

Create Helmholtz fsii flux elements.

622 {
623 
624  // Attach to inner boundary of Helmholtz mesh (0)
625  unsigned b=0;
626 
627  // How many bulk elements are adjacent to boundary b?
628  unsigned n_element = Helmholtz_mesh_pt->nboundary_element(b);
629 
630  // Loop over the bulk elements adjacent to boundary b
631  for(unsigned e=0;e<n_element;e++)
632  {
633  // Get pointer to the bulk element that is adjacent to boundary b
634  HELMHOLTZ_ELEMENT* bulk_elem_pt = dynamic_cast<HELMHOLTZ_ELEMENT*>(
636 
637  //Find the index of the face of element e along boundary b
638  int face_index = Helmholtz_mesh_pt->face_index_at_boundary(b,e);
639 
640  // Create element
642  <HELMHOLTZ_ELEMENT,ELASTICITY_ELEMENT>* el_pt=
644  <HELMHOLTZ_ELEMENT,ELASTICITY_ELEMENT>(bulk_elem_pt,
645  face_index);
646 
647  // Add to mesh
649 
650  // Associate element with bulk boundary (to allow it to access
651  // the boundary coordinates in the bulk mesh)
652  el_pt->set_boundary_number_in_bulk_mesh(b);
653  }
654 
655 } // end of create_helmholtz_flux_elements
Definition: helmholtz_time_harmonic_linear_elasticity_interaction.h:453

References b, e(), and oomph::FaceElement::set_boundary_number_in_bulk_mesh().

◆ create_helmholtz_fsi_flux_elements() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_helmholtz_fsi_flux_elements ( )
private

Create Helmholtz FSI flux elements.

◆ create_solid_traction_elements()

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::create_solid_traction_elements
private

Create solid traction elements.

868 {
869  // Loop over pressure loaded boundaries
870  unsigned b=0;
871  unsigned nb=3;
872  for (unsigned i=0;i<nb;i++)
873  {
874  switch(i)
875  {
876  case 0:
878  break;
879 
880  case 1:
882  break;
883 
884  case 2:
886  break;
887  }
888 
889  // We're attaching face elements to region 0
890  unsigned r=0;
891 
892  // How many bulk elements are adjacent to boundary b?
893  unsigned n_element = Solid_mesh_pt->nboundary_element_in_region(b,r);
894 
895  // Loop over the bulk elements adjacent to boundary b
896  for(unsigned e=0;e<n_element;e++)
897  {
898  // Get pointer to the bulk element that is adjacent to boundary b
899  ELASTICITY_ELEMENT* bulk_elem_pt = dynamic_cast<ELASTICITY_ELEMENT*>(
900  Solid_mesh_pt->boundary_element_in_region_pt(b,r,e));
901 
902  //Find the index of the face of element e along boundary b
903  int face_index = Solid_mesh_pt->face_index_at_boundary_in_region(b,r,e);
904 
905  // Create element
908  (bulk_elem_pt,face_index);
909 
910  // Add to mesh
912 
913  // Associate element with bulk boundary (to allow it to access
914  // the boundary coordinates in the bulk mesh)
916 
917  //Set the traction function
919  }
920  }
921 } // end of create_traction_elements
unsigned Rib_divider_boundary_id
Definition: unstructured_acoustic_fsi.cc:259
Mesh * Solid_traction_mesh_pt
Pointer to mesh of solid traction elements.
Definition: unstructured_acoustic_fsi.cc:229
unsigned Lower_inner_boundary_id
Boundary ID of lower inner boundary.
Definition: unstructured_acoustic_fsi.cc:253
unsigned Upper_inner_boundary_id
Boundary ID of upper inner boundary.
Definition: unstructured_acoustic_fsi.cc:250
void set_boundary_number_in_bulk_mesh(const unsigned &b)
Set function for the boundary number in bulk mesh.
Definition: elements.h:4482
Definition: time_harmonic_linear_elasticity_traction_elements.h:77
void(*&)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction) traction_fct_pt()
Reference to the traction function pointer.
Definition: time_harmonic_linear_elasticity_traction_elements.h:169
int nb
Definition: level2_impl.h:286
void pressure_load(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double > > &traction)
Pressure load (real and imag part)
Definition: unstructured_acoustic_fsi.cc:151

References b, e(), i, nb, Global_Parameters::pressure_load(), UniformPSDSelfTest::r, oomph::FaceElement::set_boundary_number_in_bulk_mesh(), and oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::traction_fct_pt.

◆ delete_face_elements() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::delete_face_elements ( Mesh *const &  boundary_mesh_pt)
private

Delete (face) elements in specified mesh.

Delete face elements and wipe the mesh.

552 {
553  // How many surface elements are in the surface mesh
554  unsigned n_element = boundary_mesh_pt->nelement();
555 
556  // Loop over the surface elements
557  for(unsigned e=0;e<n_element;e++)
558  {
559  // Kill surface element
560  delete boundary_mesh_pt->element_pt(e);
561  }
562 
563  // Wipe the mesh
564  boundary_mesh_pt->flush_element_and_node_storage();
565 
566 } // end of delete_face_elements
void flush_element_and_node_storage()
Definition: mesh.h:407

References e(), oomph::Mesh::element_pt(), oomph::Mesh::flush_element_and_node_storage(), and oomph::Mesh::nelement().

◆ delete_face_elements() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::delete_face_elements ( Mesh *const &  boundary_mesh_pt)
private

Delete (face) elements in specified mesh.

◆ doc_solution() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::doc_solution

Doc the solution.

728 {
729 
730  ofstream some_file,some_file2;
731  char filename[100];
732 
733  // Number of plot points
734  unsigned n_plot=5;
735 
736  // Compute/output the radiated power
737  //----------------------------------
738  sprintf(filename,"%s/power%i.dat",Doc_info.directory().c_str(),
739  Doc_info.number());
740  some_file.open(filename);
741 
742  // Accumulate contribution from elements
743  double power=0.0;
744  unsigned nn_element=Helmholtz_outer_boundary_mesh_pt->nelement();
745  for(unsigned e=0;e<nn_element;e++)
746  {
750  power += el_pt->global_power_contribution(some_file);
751  }
752  some_file.close();
753  oomph_info << "Step: " << Doc_info.number()
754  << ": Q=" << Global_Parameters::Q << "\n"
755  << " k_squared=" << Global_Parameters::K_squared << "\n"
756  << " density ratio=" << Global_Parameters::Density_ratio << "\n"
757  << " omega_sq=" << Global_Parameters::Omega_sq << "\n"
758  << " Total radiated power " << power << "\n"
759  << " Axisymmetric radiated power " << "\n"
761  << std::endl;
762 
763 
764  // Write trace file
769  << power << " "
771  << std::endl;
772 
773  std::ostringstream case_string;
774  case_string << "TEXT X=10,Y=90, T=\"Q="
776  << ", k<sup>2</sup>="
778  << ", density ratio="
780  << ", omega_sq="
782  << "\"\n";
783 
784 
785  // Output displacement field
786  //--------------------------
787  sprintf(filename,"%s/elast_soln%i.dat",Doc_info.directory().c_str(),
788  Doc_info.number());
789  some_file.open(filename);
790  Solid_mesh_pt->output(some_file,n_plot);
791  some_file.close();
792 
793 
794  // Output fsi traction elements
795  //-----------------------------
796  sprintf(filename,"%s/traction_soln%i.dat",Doc_info.directory().c_str(),
797  Doc_info.number());
798  some_file.open(filename);
799  FSI_traction_mesh_pt->output(some_file,n_plot);
800  some_file.close();
801 
802 
803  // Output Helmholtz fsi flux elements
804  //-----------------------------------
805  sprintf(filename,"%s/flux_bc_soln%i.dat",Doc_info.directory().c_str(),
806  Doc_info.number());
807  some_file.open(filename);
808  Helmholtz_fsi_flux_mesh_pt->output(some_file,n_plot);
809  some_file.close();
810 
811 
812  // Output Helmholtz
813  //-----------------
814  sprintf(filename,"%s/helmholtz_soln%i.dat",Doc_info.directory().c_str(),
815  Doc_info.number());
816  some_file.open(filename);
817  Helmholtz_mesh_pt->output(some_file,n_plot);
818  some_file << case_string.str();
819  some_file.close();
820 
821 
822  // Output exact solution for Helmholtz
823  //------------------------------------
824  sprintf(filename,"%s/exact_helmholtz_soln%i.dat",Doc_info.directory().c_str(),
825  Doc_info.number());
826  some_file.open(filename);
827  Helmholtz_mesh_pt->output_fct(some_file,n_plot,
829  some_file.close();
830 
831 
832  cout << "Doced for Q=" << Global_Parameters::Q << " (step "
833  << Doc_info.number() << ")\n";
834 
835  // Increment label for output files
836  Doc_info.number()++;
837 
838 } //end doc
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
Definition: helmholtz_bc_elements.h:168
double global_power_contribution()
Definition: helmholtz_bc_elements.h:248
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
Output a given Vector function at f(n_plot) points in each element.
Definition: mesh.cc:2199
double Density_ratio
Density ratio: solid to fluid.
Definition: acoustic_fsi.cc:70
double exact_axisym_radiated_power()
Exact radiated power for axisymmetric solution.
Definition: acoustic_fsi.cc:205
void exact_axisym_potential(const Vector< double > &x, Vector< double > &soln)
Exact solution for Helmholtz potential for axisymmetric solution.
Definition: acoustic_fsi.cc:195

References Global_Parameters::Density_ratio, oomph::DocInfo::directory(), GlobalParameters::Doc_info, e(), Global_Parameters::exact_axisym_potential(), Global_Parameters::exact_axisym_radiated_power(), MergeRestartFiles::filename, oomph::HelmholtzBCElementBase< ELEMENT >::global_power_contribution(), Global_Parameters::K_squared, oomph::DocInfo::number(), Global_Parameters::Omega_sq, oomph::oomph_info, Global_Parameters::Q, and oomph::Problem_Parameter::Trace_file.

◆ doc_solution() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::doc_solution ( )

Doc the solution.

◆ setup_interaction() [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::setup_interaction
private

Setup interaction.

Setup interaction between two fields.

707 {
708  // Setup Helmholtz "pressure" load on traction elements
709  unsigned boundary_in_helmholtz_mesh=0;
711  <HELMHOLTZ_ELEMENT,2>
712  (this,boundary_in_helmholtz_mesh,Helmholtz_mesh_pt,FSI_traction_mesh_pt);
713 
714  // Setup Helmholtz flux from normal displacement interaction
715  unsigned boundary_in_solid_mesh=2;
717  <ELASTICITY_ELEMENT,2>(
718  this,boundary_in_solid_mesh,Solid_mesh_pt,Helmholtz_fsi_flux_mesh_pt);
719 }
void setup_bulk_elements_adjacent_to_face_mesh(Problem *problem_pt, Vector< unsigned > &boundary_in_bulk_mesh, Mesh *const &bulk_mesh_pt, Vector< Mesh * > &face_mesh_pt, const unsigned &interaction=0)
/ Templated helper functions for multi-domain methods using locate_zeta
Definition: multi_domain.template.cc:72

References oomph::Multi_domain_functions::setup_bulk_elements_adjacent_to_face_mesh().

◆ setup_interaction() [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
void CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::setup_interaction ( )
private

Setup interaction.

Member Data Documentation

◆ Doc_info

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
DocInfo CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Doc_info
private

DocInfo object for output.

◆ FSI_traction_mesh_pt

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
Mesh * CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::FSI_traction_mesh_pt
private

Pointer to mesh of FSI traction elements.

◆ Helmholtz_fsi_flux_mesh_pt

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
Mesh * CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Helmholtz_fsi_flux_mesh_pt
private

Pointer to mesh of Helmholtz FSI flux elements.

◆ Helmholtz_mesh_pt

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
TreeBasedRefineableMeshBase * CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Helmholtz_mesh_pt
private

Pointer to Helmholtz mesh.

◆ Helmholtz_outer_boundary_mesh_pt [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
HelmholtzDtNMesh<HELMHOLTZ_ELEMENT>* CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Helmholtz_outer_boundary_mesh_pt
private

Pointer to mesh containing the DtN elements.

◆ Helmholtz_outer_boundary_mesh_pt [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
Mesh* CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Helmholtz_outer_boundary_mesh_pt
private

Pointer to mesh containing the ABC elements.

◆ Lower_inner_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Lower_inner_boundary_id
private

Boundary ID of lower inner boundary.

◆ Lower_symmetry_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Lower_symmetry_boundary_id
private

Boundary ID of lower symmetry boundary.

◆ Outer_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Outer_boundary_id
private

Boundary ID of outer boundary.

◆ Rib_divider_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Rib_divider_boundary_id
private

◆ Solid_mesh_pt [1/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
TreeBasedRefineableMeshBase* CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Solid_mesh_pt
private

Pointer to solid mesh.

◆ Solid_mesh_pt [2/2]

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
RefineableTriangleMesh<ELASTICITY_ELEMENT>* CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Solid_mesh_pt
private

Pointer to refineable solid mesh.

◆ Solid_traction_mesh_pt

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
Mesh* CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Solid_traction_mesh_pt
private

Pointer to mesh of solid traction elements.

◆ Trace_file

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
ofstream CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Trace_file
private

Trace file.

◆ Upper_inner_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Upper_inner_boundary_id
private

Boundary ID of upper inner boundary.

◆ Upper_symmetry_boundary_id

template<class ELASTICITY_ELEMENT , class HELMHOLTZ_ELEMENT >
unsigned CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::Upper_symmetry_boundary_id
private

Boundary ID of upper symmetry boundary.


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