TorusProblem< ELEMENT > Class Template Reference

Solve the Axisymmetric Navier Stokes equations in a torus. More...

+ Inheritance diagram for TorusProblem< ELEMENT >:

Public Member Functions

 TorusProblem (const unsigned &max_refinement_level, const double &min_error_target, const double &max_error_target)
 
void set_initial_condition ()
 Set the initial conditions: all nodes have zero velocity. More...
 
void set_boundary_conditions (const double &time)
 Set boundary conditions on the walls. More...
 
void solve_system (const double &dt, const unsigned &nstep)
 Function that is used to run the parameter study. More...
 
RefineableFullCircleMesh< ELEMENT > * mesh_pt ()
 Return a pointer to the specific mesh used. More...
 
void actions_before_implicit_timestep ()
 Update the problem specs before next timestep: More...
 
void actions_after_adapt ()
 
- 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 ()
 
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 ()
 

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_solve ()
 
virtual void actions_after_newton_solve ()
 
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 ()
 
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 TorusProblem< ELEMENT >

Solve the Axisymmetric Navier Stokes equations in a torus.

Constructor & Destructor Documentation

◆ TorusProblem()

template<class ELEMENT >
TorusProblem< ELEMENT >::TorusProblem ( const unsigned max_refinement_level,
const double min_error_target,
const double max_error_target 
)

Constructor taking the maximum refinement level and the minimum and maximum error targets.

Constructor: specify the maximum refinement level, the minimum and maximum error targets

