ScatteringProblem< ELEMENT > Class Template Reference

Problem class to compute scattering of planar wave from unit disk. More...

+ Inheritance diagram for ScatteringProblem< ELEMENT >:

Public Member Functions

 ScatteringProblem ()
 Constructor. More...
 
 ~ScatteringProblem ()
 Destructor (empty) More...
 
void doc_solution (DocInfo &doc_info)
 Doc the solution: doc_info contains labels/output directory etc. More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (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 prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void create_outer_bc_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete boundary face elements and wipe the surface mesh. More...
 
void set_prescribed_incoming_flux_pt ()
 
void setup_outer_boundary ()
 Set up boundary condition elements on outer boundary. More...
 
 ScatteringProblem ()
 Constructor. More...
 
 ~ScatteringProblem ()
 Destructor (empty) More...
 
void doc_solution (DocInfo &doc_info)
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (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 prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void create_outer_bc_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete boundary face elements and wipe the surface mesh. More...
 
void set_prescribed_incoming_flux_pt ()
 
void setup_outer_boundary ()
 Set up boundary condition elements on outer boundary. More...
 
 ScatteringProblem ()
 Constructor. More...
 
 ~ScatteringProblem ()
 Destructor (empty) More...
 
void doc_solution (DocInfo &doc_info)
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (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 prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void create_outer_bc_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete boundary face elements and wipe the surface mesh. More...
 
void set_prescribed_incoming_flux_pt ()
 
void setup_outer_boundary ()
 Set up boundary condition elements on outer boundary. 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 ()
 

Public Attributes

TwoDAnnularMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
HelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
 
MeshHelmholtz_inner_boundary_mesh_pt
 
TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
ofstream Trace_file
 Trace file. More...
 
MESH * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
- Public Attributes inherited from oomph::Problem
bool Shut_up_in_newton_solve
 

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...
 
- 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 ELEMENT>
class ScatteringProblem< ELEMENT >

Problem class to compute scattering of planar wave from unit disk.

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

Constructor & Destructor Documentation

◆ ScatteringProblem() [1/3]

template<class ELEMENT , class MESH >
ScatteringProblem< ELEMENT, MESH >::ScatteringProblem

Constructor.

Constructor for Helmholtz problem.

309 {
310 
311  // Setup "bulk" mesh
312 
313  // # of elements in theta
314  unsigned n_theta=15;
315 
316  // # of elements in radius
317  unsigned n_r=5;
318 
319  // Inner radius
320  double a=1.0;
321 
322  // Thickness of annular computational domain
323  double h=0.5;
324 
325  // Set outer radius
327 
328  // Mesh is periodic
329  bool periodic=true;
330 
331  // Full circle
332  double azimuthal_fraction=1.0;
333 
334 #ifdef ADAPTIVE
335 
336  // Build "bulk" mesh
337  Bulk_mesh_pt=
339  azimuthal_fraction,n_theta,n_r,a,h);
340 
341  // Create/set error estimator
342  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
343 
344  // Choose error tolerances to force some uniform refinement
345  Bulk_mesh_pt->min_permitted_error()=0.004;
346  Bulk_mesh_pt->max_permitted_error()=0.01;
347 
348 #else
349 
350  // Build "bulk" mesh
351  Bulk_mesh_pt=
352  new TwoDAnnularMesh<ELEMENT>(periodic,
353  azimuthal_fraction,n_theta,n_r,a,h);
354 
355 #endif
356 
357  // Pointer to mesh containing the Helmholtz outer boundary condition
358  // elements. Specify outer radius and number of Fourier terms to be
359  // used in gamma integral
362 
363  // Pointer to mesh containing the Helmholtz inner boundary condition
364  // elements. Specify outer radius
366 
367  // Create prescribed-flux elements from all elements that are
368  // adjacent to the inner boundary , but add them to a separate mesh.
370 
371  // Create outer boundary elements from all elements that are
372  // adjacent to the outer boundary , but add them to a separate mesh.
374 
375  // Add the several sub meshes to the problem
379 
380  // Build the Problem's global mesh from its various sub-meshes
382 
383  // Complete the build of all elements so they are fully functional
384 
385  // Loop over the Helmholtz bulk elements to set up element-specific
386  // things that cannot be handled by constructor: Pass pointer to
387  // wave number squared
388  unsigned n_element = Bulk_mesh_pt->nelement();
389  for(unsigned e=0;e<n_element;e++)
390  {
391  // Upcast from GeneralisedElement to Helmholtz bulk element
392  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
393 
394  //Set the k_squared pointer
395  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
396  }
397 
398  // Set up elements on outer boundary
400 
401  // Set pointer to prescribed flux function for flux elements
403 
404  // Setup equation numbering scheme
405  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
406 
407 } // end of constructor
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void create_flux_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
Definition: scattering.cc:611
Mesh * Helmholtz_inner_boundary_mesh_pt
Definition: scattering.cc:297
HelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
Definition: scattering.cc:293
TwoDAnnularMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: scattering.cc:287
void create_outer_bc_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
Definition: scattering.cc:645
void set_prescribed_incoming_flux_pt()
Definition: scattering.cc:493
void setup_outer_boundary()
Set up boundary condition elements on outer boundary.
Definition: scattering.cc:451
Definition: helmholtz_bc_elements.h:562
Definition: mesh.h:67
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
Definition: annular_mesh.template.h:108
Definition: annular_mesh.template.h:46
Definition: error_estimator.h:266
const Scalar * a
Definition: level2_cplx_impl.h:32
double Outer_radius
Radius of outer boundary (must be a circle!)
Definition: helmholtz_point_source.cc:71
double K_squared
Square of the wavenumber.
Definition: helmholtz_point_source.cc:60
unsigned N_fourier
Definition: scattering.cc:68

References a, e(), GlobalParameters::K_squared, GlobalParameters::N_fourier, and GlobalParameters::Outer_radius.

◆ ~ScatteringProblem() [1/3]

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

Destructor (empty)

231 {}

◆ ScatteringProblem() [2/3]

template<class ELEMENT >
ScatteringProblem< ELEMENT >::ScatteringProblem ( )

Constructor.

◆ ~ScatteringProblem() [2/3]

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

Destructor (empty)

230 {}

◆ ScatteringProblem() [3/3]

template<class ELEMENT >
ScatteringProblem< ELEMENT >::ScatteringProblem ( )

Constructor.

◆ ~ScatteringProblem() [3/3]

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

Destructor (empty)

231 {}

Member Function Documentation

◆ actions_after_adapt() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::actions_after_adapt
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Actions after adapt: Rebuild the face element meshes.

Reimplemented from oomph::Problem.

430 {
431  // Create prescribed-flux elements and BC elements
432  // from all elements that are adjacent to the boundaries and add them to
433  // Helmholtz_boundary_meshes
436 
437  // Rebuild the Problem's global mesh from its various sub-meshes
439 
440  // Set pointer to prescribed flux function and DtN mesh
443 
444 }// end of actions_after_adapt
void rebuild_global_mesh()
Definition: problem.cc:1533

◆ actions_after_adapt() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_adapt() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::actions_after_adapt ( )
virtual

Actions after adapt: Rebuild the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_after_newton_solve() [1/3]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

241 {}

◆ actions_after_newton_solve() [2/3]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

240 {}

◆ actions_after_newton_solve() [3/3]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

241 {}

◆ actions_before_adapt() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::actions_before_adapt
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Actions before adapt: Wipe the mesh of face elements.

Reimplemented from oomph::Problem.

414 {
415  // Kill the flux elements and wipe the boundary meshs
418 
419  // Rebuild the Problem's global mesh from its various sub-meshes
421 
422 }// end of actions_before_adapt
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete boundary face elements and wipe the surface mesh.
Definition: scattering.cc:688

◆ actions_before_adapt() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_adapt() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::actions_before_adapt ( )
virtual

Actions before adapt: Wipe the mesh of prescribed flux elements.

Reimplemented from oomph::Problem.

◆ actions_before_newton_convergence_check() [1/3]

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

245  {
247  {
248  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
249  }
250  }
bool DtN_BC
Definition: helmholtz_point_source.cc:64

References GlobalParameters::DtN_BC.

◆ actions_before_newton_convergence_check() [2/3]

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

244  {
246  {
247  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
248  }
249  }

References GlobalParameters::DtN_BC.

◆ actions_before_newton_convergence_check() [3/3]

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

245  {
247  {
248  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
249  }
250  }

References GlobalParameters::DtN_BC.

◆ actions_before_newton_solve() [1/3]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

238 {}

◆ actions_before_newton_solve() [2/3]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

237 {}

◆ actions_before_newton_solve() [3/3]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

238 {}

◆ create_flux_elements() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_inner_boundary_mesh_pt 
)

Create Helmholtz flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified surface Mesh

Create Helmholtz inner Flux Elements on the b-th boundary of the Mesh object pointed to by bulk_mesh_pt and add the elements to the Mesh object pointed to by helmholtz_inner_boundary_mesh_pt

613 {
614  // Loop over the bulk elements adjacent to boundary b
615  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
616  for(unsigned e=0;e<n_element;e++)
617  {
618  // Get pointer to the bulk element that is adjacent to boundary b
619  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
620  bulk_mesh_pt->boundary_element_pt(b,e));
621 
622  //Find the index of the face of element e along boundary b
623  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
624 
625  // Build the corresponding prescribed incoming-flux element
626  HelmholtzFluxElement<ELEMENT>* flux_element_pt = new
627  HelmholtzFluxElement<ELEMENT>(bulk_elem_pt,face_index);
628 
629  //Add the prescribed incoming-flux element to the surface mesh
630  helmholtz_inner_boundary_mesh_pt->add_element_pt(flux_element_pt);
631 
632  } //end of loop over bulk elements adjacent to boundary b
633 
634 } // end of create_flux_elements
Scalar * b
Definition: benchVecAdd.cpp:17
Definition: helmholtz_flux_elements.h:53
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

