PMLProblem< ELEMENT > Class Template Reference
+ Inheritance diagram for PMLProblem< ELEMENT >:

Public Member Functions

 PMLProblem ()
 Constructor. More...
 
 ~PMLProblem ()
 Destructor (empty) 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 doc_solution (DocInfo &doc_info)
 Doc the solution: doc_info contains labels/output directory etc. More...
 
void create_pml_meshes ()
 Create PML meshes. More...
 
void apply_boundary_conditions ()
 Apply boundary conditions. More...
 
 PMLProblem ()
 Constructor. More...
 
 ~PMLProblem ()
 Destructor (empty) 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 doc_solution (DocInfo &doc_info)
 
void create_pml_meshes ()
 Create PML meshes. More...
 
void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
 
void create_power_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_power_boundary_mesh_pt)
 
- Public Member Functions inherited from oomph::Problem
virtual void debug_hook_fct (const unsigned &i)
 
void set_analytic_dparameter (double *const &parameter_pt)
 
void unset_analytic_dparameter (double *const &parameter_pt)
 
bool is_dparameter_calculated_analytically (double *const &parameter_pt)
 
void set_analytic_hessian_products ()
 
void unset_analytic_hessian_products ()
 
bool are_hessian_products_calculated_analytically ()
 
void set_pinned_values_to_zero ()
 
bool distributed () const
 
virtual void actions_before_adapt ()
 
virtual void actions_after_adapt ()
 Actions that are to be performed after a mesh adaptation. More...
 
OomphCommunicatorcommunicator_pt ()
 access function to the oomph-lib communicator More...
 
const OomphCommunicatorcommunicator_pt () const
 access function to the oomph-lib communicator, const version More...
 
 Problem ()
 
 Problem (const Problem &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Problem &)=delete
 Broken assignment operator. More...
 
virtual ~Problem ()
 Virtual destructor to clean up memory. More...
 
Mesh *& mesh_pt ()
 Return a pointer to the global mesh. More...
 
Mesh *const & mesh_pt () const
 Return a pointer to the global mesh (const version) More...
 
Mesh *& mesh_pt (const unsigned &imesh)
 
Mesh *const & mesh_pt (const unsigned &imesh) const
 Return a pointer to the i-th submesh (const version) More...
 
unsigned nsub_mesh () const
 Return number of submeshes. More...
 
unsigned add_sub_mesh (Mesh *const &mesh_pt)
 
void flush_sub_meshes ()
 
void build_global_mesh ()
 
void rebuild_global_mesh ()
 
LinearSolver *& linear_solver_pt ()
 Return a pointer to the linear solver object. More...
 
LinearSolver *const & linear_solver_pt () const
 Return a pointer to the linear solver object (const version) More...
 
LinearSolver *& mass_matrix_solver_for_explicit_timestepper_pt ()
 
LinearSolvermass_matrix_solver_for_explicit_timestepper_pt () const
 
EigenSolver *& eigen_solver_pt ()
 Return a pointer to the eigen solver object. More...
 
EigenSolver *const & eigen_solver_pt () const
 Return a pointer to the eigen solver object (const version) More...
 
Time *& time_pt ()
 Return a pointer to the global time object. More...
 
Timetime_pt () const
 Return a pointer to the global time object (const version). More...
 
doubletime ()
 Return the current value of continuous time. More...
 
double time () const
 Return the current value of continuous time (const version) More...
 
TimeStepper *& time_stepper_pt ()
 
const TimeSteppertime_stepper_pt () const
 
TimeStepper *& time_stepper_pt (const unsigned &i)
 Return a pointer to the i-th timestepper. More...
 
ExplicitTimeStepper *& explicit_time_stepper_pt ()
 Return a pointer to the explicit timestepper. More...
 
unsigned long set_timestepper_for_all_data (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data=false)
 
virtual void shift_time_values ()
 Shift all values along to prepare for next timestep. More...
 
AssemblyHandler *& assembly_handler_pt ()
 Return a pointer to the assembly handler object. More...
 
AssemblyHandler *const & assembly_handler_pt () const
 Return a pointer to the assembly handler object (const version) More...
 
doubleminimum_dt ()
 Access function to min timestep in adaptive timestepping. More...
 
doublemaximum_dt ()
 Access function to max timestep in adaptive timestepping. More...
 
unsignedmax_newton_iterations ()
 Access function to max Newton iterations before giving up. More...
 
void problem_is_nonlinear (const bool &prob_lin)
 Access function to Problem_is_nonlinear. More...
 
doublemax_residuals ()
 
booltime_adaptive_newton_crash_on_solve_fail ()
 Access function for Time_adaptive_newton_crash_on_solve_fail. More...
 
doublenewton_solver_tolerance ()
 
void add_time_stepper_pt (TimeStepper *const &time_stepper_pt)
 
void set_explicit_time_stepper_pt (ExplicitTimeStepper *const &explicit_time_stepper_pt)
 
void initialise_dt (const double &dt)
 
void initialise_dt (const Vector< double > &dt)
 
Data *& global_data_pt (const unsigned &i)
 Return a pointer to the the i-th global data object. More...
 
void add_global_data (Data *const &global_data_pt)
 
void flush_global_data ()
 
LinearAlgebraDistribution *const & dof_distribution_pt () const
 Return the pointer to the dof distribution (read-only) More...
 
unsigned long ndof () const
 Return the number of dofs. More...
 