185 {
186 
187  using namespace Global_Physical_Variables;
188 
189  //Create a timestepper
191 
192  //Create the domain for the mesh, which consists of a circle of
193  //radius Radius and centred at (1/Delta, 0)
194  GeomObject* area_pt = new FilledCircle(1.0/Delta,0.0,Radius);
195 
196  //Define pi
197  const double pi = MathematicalConstants::Pi;
198 
199  //Set the positions of the angles that divide the outer ring
200  //These must be in the range -pi,pi, ordered from smallest to
201  //largest
202  Vector<double> theta_positions(4);
203  theta_positions[0] = -0.75*pi;
204  theta_positions[1] = -0.25*pi;
205  theta_positions[2] = 0.25*pi;
206  theta_positions[3] = 0.75*pi;
207 
208  //Define the radial fraction of the central box (always halfway
209  //along the radius)
210  Vector<double> radial_frac(4,0.5);
211 
212  //Now create the mesh
213  Problem::mesh_pt() = new RefineableFullCircleMesh<ELEMENT>(
214  area_pt,theta_positions,radial_frac,time_stepper_pt());
215 
216  // Set error estimator
218  mesh_pt()->spatial_error_estimator_pt()=error_estimator_pt;
219 
220  // Maximum number of refinements (increase this if you want a finer mesh)
221  mesh_pt()->max_refinement_level() = max_refinement_level;
222 
223  // Error targets for adaptive refinement
224  mesh_pt()->max_permitted_error() = max_error_target;
225  mesh_pt()->min_permitted_error() = min_error_target;
226 
227  //Build iterative linear solver
228  GMRES<CRDoubleMatrix>* iterative_linear_solver_pt = new
230 
231  // Set maximum number of iterations
232  iterative_linear_solver_pt->max_iter() = 100;
233 
234  // Set tolerance
235  iterative_linear_solver_pt->tolerance() = 1.0e-8;
236 
239  //Set the mesh
240  prec_pt->set_navier_stokes_mesh(this->mesh_pt());
241 
242  //Set the preconditioner
243  iterative_linear_solver_pt->preconditioner_pt() = prec_pt;
244 
245  //Set the linear solver
246  this->linear_solver_pt() = iterative_linear_solver_pt;
247 
248  //Precondition the pressure block with AMG if we have Hypre
249 #ifdef OOMPH_HAS_HYPRE
250 //Trap because HYPRE can't handle the case when OOMPH_HAS_MPI, but we
251 //run in serial
252 #ifndef OOMPH_HAS_MPI
253 
254  //Set up the internal preconditioners
255  Preconditioner* P_matrix_preconditioner_pt = new HyprePreconditioner;
256 
257  // Set parameters for use as preconditioner on Poisson-type problem
259  static_cast<HyprePreconditioner*>(P_matrix_preconditioner_pt));
260 
261  // Use Hypre for the Schur complement block
262  prec_pt->set_p_preconditioner(P_matrix_preconditioner_pt);
263 
264  // Shut up!
265  static_cast<HyprePreconditioner*>(P_matrix_preconditioner_pt)->
266  disable_doc_time();
267 #endif
268 #endif
269 
270 //Set Block diagonal preconditioner for the momentum block
271 BlockDiagonalPreconditioner<CRDoubleMatrix>* F_matrix_preconditioner_pt =
273 
274 // Set mesh
275 F_matrix_preconditioner_pt->add_mesh(this->mesh_pt());
276 
277 //Use AMG on the momentum blocks if we have Hypre
278 #ifdef OOMPH_HAS_HYPRE
279 //Trap because HYPRE can't handle the case when
280 //OOMPH_HAS_MPI, but we run in serial
281 #ifndef OOMPH_HAS_MPI
283  (F_matrix_preconditioner_pt)->set_subsidiary_preconditioner_function
284  (Hypre_Subsidiary_Preconditioner_Helper::set_hypre_preconditioner);
285 #endif
286 #endif
287 
288  // Use Hypre for momentum block
289  prec_pt->set_f_preconditioner(F_matrix_preconditioner_pt);
290 
291  //Loop over all the (fluid) elements
292  unsigned n_element = mesh_pt()->nelement();
293  for(unsigned e=0;e<n_element;e++)
294  {
295  //Cast to the particular element type, this is necessary because
296  //the base elements don't have the member functions that we're about
297  //to call!
298  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(e));
299 
300  //There is no need for ALE
301  el_pt->disable_ALE();
302 
303  //Set the Reynolds number for each element
304  //(yes we could have different Reynolds number in each element!!)
305  el_pt->re_pt() = &Re;
306  //Set the product of Reynolds and Strouhal numbers
307  el_pt->re_st_pt() = &Re;
308  }
309 
310  //Let this problem be conventional form by setting gamma to zero
311  ELEMENT::Gamma[0] = 0.0; //r-momentum
312  ELEMENT::Gamma[1] = 0.0; //z-momentum
313 
314  //Set the boundary conditions (no slip on the torus walls)
315  //Loop over the nodes on the (only) mesh boundary
316  unsigned n_boundary_node = mesh_pt()->nboundary_node(0);
317  for(unsigned n=0;n<n_boundary_node;n++)
318  {
319  //Pin all three velocity components on the wall
320  for(unsigned i=0;i<3;i++)
321  {
322  mesh_pt()->boundary_node_pt(0,n)->pin(i);
323  }
324  }
325 
326  // Pin redudant pressure dofs
327  RefineableAxisymmetricNavierStokesEquations::
328  pin_redundant_nodal_pressures(mesh_pt()->element_pt());
329 
330  //Pin a single pressure value
331  dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(0))->fix_pressure(0,0.0);
332 
333  //Setup all the equation numbering and look-up schemes
334  std::cout << assign_eqn_numbers() << std::endl;
335 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Definition: torus.cc:63
RefineableFullCircleMesh< ELEMENT > * mesh_pt()
Return a pointer to the specific mesh used.
Definition: torus.cc:153
Definition: general_purpose_block_preconditioners.h:321
The GMRES method.
Definition: iterative_linear_solver.h:1227
void add_mesh(const Mesh *mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
Definition: general_purpose_block_preconditioners.h:191
Definition: geom_objects.h:101
Definition: hypre_solver.h:826
Preconditioner *& preconditioner_pt()
Access function to preconditioner.
Definition: iterative_linear_solver.h:95
double & tolerance()
Access to convergence tolerance.
Definition: iterative_linear_solver.h:107
unsigned & max_iter()
Access to max. number of iterations.
Definition: iterative_linear_solver.h:113
Definition: navier_stokes_preconditioners.h:607
void set_p_preconditioner(Preconditioner *new_p_preconditioner_pt)
Function to set a new pressure matrix preconditioner (inexact solver)
Definition: navier_stokes_preconditioners.h:745
void set_f_preconditioner(Preconditioner *new_f_preconditioner_pt)
Function to set a new momentum matrix preconditioner (inexact solver)
Definition: navier_stokes_preconditioners.h:769
void set_navier_stokes_mesh(Mesh *mesh_pt, const bool &allow_multiple_element_type_in_navier_stokes_mesh=false)
Definition: navier_stokes_preconditioners.h:732
Definition: preconditioner.h:54
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Definition: problem.h:1466
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: full_circle_mesh.template.h:114
Definition: error_estimator.h:266
double Pi
Definition: two_d_biharmonic.cc:235
double Re
Reynolds number.
Definition: axisym_navier_stokes/counter_rotating_disks/counter_rotating_disks.cc:61
double Gamma
Aspect ratio (cylinder height / cylinder radius)
Definition: axisym_navier_stokes/counter_rotating_disks/counter_rotating_disks.cc:70
Global variables.
Definition: TwenteMeshGluing.cpp:60
double Radius
Radius of the pipe.
Definition: pipe.cc:55
double Delta
The curvature of the torus.
Definition: torus.cc:51
Z2ErrorEstimator * error_estimator_pt
Definition: MortaringCantileverCompareToNonMortaring.cpp:190
const Mdouble pi
Definition: ExtendedMath.h:23
void set_defaults_for_2D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Definition: hypre_solver.cc:45

References oomph::GeneralPurposeBlockPreconditioner< MATRIX >::add_mesh(), Global_Physical_Variables::Delta, e(), MeshRefinement::error_estimator_pt, GlobalPhysicalVariables::Gamma, i, oomph::IterativeLinearSolver::max_iter(), n, constants::pi, BiharmonicTestFunctions2::Pi, oomph::IterativeLinearSolver::preconditioner_pt(), Global_Physical_Variables::Radius, GlobalPhysicalVariables::Re, oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem(), oomph::NavierStokesSchurComplementPreconditioner::set_f_preconditioner(), oomph::NavierStokesSchurComplementPreconditioner::set_navier_stokes_mesh(), oomph::NavierStokesSchurComplementPreconditioner::set_p_preconditioner(), and oomph::IterativeLinearSolver::tolerance().

Member Function Documentation

◆ actions_after_adapt()

template<class ELEMENT >
void TorusProblem< ELEMENT >::actions_after_adapt ( )
inlinevirtual

After adaptation: Pin pressure again (the previously pinned value might have disappeared) and pin redudant pressure dofs.

Reimplemented from oomph::Problem.

162  {
163  // Unpin all pressure dofs
164  RefineableAxisymmetricNavierStokesEquations::
165  unpin_all_pressure_dofs(mesh_pt()->element_pt());
166 
167  // Pin redudant pressure dofs
168  RefineableAxisymmetricNavierStokesEquations::
169  pin_redundant_nodal_pressures(mesh_pt()->element_pt());
170 
171  //Pin a single pressure value
172  dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(0))->fix_pressure(0,0.0);
173  }

