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

Public Member Functions

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

Private Attributes

MeshBulk_mesh_pt
 Pointer to the "bulk" mesh. More...
 
MeshSurface_mesh_P_pt
 Pointer to the "surface" mesh. More...
 
PreconditionerPrec_pt
 Preconditioner. More...
 
PreconditionerNavier_stokes_prec_pt
 Preconditioner for the Navier-Stokes block. More...
 
PreconditionerF_preconditioner_pt
 Preconditioner for the momentum block. More...
 
PreconditionerP_preconditioner_pt
 Preconditioner for the pressure block. More...
 
IterativeLinearSolverSolver_pt
 Iterative linear solver. More...
 
unsigned Bottom_bound
 Bottom boundary. More...
 
unsigned Right_bound
 Right boundary. More...
 
unsigned Top_bound
 Top boundary. More...
 
unsigned Left_bound
 Left boundary. 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_after_newton_solve ()
 
virtual void actions_before_newton_convergence_check ()
 
virtual void actions_before_newton_step ()
 
virtual void actions_before_implicit_timestep ()
 
virtual void actions_after_implicit_timestep ()
 
virtual void actions_after_implicit_timestep_and_error_estimation ()
 
virtual void actions_before_explicit_timestep ()
 Actions that should be performed before each explicit time step. More...
 
virtual void actions_after_explicit_timestep ()
 Actions that should be performed after each explicit time step. More...
 
virtual void actions_before_read_unstructured_meshes ()
 
virtual void actions_after_read_unstructured_meshes ()
 
virtual void actions_after_change_in_global_parameter (double *const &parameter_pt)
 
virtual void actions_after_change_in_bifurcation_parameter ()
 
virtual void actions_after_parameter_increase (double *const &parameter_pt)
 
doubledof_derivative (const unsigned &i)
 
doubledof_current (const unsigned &i)
 
virtual void set_initial_condition ()
 
virtual double global_temporal_error_norm ()
 
unsigned newton_solve_continuation (double *const &parameter_pt)
 
unsigned newton_solve_continuation (double *const &parameter_pt, DoubleVector &z)
 
void calculate_continuation_derivatives (double *const &parameter_pt)
 
void calculate_continuation_derivatives (const DoubleVector &z)
 
void calculate_continuation_derivatives_fd (double *const &parameter_pt)
 
bool does_pointer_correspond_to_problem_data (double *const &parameter_pt)
 
void set_consistent_pinned_values_for_continuation ()
 
- Protected Attributes inherited from oomph::Problem
Vector< Problem * > Copy_of_problem_pt
 
std::map< double *, boolCalculate_dparameter_analytic
 
bool Calculate_hessian_products_analytic
 
LinearAlgebraDistributionDof_distribution_pt
 
Vector< double * > Dof_pt
 Vector of pointers to dofs. More...
 
DoubleVectorWithHaloEntries Element_count_per_dof
 
double Relaxation_factor
 
double Newton_solver_tolerance
 
unsigned Max_newton_iterations
 Maximum number of Newton iterations. More...
 
unsigned Nnewton_iter_taken
 
Vector< doubleMax_res
 Maximum residuals at start and after each newton iteration. More...
 
double Max_residuals
 
bool Time_adaptive_newton_crash_on_solve_fail
 
bool Jacobian_reuse_is_enabled
 Is re-use of Jacobian in Newton iteration enabled? Default: false. More...
 
bool Jacobian_has_been_computed
 
bool Problem_is_nonlinear
 
bool Pause_at_end_of_sparse_assembly
 
bool Doc_time_in_distribute
 
unsigned Sparse_assembly_method
 
unsigned Sparse_assemble_with_arrays_initial_allocation
 
unsigned Sparse_assemble_with_arrays_allocation_increment
 
Vector< Vector< unsigned > > Sparse_assemble_with_arrays_previous_allocation
 
double Numerical_zero_for_sparse_assembly
 
