FourierDecomposedHelmholtzProblem< ELEMENT > Class Template Reference

Problem class. More...

+ Inheritance diagram for FourierDecomposedHelmholtzProblem< ELEMENT >:

Public Member Functions

 FourierDecomposedHelmholtzProblem ()
 Constructor. More...
 
 ~FourierDecomposedHelmholtzProblem ()
 Destructor (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem after solve (empty) More...
 
void doc_solution (DocInfo &doc_info)
 Doc the solution: doc_info contains labels/output directory etc. More...
 
void actions_before_newton_convergence_check ()
 Recompute gamma integral before checking Newton residuals. More...
 
void check_gamma (DocInfo &doc_info)
 Check gamma computation. More...
 
 FourierDecomposedHelmholtzProblem ()
 Constructor. More...
 
 ~FourierDecomposedHelmholtzProblem ()
 Destructor (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem after solve (empty) More...
 
void doc_solution (DocInfo &doc_info)
 
void actions_before_newton_convergence_check ()
 Recompute gamma integral before checking Newton residuals. More...
 
void actions_before_adapt ()
 Actions before adapt: Wipe the mesh of prescribed flux elements. More...
 
void actions_after_adapt ()
 Actions after adapt: Rebuild the mesh of prescribed flux elements. More...
 
void check_gamma (DocInfo &doc_info)
 Check gamma computation. 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 ()
 Create BC elements on outer boundary. More...
 
void create_flux_elements_on_inner_boundary ()
 Create flux elements on inner boundary. More...
 
void create_outer_bc_elements ()
 Create BC elements on outer boundary. More...
 
void create_flux_elements_on_inner_boundary ()
 Create flux elements on inner boundary. More...
 
void delete_face_elements (Mesh *const &boundary_mesh_pt)
 Delete boundary face elements and wipe the surface mesh. More...
 

Private Attributes

AnnularQuadMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to bulk mesh. More...
 
FourierDecomposedHelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
 
MeshHelmholtz_inner_boundary_mesh_pt
 on the inner boundary More...
 
TriangleMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
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 FourierDecomposedHelmholtzProblem< ELEMENT >

Problem class.

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

Constructor & Destructor Documentation

◆ FourierDecomposedHelmholtzProblem() [1/2]

Constructor.

Constructor for Fourier-decomposed Helmholtz problem.

438 {
439 
440 // Build annular mesh
441 
442 // # of elements in r-direction
443  unsigned n_r=10*ProblemParameters::El_multiplier;
444 
445  // # of elements in theta-direction
446  unsigned n_theta=10*ProblemParameters::El_multiplier;
447 
448  // Domain boundaries in theta-direction
449  double theta_min=-90.0;
450  double theta_max=90.0;
451 
452  // Domain boundaries in r-direction
453  double r_min=1.0;
454  double r_max=3.0;
455 
456  // Build and assign mesh
457  Bulk_mesh_pt =
458  new AnnularQuadMesh<ELEMENT>(n_r,n_theta,r_min,r_max,theta_min,theta_max);
459 
460  // Create mesh for DtN elements on outer boundary
464 
465  // Populate it with elements
467 
468  // Create flux elements on inner boundary
471 
472  // Add the several sub meshes to the problem
476 
477  // Build the Problem's global mesh from its various sub-meshes
479 
480  // Complete the build of all elements so they are fully functional
481  unsigned n_element = Bulk_mesh_pt->nelement();
482  for(unsigned i=0;i<n_element;i++)
483  {
484  // Upcast from GeneralsedElement to the present element
485  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(i));
486 
487  //Set the k_squared pointer
488  el_pt->k_squared_pt()=&ProblemParameters::K_squared;
489 
490  // Set pointer to Fourier wave number
491  el_pt->fourier_wavenumber_pt()=&ProblemParameters::N_fourier;
492  }
493 
494  // Setup equation numbering scheme
495  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
496 
497 } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: sphere_scattering.cc:56
void create_outer_bc_elements()
Create BC elements on outer boundary.
Definition: sphere_scattering.cc:658
Mesh * Helmholtz_inner_boundary_mesh_pt
on the inner boundary
Definition: sphere_scattering.cc:426
void create_flux_elements_on_inner_boundary()
Create flux elements on inner boundary.
Definition: sphere_scattering.cc:697
FourierDecomposedHelmholtzDtNMesh< ELEMENT > * Helmholtz_outer_boundary_mesh_pt
Definition: sphere_scattering.cc:422
AnnularQuadMesh< ELEMENT > * Bulk_mesh_pt
Pointer to bulk mesh.
Definition: sphere_scattering.cc:418
Definition: fourier_decomposed_helmholtz_bc_elements.h:351
Definition: mesh.h:67
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
double r_min
Definition: two_d_biharmonic.cc:237
double r_max
Definition: two_d_biharmonic.cc:238
unsigned El_multiplier
Multiplier for number of elements.
Definition: sphere_scattering.cc:363
unsigned Nterms_for_DtN
Number of terms in computation of DtN boundary condition.
Definition: sphere_scattering.cc:231
double K_squared
Square of the wavenumber.
Definition: sphere_scattering.cc:225
int N_fourier
Fourier wave number.
Definition: sphere_scattering.cc:228

References ProblemParameters::El_multiplier, i, ProblemParameters::K_squared, ProblemParameters::N_fourier, ProblemParameters::Nterms_for_DtN, BiharmonicTestFunctions2::r_max, and BiharmonicTestFunctions2::r_min.

◆ ~FourierDecomposedHelmholtzProblem() [1/2]

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

Destructor (empty)

388 {}

◆ FourierDecomposedHelmholtzProblem() [2/2]

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

Constructor.

◆ ~FourierDecomposedHelmholtzProblem() [2/2]

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

Destructor (empty)

333 {}

Member Function Documentation

◆ actions_after_adapt()

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

438 {
439 
440 
441  // Complete the build of all elements so they are fully functional
442 
443  // Loop over the Helmholtz bulk elements to set up element-specific
444  // things that cannot be handled by constructor: Pass pointer to
445  // wave number squared
446  unsigned n_element = Bulk_mesh_pt->nelement();
447  for(unsigned e=0;e<n_element;e++)
448  {
449  // Upcast from GeneralisedElement to Helmholtz bulk element
450  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
451 
452  //Set the k_squared pointer
453  el_pt->k_squared_pt() = &ProblemParameters::K_squared;
454 
455  // Set pointer to Fourier wave number
456  el_pt->fourier_wavenumber_pt()=&ProblemParameters::N_fourier;
457  }
458 
459  // Create prescribed-flux elements and BC elements
460  // from all elements that are adjacent to the boundaries and add them to
461  // Helmholtz_boundary_meshes
463  if (!CommandLineArgs::command_line_flag_has_been_set("--square_domain"))
464  {
466  }
467 
468  // Rebuild the Problem's global mesh from its various sub-meshes
470 
471 }// end of actions_after_adapt
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void rebuild_global_mesh()
Definition: problem.cc:1533
bool command_line_flag_has_been_set(const std::string &flag)
Definition: oomph_utilities.cc:501

References oomph::CommandLineArgs::command_line_flag_has_been_set(), e(), ProblemParameters::K_squared, and ProblemParameters::N_fourier.

◆ actions_after_newton_solve() [1/2]

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

Update the problem after solve (empty)

Reimplemented from oomph::Problem.

394 {}

◆ actions_after_newton_solve() [2/2]

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

Update the problem after solve (empty)

Reimplemented from oomph::Problem.

339 {}

◆ actions_before_adapt()

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

419 {
420  // Kill the flux elements and wipe the boundary meshs
421  if (!CommandLineArgs::command_line_flag_has_been_set("--square_domain"))
422  {
424  }
426 
427  // Rebuild the Problem's global mesh from its various sub-meshes
429 
430 }// end of actions_before_adapt
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete boundary face elements and wipe the surface mesh.
Definition: unstructured_sphere_scattering.cc:373

References oomph::CommandLineArgs::command_line_flag_has_been_set().

◆ actions_before_newton_convergence_check() [1/2]

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

402  {
403  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
404  }

◆ actions_before_newton_convergence_check() [2/2]

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

Recompute gamma integral before checking Newton residuals.

Reimplemented from oomph::Problem.

347  {
348  if (!CommandLineArgs::command_line_flag_has_been_set("--square_domain"))
349  {
350  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
351  }
352  }

References oomph::CommandLineArgs::command_line_flag_has_been_set().

◆ actions_before_newton_solve() [1/2]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

391 {}

◆ actions_before_newton_solve() [2/2]

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

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

336 {}

◆ check_gamma() [1/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::check_gamma ( DocInfo doc_info)

Check gamma computation.

Check gamma computation: \( \gamma = -du/dn \).

506 {
507 
508 
509  // Compute gamma stuff
510  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
511 
512  ofstream some_file;
513  char filename[100];
514 
515  sprintf(filename,"%s/gamma_test%i.dat",doc_info.directory().c_str(),
516  doc_info.number());
517  some_file.open(filename);
518 
519  //first loop over elements e
520  unsigned nel=Helmholtz_outer_boundary_mesh_pt->nelement();
521  for (unsigned e=0;e<nel;e++)
522  {
523  // Get a pointer to element
526  (Helmholtz_outer_boundary_mesh_pt->element_pt(e));
527 
528  //Set the value of n_intpt
529  const unsigned n_intpt =el_pt->integral_pt()->nweight();
530 
531  // Get gamma at all gauss points in element
533  Helmholtz_outer_boundary_mesh_pt->gamma_at_gauss_point(el_pt));
534 
535  //Loop over the integration points
536  for(unsigned ipt=0;ipt<n_intpt;ipt++)
537  {
538  //Allocate and initialise coordiante
539  Vector<double> x(el_pt->dim()+1,0.0);
540 
541  //Set the Vector to hold local coordinates
542  unsigned n=el_pt->dim();
543  Vector<double> s(n,0.0);
544  for(unsigned i=0;i<n;i++)
545  {
546  s[i]=el_pt->integral_pt()->knot(ipt,i);
547  }
548 
549  //Get the coordinates of the integration point
550  el_pt->interpolated_x(s,x);
551 
552  complex<double> flux;
554  some_file << atan2(x[0],x[1]) << " "
555  << gamma[ipt].real() << " "
556  << gamma[ipt].imag() << " "
557  << flux.real() << " "
558  << flux.imag() << " "
559  << std::endl;
560 
561  }// end of loop over integration points
562 
563  }// end of loop over elements
564 
565  some_file.close();
566 
567 }//end of output_gamma
AnnoyingScalar atan2(const AnnoyingScalar &y, const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:139
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
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
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
unsigned dim() const
Definition: elements.h:2611
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
Definition: fourier_decomposed_helmholtz_bc_elements.h:422
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
Definition: oomph-lib/src/generic/Vector.h:58
RealScalar s
Definition: level1_cplx_impl.h:130
string filename
Definition: MergeRestartFiles.py:39
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void exact_minus_dudr(const Vector< double > &x, std::complex< double > &flux)
Definition: sphere_scattering.cc:302
Mdouble gamma(Mdouble gamma_in)
This is the gamma function returns the true value for the half integer value.
Definition: ExtendedMath.cc:116
list x
Definition: plotDoE.py:28

References atan2(), oomph::FiniteElement::dim(), oomph::DocInfo::directory(), e(), ProblemParameters::exact_minus_dudr(), MergeRestartFiles::filename, ProblemParameters::flux(), mathsFunc::gamma(), i, oomph::FiniteElement::integral_pt(), oomph::FaceElement::interpolated_x(), oomph::Integral::knot(), n, oomph::DocInfo::number(), oomph::Integral::nweight(), s, and plotDoE::x.

◆ check_gamma() [2/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::check_gamma ( DocInfo doc_info)

Check gamma computation.

◆ create_flux_elements_on_inner_boundary() [1/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::create_flux_elements_on_inner_boundary
private

Create flux elements on inner boundary.

698 {
699  // Apply flux bc on inner boundary (boundary 3)
700  unsigned b=3;
701 
702 // Loop over the bulk elements adjacent to boundary b
703  unsigned n_element = Bulk_mesh_pt->nboundary_element(b);
704  for(unsigned e=0;e<n_element;e++)
705  {
706  // Get pointer to the bulk element that is adjacent to boundary b
707  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
708  Bulk_mesh_pt->boundary_element_pt(b,e));
709 
710  //Find the index of the face of element e along boundary b
711  int face_index = Bulk_mesh_pt->face_index_at_boundary(b,e);
712 
713  // Build the corresponding prescribed incoming-flux element
714  FourierDecomposedHelmholtzFluxElement<ELEMENT>* flux_element_pt = new
715  FourierDecomposedHelmholtzFluxElement<ELEMENT>(bulk_elem_pt,face_index);
716 
717  //Add the prescribed incoming-flux element to the surface mesh
719 
720  // Set the pointer to the prescribed flux function
721  flux_element_pt->flux_fct_pt() = &ProblemParameters::exact_minus_dudr;
722 
723  } //end of loop over bulk elements adjacent to boundary b
724 
725 } // end of create flux elements on inner boundary
Scalar * b
Definition: benchVecAdd.cpp:17
Definition: fourier_decomposed_helmholtz_flux_elements.h:56
FourierDecomposedHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
Definition: fourier_decomposed_helmholtz_flux_elements.h:92
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617

References b, e(), ProblemParameters::exact_minus_dudr(), and oomph::FourierDecomposedHelmholtzFluxElement< ELEMENT >::flux_fct_pt().

◆ create_flux_elements_on_inner_boundary() [2/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::create_flux_elements_on_inner_boundary ( )
private

Create flux elements on inner boundary.

◆ create_outer_bc_elements() [1/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::create_outer_bc_elements
private

Create BC elements on outer boundary.

659 {
660  // Outer boundary is boundary 1:
661  unsigned b=1;
662 
663  // Loop over the bulk elements adjacent to boundary b?
664  unsigned n_element = Bulk_mesh_pt->nboundary_element(b);
665  for(unsigned e=0;e<n_element;e++)
666  {
667  // Get pointer to the bulk element that is adjacent to boundary b
668  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
669  Bulk_mesh_pt->boundary_element_pt(b,e));
670 
671  //Find the index of the face of element e along boundary b
672  int face_index = Bulk_mesh_pt->face_index_at_boundary(b,e);
673 
674  // Build the corresponding DtN element
677  face_index);
678 
679  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
680  Helmholtz_outer_boundary_mesh_pt->add_element_pt(flux_element_pt);
681 
682  // Set pointer to the mesh that contains all the boundary condition
683  // elements on this boundary
684  flux_element_pt->
685  set_outer_boundary_mesh_pt(Helmholtz_outer_boundary_mesh_pt);
686  }
687 
688 } // end of create_outer_bc_elements

References b, and e().

◆ create_outer_bc_elements() [2/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< ELEMENT >::create_outer_bc_elements ( )
private

Create BC elements on outer boundary.

◆ delete_face_elements()

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

Delete boundary face elements and wipe the surface mesh.

374  {
375  // Loop over the surface elements
376  unsigned n_element = boundary_mesh_pt->nelement();
377  for(unsigned e=0;e<n_element;e++)
378  {
379  // Kill surface element
380  delete boundary_mesh_pt->element_pt(e);
381  }
382 
383  // Wipe the mesh
384  boundary_mesh_pt->flush_element_and_node_storage();
385 
386  }
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() [1/2]

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

575 {
576 
577  ofstream some_file;
578  char filename[100];
579 
580  // Number of plot points: npts x npts
581  unsigned npts=5;
582 
583  // Output solution
584  //-----------------
585  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
586  doc_info.number());
587  some_file.open(filename);
588  Bulk_mesh_pt->output(some_file,npts);
589  some_file.close();
590 
591 
592  // Output exact solution
593  //----------------------
594  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
595  doc_info.number());
596  some_file.open(filename);
597  Bulk_mesh_pt->output_fct(some_file,npts,ProblemParameters::get_exact_u);
598  some_file.close();
599 
600 
601  // Doc error and return of the square of the L2 error
602  //---------------------------------------------------
603  double error,norm;
604  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
605  doc_info.number());
606  some_file.open(filename);
607  Bulk_mesh_pt->compute_error(some_file,ProblemParameters::get_exact_u,
608  error,norm);
609  some_file.close();
610 
611  // Doc L2 error and norm of solution
612  cout << "\nNorm of error : " << sqrt(error) << std::endl;
613  cout << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
614 
615  // Check gamma computation
616  check_gamma(doc_info);
617 
618 
619  // Compute/output the radiated power
620  //----------------------------------
621  sprintf(filename,"%s/power%i.dat",doc_info.directory().c_str(),
622  doc_info.number());
623  some_file.open(filename);
624  sprintf(filename,"%s/total_power%i.dat",doc_info.directory().c_str(),
625  doc_info.number());
626 
627  // Accumulate contribution from elements
628  double power=0.0;
629  unsigned nn_element=Helmholtz_outer_boundary_mesh_pt->nelement();
630  for(unsigned e=0;e<nn_element;e++)
631  {
634  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
635  power += el_pt->global_power_contribution(some_file);
636  }
637  some_file.close();
638 
639  // Output total power
640  oomph_info << "Radiated power: "
643  << power << std::endl;
644  some_file.open(filename);
645  some_file << ProblemParameters::N_fourier << " "
647  << power << std::endl;
648  some_file.close();
649 
650 } // end of doc
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
void check_gamma(DocInfo &doc_info)
Check gamma computation.
Definition: sphere_scattering.cc:505
Definition: fourier_decomposed_helmholtz_bc_elements.h:62
double global_power_contribution()
Definition: fourier_decomposed_helmholtz_bc_elements.h:145
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector of size 2, containing real and imag parts.
Definition: sphere_scattering.cc:243
int error
Definition: calibrate.py:297
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::DocInfo::directory(), e(), calibrate::error, MergeRestartFiles::filename, ProblemParameters::get_exact_u(), oomph::FourierDecomposedHelmholtzBCElementBase< ELEMENT >::global_power_contribution(), ProblemParameters::K_squared, ProblemParameters::N_fourier, oomph::DocInfo::number(), oomph::oomph_info, and sqrt().

◆ doc_solution() [2/2]

template<class ELEMENT >
void FourierDecomposedHelmholtzProblem< 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 [1/2]

template<class ELEMENT >
AnnularQuadMesh<ELEMENT>* FourierDecomposedHelmholtzProblem< ELEMENT >::Bulk_mesh_pt
private

Pointer to bulk mesh.

◆ Bulk_mesh_pt [2/2]

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

Pointer to the "bulk" mesh.

◆ Helmholtz_inner_boundary_mesh_pt

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

on the inner boundary

Mesh of face elements that apply the prescribed flux on the inner boundary

◆ Helmholtz_outer_boundary_mesh_pt

template<class ELEMENT >
FourierDecomposedHelmholtzDtNMesh< ELEMENT > * FourierDecomposedHelmholtzProblem< ELEMENT >::Helmholtz_outer_boundary_mesh_pt
private

Pointer to mesh containing the DtN boundary condition elements

◆ Trace_file

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

Trace file.


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