◆ actions_before_implicit_timestep()

template<class ELEMENT >
void TorusProblem< ELEMENT >::actions_before_implicit_timestep ( )
inlinevirtual

Update the problem specs before next timestep:

Reimplemented from oomph::Problem.

void set_boundary_conditions(const double &time)
Set boundary conditions on the walls.
Definition: torus.cc:342
double & time()
Return the current value of continuous time.
Definition: problem.cc:11531

◆ mesh_pt()

template<class ELEMENT >
RefineableFullCircleMesh<ELEMENT>* TorusProblem< ELEMENT >::mesh_pt ( )
inline

Return a pointer to the specific mesh used.

154  {return dynamic_cast<RefineableFullCircleMesh<ELEMENT>*>(Problem::mesh_pt());}

◆ set_boundary_conditions()

template<class ELEMENT >
void TorusProblem< ELEMENT >::set_boundary_conditions ( const double time)

Set boundary conditions on the walls.

Set the boundary conditions as a function of time we are going to spin up the torus

343 {
344  //NOTE: The default value of all parameters is zero, so we need only
345  //set the values that are non-zero on the boundary, i.e. the swirl
346 
347  //Loop over the nodes on the boundary
348  unsigned n_boundary_node = mesh_pt()->nboundary_node(0);
349  //Loop over the nodes on the boundary
350  for(unsigned n=0;n<n_boundary_node;n++)
351  {
352  //Get the radial values
353  double r = mesh_pt()->boundary_node_pt(0,n)->x(0);
354  //Set the value of the w-velocity
355  //Fast(ish) spin-up
356  mesh_pt()->boundary_node_pt(0,n)->set_value(2,r*(1.0 - exp(-20.0*time)));
357  }
358 }
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 exp(const bfloat16 &a)
Definition: BFloat16.h:615
r
Definition: UniformPSDSelfTest.py:20