unsigned ntime_stepper () const
 Return the number of time steppers. More...
 
unsigned nglobal_data () const
 Return the number of global data values. More...
 
unsigned self_test ()
 Self-test: Check meshes and global data. Return 0 for OK. More...
 
void enable_store_local_dof_pt_in_elements ()
 
void disable_store_local_dof_pt_in_elements ()
 
unsigned long assign_eqn_numbers (const bool &assign_local_eqn_numbers=true)
 
void describe_dofs (std::ostream &out= *(oomph_info.stream_pt())) const
 
void enable_discontinuous_formulation ()
 
void disable_discontinuous_formulation ()
 
void get_dofs (DoubleVector &dofs) const
 
void get_dofs (const unsigned &t, DoubleVector &dofs) const
 Return vector of the t'th history value of all dofs. More...
 
void set_dofs (const DoubleVector &dofs)
 Set the values of the dofs. More...
 
void set_dofs (const unsigned &t, DoubleVector &dofs)
 Set the history values of the dofs. More...
 
void set_dofs (const unsigned &t, Vector< double * > &dof_pt)
 
void add_to_dofs (const double &lambda, const DoubleVector &increment_dofs)
 Add lambda x incremenet_dofs[l] to the l-th dof. More...
 
doubleglobal_dof_pt (const unsigned &i)
 
doubledof (const unsigned &i)
 i-th dof in the problem More...
 
double dof (const unsigned &i) const
 i-th dof in the problem (const version) More...
 
double *& dof_pt (const unsigned &i)
 Pointer to i-th dof in the problem. More...
 
doubledof_pt (const unsigned &i) const
 Pointer to i-th dof in the problem (const version) More...
 
virtual void get_inverse_mass_matrix_times_residuals (DoubleVector &Mres)
 
virtual void get_dvaluesdt (DoubleVector &f)
 
virtual void get_residuals (DoubleVector &residuals)
 Get the total residuals Vector for the problem. More...
 
virtual void get_jacobian (DoubleVector &residuals, DenseDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CRDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, CCDoubleMatrix &jacobian)
 
virtual void get_jacobian (DoubleVector &residuals, SumOfMatrices &jacobian)
 
void get_fd_jacobian (DoubleVector &residuals, DenseMatrix< double > &jacobian)
 Get the full Jacobian by finite differencing. More...
 
void get_derivative_wrt_global_parameter (double *const &parameter_pt, DoubleVector &result)
 
void get_hessian_vector_products (DoubleVectorWithHaloEntries const &Y, Vector< DoubleVectorWithHaloEntries > const &C, Vector< DoubleVectorWithHaloEntries > &product)
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &steady=true)
 Solve the eigenproblem. More...
 
void solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true)
 
virtual void get_eigenproblem_matrices (CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0)
 
void assign_eigenvector_to_dofs (DoubleVector &eigenvector)
 Assign the eigenvector passed to the function to the dofs. More...
 
void add_eigenvector_to_dofs (const double &epsilon, const DoubleVector &eigenvector)
 
void store_current_dof_values ()
 Store the current values of the degrees of freedom. More...
 
void restore_dof_values ()
 Restore the stored values of the degrees of freedom. More...
 
void enable_jacobian_reuse ()
 
void disable_jacobian_reuse ()
 Disable recycling of Jacobian in Newton iteration. More...
 
bool jacobian_reuse_is_enabled ()
 Is recycling of Jacobian in Newton iteration enabled? More...
 
booluse_predictor_values_as_initial_guess ()
 
void newton_solve ()
 Use Newton method to solve the problem. More...
 
void enable_globally_convergent_newton_method ()
 enable globally convergent Newton method More...
 
void disable_globally_convergent_newton_method ()
 disable globally convergent Newton method More...
 
void newton_solve (unsigned const &max_adapt)
 
void steady_newton_solve (unsigned const &max_adapt=0)
 
void copy (Problem *orig_problem_pt)
 
virtual Problemmake_copy ()
 
virtual void read (std::ifstream &restart_file, bool &unsteady_restart)
 
virtual void read (std::ifstream &restart_file)
 
virtual void dump (std::ofstream &dump_file) const
 
void dump (const std::string &dump_file_name) const
 
void delete_all_external_storage ()
 
virtual void symmetrise_eigenfunction_for_adaptive_pitchfork_tracking ()
 
doublebifurcation_parameter_pt () const
 
void get_bifurcation_eigenfunction (Vector< DoubleVector > &eigenfunction)
 
void activate_fold_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const bool &block_solve=true)
 
void activate_bifurcation_tracking (double *const &parameter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation, const bool &block_solve=true)
 
void activate_pitchfork_tracking (double *const &parameter_pt, const DoubleVector &symmetry_vector, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const bool &block_solve=true)
 
void activate_hopf_tracking (double *const &parameter_pt, const double &omega, const DoubleVector &null_real, const DoubleVector &null_imag, const bool &block_solve=true)
 
void deactivate_bifurcation_tracking ()
 
void reset_assembly_handler_to_default ()
 Reset the system to the standard non-augemented state. More...
 
double arc_length_step_solve (double *const &parameter_pt, const double &ds, const unsigned &max_adapt=0)
 
double arc_length_step_solve (Data *const &data_pt, const unsigned &data_index, const double &ds, const unsigned &max_adapt=0)
 
void reset_arc_length_parameters ()
 