References oomph::Mesh::add_element_pt(), b, oomph::Mesh::boundary_element_pt(), e(), oomph::Mesh::face_index_at_boundary(), and oomph::Mesh::nboundary_element().

◆ create_flux_elements() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_inner_boundary_mesh_pt 
)

Create Helmholtz flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified surface Mesh

◆ create_flux_elements() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_inner_boundary_mesh_pt 
)

Create Helmholtz flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified surface Mesh

◆ create_outer_bc_elements() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::create_outer_bc_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_outer_boundary_mesh_pt 
)

Create BC elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified survace Mesh

Create outer BC elements on the b-th boundary of the Mesh object pointed to by bulk_mesh_pt and add the elements to the Mesh object pointed to by helmholtz_outer_boundary_mesh_pt.

647 {
648  // Loop over the bulk elements adjacent to boundary b?
649  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
650  for(unsigned e=0;e<n_element;e++)
651  {
652  // Get pointer to the bulk element that is adjacent to boundary b
653  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
654  bulk_mesh_pt->boundary_element_pt(b,e));
655 
656  //Find the index of the face of element e along boundary b
657  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
658 
659  // Build the corresponding outer flux element
660 
661  // Dirichlet to Neumann boundary conditon
663  {
664  HelmholtzDtNBoundaryElement<ELEMENT>* flux_element_pt = new
665  HelmholtzDtNBoundaryElement<ELEMENT>(bulk_elem_pt,face_index);
666 
667  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
668  helmholtz_outer_boundary_mesh_pt->add_element_pt(flux_element_pt);
669  }
670  // ABCs BC
671  else
672  {
673  HelmholtzAbsorbingBCElement<ELEMENT>* flux_element_pt = new
674  HelmholtzAbsorbingBCElement<ELEMENT>(bulk_elem_pt,face_index);
675 
676  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
677  helmholtz_outer_boundary_mesh_pt->add_element_pt(flux_element_pt);
678  }
679  } //end of loop over bulk elements adjacent to boundary b
680 } // end of create_outer_bc_elements
Definition: helmholtz_bc_elements.h:639
Definition: helmholtz_bc_elements.h:1121

