LinearWaveProblem< ELEMENT, TIMESTEPPER > Class Template Reference

LinearWave problem in rectanglular domain. More...

+ Inheritance diagram for LinearWaveProblem< ELEMENT, TIMESTEPPER >:

Public Member Functions

 LinearWaveProblem (const unsigned &nx, const unsigned &ny, const bool &impulsive_start, LinearWaveEquations< 2 >::LinearWaveSourceFctPt source_fct_pt)
 Constructor for LinearWave problem. More...
 
 ~LinearWaveProblem ()
 Destructor (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after solve (empty) More...
 
void actions_before_implicit_timestep ()
 
void set_initial_condition ()
 Set initial condition (incl history values) More...
 
void doc_solution (DocInfo &doc_info)
 Doc the solution. More...
 
void unsteady_run ()
 Do unsteady run. More...
 
 LinearWaveProblem (const unsigned &nx, const unsigned &ny, LinearWaveEquations< 2 >::LinearWaveSourceFctPt source_fct_pt)
 Constructor for LinearWave problem. More...
 
 ~LinearWaveProblem ()
 Destructor (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after solve (empty) More...
 
void actions_before_implicit_timestep ()
 
void set_initial_condition ()
 
void doc_solution (DocInfo &doc_info)
 Doc the solution. More...
 
void unsteady_run ()
 Do unsteady run. More...
 
 LinearWaveProblem (const unsigned &nx, const unsigned &ny, const bool &impulsive_start, LinearWaveEquations< 2 >::LinearWaveSourceFctPt source_fct_pt)
 
 ~LinearWaveProblem ()
 Destructor (empty) More...
 
void actions_after_newton_solve ()
 Update the problem specs after solve (empty) More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_implicit_timestep ()
 Update the problem specs after solve (empty) More...
 
void actions_before_implicit_timestep ()
 
RefineableRectangularQuadMesh< ELEMENT > * mesh_pt ()
 Access function for the specific mesh. More...
 
void set_initial_condition ()
 Set initial condition (incl history values) More...
 
void doc_solution (DocInfo &doc_info)
 Doc the solution. More...
 
void unsteady_run ()
 Do unsteady run. 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
 
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 Member Functions

void create_flux_elements (const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
 

Private Attributes

ofstream Trace_file
 
bool Impulsive_start
 
RectangularQuadMesh< ELEMENT > * Bulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshSurface_mesh_pt
 Pointer to the "surface" mesh. 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_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 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 TIMESTEPPER>
class LinearWaveProblem< ELEMENT, TIMESTEPPER >

LinearWave problem in rectanglular domain.

LinearWave problem with flux boundary conditions in rectangle.

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

Constructor & Destructor Documentation

◆ LinearWaveProblem() [1/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::LinearWaveProblem ( const unsigned nx,
const unsigned ny,
const bool impulsive_start,
LinearWaveEquations< 2 >::LinearWaveSourceFctPt  source_fct_pt 
)

Constructor for LinearWave problem.

Constructor: pass number of elements in x and y directions, bool indicating impulsive or "smooth" start, and pointer to source function

214  :
215  Impulsive_start(impulsive_start)
216 {
217 
218  //Allocate the timestepper -- this constructs the time object as well
219  add_time_stepper_pt(new TIMESTEPPER());
220 
221  // Set up parameters for exact solution
222  //-------------------------------------
223 
224  // Steepness of tanh profile
226 
227  // Orientation of step wave
229 
230  // Set up mesh
231  //------------
232 
233  // # of elements in x-direction
234  unsigned Nx=nx;
235 
236  // # of elements in y-direction
237  unsigned Ny=ny;
238 
239  // Domain length in x-direction
240  double Lx=1.0;
241 
242  // Domain length in y-direction
243  double Ly=2.0;
244 
245  // Build and assign mesh
246  Problem::mesh_pt()=new RectangularQuadMesh<ELEMENT>(
248 
249  // Set the boundary conditions for this problem: All nodes are
250  // free by default -- just pin the ones that have Dirichlet conditions
251  // here.
252  unsigned num_bound = mesh_pt()->nboundary();
253  for(unsigned ibound=0;ibound<num_bound;ibound++)
254  {
255  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
256  for (unsigned inod=0;inod<num_nod;inod++)
257  {
258  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
259  }
260  } //end of boundary conditions
261 
262  // Complete build of elements
263  // --------------------------
264 
265  // Loop over the elements to set up element-specific
266  // things that cannot be handled by constructor
267 unsigned n_element = mesh_pt()->nelement();
268 for(unsigned i=0;i<n_element;i++)
269  {
270  // Upcast from GeneralisedElement to the present element
271  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
272 
273  //Set the source function pointer
274  el_pt->source_fct_pt() = source_fct_pt;
275  }
276 
277  // Setup equation numbering scheme
278  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
279 
280 } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
RefineableRectangularQuadMesh< ELEMENT > * mesh_pt()
Access function for the specific mesh.
Definition: two_d_linear_wave_adapt.cc:190
bool Impulsive_start
Definition: two_d_linear_wave.cc:202
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
TimeStepper *& time_stepper_pt()
Definition: problem.h:1524
Definition: rectangular_quadmesh.template.h:59
double Pi
Definition: two_d_biharmonic.cc:235
unsigned Nx
Number of elements in each direction (used by SimpleCubicMesh)
Definition: structured_cubic_point_source.cc:114
unsigned Ny
Definition: structured_cubic_point_source.cc:115
double Ly
Length of domain in y direction.
Definition: periodic_load.cc:58
double Lx
Length of domain in x direction.
Definition: periodic_load.cc:55
const unsigned nx
Definition: ConstraintElementsUnitTest.cpp:30
const unsigned ny
Definition: ConstraintElementsUnitTest.cpp:31
double Alpha
Parameter for steepness of step.
Definition: two_d_linear_wave.cc:52
double Phi
Orientation of step wave.
Definition: two_d_linear_wave.cc:55
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: gen_axisym_advection_diffusion_elements.h:229

References oomph::Problem::add_time_stepper_pt(), TanhSolnForLinearWave::Alpha, oomph::Problem::assign_eqn_numbers(), i, Global_Parameters::Lx, Global_Parameters::Ly, LinearWaveProblem< ELEMENT, TIMESTEPPER >::mesh_pt(), Mesh_Parameters::nx, GlobalParameters::Nx, Mesh_Parameters::ny, GlobalParameters::Ny, TanhSolnForLinearWave::Phi, BiharmonicTestFunctions2::Pi, oomph::source_fct_pt(), and oomph::Problem::time_stepper_pt().

◆ ~LinearWaveProblem() [1/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::~LinearWaveProblem ( )
inline

Destructor (empty)

121 {}

◆ LinearWaveProblem() [2/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::LinearWaveProblem ( const unsigned nx,
const unsigned ny,
LinearWaveEquations< 2 >::LinearWaveSourceFctPt  source_fct_pt 
)

Constructor for LinearWave problem.

Constructor: pass number of elements in x and y directions and pointer to source function.

260 {
261 
262  //Create the timestepper -- this constructs the time object as well
263  add_time_stepper_pt(new TIMESTEPPER());
264 
265 
266  // Set up parameters for exact solution
267  //-------------------------------------
268 
269  // Steepness of tanh profile
271 
272  // Orientation of step wave
274 
275 
276 
277  // Set up mesh
278  //------------
279 
280  // # of elements in x-direction
281  unsigned Nx=nx;
282 
283  // # of elements in y-direction
284  unsigned Ny=ny;
285 
286  // Domain length in x-direction
287  double Lx=1.0;
288 
289  // Domain length in y-direction
290  double Ly=2.0;
291 
292  // Build bulk mesh
294 
295  // Create "surface mesh" that will contain only the prescribed-flux
296  // elements. The constructor just creates the mesh without
297  // giving it any elements, nodes, etc.
298  Surface_mesh_pt = new Mesh;
299 
300  // Create prescribed-flux elements from all elements that are
301  // adjacent to boundary 1, but add them to a separate mesh.
302  // Note that this is exactly the same function as used in the
303  // single mesh version of the problem, we merely pass different Mesh pointers.
305 
306  // Add the two sub meshes to the problem
309 
310  // Combine all submeshes into a single Mesh
312 
313 
314 
315  // Set the boundary conditions for this problem: All nodes are
316  // free by default -- just pin the ones that have Dirichlet conditions
317  // here.
318  unsigned num_bound = Bulk_mesh_pt->nboundary();
319  for(unsigned ibound=0;ibound<num_bound;ibound++)
320  {
321  if (ibound!=1) // Not on the flux boundary
322  {
323  unsigned num_nod= Bulk_mesh_pt->nboundary_node(ibound);
324  for (unsigned inod=0;inod<num_nod;inod++)
325  {
326  Bulk_mesh_pt->boundary_node_pt(ibound,inod)->pin(0);
327  }
328  }
329  }
330 
331  // Complete build of bulk elements
332  // -------------------------------
333 
334  // Loop over the elements to set up element-specific
335  // things that cannot be handled by constructor
336  unsigned n_element = Bulk_mesh_pt->nelement();
337  for(unsigned i=0;i<n_element;i++)
338  {
339  // Upcast from FiniteElement to the present element
340  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(i));
341 
342  //Set the source function pointer
343  el_pt->source_fct_pt() = source_fct_pt;
344  }
345 
346 
347  // Complete build of flux elements
348  // -------------------------------
349 
350  // Loop over surface elements to set up element-specific
351  // things that cannot be handled by constructor
352  n_element = Surface_mesh_pt->nelement();
353  for (unsigned e=0;e<n_element;e++)
354  {
355  // Upcast from GeneralisedElement to LinearWave flux element
357  dynamic_cast< LinearWaveFluxElement<ELEMENT>*>(
359 
360  el_pt->flux_fct_pt() =
362  }
363 
364  //Setup equation numbering scheme
365  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
366 
367 } // end of constructor
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Mesh * Surface_mesh_pt
Pointer to the "surface" mesh.
Definition: two_d_linear_wave_flux.cc:244
void create_flux_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
Definition: two_d_linear_wave_flux.cc:423
RectangularQuadMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: two_d_linear_wave_flux.cc:241
Definition: linear_wave_flux_elements.h:53
LinearWavePrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
Definition: linear_wave_flux_elements.h:84
Definition: mesh.h:67
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
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
void prescribed_flux_on_fixed_y_boundary(const double &time, const Vector< double > &x, double &flux)
Prescribed flux on a fixed y max boundary.
Definition: two_d_linear_wave_flux.cc:108

References TanhSolnForLinearWave::Alpha, e(), oomph::LinearWaveFluxElement< ELEMENT >::flux_fct_pt(), i, Global_Parameters::Lx, Global_Parameters::Ly, Mesh_Parameters::nx, GlobalParameters::Nx, Mesh_Parameters::ny, GlobalParameters::Ny, TanhSolnForLinearWave::Phi, BiharmonicTestFunctions2::Pi, TanhSolnForLinearWave::prescribed_flux_on_fixed_y_boundary(), and oomph::source_fct_pt().

◆ ~LinearWaveProblem() [2/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::~LinearWaveProblem ( )
inline

Destructor (empty)

151 {}

◆ LinearWaveProblem() [3/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::LinearWaveProblem ( const unsigned nx,
const unsigned ny,
const bool impulsive_start,
LinearWaveEquations< 2 >::LinearWaveSourceFctPt  source_fct_pt 
)

Constructor: pass number of elements in x and y directions bool indicating impulsive or "smooth" start and pointer to source function

◆ ~LinearWaveProblem() [3/3]

template<class ELEMENT , class TIMESTEPPER >
LinearWaveProblem< ELEMENT, TIMESTEPPER >::~LinearWaveProblem ( )
inline

Destructor (empty)

121 {}

Member Function Documentation

◆ actions_after_implicit_timestep() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

130 {}

◆ actions_after_implicit_timestep() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

160 {}

◆ actions_after_implicit_timestep() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_implicit_timestep ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

130 {}

◆ actions_after_newton_solve() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_newton_solve ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

124 {}

◆ actions_after_newton_solve() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_newton_solve ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

154 {}

◆ actions_after_newton_solve() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_after_newton_solve ( )
inlinevirtual

Update the problem specs after solve (empty)

Reimplemented from oomph::Problem.

124 {}

◆ actions_before_implicit_timestep() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_implicit_timestep ( )
inlinevirtual

Update the problem specs before next timestep: Set time-dependent Dirchlet boundary from exact solution.

Reimplemented from oomph::Problem.

135  {
137  initial_value_fct(1);
139  initial_veloc_fct(1);
141  initial_accel_fct(1);
142 
143  // Assign values for analytical value, veloc and accel:
144  initial_value_fct[0]=&TanhSolnForLinearWave::exact_u;
145  initial_veloc_fct[0]=&TanhSolnForLinearWave::exact_dudt;
146  initial_accel_fct[0]=&TanhSolnForLinearWave::exact_d2udt2;
147 
148  // Loop over boundaries
149  unsigned num_bound=mesh_pt()->nboundary();
150  for (unsigned ibound=0;ibound<num_bound;ibound++)
151  {
152  // Loop over boundary nodes
153  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
154  for (unsigned inod=0;inod<num_nod;inod++)
155  {
156  // Set the boundary condition from the exact solution
157  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
158 
159  bool use_direct_assignment=false;
160  if (use_direct_assignment)
161  {
162  // Set nodal coordinates for evaluation of BC:
163  Vector<double> x(2);
164  x[0]=nod_pt->x(0);
165  x[1]=nod_pt->x(1);
166 
167  // Set exact solution at current time
168  nod_pt->set_value(0,
170  }
171  else
172  {
173  // Get timestepper
174  TIMESTEPPER* timestepper_pt=dynamic_cast<TIMESTEPPER*>
175  (time_stepper_pt());
176 
177  // Assign the history values
178  timestepper_pt->assign_initial_data_values(nod_pt,
179  initial_value_fct,
180  initial_veloc_fct,
181  initial_accel_fct);
182  }
183  }
184  }
185  } // end of actions before timestep
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Time *& time_pt()
Return a pointer to the global time object.
Definition: problem.h:1504
double & time()
Return the current value of continuous time.
Definition: problem.cc:11531
Definition: oomph-lib/src/generic/Vector.h:58
double exact_d2udt2(const double &time, const Vector< double > &x)
2nd time-deriv of exact solution
Definition: two_d_linear_wave.cc:73
double exact_dudt(const double &time, const Vector< double > &x)
1st time-deriv of exact solution
Definition: two_d_linear_wave.cc:65
double exact_u(const double &time, const Vector< double > &x)
Exact solution.
Definition: two_d_linear_wave.cc:58
list x
Definition: plotDoE.py:28

References TanhSolnForLinearWave::exact_d2udt2(), TanhSolnForLinearWave::exact_dudt(), TanhSolnForLinearWave::exact_u(), oomph::Data::set_value(), plotDoE::x, and oomph::Node::x().

◆ actions_before_implicit_timestep() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_implicit_timestep ( )
inlinevirtual

Update the problem specs before next timestep: Set Dirchlet boundary conditions from exact solution.

Reimplemented from oomph::Problem.

165  {
167  initial_value_fct(1);
169  initial_veloc_fct(1);
171  initial_accel_fct(1);
172 
173  // Assign values for analytical value, veloc and accel:
174  initial_value_fct[0]=&TanhSolnForLinearWave::exact_u;
175  initial_veloc_fct[0]=&TanhSolnForLinearWave::exact_dudt;
176  initial_accel_fct[0]=&TanhSolnForLinearWave::exact_d2udt2;
177 
178  // Loop over boundaries
179  unsigned num_bound=Bulk_mesh_pt->nboundary();
180  for (unsigned ibound=0;ibound<num_bound;ibound++)
181  {
182  // Ignore boundary 1 where flux condition is applied
183  if (ibound!=1)
184  {
185  // Loop over boundary nodes
186  unsigned num_nod=Bulk_mesh_pt->nboundary_node(ibound);
187  for (unsigned inod=0;inod<num_nod;inod++)
188  {
189  // Set the boundary condition from the exact solution
190  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(ibound,inod);
191 
192  bool use_direct_assignment=false;
193  if (use_direct_assignment)
194  {
195  // Set nodal coordinates for evaluation of BC:
196  Vector<double> x(2);
197  x[0]=nod_pt->x(0);
198  x[1]=nod_pt->x(1);
199 
200  // Set exact solution at current time
201  nod_pt->
202  set_value(0,
204  }
205  else
206  {
207  // Get timestepper
208  TIMESTEPPER* timestepper_pt=dynamic_cast<TIMESTEPPER*>
209  (time_stepper_pt());
210 
211  // Assign the history values
212  timestepper_pt->assign_initial_data_values(nod_pt,
213  initial_value_fct,
214  initial_veloc_fct,
215  initial_accel_fct);
216  }
217  }
218  }
219  }
220  } // end of actions before timestep

References TanhSolnForLinearWave::exact_d2udt2(), TanhSolnForLinearWave::exact_dudt(), TanhSolnForLinearWave::exact_u(), plotDoE::x, and oomph::Node::x().

◆ actions_before_implicit_timestep() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_implicit_timestep ( )
inlinevirtual

Update the problem specs before next timestep: Set time-dependent Dirchlet boundary from exact solution.

Reimplemented from oomph::Problem.

135  {
137  initial_value_fct(1);
139  initial_veloc_fct(1);
141  initial_accel_fct(1);
142 
143  // Assign values for analytical value, veloc and accel:
144  initial_value_fct[0]=&TanhSolnForLinearWave::exact_u;
145  initial_veloc_fct[0]=&TanhSolnForLinearWave::exact_dudt;
146  initial_accel_fct[0]=&TanhSolnForLinearWave::exact_d2udt2;
147 
148  // Loop over boundaries
149  unsigned num_bound=mesh_pt()->nboundary();
150  for (unsigned ibound=0;ibound<num_bound;ibound++)
151  {
152  // Loop over boundary nodes
153  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
154  for (unsigned inod=0;inod<num_nod;inod++)
155  {
156  // Set the boundary condition from the exact solution
157  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
158 
159  bool use_direct_assignment=false;
160  if (use_direct_assignment)
161  {
162  // Set nodal coordinates for evaluation of BC:
163  Vector<double> x(2);
164  x[0]=nod_pt->x(0);
165  x[1]=nod_pt->x(1);
166 
167  // Set exact solution at current time
168  nod_pt->
169  set_value(0,
171  }
172  else
173  {
174  // Get timestepper
175  TIMESTEPPER* timestepper_pt=dynamic_cast<TIMESTEPPER*>
176  (time_stepper_pt());
177 
178  // Assign the history values
179  timestepper_pt->assign_initial_data_values(nod_pt,
180  initial_value_fct,
181  initial_veloc_fct,
182  initial_accel_fct);
183  }
184  }
185  }
186  } // end of actions before timestep

References TanhSolnForLinearWave::exact_d2udt2(), TanhSolnForLinearWave::exact_dudt(), TanhSolnForLinearWave::exact_u(), plotDoE::x, and oomph::Node::x().

◆ actions_before_newton_solve() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

127 {}

◆ actions_before_newton_solve() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

157 {}

◆ actions_before_newton_solve() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

127 {}

◆ create_flux_elements()

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::create_flux_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)
private

Create LinearWave flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the Mesh object pointed to by surface_mesh_pt

Create LinearWave 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 surface_mesh_pt.

425 {
426  // How many bulk elements are adjacent to boundary b?
427  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
428 
429  // Loop over the bulk elements adjacent to boundary b?
430  for(unsigned e=0;e<n_element;e++)
431  {
432  // Get pointer to the bulk element that is adjacent to boundary b
433  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
434  bulk_mesh_pt->boundary_element_pt(b,e));
435 
436  //Find the index of the face of element e along boundary b
437  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
438 
439  // Build the corresponding prescribed-flux element
440  LinearWaveFluxElement<ELEMENT>* flux_element_pt = new
441  LinearWaveFluxElement<ELEMENT>(bulk_elem_pt,face_index);
442 
443  //Add the prescribed-flux element to the surface mesh
444  surface_mesh_pt->add_element_pt(flux_element_pt);
445 
446  } //end of loop over bulk elements adjacent to boundary b
447 
448 } // end of create flux elements
Scalar * b
Definition: benchVecAdd.cpp:17
int face_index_at_boundary(const unsigned &b, const unsigned &e) const
Definition: mesh.h:896
unsigned nboundary_element(const unsigned &b) const
Return number of finite elements that are adjacent to boundary b.
Definition: mesh.h:878
FiniteElement * boundary_element_pt(const unsigned &b, const unsigned &e) const
Return pointer to e-th finite element on boundary b.
Definition: mesh.h:840
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617

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/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::doc_solution ( DocInfo doc_info)

Doc the solution.

405 {
406 
407  ofstream some_file;
408  char filename[100];
409 
410  // Number of plot points
411  unsigned npts;
412  npts=5;
413 
414  cout << std::endl;
415  cout << "=================================================" << std::endl;
416  cout << "Docing solution for t=" << time_pt()->time() << std::endl;
417  cout << "=================================================" << std::endl;
418 
419  // Output solution
420  //-----------------
421  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
422  doc_info.number());
423  some_file.open(filename);
424  mesh_pt()->output(some_file,npts);
425  some_file << "TEXT X=2.5,Y=93.6,F=HELV,HU=POINT,C=BLUE,H=26,T=\"time = "
426  << time_pt()->time() << "\"";
427  some_file << "GEOMETRY X=2.5,Y=98,T=LINE,C=BLUE,LT=0.4" << std::endl;
428  some_file << "1" << std::endl;
429  some_file << "2" << std::endl;
430  some_file << " 0 0" << std::endl;
431  some_file << time_pt()->time()*20.0 << " 0" << std::endl;
432  some_file.close();
433 
434  // Output exact solution
435  //----------------------
436  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
437  doc_info.number());
438  oomph_info << " FILENAME: " << filename << std::endl;
439  some_file.open(filename);
440  mesh_pt()->output_fct(some_file,npts,time_pt()->time(),
442  some_file.close();
443 
444  // Doc error
445  //----------
446  double error,norm;
447  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
448  doc_info.number());
449  some_file.open(filename);
450  mesh_pt()->compute_error(some_file,
452  time_pt()->time(),
453  error,norm);
454  some_file.close();
455  cout << "error: " << error << std::endl;
456  cout << "norm : " << norm << std::endl << std::endl;
457 
458  // Write trace file
459  Trace_file << time_pt()->time() << " " << time_pt()->dt()
460  << " " << mesh_pt()->nelement() << " "
461  << error << " " << norm << std::endl;
462 
463 } // end of doc solution
ofstream Trace_file
Definition: two_d_linear_wave.cc:199
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 & dt(const unsigned &t=0)
Definition: timesteppers.h:136
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
string filename
Definition: MergeRestartFiles.py:39
void get_exact_u(const double &time, const Vector< double > &x, Vector< double > &u)
Exact solution as a vector.
Definition: two_d_linear_wave.cc:82
int error
Definition: calibrate.py:297
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::DocInfo::directory(), calibrate::error, MergeRestartFiles::filename, TanhSolnForLinearWave::get_exact_u(), oomph::DocInfo::number(), oomph::oomph_info, and oomph::Problem_Parameter::Trace_file.

◆ doc_solution() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::doc_solution ( DocInfo doc_info)

Doc the solution.

◆ doc_solution() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::doc_solution ( DocInfo doc_info)

Doc the solution.

◆ mesh_pt()

template<class ELEMENT , class TIMESTEPPER >
RefineableRectangularQuadMesh<ELEMENT>* LinearWaveProblem< ELEMENT, TIMESTEPPER >::mesh_pt ( )
inline

Access function for the specific mesh.

191  {
192  return dynamic_cast<RefineableRectangularQuadMesh<ELEMENT>*>(
193  Problem::mesh_pt());
194  }
Definition: rectangular_quadmesh.template.h:326

Referenced by LinearWaveProblem< ELEMENT, TIMESTEPPER >::LinearWaveProblem().

◆ set_initial_condition() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::set_initial_condition
virtual

Set initial condition (incl history values)

Set initial condition.

Reimplemented from oomph::Problem.

290 {
291 
292  // Get timestepper
293  TIMESTEPPER* timestepper_pt=dynamic_cast<TIMESTEPPER*>(time_stepper_pt());
294 
295 
296  // Impulsive start
297  //----------------
298  if (Impulsive_start)
299  {
300  // Loop over the nodes to set initial conditions everywhere
301  unsigned num_nod=mesh_pt()->nnode();
302  for (unsigned jnod=0;jnod<num_nod;jnod++)
303  {
304  // Pointer to node
305  Node* nod_pt=mesh_pt()->node_pt(jnod);
306 
307  // Get nodal coordinates
308  Vector<double> x(2);
309  x[0]=nod_pt->x(0);
310  x[1]=nod_pt->x(1);
311 
312  // Assign initial value from exact solution
314 
315  // Set history values so that they are consistent with an impulsive
316  // start from this value
317  timestepper_pt->assign_initial_values_impulsive(nod_pt);
318  }
319  } // end impulsive start
320 
321  // "Smooth" start from analytical time history
322  //--------------------------------------------
323  else
324  {
325 
326  // Vector of function pointers to functions that specify the
327  // value, and the first and second time-derivatives of the
328  // function used as the initial condition
330  initial_value_fct(1);
332  initial_veloc_fct(1);
334  initial_accel_fct(1);
335 
336  // Assign values for analytical value, veloc and accel:
337  initial_value_fct[0]=&TanhSolnForLinearWave::exact_u;
338  initial_veloc_fct[0]=&TanhSolnForLinearWave::exact_dudt;
339  initial_accel_fct[0]=&TanhSolnForLinearWave::exact_d2udt2;
340 
341  // Assign Newmark history values so that Newmark approximations
342  // for velocity and accel are correct at initial time:
343 
344  // Loop over the nodes to set initial conditions everywhere
345  unsigned num_nod=mesh_pt()->nnode();
346  for (unsigned jnod=0;jnod<num_nod;jnod++)
347  {
348  // Pointer to node
349  Node* nod_pt=mesh_pt()->node_pt(jnod);
350 
351  // Assign the history values
352  timestepper_pt->assign_initial_data_values(nod_pt,
353  initial_value_fct,
354  initial_veloc_fct,
355  initial_accel_fct);
356  } // end of smooth start
357 
358 
359  // Paranoia: Check that the initial values were assigned correctly
360  double err_max=0.0;
361  for (unsigned jnod=0;jnod<num_nod;jnod++)
362  {
363  // Pointer to node
364  Node* nod_pt=mesh_pt()->node_pt(jnod);
365 
366  // Get nodal coordinates
367  Vector<double> x(2);
368  x[0]=nod_pt->x(0);
369  x[1]=nod_pt->x(1);
370 
371  // Get exact value and first and second time-derivatives
372  double u_exact=
374  double dudt_exact=
376  double d2udt2_exact=
378 
379  // Get Newmark approximations for zero-th, first and second
380  // time-derivatives of the nodal values.
381  double u_fe=timestepper_pt->time_derivative(0,nod_pt,0);
382  double dudt_fe=timestepper_pt->time_derivative(1,nod_pt,0);
383  double d2udt2_fe=timestepper_pt->time_derivative(2,nod_pt,0);
384 
385  // Error
386  double error=sqrt(pow(u_exact-u_fe,2)+
387  pow(dudt_exact-dudt_fe,2)+
388  pow(d2udt2_exact-d2udt2_fe,2));
389  if (error>err_max) err_max=error;
390  }
391  cout << "Max. error in assignment of initial condition "
392  << err_max << std::endl;
393  }
394 
395 
396 } // end of set initial condition
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625

References calibrate::error, TanhSolnForLinearWave::exact_d2udt2(), TanhSolnForLinearWave::exact_dudt(), TanhSolnForLinearWave::exact_u(), Eigen::ArrayBase< Derived >::pow(), oomph::Data::set_value(), sqrt(), plotDoE::x, and oomph::Node::x().

◆ set_initial_condition() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::set_initial_condition ( )
virtual

Set initial condition (incl history values) according to specified function.

Reimplemented from oomph::Problem.

◆ set_initial_condition() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::set_initial_condition ( )
virtual

Set initial condition (incl history values)

Reimplemented from oomph::Problem.

◆ unsteady_run() [1/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::unsteady_run

Do unsteady run.

Perform run up to specified time.

Unsteady run.

473 {
474 
475  // Setup labels for output
476  DocInfo doc_info;
477 
478  // Output directory
479  if (Impulsive_start)
480  {
481  doc_info.set_directory("RESLT_impulsive");
482  }
483  else
484  {
485  doc_info.set_directory("RESLT_smooth");
486  }
487 
488  // Open trace file
489  char filename[100];
490  sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
491  Trace_file.open(filename);
492 
493  // Initialise time
494  double time0=0.0;
495  time_pt()->time()=time0;
496 
497  // Set initial timestep
498  double dt=0.005;
499  time_pt()->initialise_dt(dt);
500 
501  // Set IC
503 
504  //Output initial condition
505  doc_solution(doc_info);
506 
507  //Increment counter for solutions
508  doc_info.number()++;
509 
510  // Maximum time
511  double t_max=4.0;
512 
513  // Number of steps
514  unsigned nstep=unsigned(t_max/dt);
515 
516  // If validation run only do 2 timesteps
517  if (CommandLineArgs::Argc>1)
518  {
519  nstep=2;
520  cout << "Validation run -- only doing two timesteps." << std::endl;
521  }
522 
523  // Timestepping loop
524  for (unsigned istep=0;istep<nstep;istep++)
525  {
526  //Take fixed timestep without spatial adaptivity
528 
529  //Output solution
530  doc_solution(doc_info);
531 
532  //Increment counter for solutions
533  doc_info.number()++;
534  }
535 
536  // Close trace file
537  Trace_file.close();
538 
539 } // end of unsteady run
void doc_solution(DocInfo &doc_info)
Doc the solution.
Definition: two_d_linear_wave.cc:404
void set_initial_condition()
Set initial condition (incl history values)
Definition: two_d_linear_wave.cc:289
Definition: oomph_utilities.h:499
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
void unsteady_newton_solve(const double &dt)
Definition: problem.cc:10953
void initialise_dt(const double &dt_)
Set all timesteps to the same value, dt.
Definition: timesteppers.h:99
int Argc
Number of arguments + 1.
Definition: oomph_utilities.cc:407

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

◆ unsteady_run() [2/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::unsteady_run ( )

Do unsteady run.

◆ unsteady_run() [3/3]

template<class ELEMENT , class TIMESTEPPER >
void LinearWaveProblem< ELEMENT, TIMESTEPPER >::unsteady_run ( )

Do unsteady run.

Member Data Documentation

◆ Bulk_mesh_pt

template<class ELEMENT , class TIMESTEPPER >
RectangularQuadMesh<ELEMENT>* LinearWaveProblem< ELEMENT, TIMESTEPPER >::Bulk_mesh_pt
private

Pointer to the "bulk" mesh.

◆ Impulsive_start

template<class ELEMENT , class TIMESTEPPER >
bool LinearWaveProblem< ELEMENT, TIMESTEPPER >::Impulsive_start
private

◆ Surface_mesh_pt

template<class ELEMENT , class TIMESTEPPER >
Mesh* LinearWaveProblem< ELEMENT, TIMESTEPPER >::Surface_mesh_pt
private

Pointer to the "surface" mesh.

◆ Trace_file

template<class ELEMENT , class TIMESTEPPER >
ofstream LinearWaveProblem< ELEMENT, TIMESTEPPER >::Trace_file
private

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