intsign_of_jacobian ()
 
void explicit_timestep (const double &dt, const bool &shift_values=true)
 Take an explicit timestep of size dt. More...
 
void unsteady_newton_solve (const double &dt)
 
void unsteady_newton_solve (const double &dt, const bool &shift_values)
 
void unsteady_newton_solve (const double &dt, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const bool &first, const bool &shift=true)
 
double doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const unsigned &suppress_resolve_after_spatial_adapt_flag, const bool &first, const bool &shift=true)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon)
 
double adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon, const bool &shift_values)
 
void assign_initial_values_impulsive ()
 
void assign_initial_values_impulsive (const double &dt)
 
void calculate_predictions ()
 Calculate predictions. More...
 
void enable_mass_matrix_reuse ()
 
void disable_mass_matrix_reuse ()
 
bool mass_matrix_reuse_is_enabled ()
 Return whether the mass matrix is being reused. More...
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void refine_uniformly (DocInfo &doc_info)
 
void refine_uniformly_and_prune (DocInfo &doc_info)
 
void refine_uniformly ()
 
void refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform refinement for submesh i_mesh with documentation. More...
 
void refine_uniformly (const unsigned &i_mesh)
 Do uniform refinement for submesh i_mesh without documentation. More...
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh)
 
void p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info)
 
void p_refine_uniformly (DocInfo &doc_info)
 
void p_refine_uniformly_and_prune (DocInfo &doc_info)
 
void p_refine_uniformly ()
 
void p_refine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-refinement for submesh i_mesh with documentation. More...
 
void p_refine_uniformly (const unsigned &i_mesh)
 Do uniform p-refinement for submesh i_mesh without documentation. More...
 
void refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void refine_selected_elements (const unsigned &i_mesh, const Vector< RefineableElement * > &elements_to_be_refined_pt)
 
void refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void refine_selected_elements (const Vector< Vector< RefineableElement * >> &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined)
 
void p_refine_selected_elements (const unsigned &i_mesh, const Vector< PRefineableElement * > &elements_to_be_refined_pt)
 
void p_refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined)
 
void p_refine_selected_elements (const Vector< Vector< PRefineableElement * >> &elements_to_be_refined_pt)
 
unsigned unrefine_uniformly ()
 
unsigned unrefine_uniformly (const unsigned &i_mesh)
 
void p_unrefine_uniformly (DocInfo &doc_info)
 
void p_unrefine_uniformly (const unsigned &i_mesh, DocInfo &doc_info)
 Do uniform p-unrefinement for submesh i_mesh without documentation. More...
 
void adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void adapt ()
 
void p_adapt (unsigned &n_refined, unsigned &n_unrefined)
 
void p_adapt ()
 
void adapt_based_on_error_estimates (unsigned &n_refined, unsigned &n_unrefined, Vector< Vector< double >> &elemental_error)
 
void adapt_based_on_error_estimates (Vector< Vector< double >> &elemental_error)
 
void get_all_error_estimates (Vector< Vector< double >> &elemental_error)
 
void doc_errors (DocInfo &doc_info)
 Get max and min error for all elements in submeshes. More...
 
void doc_errors ()
 Get max and min error for all elements in submeshes. More...
 
void enable_info_in_newton_solve ()
 
void disable_info_in_newton_solve ()
 Disable the output of information when in the newton solver. More...
 
- Public Member Functions inherited from oomph::ExplicitTimeSteppableObject
 ExplicitTimeSteppableObject ()
 Empty constructor. More...
 
 ExplicitTimeSteppableObject (const ExplicitTimeSteppableObject &)=delete
 Broken copy constructor. More...
 
void operator= (const ExplicitTimeSteppableObject &)=delete
 Broken assignment operator. More...
 
virtual ~ExplicitTimeSteppableObject ()
 Empty destructor. More...
 
virtual void actions_before_explicit_stage ()
 
virtual void actions_after_explicit_stage ()
 

Private Attributes

TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshPML_right_mesh_pt
 Pointer to the right PML mesh. More...
 
MeshPML_top_mesh_pt
 Pointer to the top PML mesh. More...
 
MeshPML_left_mesh_pt
 Pointer to the left PML mesh. More...
 
MeshPML_bottom_mesh_pt
 Pointer to the bottom PML mesh. More...
 
MeshPML_top_right_mesh_pt
 Pointer to the top right corner PML mesh. More...
 
MeshPML_top_left_mesh_pt
 Pointer to the top left corner PML mesh. More...
 
MeshPML_bottom_right_mesh_pt
 Pointer to the bottom right corner PML mesh. More...
 
MeshPML_bottom_left_mesh_pt
 Pointer to the bottom left corner PML mesh. More...
 
ofstream Trace_file
 Trace file. More...
 
MeshHelmholtz_inner_boundary_mesh_pt
 
MeshHelmholtz_power_boundary_mesh_pt
 Pointer to mesh of elements that compute the radiated power. 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_convergence_check ()
 
virtual void actions_before_newton_step ()
 
virtual void actions_after_newton_step ()
 
virtual void actions_before_implicit_timestep ()
 
virtual void actions_after_implicit_timestep ()
 
virtual void actions_after_implicit_timestep_and_error_estimation ()
 
virtual void actions_before_explicit_timestep ()
 Actions that should be performed before each explicit time step. More...
 
virtual void actions_after_explicit_timestep ()
 Actions that should be performed after each explicit time step. More...
 