References oomph::Mesh::add_element_pt(), b, oomph::Mesh::boundary_element_pt(), GlobalParameters::DtN_BC, e(), oomph::Mesh::face_index_at_boundary(), and oomph::Mesh::nboundary_element().

◆ create_outer_bc_elements() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::create_outer_bc_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_outer_boundary_mesh_pt 
)

Create BC elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified survace Mesh

◆ create_outer_bc_elements() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::create_outer_bc_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_outer_boundary_mesh_pt 
)

Create BC elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified survace Mesh

◆ delete_face_elements() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::delete_face_elements ( Mesh *const &  boundary_mesh_pt)

Delete boundary face elements and wipe the surface mesh.

Delete face elements and wipe the boundary mesh.

689 {
690  // Loop over the surface elements
691  unsigned n_element = boundary_mesh_pt->nelement();
692  for(unsigned e=0;e<n_element;e++)
693  {
694  // Kill surface element
695  delete boundary_mesh_pt->element_pt(e);
696  }
697 
698  // Wipe the mesh
699  boundary_mesh_pt->flush_element_and_node_storage();
700 
701 } // end of delete_outer_face_elements
void flush_element_and_node_storage()
Definition: mesh.h:407
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590

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

◆ delete_face_elements() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::delete_face_elements ( Mesh *const &  boundary_mesh_pt)

Delete boundary face elements and wipe the surface mesh.

◆ delete_face_elements() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::delete_face_elements ( Mesh *const &  boundary_mesh_pt)