References Eigen::bfloat16_impl::exp(), n, and UniformPSDSelfTest::r.

◆ set_initial_condition()

template<class ELEMENT >
void TorusProblem< ELEMENT >::set_initial_condition ( )
inlinevirtual

Set the initial conditions: all nodes have zero velocity.

Reimplemented from oomph::Problem.

135  {
136  const unsigned n_node = mesh_pt()->nnode();
137  for(unsigned n=0;n<n_node;n++)
138  {
139  for(unsigned i=0;i<3;i++)
140  {
141  mesh_pt()->node_pt(n)->set_value(i,0.0);
142  }
143  }
144  }

References i, and n.

◆ solve_system()

template<class ELEMENT >
void TorusProblem< ELEMENT >::solve_system ( const double dt,
const unsigned nstep 
)

Function that is used to run the parameter study.

Solve the system for a number of different values of the Reynolds number.

367 {
368  using namespace Global_Physical_Variables;
369 
370  //Open a trace file
371  ofstream trace("time_trace.dat");
372 
373  //Define a string that we can set to be the name of the output file
374  char filename[100];
375  //Define an output filestream
376  ofstream file;
377 
378  //Set the Reynolds number
379  Re = 1000.0;
380 
381  //Set an impulsive start from rest
383 
384  //Output intital data
385  trace << time() << " " << mesh_pt()->boundary_node_pt(0,0)->value(2)
386  << " " << mesh_pt()->node_pt(0)->value(0)
387  << " " << mesh_pt()->node_pt(0)->value(1)
388  << " " << mesh_pt()->node_pt(0)->value(2) << std::endl;
389 
390  //Increase the maximum value of the residuals to get
391  //past the first few steps
392  Max_residuals = 50.0;
393 
394  //Now perform the first timestep with 2 steps of refinement
395  unsteady_newton_solve(dt,2,true);
396 
397  trace << time() << " " << mesh_pt()->boundary_node_pt(0,0)->value(2)
398  << " " << mesh_pt()->node_pt(0)->value(0)
399  << " " << mesh_pt()->node_pt(0)->value(1)
400  << " " << mesh_pt()->node_pt(0)->value(2) << std::endl;
401 
402  //Output data after the first timestep
403  //Create the filename, including the array index
404  sprintf(filename,"soln_Re%g_t%g.dat",Re,time());
405  //Actually, write the data
406  file.open(filename);
407  mesh_pt()->output(file,5);
408  file.close();
409 
410  //Now do the other steps with only one adaptation per step
411  for(unsigned n=1;n<nstep;n++)
412  {
413  //Solve the problem
414  unsteady_newton_solve(dt,1,false);
415 
416  trace << time() << " " << mesh_pt()->boundary_node_pt(0,0)->value(2)
417  << " " << mesh_pt()->node_pt(0)->value(0)
418  << " " << mesh_pt()->node_pt(0)->value(1)
419  << " " << mesh_pt()->node_pt(0)->value(2) << std::endl;
420 
421  //Output data at each step
422  //Create the filename, including the array index
423  sprintf(filename,"soln_Re%g_t%g.dat",Re,time());
424  //Actually, write the data
425  file.open(filename);
426  mesh_pt()->output(file,5);
427  file.close();
428  }
429 
430  //Close the trace file
431  trace.close();
432 }
void assign_initial_values_impulsive()
Definition: problem.cc:11499
double Max_residuals
Definition: problem.h:610
void unsteady_newton_solve(const double &dt)
Definition: problem.cc:10953
string filename
Definition: MergeRestartFiles.py:39

References MergeRestartFiles::filename, n, and GlobalPhysicalVariables::Re.


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