virtual void actions_before_read_unstructured_meshes ()
 
virtual void actions_after_read_unstructured_meshes ()
 
virtual void actions_after_change_in_global_parameter (double *const &parameter_pt)
 
virtual void actions_after_change_in_bifurcation_parameter ()
 
virtual void actions_after_parameter_increase (double *const &parameter_pt)
 
doubledof_derivative (const unsigned &i)
 
doubledof_current (const unsigned &i)
 
virtual void set_initial_condition ()
 
virtual double global_temporal_error_norm ()
 
unsigned newton_solve_continuation (double *const &parameter_pt)
 
unsigned newton_solve_continuation (double *const &parameter_pt, DoubleVector &z)
 
void calculate_continuation_derivatives (double *const &parameter_pt)
 
void calculate_continuation_derivatives (const DoubleVector &z)
 
void calculate_continuation_derivatives_fd (double *const &parameter_pt)
 
bool does_pointer_correspond_to_problem_data (double *const &parameter_pt)
 
void set_consistent_pinned_values_for_continuation ()
 
- Protected Attributes inherited from oomph::Problem
Vector< Problem * > Copy_of_problem_pt
 
std::map< double *, boolCalculate_dparameter_analytic
 
bool Calculate_hessian_products_analytic
 
LinearAlgebraDistributionDof_distribution_pt
 
Vector< double * > Dof_pt
 Vector of pointers to dofs. More...
 
DoubleVectorWithHaloEntries Element_count_per_dof
 
double Relaxation_factor
 
double Newton_solver_tolerance
 
unsigned Max_newton_iterations
 Maximum number of Newton iterations. More...
 
unsigned Nnewton_iter_taken
 
Vector< doubleMax_res
 Maximum residuals at start and after each newton iteration. More...
 
double Max_residuals
 
bool Time_adaptive_newton_crash_on_solve_fail
 
bool Jacobian_reuse_is_enabled
 Is re-use of Jacobian in Newton iteration enabled? Default: false. More...
 
bool Jacobian_has_been_computed
 
bool Problem_is_nonlinear
 
bool Pause_at_end_of_sparse_assembly
 
bool Doc_time_in_distribute
 
unsigned Sparse_assembly_method
 
unsigned Sparse_assemble_with_arrays_initial_allocation
 
unsigned Sparse_assemble_with_arrays_allocation_increment
 
Vector< Vector< unsigned > > Sparse_assemble_with_arrays_previous_allocation
 
double Numerical_zero_for_sparse_assembly
 
double FD_step_used_in_get_hessian_vector_products
 
bool Mass_matrix_reuse_is_enabled
 
bool Mass_matrix_has_been_computed
 
bool Discontinuous_element_formulation
 
double Minimum_dt
 Minimum desired dt: if dt falls below this value, exit. More...
 
double Maximum_dt
 Maximum desired dt. More...
 
double DTSF_max_increase
 
double DTSF_min_decrease
 
double Minimum_dt_but_still_proceed
 
bool Scale_arc_length
 Boolean to control whether arc-length should be scaled. More...
 
double Desired_proportion_of_arc_length
 Proportion of the arc-length to taken by the parameter. More...
 
double Theta_squared
 
int Sign_of_jacobian
 Storage for the sign of the global Jacobian. More...
 
double Continuation_direction
 
double Parameter_derivative
 Storage for the derivative of the global parameter wrt arc-length. More...
 
double Parameter_current
 Storage for the present value of the global parameter. More...
 
bool Use_continuation_timestepper
 Boolean to control original or new storage of dof stuff. More...
 
unsigned Dof_derivative_offset
 
unsigned Dof_current_offset
 
Vector< doubleDof_derivative
 Storage for the derivative of the problem variables wrt arc-length. More...
 
Vector< doubleDof_current
 Storage for the present values of the variables. More...
 
double Ds_current
 Storage for the current step value. More...
 
unsigned Desired_newton_iterations_ds
 
double Minimum_ds
 Minimum desired value of arc-length. More...
 
bool Bifurcation_detection
 Boolean to control bifurcation detection via determinant of Jacobian. More...
 
bool Bisect_to_find_bifurcation
 Boolean to control wheter bisection is used to located bifurcation. More...
 
bool First_jacobian_sign_change
 Boolean to indicate whether a sign change has occured in the Jacobian. More...
 
bool Arc_length_step_taken
 Boolean to indicate whether an arc-length step has been taken. More...
 
bool Use_finite_differences_for_continuation_derivatives
 
OomphCommunicatorCommunicator_pt
 The communicator for this problem. More...
 
bool Always_take_one_newton_step
 
double Timestep_reduction_factor_after_nonconvergence
 
bool Keep_temporal_error_below_tolerance
 
- Static Protected Attributes inherited from oomph::Problem
static ContinuationStorageScheme Continuation_time_stepper
 Storage for the single static continuation timestorage object. More...
 

Detailed Description

template<class ELEMENT>
class PMLProblem< ELEMENT >

////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// Problem class to demonstrate use of perfectly matched layers for Helmholtz problems.

Constructor & Destructor Documentation

◆ PMLProblem() [1/2]

template<class ELEMENT >
PMLProblem< ELEMENT >::PMLProblem

Constructor.

Constructor for Helmholtz problem.

