HelmholtzPointSourceProblem< ELEMENT > Class Template Reference

Problem class. More...

+ Inheritance diagram for HelmholtzPointSourceProblem< ELEMENT >:

Public Member Functions

 HelmholtzPointSourceProblem ()
 Constructor. More...
 
 ~HelmholtzPointSourceProblem ()
 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...
 
- 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_outer_bc_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete boundary face elements and wipe the surface mesh. More...
 
void setup_outer_boundary ()
 Set up boundary condition elements on outer boundary. More...
 
void setup_point_source ()
 Set up point source. More...
 

Private Attributes

TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshAsGeomObjectMesh_as_geom_obj_pt
 Mesh as geometric object representation of bulk mesh. More...
 
HelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
 
ofstream Trace_file
 Trace file. More...
 

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 ELEMENT>
class HelmholtzPointSourceProblem< ELEMENT >

Problem class.

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

Constructor & Destructor Documentation

◆ HelmholtzPointSourceProblem()

template<class ELEMENT >
HelmholtzPointSourceProblem< ELEMENT >::HelmholtzPointSourceProblem

Constructor.

Constructor for Helmholtz problem.

Initialise

331 {
332 
333  // Open trace file
334  Trace_file.open("RESLT/trace.dat");
335 
338 
339  // Setup "bulk" mesh
340 
341  // Create circle representing outer boundary
342  double x_c=0.0;
343  double y_c=0.0;
344  Circle* outer_circle_pt=new Circle(x_c,y_c,
346 
347  // Outer boundary
348  //---------------
349  TriangleMeshClosedCurve* outer_boundary_pt=0;
350 
351  unsigned n_segments=40;
352  Vector<TriangleMeshCurveSection*> outer_boundary_line_pt(2);
353 
354  // The intrinsic coordinates for the beginning and end of the curve
355  double s_start = 0.0;
356  double s_end = MathematicalConstants::Pi;
357  unsigned boundary_id = 0;
358  outer_boundary_line_pt[0]=
359  new TriangleMeshCurviLine(outer_circle_pt,
360  s_start,
361  s_end,
362  n_segments,
363  boundary_id);
364 
365  // The intrinsic coordinates for the beginning and end of the curve
366  s_start = MathematicalConstants::Pi;
367  s_end = 2.0*MathematicalConstants::Pi;
368  boundary_id = 1;
369  outer_boundary_line_pt[1]=
370  new TriangleMeshCurviLine(outer_circle_pt,
371  s_start,
372  s_end,
373  n_segments,
374  boundary_id);
375 
376  // Create closed curve for outer boundary
377  outer_boundary_pt=new TriangleMeshClosedCurve(outer_boundary_line_pt);
378 
379 
380  // Use the TriangleMeshParameters object for helping on the manage of the
381  // TriangleMesh parameters. The only parameter that needs to take is the
382  // outer boundary.
383  TriangleMeshParameters triangle_mesh_parameters(outer_boundary_pt);
384 
385 #ifdef ADAPTIVE
386 
387  // Build "bulk" mesh
388  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
389 
390  // Create/set error estimator
391  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
392 
393  // Choose error tolerances to force some uniform refinement
394  Bulk_mesh_pt->min_permitted_error()=0.00004;
395  Bulk_mesh_pt->max_permitted_error()=0.0001;
396 
397 #else
398 
399  // Build "bulk" mesh
400  Bulk_mesh_pt=new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
401 
402 #endif
403 
404  // Pointer to mesh containing the Helmholtz outer boundary condition
405  // elements. Specify outer radius and number of Fourier terms to be
406  // used in gamma integral
407  unsigned n_fourier=10;
410  n_fourier);
411 
412  // Create outer boundary elements from all elements that are
413  // adjacent to the outer boundary , but add them to a separate mesh.
416 
417  // Add the several sub meshes to the problem
420 
421  // Build the Problem's global mesh from its various sub-meshes
423 
424  // Complete the build of all elements so they are fully functional
425 
426  // Loop over the Helmholtz bulk elements to set up element-specific
427  // things that cannot be handled by constructor: Pass pointer to
428  // wave number squared
429  unsigned n_element = Bulk_mesh_pt->nelement();
430  for(unsigned e=0;e<n_element;e++)
431  {
432  // Upcast from GeneralisedElement to Helmholtz bulk element
433  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
434 
435  //Set the k_squared pointer
436  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
437  }
438 
439  // Set up elements on outer boundary
441 
442  // Setup point source
444 
445  // Setup equation numbering scheme
446  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
447 
448 } // end of constructor
Array< double, 1, 3 > e(1./3., 0.5, 2.)
HelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
Definition: helmholtz_point_source.cc:316
ofstream Trace_file
Trace file.
Definition: helmholtz_point_source.cc:319
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: helmholtz_point_source.cc:307
MeshAsGeomObject * Mesh_as_geom_obj_pt
Mesh as geometric object representation of bulk mesh.
Definition: helmholtz_point_source.cc:312
void setup_point_source()
Set up point source.
Definition: helmholtz_point_source.cc:511
void setup_outer_boundary()
Set up boundary condition elements on outer boundary.
Definition: helmholtz_point_source.cc:540
void create_outer_bc_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
Definition: helmholtz_point_source.cc:643
Definition: geom_objects.h:873
Definition: helmholtz_bc_elements.h:562
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
Unstructured refineable Triangle Mesh.
Definition: triangle_mesh.template.h:2249
Base class defining a closed curve for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1339
Definition: unstructured_two_d_mesh_geometry_base.h:662
Definition: triangle_mesh.template.h:94
Definition: triangle_mesh.template.h:424
Definition: oomph-lib/src/generic/Vector.h:58
Definition: error_estimator.h:266
double Pi
Definition: two_d_biharmonic.cc:235
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