double FD_step_used_in_get_hessian_vector_products
 
bool Mass_matrix_reuse_is_enabled
 
bool Mass_matrix_has_been_computed
 
bool Discontinuous_element_formulation
 
double Minimum_dt
 Minimum desired dt: if dt falls below this value, exit. More...
 
double Maximum_dt
 Maximum desired dt. More...
 
double DTSF_max_increase
 
double DTSF_min_decrease
 
double Minimum_dt_but_still_proceed
 
bool Scale_arc_length
 Boolean to control whether arc-length should be scaled. More...
 
double Desired_proportion_of_arc_length
 Proportion of the arc-length to taken by the parameter. More...
 
double Theta_squared
 
int Sign_of_jacobian
 Storage for the sign of the global Jacobian. More...
 
double Continuation_direction
 
double Parameter_derivative
 Storage for the derivative of the global parameter wrt arc-length. More...
 
double Parameter_current
 Storage for the present value of the global parameter. More...
 
bool Use_continuation_timestepper
 Boolean to control original or new storage of dof stuff. More...
 
unsigned Dof_derivative_offset
 
unsigned Dof_current_offset
 
Vector< doubleDof_derivative
 Storage for the derivative of the problem variables wrt arc-length. More...
 
Vector< doubleDof_current
 Storage for the present values of the variables. More...
 
double Ds_current
 Storage for the current step value. More...
 
unsigned Desired_newton_iterations_ds
 
double Minimum_ds
 Minimum desired value of arc-length. More...
 
bool Bifurcation_detection
 Boolean to control bifurcation detection via determinant of Jacobian. More...
 
bool Bisect_to_find_bifurcation
 Boolean to control wheter bisection is used to located bifurcation. More...
 
bool First_jacobian_sign_change
 Boolean to indicate whether a sign change has occured in the Jacobian. More...
 
bool Arc_length_step_taken
 Boolean to indicate whether an arc-length step has been taken. More...
 
bool Use_finite_differences_for_continuation_derivatives
 
OomphCommunicatorCommunicator_pt
 The communicator for this problem. More...
 
bool Always_take_one_newton_step
 
double Timestep_reduction_factor_after_nonconvergence
 
bool Keep_temporal_error_below_tolerance
 
- Static Protected Attributes inherited from oomph::Problem
static ContinuationStorageScheme Continuation_time_stepper
 Storage for the single static continuation timestorage object. More...
 

Constructor & Destructor Documentation

◆ TiltedCavityProblem()

template<class ELEMENT >
TiltedCavityProblem< ELEMENT >::TiltedCavityProblem

Problem constructor.

Constructor: Pass number of elements in x and y directions and lengths