166 {
167 
168  // Open trace file
169  Trace_file.open("RESLT/trace.dat");
170 
171  // Create circle representing inner boundary
172  double a=0.2;
173  double x_c=0.0;
174  double y_c=0.0;
175  Circle* inner_circle_pt=new Circle(x_c,y_c,a);
176 
177  // Outer boundary
178  //---------------
179  TriangleMeshClosedCurve* outer_boundary_pt=0;
180 
181  unsigned n_segments = 20;
182  Vector<TriangleMeshCurveSection*> outer_boundary_line_pt(4);
183 
184  // Each polyline only has three vertices, provide storage for their
185  // coordinates
186  Vector<Vector<double> > vertex_coord(2);
187  for(unsigned i=0;i<2;i++)
188  {
189  vertex_coord[i].resize(2);
190  }
191 
192  // First polyline
193  vertex_coord[0][0]=-2.0;
194  vertex_coord[0][1]=-2.0;
195  vertex_coord[1][0]=-2.0;
196  vertex_coord[1][1]=2.0;
197 
198  // Build the 1st boundary polyline
199  unsigned boundary_id=2;
200  outer_boundary_line_pt[0] = new TriangleMeshPolyLine(vertex_coord,
201  boundary_id);
202 
203  // Second boundary polyline
204  vertex_coord[0][0]=-2.0;
205  vertex_coord[0][1]=2.0;
206  vertex_coord[1][0]=2.0;
207  vertex_coord[1][1]=2.0;
208 
209  // Build the 2nd boundary polyline
210  boundary_id=3;
211  outer_boundary_line_pt[1] = new TriangleMeshPolyLine(vertex_coord,
212  boundary_id);
213 
214  // Third boundary polyline
215  vertex_coord[0][0]=2.0;
216  vertex_coord[0][1]=2.0;
217  vertex_coord[1][0]=2.0;
218  vertex_coord[1][1]=-2.0;
219 
220  // Build the 3rd boundary polyline
221  boundary_id=4;
222  outer_boundary_line_pt[2] = new TriangleMeshPolyLine(vertex_coord,
223  boundary_id);
224 
225  // Fourth boundary polyline
226  vertex_coord[0][0]=2.0;
227  vertex_coord[0][1]=-2.0;
228  vertex_coord[1][0]=-2.0;
229  vertex_coord[1][1]=-2.0;
230 
231  // Build the 4th boundary polyline
232  boundary_id=5;
233  outer_boundary_line_pt[3] = new TriangleMeshPolyLine(vertex_coord,
234  boundary_id);
235 
236  // Create the triangle mesh polygon for outer boundary
237  outer_boundary_pt = new TriangleMeshPolygon(outer_boundary_line_pt);
238 
239  // Inner circular boundary
240  //------------------------
241  Vector<TriangleMeshCurveSection*> inner_boundary_line_pt(2);
242 
243  // The intrinsic coordinates for the beginning and end of the curve
244  double s_start = 0.0;
245  double s_end = MathematicalConstants::Pi;
246  boundary_id = 0;
247  inner_boundary_line_pt[0]=
248  new TriangleMeshCurviLine(inner_circle_pt,
249  s_start,
250  s_end,
251  n_segments,
252  boundary_id);
253 
254  // The intrinsic coordinates for the beginning and end of the curve
255  s_start = MathematicalConstants::Pi;
256  s_end = 2.0*MathematicalConstants::Pi;
257  boundary_id = 1;
258  inner_boundary_line_pt[1]=
259  new TriangleMeshCurviLine(inner_circle_pt,
260  s_start,
261  s_end,
262  n_segments,
263  boundary_id);
264 
265 
266  // Combine to hole
267  //----------------
269  Vector<double> hole_coords(2);
270  hole_coords[0]=0.0;
271  hole_coords[1]=0.0;
272  hole_pt[0]=new TriangleMeshClosedCurve(inner_boundary_line_pt,
273  hole_coords);
274 
275 
276  // Use the TriangleMeshParameters object for helping on the manage
277  // of the TriangleMesh parameters. The only parameter that needs to take
278  // is the outer boundary.
279  TriangleMeshParameters triangle_mesh_parameters(outer_boundary_pt);
280 
281  // Specify the closed curve using the TriangleMeshParameters object
282  triangle_mesh_parameters.internal_closed_curve_pt() = hole_pt;
283 
284  // Target element size in bulk mesh
285  triangle_mesh_parameters.element_area() = 0.1;
286 
287 #ifdef ADAPTIVE
288 
289  // Build adaptive "bulk" mesh
290  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
291 
292  // Create/set error estimator
293  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
294 
295  // Choose error tolerances to force some uniform refinement
296  Bulk_mesh_pt->min_permitted_error()=0.00004;
297  Bulk_mesh_pt->max_permitted_error()=0.0001;
298 
299 #else
300 
301  // Build "bulk" mesh
302  Bulk_mesh_pt=new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
303 
304 #endif
305 
306  // Create the main triangular mesh
308 
309  // Create PML meshes and add them to the global mesh
311 
312  // Build the entire mesh from its submeshes
314 
315  // Let's have a look where the boundaries are
316  this->mesh_pt()->output("global_mesh.dat");
317  this->mesh_pt()->output_boundaries("global_mesh_boundary.dat");
318 
319  // Complete the build of all elements so they are fully functional
320  unsigned n_element = this->mesh_pt()->nelement();
321  for(unsigned e=0;e<n_element;e++)
322  {
323  // Upcast from GeneralisedElement to Helmholtz bulk element
324  PMLHelmholtzEquations<2> *el_pt =
325  dynamic_cast<PMLHelmholtzEquations<2>*>(mesh_pt()->element_pt(e));
326 
327  //Set the k_squared double pointer
329  }
330 
331  // Apply boundary conditions
333 
334  // Setup equation numbering scheme
335  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
336 
337 } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void apply_boundary_conditions()
Apply boundary conditions.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:422
ofstream Trace_file
Trace file.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:155
void create_pml_meshes()
Create PML meshes.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:530
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:125
Definition: geom_objects.h:873
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
void output_boundaries(std::ostream &outfile)
Output the nodes on the boundaries (into separate tecplot zones)
Definition: mesh.cc:1064
Definition: pml_helmholtz_elements.h:62
double *& k_squared_pt()
Get pointer to k_squared.
Definition: pml_helmholtz_elements.h:97
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
Mesh *& mesh_pt()
Return a pointer to the global mesh.
Definition: problem.h:1280
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
Class defining a polyline for use in Triangle Mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:868
Class defining a closed polygon for the Triangle mesh generation.
Definition: unstructured_two_d_mesh_geometry_base.h:1451
Definition: triangle_mesh.template.h:424
Definition: oomph-lib/src/generic/Vector.h:58
Definition: error_estimator.h:266
const Scalar * a
Definition: level2_cplx_impl.h:32
double Pi
Definition: two_d_biharmonic.cc:235
double K_squared
Square of the wavenumber.
Definition: helmholtz_point_source.cc:60