References e(), GlobalParameters::K_squared, GlobalParameters::Outer_radius, oomph::MathematicalConstants::Pi, and oomph::Problem_Parameter::Trace_file.

◆ ~HelmholtzPointSourceProblem()

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

Destructor (empty)

255 {}

Member Function Documentation

◆ actions_after_adapt()

template<class ELEMENT >
void HelmholtzPointSourceProblem< ELEMENT >::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.

471 {
472 
473  // Complete the build of all elements so they are fully functional
474 
475  // Loop over the Helmholtz bulk elements to set up element-specific
476  // things that cannot be handled by constructor: Pass pointer to
477  // wave number squared
478  unsigned n_element = Bulk_mesh_pt->nelement();
479  for(unsigned e=0;e<n_element;e++)
480  {
481  // Upcast from GeneralisedElement to Helmholtz bulk element
482  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
483 
484  //Set the k_squared pointer
485  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
486  }
487 
488  // Create BC elements
489  // from all elements that are adjacent to the boundaries and add them to
490  // Helmholtz_boundary_meshes
493 
494  // Rebuild the Problem's global mesh from its various sub-meshes
496 
497  // Set up DtN mesh
499 
500  // Setup point source
502 
503 }// end of actions_after_adapt
void rebuild_global_mesh()
Definition: problem.cc:1533

References e(), and GlobalParameters::K_squared.

◆ actions_after_newton_solve()

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

265 {}

◆ actions_before_adapt()

template<class ELEMENT >
void HelmholtzPointSourceProblem< ELEMENT >::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.

455 {
456  // Wipe the boundary mesh
458 
459 
460  // Rebuild the Problem's global mesh from its various sub-meshes
462 
463 }// end of actions_before_adapt
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete boundary face elements and wipe the surface mesh.
Definition: helmholtz_point_source.cc:686

◆ actions_before_newton_convergence_check()

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

269  {
271  {
272  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
273  }
274  }
bool DtN_BC
Definition: helmholtz_point_source.cc:64

References GlobalParameters::DtN_BC.

◆ actions_before_newton_solve()

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

262 {}

◆ create_outer_bc_elements()

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

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.

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

◆ delete_face_elements()

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

Delete boundary face elements and wipe the surface mesh.

Delete face elements and wipe the boundary mesh.

687 {
688  // Loop over the surface elements
689  unsigned n_element = boundary_mesh_pt->nelement();
690  for(unsigned e=0;e<n_element;e++)
691  {
692  // Kill surface element
693  delete boundary_mesh_pt->element_pt(e);
694  }
695 
696  // Wipe the mesh
697  boundary_mesh_pt->flush_element_and_node_storage();
698 
699 } // 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().

◆ doc_solution()

template<class ELEMENT >
void HelmholtzPointSourceProblem< ELEMENT >::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