313 {
314  // Alias the namespace for convenience
315  namespace GV = Global_Variables;
316 
317  // Assign the enumerations for the boundaries.
318  Bottom_bound = 0;
319  Right_bound = 1;
320  Top_bound = 2;
321  Left_bound = 3;
322 
323  // First we create the tilted cavity mesh.
325  GV::Lx,GV::Ly,
326  GV::Ang_rad);
327 
328  // Set the boundary conditions, recall that the boundaries are:
329  //
330  // 2 non slip
331  // ----------
332  // | |
333  // 3 Inflow| |1 P.O. (parallel outflow)
334  // | |
335  // ----------
336  // 0 non slip
337 
338  // Parallel outflow boundary is located at 1.
339  const unsigned po_bound = Right_bound;
340 
341  // Create a "surface mesh" that will contain only
342  // ImposeParallelOutflowElements in boundary 1
343  // The constructor just creates the mesh without
344  // giving it any elements, nodes, etc.
345  Surface_mesh_P_pt = new Mesh;
346 
347  // Create ImposeParallelOutflowElement from all elements that are
348  // adjacent to the Neumann boundary.
350  Bulk_mesh_pt,
352 
353  // Add the two meshes to the problem.
356 
357  // combine all sub-meshes into a single mesh.
359 
360  // All nodes are free by default
361  // just pin the ones that have Dirichlet conditions
362  const unsigned num_bound = mesh_pt()->nboundary();
363  for(unsigned ibound=0;ibound<num_bound;ibound++)
364  {
365  if(ibound != po_bound)
366  {
367  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
368  for (unsigned inod=0;inod<num_nod;inod++)
369  {
370  // Get node
371  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
372 
373  nod_pt->pin(0);
374  nod_pt->pin(1);
375 
376  nod_pt->set_value(0,0);
377  nod_pt->set_value(1,0);
378  }
379  }
380  }
381 
382  //Complete the problem setup to make the elements fully functional
383 
384  //Loop over the elements
385  unsigned n_el = Bulk_mesh_pt->nelement();
386  for(unsigned e=0;e<n_el;e++)
387  {
388  //Cast to a fluid element
389  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
390 
391  //Set the Reynolds number, etc
392  el_pt->re_pt() = &GV::Re;
393  } // for(unsigned e=0;e<n_el;e++)
394 
395  //Assign equation numbers
396  oomph_info << "\n equation numbers : "
397  << assign_eqn_numbers() << std::endl;
398 
399 
400  // We choose either trilinos or oomph-lib's GMRES as our linear solver
401  if(GV::Use_trilinos)
402  {
403 #ifdef OOMPH_HAS_TRILINOS
404  // Create the trilinos solver.
405  TrilinosAztecOOSolver* trilinos_solver_pt = new TrilinosAztecOOSolver;
406  trilinos_solver_pt->solver_type() = TrilinosAztecOOSolver::GMRES;
407 
408  // Store the solver pointer.
409  Solver_pt = trilinos_solver_pt;
410 #endif
411  }
412  else
413  {
414  // Create oomph-lib iterative linear solver.
416 
417  // We use RHS preconditioning. Note that by default,
418  // left hand preconditioning is used.
419  static_cast<GMRES<CRDoubleMatrix>*>(solver_pt)
420  ->set_preconditioner_RHS();
421 
422  // Store the solver pointer.
423  Solver_pt = solver_pt;
424  }
425 
426  // Set up other solver parameters.
427  Solver_pt->tolerance() = 1.0e-6;
428 
429  // Overview of solvers:
430  //
431  // The Jacobian takes the block form:
432  //
433  // | F_ns | L^T |
434  // |------------|
435  // | L | 0 |
436  //
437  // where L correspond to the constrained block,
438  // F_ns is the Navier-Stokes block with the following block structure
439  //
440  // | F | B^T |
441  // |----------|
442  // | B | 0 |
443  //
444  //
445  // The lagrange enforced flow preconditioner takes the form:
446  // F_ns + L^T inv(W) L |
447  // --------------------------
448  // | W
449  //
450  // where W=LL^T
451  //
452  // We use SuperLU to solve the W block (2,2)
453  //
454  // For the (1,1) block, we can use SuperLU or the LSC preconditioner.
455 
456 
457  // Create the preconditioner
460 
461  // Create the vector of mesh pointers!
463  mesh_pt.resize(2,0);
464  mesh_pt[0] = Bulk_mesh_pt;
466 
467  lgr_prec_pt->set_meshes(mesh_pt);
468 
470  if(GV::Use_lsc)
471  {
472  // Create the NS LSC preconditioner.
473  lsc_prec_pt = new NavierStokesSchurComplementPreconditioner(this);
474  lsc_prec_pt->set_navier_stokes_mesh(Bulk_mesh_pt);
475  lsc_prec_pt->use_lsc();
476  lgr_prec_pt->set_navier_stokes_preconditioner(lsc_prec_pt);
477 
479  {
480  if(GV::Visc == 0)
481  {
486  }
487  else
488  {
493  }
494  }
495 
497  {
502  }
503  }
504  else
505  {
507  }
508 
509  // Store the preconditioner pointers.
510  Navier_stokes_prec_pt = lsc_prec_pt;
511  Prec_pt = lgr_prec_pt;
512 
513  // Pass the preconditioner to the solver.
514  Solver_pt->preconditioner_pt() = lgr_prec_pt;
515 
516  // Pass the solver to the problem.
517  this->linear_solver_pt() = Solver_pt;
518 
519  // Set the Newton solver tolerance.
520  this->newton_solver_tolerance() = 1.0e-6;
521 }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
unsigned Bottom_bound
Bottom boundary.
Definition: two_d_tilted_square.cc:294
void create_parall_outflow_lagrange_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &surface_mesh_pt)
Definition: two_d_tilted_square.cc:553
Mesh * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Definition: two_d_tilted_square.cc:262
unsigned Left_bound
Left boundary.
Definition: two_d_tilted_square.cc:303
Mesh * Surface_mesh_P_pt
Pointer to the "surface" mesh.
Definition: two_d_tilted_square.cc:265
Preconditioner * F_preconditioner_pt
Preconditioner for the momentum block.
Definition: two_d_tilted_square.cc:274
Preconditioner * Prec_pt
Preconditioner.
Definition: two_d_tilted_square.cc:268
Preconditioner * P_preconditioner_pt
Preconditioner for the pressure block.
Definition: two_d_tilted_square.cc:277
Preconditioner * Navier_stokes_prec_pt
Preconditioner for the Navier-Stokes block.
Definition: two_d_tilted_square.cc:271
unsigned Right_bound
Right boundary.
Definition: two_d_tilted_square.cc:297
unsigned Top_bound
Top boundary.
Definition: two_d_tilted_square.cc:300
IterativeLinearSolver * Solver_pt
Iterative linear solver.
Definition: two_d_tilted_square.cc:280
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
The GMRES method.
Definition: iterative_linear_solver.h:1227
Definition: iterative_linear_solver.h:54
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
Definition: lagrange_enforced_flow_preconditioner.h:168
void set_navier_stokes_preconditioner(Preconditioner *new_ns_preconditioner_pt=0)
Definition: lagrange_enforced_flow_preconditioner.cc:1214
void set_meshes(const Vector< Mesh * > &mesh_pt)
Definition: lagrange_enforced_flow_preconditioner.cc:419
void set_superlu_for_navier_stokes_preconditioner()
Definition: lagrange_enforced_flow_preconditioner.h:289
Definition: mesh.h:67
unsigned long nboundary_node(const unsigned &ibound) const
Return number of nodes on a particular boundary.
Definition: mesh.h:833
unsigned nboundary() const
Return number of boundaries.
Definition: mesh.h:827
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
Node *& boundary_node_pt(const unsigned &b, const unsigned &n)
Return pointer to node n on boundary b.
Definition: mesh.h:493
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
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
void use_lsc()
Use LSC version of the preconditioner.
Definition: navier_stokes_preconditioners.h:782
Definition: nodes.h:906
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Definition: problem.h:1466
void build_global_mesh()
Definition: problem.cc:1493
double & newton_solver_tolerance()
Definition: problem.h:1621
Mesh *& mesh_pt()
Return a pointer to the global mesh.
Definition: problem.h:1280
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
Definition: two_d_tilted_square.cc:131
Definition: trilinos_solver.h:267
unsigned & solver_type()
Access function to Solver_type.
Definition: trilinos_solver.h:442
Namespace for physical parameters.
Definition: two_d_tilted_square.cc:59
bool Use_trilinos
Use trilinos?
Definition: two_d_tilted_square.cc:90
bool Use_lsc
Use LSC preconditioner for the Navier-Stokes block?
Definition: two_d_tilted_square.cc:93
unsigned Visc
The viscous term.
Definition: two_d_tilted_square.cc:77
bool Use_amg_for_p
Use Boomer AMG for the pressure block?
Definition: two_d_tilted_square.cc:99
bool Use_amg_for_f
Use Boomer AMG for the momentum block?
Definition: two_d_tilted_square.cc:96
static const double Ly
Definition: two_d_tilted_square.cc:74
unsigned Noel
Number of elements in 1D.
Definition: two_d_tilted_square.cc:87
double Ang_rad
Definition: two_d_tilted_square.cc:84
double Re
Reynolds number.
Definition: two_d_tilted_square.cc:80
static const double Lx
The domain length in the x and y direction respectively.
Definition: two_d_tilted_square.cc:73
Preconditioner * boomer_amg_for_2D_momentum_stressdiv_visc()
Definition: lagrange_enforced_flow_preconditioner.cc:107
Preconditioner * boomer_amg_for_2D_poisson_problem()
Definition: lagrange_enforced_flow_preconditioner.cc:225
Preconditioner * boomer_amg_for_2D_momentum_simple_visc()
Definition: lagrange_enforced_flow_preconditioner.cc:67
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References Global_Variables::Ang_rad, oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_simple_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_momentum_stressdiv_visc(), oomph::Lagrange_Enforced_Flow_Preconditioner_Subsidiary_Operator_Helper::boomer_amg_for_2D_poisson_problem(), e(), Global_Variables::Lx, Global_Variables::Ly, Global_Variables::Noel, oomph::oomph_info, oomph::Data::pin(), Global_Variables::Re, oomph::NavierStokesSchurComplementPreconditioner::set_f_preconditioner(), oomph::LagrangeEnforcedFlowPreconditioner::set_meshes(), oomph::NavierStokesSchurComplementPreconditioner::set_navier_stokes_mesh(), oomph::LagrangeEnforcedFlowPreconditioner::set_navier_stokes_preconditioner(), oomph::NavierStokesSchurComplementPreconditioner::set_p_preconditioner(), oomph::LagrangeEnforcedFlowPreconditioner::set_superlu_for_navier_stokes_preconditioner(), oomph::Data::set_value(), oomph::TrilinosAztecOOSolver::solver_type(), oomph::IterativeLinearSolver::tolerance(), Global_Variables::Use_amg_for_f, Global_Variables::Use_amg_for_p, Global_Variables::Use_lsc, oomph::NavierStokesSchurComplementPreconditioner::use_lsc(), Global_Variables::Use_trilinos, and Global_Variables::Visc.