References a, e(), oomph::TriangleMeshParameters::element_area(), i, oomph::TriangleMeshParameters::internal_closed_curve_pt(), GlobalParameters::K_squared, oomph::PMLHelmholtzEquations< DIM >::k_squared_pt(), oomph::MathematicalConstants::Pi, and oomph::Problem_Parameter::Trace_file.

◆ ~PMLProblem() [1/2]

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

Destructor (empty)

84 {}

◆ PMLProblem() [2/2]

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

Constructor.

◆ ~PMLProblem() [2/2]

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

Destructor (empty)

251 {};

Member Function Documentation

◆ actions_after_newton_solve() [1/2]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

90 {}

◆ actions_after_newton_solve() [2/2]

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

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

257 {};

◆ actions_before_newton_solve() [1/2]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

87 {}

◆ actions_before_newton_solve() [2/2]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

254 {} ;

◆ apply_boundary_conditions()

template<class ELEMENT >
void PMLProblem< ELEMENT >::apply_boundary_conditions

Apply boundary conditions.

423 {
424 
425  // Boundary conditions are set on the surface of the circle
426  // as a constant nonzero Dirichlet boundary condition
427  unsigned n_bound = Bulk_mesh_pt->nboundary();
428 
429  for(unsigned b=0;b<n_bound;b++)
430  {
431  unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
432  for (unsigned n=0;n<n_node;n++)
433  {
434  if ((b==0) || (b==1))
435  {
436  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(b,n);
437  nod_pt->pin(0);
438  nod_pt->pin(1);
439 
440  nod_pt->set_value(0,0.1);
441  nod_pt->set_value(1,0.0);
442  }
443  }
444  }
445 
446 }// end of apply_boundary_conditions
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Definition: nodes.h:906

References b, n, oomph::Data::pin(), and oomph::Data::set_value().

◆ create_flux_elements()

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

750 {
751  // Loop over the bulk elements adjacent to boundary b
752  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
753  for(unsigned e=0;e<n_element;e++)
754  {
755  // Get pointer to the bulk element that is adjacent to boundary b
756  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
757  bulk_mesh_pt->boundary_element_pt(b,e));
758 
759  //Find the index of the face of element e along boundary b
760  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
761 
762  // Build the corresponding prescribed incoming-flux element
763  PMLHelmholtzFluxElement<ELEMENT>* flux_element_pt = new
764  PMLHelmholtzFluxElement<ELEMENT>(bulk_elem_pt,face_index);
765 
766  //Add the prescribed incoming-flux element to the surface mesh
767  helmholtz_inner_boundary_mesh_pt->add_element_pt(flux_element_pt);
768 
769  // Set the pointer to the prescribed flux function
770  flux_element_pt->flux_fct_pt() =
772 
773  } //end of loop over bulk elements adjacent to boundary b
774 
775 } // end of create_flux_elements
int face_index_at_boundary(const unsigned &b, const unsigned &e) const
Definition: mesh.h:896
unsigned nboundary_element(const unsigned &b) const
Return number of finite elements that are adjacent to boundary b.
Definition: mesh.h:878
FiniteElement * boundary_element_pt(const unsigned &b, const unsigned &e) const
Return pointer to e-th finite element on boundary b.
Definition: mesh.h:840
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
Definition: pml_helmholtz_flux_elements.h:552
PMLHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
Definition: pml_helmholtz_flux_elements.h:582
void prescribed_incoming_flux(const Vector< double > &x, complex< double > &flux)
Definition: scattering.cc:156

References oomph::Mesh::add_element_pt(), b, oomph::Mesh::boundary_element_pt(), e(), oomph::Mesh::face_index_at_boundary(), oomph::PMLHelmholtzFluxElement< ELEMENT >::flux_fct_pt(), oomph::Mesh::nboundary_element(), and GlobalParameters::prescribed_incoming_flux().