584 {
585 
586  ofstream some_file,some_file2;
587  char filename[100];
588 
589  // Number of plot points
590  unsigned npts;
591  npts=5;
592 
593  // Total radiated power
594  double power=0.0;
595 
596  // Compute/output the radiated power
597  //----------------------------------
598  sprintf(filename,"%s/power%i.dat",doc_info.directory().c_str(),
599  doc_info.number());
600  some_file.open(filename);
601 
602  // Accumulate contribution from elements
603  unsigned nn_element=Helmholtz_outer_boundary_mesh_pt->nelement();
604  for(unsigned e=0;e<nn_element;e++)
605  {
607  dynamic_cast< HelmholtzBCElementBase<ELEMENT>*>(
608  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
609  power += el_pt->global_power_contribution(some_file);
610  }
611  some_file.close();
612  oomph_info << "Total radiated power: " << power << std::endl;
613 
614  // Output solution
615  //-----------------
616  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
617  doc_info.number());
618  some_file.open(filename);
619  Bulk_mesh_pt->output(some_file,npts);
620  some_file.close();
621 
622  // Output solution in Paraview
623  //---------------------------
624  sprintf(filename,"%s/soln%i.vtu",doc_info.directory().c_str(),
625  doc_info.number());
626  some_file.open(filename);
627  Bulk_mesh_pt->output_paraview(some_file,npts);
628  some_file.close();
629 
630  // Write power to trace file
631  Trace_file << power << std::endl;
632 
633 } // end of doc
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
string filename
Definition: MergeRestartFiles.py:39
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::DocInfo::directory(), e(), MergeRestartFiles::filename, oomph::HelmholtzBCElementBase< ELEMENT >::global_power_contribution(), oomph::DocInfo::number(), oomph::oomph_info, and oomph::Problem_Parameter::Trace_file.

◆ setup_outer_boundary()

template<class ELEMENT >
void HelmholtzPointSourceProblem< ELEMENT >::setup_outer_boundary
private

Set up boundary condition elements on outer boundary.

Set pointers for elements on outer boundary.

541 {
542  // Loop over the flux elements to pass pointer to DtN
543  // BC for the outer boundary
544  unsigned n_element=Helmholtz_outer_boundary_mesh_pt->nelement();
545  for(unsigned e=0;e<n_element;e++)
546  {
547  // Dirichlet to Neumann BC
549  {
550  // Upcast from GeneralisedElement to Helmholtz flux element
552  dynamic_cast< HelmholtzDtNBoundaryElement<ELEMENT>*>(
553  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
554 
555  // Set pointer to the mesh that contains all the boundary condition
556  // elements on this boundary
558  }
559  // ABCs BC
560  else
561  {
562  // Upcast from GeneralisedElement to appropriate type
564  dynamic_cast< HelmholtzAbsorbingBCElement<ELEMENT>*>(
565  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
566 
567  // Set pointer to outer radius of artificial boundary
569 
570  // Set order of absorbing boundary condition
572  }
573  }
574 }
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_point_source()

template<class ELEMENT >
void HelmholtzPointSourceProblem< ELEMENT >::setup_point_source
private

Set up point source.

Set point source.

512 {
513 
514  // Create mesh as geometric object
515  delete Mesh_as_geom_obj_pt;
517 
518  // Position of point source in Eulerian coordinates
519  Vector<double> x_point_source(2);
520  x_point_source[0]=0.2;
521  x_point_source[1]=0.3;
522 
523  GeomObject* sub_geom_object_pt=0;
524  Vector<double> s_point_source(2);
525  Mesh_as_geom_obj_pt->locate_zeta(x_point_source,sub_geom_object_pt,
526  s_point_source);
527 
528  // Set point force
529  std::complex<double> magnitude(2.0,3.0);
530  dynamic_cast<ELEMENT*>(sub_geom_object_pt)->setup(s_point_source,magnitude);
531 }
Definition: geom_objects.h:101
Definition: mesh_as_geometric_object.h:93
void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: mesh_as_geometric_object.h:373
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
double magnitude(const Vector< double > &a)
Get the magnitude of a vector.
Definition: oomph-lib/src/generic/Vector.h:303

References oomph::VectorHelpers::magnitude(), and Flag_definition::setup().

Member Data Documentation

◆ Bulk_mesh_pt

template<class ELEMENT >
TriangleMesh<ELEMENT>* HelmholtzPointSourceProblem< ELEMENT >::Bulk_mesh_pt
private

Pointer to the "bulk" mesh.

◆ Helmholtz_outer_boundary_mesh_pt

template<class ELEMENT >
HelmholtzDtNMesh<ELEMENT>* HelmholtzPointSourceProblem< ELEMENT >::Helmholtz_outer_boundary_mesh_pt
private

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

◆ Mesh_as_geom_obj_pt

template<class ELEMENT >
MeshAsGeomObject* HelmholtzPointSourceProblem< ELEMENT >::Mesh_as_geom_obj_pt
private

Mesh as geometric object representation of bulk mesh.

◆ Trace_file

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

Trace file.


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