Member Function Documentation

◆ actions_after_newton_step()

template<class ELEMENT >
void TiltedCavityProblem< ELEMENT >::actions_after_newton_step ( )
inlinevirtual

Update after Newton step - document the number of iterations required for the iterative solver to converge.

Reimplemented from oomph::Problem.

217  {
218  // Alias the namespace for convenience
219  namespace GV = Global_Variables;
220 
221  unsigned iters = 0;
222  // Get the iteration counts
223 #ifdef PARANOID
224  IterativeLinearSolver* iterative_solver_pt
225  = dynamic_cast<IterativeLinearSolver*>
226  (this->linear_solver_pt());
227  if(iterative_solver_pt == 0)
228  {
229  std::ostringstream error_message;
230  error_message << "Cannot cast the solver pointer." << std::endl;
231 
232  throw OomphLibError(error_message.str(),
235  }
236  else
237  {
238  iters = iterative_solver_pt->iterations();
239  GV::Iterations.push_back(iters);
240  }
241 #else
242  iters = static_cast<IterativeLinearSolver*>
243  (this->linear_solver_pt())->iterations();
244  GV::Iterations.push_back(iters);
245 #endif
246  }
virtual unsigned iterations() const =0
Number of iterations taken.
Definition: oomph_definitions.h:222
Vector< unsigned > Iterations
Storage for number of iterations during Newton steps.
Definition: two_d_tilted_square.cc:111
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Global_Variables::Iterations, oomph::IterativeLinearSolver::iterations(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ actions_before_newton_solve()

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

Update before Newton solve.

Reimplemented from oomph::Problem.

174  {
175  // Alias the namespace for convenience
176  namespace GV = Global_Variables;
177  GV::Iterations.clear();
178 
179  // Which is the inflow boundary?
180  const unsigned in_bound = Left_bound;
181 
182  // The number of nodes on a boundary
183  unsigned num_nod = mesh_pt()->nboundary_node(in_bound);
184  for(unsigned inod=0;inod<num_nod;inod++)
185  {
186  Node* nod_pt=mesh_pt()->boundary_node_pt(in_bound,inod);
187 
188  // Pin both velocity components
189  nod_pt->pin(0);
190  nod_pt->pin(1);
191 
192  // Get the x and y cartesian coordinates
193  double x0=nod_pt->x(0);
194  double x1=nod_pt->x(1);
195 
196  // Tilt x1 by -Ang_rad, this will give us the original coordinate.
197  double x1_old = x0*sin(-GV::Ang_rad) + x1*cos(-GV::Ang_rad);
198 
199  // Now calculate the parabolic inflow at this point.
200  double u0_old = (x1_old - GV::Y_min)*(GV::Y_max - x1_old);
201 
202  // Now apply the rotation to u0_old, using rotation matrices.
203  // with x = u0_old and y = 0, i.e. R*[u;0] since we have the
204  // velocity in the x direction only. There is no velocity
205  // in the y direction.
206  double u0=u0_old*cos(GV::Ang_rad);
207  double u1=u0_old*sin(GV::Ang_rad);
208 
209  nod_pt->set_value(0,u0);
210  nod_pt->set_value(1,u1);
211  } // for
212  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
static const double Y_max
Definition: two_d_tilted_square.cc:70
static const double Y_min
Min and max y value respectively.
Definition: two_d_tilted_square.cc:69
Vector< double > x1(const Vector< double > &coord)
Cartesian coordinates centered at the point (0.5,1)
Definition: poisson/poisson_with_singularity/two_d_poisson.cc:86
Vector< double > x0(2, 0.0)

References Global_Variables::Ang_rad, cos(), Global_Variables::Iterations, oomph::Data::pin(), oomph::Data::set_value(), sin(), oomph::Node::x(), Global::x0, Global_parameters::x1(), Global_Variables::Y_max, and Global_Variables::Y_min.

◆ create_parall_outflow_lagrange_elements()

template<class ELEMENT >
void TiltedCavityProblem< ELEMENT >::create_parall_outflow_lagrange_elements ( const unsigned b,
Mesh *const &  bulk_mesh_pt,
Mesh *const &  surface_mesh_pt 
)

Create lagrange 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 ImposeParallelOutflowElement 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.

556 {
557  // How many bulk elements are adjacent to boundary b?
558  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
559 
560  // Loop over the bulk elements adjacent to boundary b
561  for(unsigned e=0;e<n_element;e++)
562  {
563  // Get pointer to the bulk element that is adjacent to boundary b
564  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
565  bulk_mesh_pt->boundary_element_pt(b,e));
566 
567  // What is the index of the face of element e along boundary b?
568  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
569 
570  // Build the corresponding impose_impenetrability_element
571  ImposeParallelOutflowElement<ELEMENT>* flux_element_pt = new
573  face_index,0);
574 
575 
576  // Add the prescribed-flux element to the surface mesh
577  surface_mesh_pt->add_element_pt(flux_element_pt);
578 
579  // Loop over the nodes
580  unsigned nnod=flux_element_pt->nnode();
581  for (unsigned j=0;j<nnod;j++)
582  {
583  Node* nod_pt = flux_element_pt->node_pt(j);
584 
585  // Is the node also on boundary 0 or 2?
586  if ((nod_pt->is_on_boundary(0))||(nod_pt->is_on_boundary(2)))
587  {
588  // How many nodal values were used by the "bulk" element
589  // that originally created this node?
590  unsigned n_bulk_value=flux_element_pt->nbulk_value(j);
591 
592  // The remaining ones are Lagrange multipliers and we pin them.
593  unsigned nval=nod_pt->nvalue();
594  for (unsigned j=n_bulk_value;j<nval;j++)
595  {
596  nod_pt->pin(j);
597  }
598  }
599  }
600  }
601 }
Scalar * b
Definition: benchVecAdd.cpp:17
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
unsigned & nbulk_value(const unsigned &n)
Definition: elements.h:4845
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Definition: impose_parallel_outflow_element.h:46
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
virtual bool is_on_boundary() const
Definition: nodes.h:1373
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::Mesh::add_element_pt(), b, oomph::Mesh::boundary_element_pt(), e(), oomph::Mesh::face_index_at_boundary(), oomph::Node::is_on_boundary(), j, oomph::Mesh::nboundary_element(), oomph::FaceElement::nbulk_value(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::Data::nvalue(), and oomph::Data::pin().

◆ doc_solution()

template<class ELEMENT >
void TiltedCavityProblem< ELEMENT >::doc_solution

Doc the solution.

528 {
529  // Alias the namespace for convenience
530  namespace GV = Global_Variables;
531 
532  // Create the soln file name string.
533  std::stringstream filename;
534  filename <<"RESLT/soln"<<GV::Soln_num<<".dat";
535 
536  // Number of plot points
537  unsigned npts=5;
538 
539  // Output solution
540  std::ofstream some_file;
541  some_file.open(filename.str().c_str());
542  Bulk_mesh_pt->output(some_file,npts);
543  some_file.close();
544 }
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
unsigned Soln_num
Soln number (for doc_solution)
Definition: two_d_tilted_square.cc:102
string filename
Definition: MergeRestartFiles.py:39

References MergeRestartFiles::filename, and Global_Variables::Soln_num.

Member Data Documentation

◆ Bottom_bound

template<class ELEMENT >
unsigned TiltedCavityProblem< ELEMENT >::Bottom_bound
private

Bottom boundary.

◆ Bulk_mesh_pt

template<class ELEMENT >
Mesh* TiltedCavityProblem< ELEMENT >::Bulk_mesh_pt
private

Pointer to the "bulk" mesh.

◆ F_preconditioner_pt

template<class ELEMENT >
Preconditioner* TiltedCavityProblem< ELEMENT >::F_preconditioner_pt
private

Preconditioner for the momentum block.

◆ Left_bound

template<class ELEMENT >
unsigned TiltedCavityProblem< ELEMENT >::Left_bound
private

Left boundary.

◆ Navier_stokes_prec_pt

template<class ELEMENT >
Preconditioner* TiltedCavityProblem< ELEMENT >::Navier_stokes_prec_pt
private

Preconditioner for the Navier-Stokes block.

◆ P_preconditioner_pt

template<class ELEMENT >
Preconditioner* TiltedCavityProblem< ELEMENT >::P_preconditioner_pt
private

Preconditioner for the pressure block.

◆ Prec_pt

template<class ELEMENT >
Preconditioner* TiltedCavityProblem< ELEMENT >::Prec_pt
private

Preconditioner.

◆ Right_bound

template<class ELEMENT >
unsigned TiltedCavityProblem< ELEMENT >::Right_bound
private

Right boundary.

◆ Solver_pt

template<class ELEMENT >
IterativeLinearSolver* TiltedCavityProblem< ELEMENT >::Solver_pt
private

Iterative linear solver.

◆ Surface_mesh_P_pt

template<class ELEMENT >
Mesh* TiltedCavityProblem< ELEMENT >::Surface_mesh_P_pt
private

Pointer to the "surface" mesh.

◆ Top_bound

template<class ELEMENT >
unsigned TiltedCavityProblem< ELEMENT >::Top_bound
private

Top boundary.


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