◆ create_pml_meshes() [1/2]

template<class ELEMENT >
void PMLProblem< ELEMENT >::create_pml_meshes

Create PML meshes.

Create PML meshes and add them to the problem's sub-meshes.

531 {
532 
533  // Bulk mesh left boundary id
534  unsigned int left_boundary_id = 2;
535 
536  // Bulk mesh top boundary id
537  unsigned int top_boundary_id = 3;
538 
539  // Bulk mesh right boundary id
540  unsigned int right_boundary_id = 4;
541 
542  // Bulk mesh bottom boundary id
543  unsigned int bottom_boundary_id = 5;
544 
545  // PML width in elements for the right layer
546  unsigned n_x_right_pml = 3;
547 
548  // PML width in elements for the top layer
549  unsigned n_y_top_pml = 3;
550 
551  // PML width in elements for the left layer
552  unsigned n_x_left_pml = 3;
553 
554  // PML width in elements for the left layer
555  unsigned n_y_bottom_pml = 3;
556 
557  // Outer physical length of the PML layers
558  // defaults to 0.2, so 10% of the size of the
559  // physical domain
560  double width_x_right_pml = 0.2;
561  double width_y_top_pml = 0.2;
562  double width_x_left_pml = 0.2;
563  double width_y_bottom_pml = 0.2;
564 
565  // Build the PML meshes based on the new adapted triangular mesh
569  (Bulk_mesh_pt,right_boundary_id,
570  n_x_right_pml, width_x_right_pml);
574  (Bulk_mesh_pt, top_boundary_id,
575  n_y_top_pml, width_y_top_pml);
579  (Bulk_mesh_pt, left_boundary_id,
580  n_x_left_pml, width_x_left_pml);
584  (Bulk_mesh_pt, bottom_boundary_id,
585  n_y_bottom_pml, width_y_bottom_pml);
586 
587  // Add submeshes to the global mesh
592 
593  // Rebuild corner PML meshes
598  Bulk_mesh_pt, right_boundary_id);
599 
604  Bulk_mesh_pt, right_boundary_id);
605 
610  Bulk_mesh_pt, left_boundary_id);
611 
616  Bulk_mesh_pt, left_boundary_id);
617 
618  // Add submeshes to the global mesh
623 
624 } // end of create_pml_meshes
Mesh * PML_bottom_left_mesh_pt
Pointer to the bottom left corner PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:152
Mesh * PML_bottom_right_mesh_pt
Pointer to the bottom right corner PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:149
Mesh * PML_bottom_mesh_pt
Pointer to the bottom PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:140
Mesh * PML_top_mesh_pt
Pointer to the top PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:134
Mesh * PML_top_right_mesh_pt
Pointer to the top right corner PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:143
Mesh * PML_right_mesh_pt
Pointer to the right PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:131
Mesh * PML_top_left_mesh_pt
Pointer to the top left corner PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:146
Mesh * PML_left_mesh_pt
Pointer to the left PML mesh.
Definition: pml_helmholtz/scattering/unstructured_two_d_helmholtz.cc:137
Definition: pml_meshes.h:48
Mesh * create_bottom_right_pml_mesh(Mesh *pml_right_mesh_pt, Mesh *pml_bottom_mesh_pt, Mesh *bulk_mesh_pt, const unsigned &right_boundary_id, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:2103
Mesh * create_top_right_pml_mesh(Mesh *pml_right_mesh_pt, Mesh *pml_top_mesh_pt, Mesh *bulk_mesh_pt, const unsigned &right_boundary_id, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:1968
Mesh * create_top_pml_mesh(Mesh *bulk_mesh_pt, const unsigned &top_boundary_id, const unsigned &n_y_top_pml, const double &width_y_top_pml, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:1646
Mesh * create_top_left_pml_mesh(Mesh *pml_left_mesh_pt, Mesh *pml_top_mesh_pt, Mesh *bulk_mesh_pt, const unsigned &left_boundary_id, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:2234
Mesh * create_bottom_left_pml_mesh(Mesh *pml_left_mesh_pt, Mesh *pml_bottom_mesh_pt, Mesh *bulk_mesh_pt, const unsigned &left_boundary_id, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:2367
Mesh * create_bottom_pml_mesh(Mesh *bulk_mesh_pt, const unsigned &bottom_boundary_id, const unsigned &n_y_bottom_pml, const double &width_y_bottom_pml, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:1860
Mesh * create_left_pml_mesh(Mesh *bulk_mesh_pt, const unsigned &left_boundary_id, const unsigned &n_x_left_pml, const double &width_x_left_pml, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:1752
Mesh * create_right_pml_mesh(Mesh *bulk_mesh_pt, const unsigned &right_boundary_id, const unsigned &n_x_right_pml, const double &width_x_right_pml, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: pml_meshes.h:1539

References oomph::TwoDimensionalPMLHelper::create_bottom_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_right_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_right_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_pml_mesh(), and oomph::TwoDimensionalPMLHelper::create_top_right_pml_mesh().

◆ create_pml_meshes() [2/2]

template<class ELEMENT >
void PMLProblem< ELEMENT >::create_pml_meshes ( )

Create PML meshes.

◆ create_power_elements()

template<class ELEMENT >
void PMLProblem< ELEMENT >::create_power_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  helmholtz_power_boundary_mesh_pt 
)

Create Helmholtz power 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_power_boundary_mesh_pt

786 {
787  // Loop over the bulk elements adjacent to boundary b
788  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
789  for(unsigned e=0;e<n_element;e++)
790  {
791  // Get pointer to the bulk element that is adjacent to boundary b
792  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
793  bulk_mesh_pt->boundary_element_pt(b,e));
794 
795  //Find the index of the face of element e along boundary b
796  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
797 
798  // Build the corresponding prescribed power element
799  PMLHelmholtzPowerElement<ELEMENT>* power_element_pt = new
800  PMLHelmholtzPowerElement<ELEMENT>(bulk_elem_pt,face_index);
801 
802  //Add the prescribed power element to the surface mesh
803  helmholtz_power_boundary_mesh_pt->add_element_pt(power_element_pt);
804 
805  } //end of loop over bulk elements adjacent to boundary b
806 
807 } // end of create_power_elements
Definition: pml_helmholtz_flux_elements.h:53

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

◆ doc_solution() [1/2]

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

454 {
455 
456  ofstream some_file,some_file2;
457  char filename[100];
458 
459  // Number of plot points
460  unsigned npts;
461  npts=5;
462 
463  // Output solution
464  //-----------------
465  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
466  doc_info.number());
467  some_file.open(filename);
468  Bulk_mesh_pt->output(some_file,npts);
469  some_file.close();
470 
471  // Output coarse solution
472  //-----------------------
473  sprintf(filename,"%s/coarse_soln%i.dat",doc_info.directory().c_str(),
474  doc_info.number());
475  some_file.open(filename);
476  unsigned npts_coarse=2;
477  Bulk_mesh_pt->output(some_file,npts_coarse);
478  some_file.close();
479 
480 
481  // Output solution within pml domains
482  //-----------------------------------
483  sprintf(filename,"%s/pml_soln%i.dat",doc_info.directory().c_str(),
484  doc_info.number());
485  some_file.open(filename);
486  PML_top_mesh_pt->output(some_file,npts);
487  PML_right_mesh_pt->output(some_file,npts);
488  PML_bottom_mesh_pt->output(some_file,npts);
489  PML_left_mesh_pt->output(some_file,npts);
490  PML_top_right_mesh_pt->output(some_file,npts);
491  PML_bottom_right_mesh_pt->output(some_file,npts);
492  PML_top_left_mesh_pt->output(some_file,npts);
493  PML_bottom_left_mesh_pt->output(some_file,npts);
494  some_file.close();
495 
496 
497 
498 
499  // Write norm of solution to trace file
500  double norm=0.0;
501  Bulk_mesh_pt->compute_norm(norm);
502  Trace_file << norm << std::endl;
503 
504  // //Do animation of Helmholtz solution
505  // //-----------------------------------
506  // unsigned nstep=40;
507  // for (unsigned i=0;i<nstep;i++)
508  // {
509  // sprintf(filename,"%s/helmholtz_animation%i_frame%i.dat",
510  // doc_info.directory().c_str(),
511  // doc_info.number(),i);
512  // some_file.open(filename);
513  // double phi=2.0*MathematicalConstants::Pi*double(i)/double(nstep-1);
514  // unsigned nelem=Bulk_mesh_pt->nelement();
515  // for (unsigned e=0;e<nelem;e++)
516  // {
517  // ELEMENT* el_pt=dynamic_cast<ELEMENT*>(
518  // Bulk_mesh_pt->element_pt(e));
519  // el_pt->output_real(some_file,phi,npts);
520  // }
521  // some_file.close();
522  // }
523 
524 } // 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
string filename
Definition: MergeRestartFiles.py:39

References oomph::DocInfo::directory(), MergeRestartFiles::filename, oomph::DocInfo::number(), and oomph::Problem_Parameter::Trace_file.

◆ doc_solution() [2/2]

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

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

Member Data Documentation

◆ Bulk_mesh_pt

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

Pointer to the "bulk" mesh.

◆ Helmholtz_inner_boundary_mesh_pt

template<class ELEMENT >
Mesh* PMLProblem< ELEMENT >::Helmholtz_inner_boundary_mesh_pt
private

Pointer to the mesh containing the Helmholtz inner boundary condition elements

◆ Helmholtz_power_boundary_mesh_pt

template<class ELEMENT >
Mesh* PMLProblem< ELEMENT >::Helmholtz_power_boundary_mesh_pt
private

Pointer to mesh of elements that compute the radiated power.

◆ PML_bottom_left_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_bottom_left_mesh_pt
private

Pointer to the bottom left corner PML mesh.

◆ PML_bottom_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_bottom_mesh_pt
private

Pointer to the bottom PML mesh.

◆ PML_bottom_right_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_bottom_right_mesh_pt
private

Pointer to the bottom right corner PML mesh.

◆ PML_left_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_left_mesh_pt
private

Pointer to the left PML mesh.

◆ PML_right_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_right_mesh_pt
private

Pointer to the right PML mesh.

◆ PML_top_left_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_top_left_mesh_pt
private

Pointer to the top left corner PML mesh.

◆ PML_top_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_top_mesh_pt
private

Pointer to the top PML mesh.

◆ PML_top_right_mesh_pt

template<class ELEMENT >
Mesh * PMLProblem< ELEMENT >::PML_top_right_mesh_pt
private

Pointer to the top right corner PML mesh.

◆ Trace_file

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

Trace file.


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