Delete boundary face elements and wipe the surface mesh.

◆ doc_solution() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::doc_solution ( DocInfo doc_info)

Doc the solution: doc_info contains labels/output directory etc.

Doc the solution. DocInfo object stores flags/labels for where the output gets written to

518 {
519 
520  ofstream some_file,some_file2;
521  char filename[100];
522 
523  // Number of plot points
524  unsigned npts;
525  npts=5;
526 
527  // Compute/output the radiated power
528  //----------------------------------
529  sprintf(filename,"%s/power%i.dat",doc_info.directory().c_str(),
530  doc_info.number());
531  some_file.open(filename);
532 
533  // Accumulate contribution from elements
534  double power=0.0;
535  unsigned nn_element=Helmholtz_outer_boundary_mesh_pt->nelement();
536  for(unsigned e=0;e<nn_element;e++)
537  {
539  dynamic_cast< HelmholtzBCElementBase<ELEMENT>*>(
540  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
541  power += el_pt->global_power_contribution(some_file);
542  }
543  some_file.close();
544  oomph_info << "Total radiated power: " << power << std::endl;
545 
546  // Output solution
547  //-----------------
548  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
549  doc_info.number());
550  some_file.open(filename);
551  Bulk_mesh_pt->output(some_file,npts);
552  some_file.close();
553 
554  // Output exact solution
555  //----------------------
556  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
557  doc_info.number());
558  some_file.open(filename);
559  Bulk_mesh_pt->output_fct(some_file,npts,GlobalParameters::get_exact_u);
560  some_file.close();
561 
562  double error,norm;
563  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
564  doc_info.number());
565  some_file.open(filename);
566  Bulk_mesh_pt->compute_error(some_file,GlobalParameters::get_exact_u,
567  error,norm);
568  some_file.close();
569 
570  // Doc L2 error and norm of solution
571  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
572  oomph_info << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
573 
574 
575  // Do animation of Helmholtz solution
576  //-----------------------------------
577  unsigned nstep=40;
578  for (unsigned i=0;i<nstep;i++)
579  {
580  sprintf(filename,"%s/helmholtz_animation%i_frame%i.dat",
581  doc_info.directory().c_str(),
582  doc_info.number(),i);
583  some_file.open(filename);
584  sprintf(filename,"%s/exact_helmholtz_animation%i_frame%i.dat",
585  doc_info.directory().c_str(),
586  doc_info.number(),i);
587  some_file2.open(filename);
588  double phi=2.0*MathematicalConstants::Pi*double(i)/double(nstep-1);
589  unsigned nelem=Bulk_mesh_pt->nelement();
590  for (unsigned e=0;e<nelem;e++)
591  {
592  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(
593  Bulk_mesh_pt->element_pt(e));
594  el_pt->output_real(some_file,phi,npts);
595  el_pt->output_real_fct(some_file2,phi,npts,
597  }
598  some_file.close();
599  some_file2.close();
600  }
601 
602 } // end of doc
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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
Definition: helmholtz_bc_elements.h:168
double global_power_contribution()
Definition: helmholtz_bc_elements.h:248
double Pi
Definition: two_d_biharmonic.cc:235
void get_exact_u(const double &time, const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition: stefan_boltzmann.cc:169
string filename
Definition: MergeRestartFiles.py:39
int error
Definition: calibrate.py:297
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::DocInfo::directory(), e(), calibrate::error, MergeRestartFiles::filename, GlobalParameters::get_exact_u(), oomph::HelmholtzBCElementBase< ELEMENT >::global_power_contribution(), i, oomph::DocInfo::number(), oomph::oomph_info, oomph::MathematicalConstants::Pi, and sqrt().

◆ doc_solution() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::doc_solution ( DocInfo doc_info)

Doc the solution. DocInfo object stores flags/labels for where the output gets written to

◆ doc_solution() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::doc_solution ( DocInfo doc_info)

Doc the solution. DocInfo object stores flags/labels for where the output gets written to

◆ set_prescribed_incoming_flux_pt() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::set_prescribed_incoming_flux_pt

Set pointer to prescribed-flux function for all elements in the surface mesh on the surface of the unit disk

Set pointer to prescribed incoming-flux function for all elements in the inner boundary

494 {
495  // Loop over the flux elements to pass pointer to prescribed
496  // flux function for the inner boundary
497  unsigned n_element=Helmholtz_inner_boundary_mesh_pt->nelement();
498  for(unsigned e=0;e<n_element;e++)
499  {
500  // Upcast from GeneralisedElement to Helmholtz flux element
502  dynamic_cast< HelmholtzFluxElement<ELEMENT>*>(
504 
505  // Set the pointer to the prescribed flux function
506  el_pt->flux_fct_pt() =
508  }
509 
510 }// end of set prescribed_incoming_flux pt
HelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
Definition: helmholtz_flux_elements.h:88
void prescribed_incoming_flux(const Vector< double > &x, complex< double > &flux)
Definition: scattering.cc:156

References e(), oomph::HelmholtzFluxElement< ELEMENT >::flux_fct_pt(), and GlobalParameters::prescribed_incoming_flux().

◆ set_prescribed_incoming_flux_pt() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::set_prescribed_incoming_flux_pt ( )

Set pointer to prescribed-flux function for all elements in the surface mesh on the surface of the unit disk

◆ set_prescribed_incoming_flux_pt() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::set_prescribed_incoming_flux_pt ( )

Set pointer to prescribed-flux function for all elements in the surface mesh on the surface of the unit disk

◆ setup_outer_boundary() [1/3]

template<class ELEMENT , class MESH >
void ScatteringProblem< ELEMENT, MESH >::setup_outer_boundary

Set up boundary condition elements on outer boundary.

Set pointers for elements on outer boundary.

452 {
453  // Loop over the flux elements to pass pointer to DtN
454  // BC for the outer boundary
455  unsigned n_element=Helmholtz_outer_boundary_mesh_pt->nelement();
456  for(unsigned e=0;e<n_element;e++)
457  {
458  // Dirichlet to Neumann BC
460  {
461  // Upcast from GeneralisedElement to Helmholtz flux element
463  dynamic_cast< HelmholtzDtNBoundaryElement<ELEMENT>*>(
464  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
465 
466  // Set pointer to the mesh that contains all the boundary condition
467  // elements on this boundary
469  }
470  // ABCs BC
471  else
472  {
473  // Upcast from GeneralisedElement to appropriate type
475  dynamic_cast< HelmholtzAbsorbingBCElement<ELEMENT>*>(
476  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
477 
478  // Set pointer to outer radius of artificial boundary
480 
481  // Set order of absorbing boundary condition
483  }
484  }
485 }
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
void set_outer_boundary_mesh_pt(HelmholtzDtNMesh< ELEMENT > *mesh_pt)
Set mesh of all DtN boundary condition elements.
Definition: helmholtz_bc_elements.h:1169
unsigned ABC_order
Flag to choose wich order to use.
Definition: helmholtz_point_source.cc:68

References GlobalParameters::ABC_order, oomph::HelmholtzAbsorbingBCElement< ELEMENT >::abc_order_pt(), GlobalParameters::DtN_BC, e(), GlobalParameters::Outer_radius, oomph::HelmholtzAbsorbingBCElement< ELEMENT >::outer_radius_pt(), and oomph::HelmholtzDtNBoundaryElement< ELEMENT >::set_outer_boundary_mesh_pt().

◆ setup_outer_boundary() [2/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::setup_outer_boundary ( )

Set up boundary condition elements on outer boundary.

◆ setup_outer_boundary() [3/3]

template<class ELEMENT >
void ScatteringProblem< ELEMENT >::setup_outer_boundary ( )

Set up boundary condition elements on outer boundary.

Member Data Documentation

◆ Bulk_mesh_pt [1/3]

template<class ELEMENT >
TwoDAnnularMesh<ELEMENT>* ScatteringProblem< ELEMENT >::Bulk_mesh_pt

Pointer to the "bulk" mesh.

◆ Bulk_mesh_pt [2/3]

template<class ELEMENT >
TriangleMesh<ELEMENT>* ScatteringProblem< ELEMENT >::Bulk_mesh_pt

Pointer to the "bulk" mesh.

◆ Bulk_mesh_pt [3/3]

template<class ELEMENT >
MESH* ScatteringProblem< ELEMENT >::Bulk_mesh_pt

Pointer to the "bulk" mesh.

◆ Helmholtz_inner_boundary_mesh_pt

template<class ELEMENT >
Mesh * ScatteringProblem< ELEMENT >::Helmholtz_inner_boundary_mesh_pt

Pointer to the mesh containing the Helmholtz inner boundary condition elements

◆ Helmholtz_outer_boundary_mesh_pt

template<class ELEMENT >
HelmholtzDtNMesh< ELEMENT > * ScatteringProblem< ELEMENT >::Helmholtz_outer_boundary_mesh_pt

Pointer to mesh containing the DtN (or ABC) boundary condition elements

◆ Trace_file

template<class ELEMENT >
ofstream ScatteringProblem< ELEMENT >::Trace_file

Trace file.


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