DDConvectionProblem< NST_ELEMENT, AD_ELEMENT > Class Template Reference
+ Inheritance diagram for DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >:

Public Member Functions

 DDConvectionProblem ()
 Constructor. More...
 
 ~DDConvectionProblem ()
 Destructor. Empty. More...
 
void actions_before_newton_solve ()
 Update the problem specs before solve (empty) More...
 
void actions_after_newton_solve ()
 Update the problem after solve (empty) More...
 
void actions_before_adapt ()
 Actions before adapt:(empty) More...
 
void actions_before_implicit_timestep ()
 
void fix_pressure (const unsigned &e, const unsigned &pdof, const double &pvalue)
 Fix pressure in element e at pressure dof pdof and set to pvalue. More...
 
void doc_solution ()
 Doc the solution. More...
 
void set_boundary_conditions (const double &time)
 Set the boundary conditions. More...
 
RectangularQuadMesh< NST_ELEMENT > * nst_mesh_pt ()
 Access function to the Navier-Stokes mesh. More...
 
RectangularQuadMesh< AD_ELEMENT > * temp_mesh_pt ()
 Access function to the Advection-Diffusion mesh. More...
 
RectangularQuadMesh< AD_ELEMENT > * conc_mesh_pt ()
 Access function to the Advection-Diffusion concentration mesh. More...
 
void get_kinetic_energy (double &E, double &Edot)
 Get kinetic energy and kinetic energy flux. 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_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 ()
 

Protected Attributes

RectangularQuadMesh< NST_ELEMENT > * Nst_mesh_pt
 
RectangularQuadMesh< AD_ELEMENT > * Temp_mesh_pt
 
RectangularQuadMesh< AD_ELEMENT > * Conc_mesh_pt
 
- 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
 

Private Attributes

DocInfo Doc_info
 DocInfo object. 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 ()
 
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 ()
 
- 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 NST_ELEMENT, class AD_ELEMENT>
class DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >

/////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// 2D double-diffusive Convection problem on three rectangular domains, discretised with Navier-Stokes and Advection-Diffusion elements. The specific type of elements is specified via the template parameters.

Constructor & Destructor Documentation

◆ DDConvectionProblem()

template<class NST_ELEMENT , class AD_ELEMENT >
DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::DDConvectionProblem

Constructor.

Constructor for convection problem.

751 {
752 
753  //Allocate a timestepper
755 
756  // Set output directory
757  Doc_info.set_directory("RESLT");
758 
759  // # of elements in x-direction
760  unsigned n_x=8;
761 
762  // # of elements in y-direction
763  unsigned n_y=8;
764 
765  // Domain length in x-direction
767 
768  // Domain length in y-direction
769  double l_y=1.0;
770 
771  // Build two standard rectangular quadmesh
772  Nst_mesh_pt =
773  new RectangularQuadMesh<NST_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
774  Temp_mesh_pt =
775  new RectangularQuadMesh<AD_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
776  Conc_mesh_pt =
777  new RectangularQuadMesh<AD_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
778 
779  // Set the boundary conditions for this problem: All nodes are
780  // free by default -- only need to pin the ones that have Dirichlet
781  // conditions here
782 
783  //Loop over the boundaries
784  unsigned num_bound = nst_mesh_pt()->nboundary();
785  for(unsigned ibound=0;ibound<num_bound;ibound++)
786  {
787  //Loop over the number of nodes on the boundry
788  unsigned num_nod= nst_mesh_pt()->nboundary_node(ibound);
789  for (unsigned inod=0;inod<num_nod;inod++)
790  {
791  //If we are on the side-walls, the v-velocity
792  //satisfies natural boundary conditions, so we only pin the
793  //first value
794  if ((ibound==1) || (ibound==3))
795  {
796  nst_mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
797  }
798  else // On the top and bottom walls, we have "stress-free" conditions
799  //which actually corresponds to transverse stress free and normal
800  //zero velocity (symmetry)
801  //Thus we pin the second value
802  {
803  nst_mesh_pt()->boundary_node_pt(ibound,inod)->pin(1);
804  }
805  }
806  }
807 
808  //Pin the zero-th pressure dof in element 0 and set its value to
809  //zero:
810  fix_pressure(0,0,0.0);
811 
812  //Loop over the boundaries of the AD mesh
813  num_bound = temp_mesh_pt()->nboundary();
814  for(unsigned ibound=0;ibound<num_bound;ibound++)
815  {
816  //Loop over the number of nodes on the boundry
817  unsigned num_nod= temp_mesh_pt()->nboundary_node(ibound);
818  for (unsigned inod=0;inod<num_nod;inod++)
819  {
820  //If we are on the side-walls, the temperature
821  //satisfies natural boundary conditions, so we don't pin anything
822  // in this mesh
823  if ((ibound==1) || (ibound==3))
824  {
825 
826  }
827  //Otherwise pin the temperature
828  else // pin all values
829  {
830  temp_mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
831  }
832  }
833  }
834 
835 
836  //Loop over the boundaries of the AD mesh
837  num_bound = conc_mesh_pt()->nboundary();
838  for(unsigned ibound=0;ibound<num_bound;ibound++)
839  {
840  //Loop over the number of nodes on the boundry
841  unsigned num_nod= conc_mesh_pt()->nboundary_node(ibound);
842  for (unsigned inod=0;inod<num_nod;inod++)
843  {
844  //If we are on the side-walls, the concentration
845  //satisfies natural boundary conditions, so we don't pin anything
846  // in this mesh
847  if ((ibound==1) || (ibound==3))
848  {
849 
850  }
851  //Otherwiwse pin the concentration
852  else // pin all values
853  {
854  conc_mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
855  }
856  }
857  }
858 
859 
860  // Complete the build of all elements so they are fully functional
861 
862  // Loop over the elements to set up element-specific
863  // things that cannot be handled by the (argument-free!) ELEMENT
864  // constructors.
865  unsigned n_nst_element = nst_mesh_pt()->nelement();
866  for(unsigned i=0;i<n_nst_element;i++)
867  {
868  // Upcast from GeneralsedElement to the present element
869  NST_ELEMENT *el_pt = dynamic_cast<NST_ELEMENT*>
870  (nst_mesh_pt()->element_pt(i));
871 
872  // Set the Reynolds number (1/Pr in our non-dimensionalisation)
874 
875  // Set ReSt (also 1/Pr in our non-dimensionalisation)
876  el_pt->re_st_pt() = &Global_Physical_Variables::Inverse_Prandtl;
877 
878  // Set the Thermal Rayleigh number
879  el_pt->ra_t_pt() = &Global_Physical_Variables::Rayleigh_T;
880 
881  // Set the Solutal Rayleigh number
882  el_pt->ra_s_pt() = &Global_Physical_Variables::Rayleigh_S;
883 
884  //Set Gravity vector
886 
887  //The mesh is fixed, so we can disable ALE
888  el_pt->disable_ALE();
889  }
890 
891  unsigned n_temp_element = temp_mesh_pt()->nelement();
892  for(unsigned i=0;i<n_temp_element;i++)
893  {
894  // Upcast from GeneralsedElement to the present element
895  AD_ELEMENT *el_pt = dynamic_cast<AD_ELEMENT*>
896  (temp_mesh_pt()->element_pt(i));
897 
898  // Set the Peclet number
899  el_pt->pe_pt() = &Global_Physical_Variables::Peclet;
900 
901  // Set the timescale to be the same as the Navier--Stokes
902  // equations (1.0)
903  el_pt->pe_st_pt() =&Global_Physical_Variables::Peclet;
904 
905  //The mesh is fixed, so we can disable ALE
906  el_pt->disable_ALE();
907  }
908 
909  unsigned n_conc_element = conc_mesh_pt()->nelement();
910  for(unsigned i=0;i<n_conc_element;i++)
911  {
912  // Upcast from GeneralsedElement to the present element
913  AD_ELEMENT *el_pt = dynamic_cast<AD_ELEMENT*>
914  (conc_mesh_pt()->element_pt(i));
915 
916  // Set the Peclet number
917  el_pt->pe_pt() = &Global_Physical_Variables::Lewis;
918 
919  // Set the Peclet number multiplied by the Strouhal number
920  el_pt->pe_st_pt() =&Global_Physical_Variables::Lewis;
921 
922  //The mesh is fixed, so we can disable ALE
923  el_pt->disable_ALE();
924  }
925 
926 
927  // Set sources for temperature
930 
931  // Set sources for concentrations
934 
935  // combine the submeshes
940 
941  // Setup equation numbering scheme
942  cout <<"Number of equations: " << assign_eqn_numbers() << endl;
943 
944 } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
RectangularQuadMesh< NST_ELEMENT > * nst_mesh_pt()
Access function to the Navier-Stokes mesh.
Definition: multimesh_dd_convection.cc:697
DocInfo Doc_info
DocInfo object.
Definition: multimesh_dd_convection.cc:736
RectangularQuadMesh< AD_ELEMENT > * Temp_mesh_pt
Definition: multimesh_dd_convection.cc:741
RectangularQuadMesh< AD_ELEMENT > * temp_mesh_pt()
Access function to the Advection-Diffusion mesh.
Definition: multimesh_dd_convection.cc:703
RectangularQuadMesh< AD_ELEMENT > * conc_mesh_pt()
Access function to the Advection-Diffusion concentration mesh.
Definition: multimesh_dd_convection.cc:709
RectangularQuadMesh< NST_ELEMENT > * Nst_mesh_pt
Definition: multimesh_dd_convection.cc:740
void fix_pressure(const unsigned &e, const unsigned &pdof, const double &pvalue)
Fix pressure in element e at pressure dof pdof and set to pvalue.
Definition: multimesh_dd_convection.cc:682
RectangularQuadMesh< AD_ELEMENT > * Conc_mesh_pt
Definition: multimesh_dd_convection.cc:742
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
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
void add_time_stepper_pt(TimeStepper *const &time_stepper_pt)
Definition: problem.cc:1545
unsigned add_sub_mesh(Mesh *const &mesh_pt)
Definition: problem.h:1330
void build_global_mesh()
Definition: problem.cc:1493
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
TimeStepper *& time_stepper_pt()
Definition: problem.h:1524
double Rayleigh_S
Solutal Rayleigh number.
Definition: dd_convection.cc:65
Vector< double > Direction_of_gravity(2)
Gravity vector.
double Inverse_Prandtl
1/Prandtl number
Definition: mpi/multi_domain/boussinesq_convection/multi_domain_boussinesq_convection.cc:53
double Lewis
The Lewis number.
Definition: dd_convection.cc:52
double Lambda
Length of domain.
Definition: dd_convection.cc:68
double Rayleigh_T
Definition: dd_convection.cc:62
double Peclet
Peclet number.
Definition: pipe.cc:49
void setup_multi_domain_interactions(Problem *problem_pt, Mesh *const &first_mesh_pt, Mesh *const &second_mesh_pt, const unsigned &first_interaction=0, const unsigned &second_interaction=0)
Definition: multi_domain.template.cc:244

References oomph::Global_Physical_Variables::Direction_of_gravity, GlobalParameters::Doc_info, i, Global_Physical_Variables::Inverse_Prandtl, Global_Physical_Variables::Lambda, Global_Physical_Variables::Lewis, Global_Physical_Variables::Peclet, Global_Physical_Variables::Rayleigh_S, Global_Physical_Variables::Rayleigh_T, oomph::DocInfo::set_directory(), and oomph::Multi_domain_functions::setup_multi_domain_interactions().

◆ ~DDConvectionProblem()

template<class NST_ELEMENT , class AD_ELEMENT >
DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::~DDConvectionProblem ( )
inline

Destructor. Empty.

656  {
657  //Delete the meshes
658  delete Conc_mesh_pt;
659  delete Temp_mesh_pt;
660  delete Nst_mesh_pt;
661  //Delete the timestepper
662  delete this->time_stepper_pt();
663  }

Member Function Documentation

◆ actions_after_newton_solve()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::actions_after_newton_solve ( )
inlinevirtual

Update the problem after solve (empty)

Reimplemented from oomph::Problem.

669 {}

◆ actions_before_adapt()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::actions_before_adapt ( )
inlinevirtual

Actions before adapt:(empty)

Reimplemented from oomph::Problem.

672 {}

◆ actions_before_implicit_timestep()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::actions_before_implicit_timestep ( )
inlinevirtual

Actions before the timestep (update the the time-dependent boundary conditions)

Reimplemented from oomph::Problem.

677  {
679  }
void set_boundary_conditions(const double &time)
Set the boundary conditions.
Definition: multimesh_dd_convection.cc:953
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

◆ actions_before_newton_solve()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::actions_before_newton_solve ( )
inlinevirtual

Update the problem specs before solve (empty)

Reimplemented from oomph::Problem.

666 {}

◆ conc_mesh_pt()

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<AD_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::conc_mesh_pt ( )
inline

Access function to the Advection-Diffusion concentration mesh.

710  {
711  return dynamic_cast<RectangularQuadMesh<AD_ELEMENT>*>(Conc_mesh_pt);
712  }

◆ doc_solution()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::doc_solution

Doc the solution.

1044 {
1045  //Declare an output stream and filename
1046  ofstream some_file;
1047  char filename[100];
1048 
1049  // Number of plot points: npts x npts
1050  unsigned npts=5;
1051 
1052  // Output whole solution (this will output elements from one mesh
1053  //---------------------- followed by the other mesh at the moment...?)
1054  sprintf(filename,"%s/soln%i.dat",Doc_info.directory().c_str(),
1055  Doc_info.number());
1056  some_file.open(filename);
1057  mesh_pt()->output(some_file,npts);
1058  some_file.close();
1059 
1060  Doc_info.number()++;
1061 } // end of doc
std::string directory() const
Output directory.
Definition: oomph_utilities.h:524
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
Mesh *& mesh_pt()
Return a pointer to the global mesh.
Definition: problem.h:1280
string filename
Definition: MergeRestartFiles.py:39

References oomph::DocInfo::directory(), GlobalParameters::Doc_info, MergeRestartFiles::filename, and oomph::DocInfo::number().

◆ fix_pressure()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::fix_pressure ( const unsigned e,
const unsigned pdof,
const double pvalue 
)
inline

Fix pressure in element e at pressure dof pdof and set to pvalue.

684  {
685  //Cast to specific element and fix pressure (only for NST element...?)
686  dynamic_cast<NST_ELEMENT*>(nst_mesh_pt()->element_pt(e))->
687  fix_pressure(pdof,pvalue);
688  } // end_of_fix_pressure
Array< double, 1, 3 > e(1./3., 0.5, 2.)

References e().

◆ get_kinetic_energy()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::get_kinetic_energy ( double E,
double Edot 
)
inline

Get kinetic energy and kinetic energy flux.

716  {
717  //Reset values to zero
718  E = 0.0; Edot=0.0;
719 
720  //Loop over the elements
721  unsigned n_element = nst_mesh_pt()->nelement();
722  for(unsigned e=0;e<n_element;e++)
723  {
724  NST_ELEMENT* elem_pt =
725  dynamic_cast<NST_ELEMENT*>(nst_mesh_pt()->element_pt(e));
726 
727  E += elem_pt->kin_energy();
728  Edot += elem_pt->d_kin_energy_dt();
729  }
730  }
double E
Elastic modulus.
Definition: TwenteMeshGluing.cpp:68

References Global_Physical_Variables::E, and e().

◆ nst_mesh_pt()

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<NST_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::nst_mesh_pt ( )
inline

Access function to the Navier-Stokes mesh.

698  {
699  return dynamic_cast<RectangularQuadMesh<NST_ELEMENT>*>(Nst_mesh_pt);
700  }

◆ set_boundary_conditions()

template<class NST_ELEMENT , class AD_ELEMENT >
void DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::set_boundary_conditions ( const double time)

Set the boundary conditions.

Set the boundary conditions as a function of continuous time

955 {
956  // Loop over all the boundaries on the NST mesh
957  unsigned num_bound=nst_mesh_pt()->nboundary();
958  for(unsigned ibound=0;ibound<num_bound;ibound++)
959  {
960  // Loop over the nodes on boundary
961  unsigned num_nod=nst_mesh_pt()->nboundary_node(ibound);
962  for(unsigned inod=0;inod<num_nod;inod++)
963  {
964  // Get pointer to node
965  Node* nod_pt=nst_mesh_pt()->boundary_node_pt(ibound,inod);
966 
967  //Set the number of velocity components
968  unsigned vel_max=2;
969 
970  //If we are on the side walls we only set the x-velocity.
971  if((ibound==1) || (ibound==3)) {vel_max = 1;}
972 
973  //Set the pinned velocities to zero on NST mesh
974  for(unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
975  }
976  }
977 
978  // Loop over all the boundaries on the AD mesh
979  num_bound=temp_mesh_pt()->nboundary();
980  for(unsigned ibound=0;ibound<num_bound;ibound++)
981  {
982  // Loop over the nodes on boundary
983  unsigned num_nod=temp_mesh_pt()->nboundary_node(ibound);
984  for(unsigned inod=0;inod<num_nod;inod++)
985  {
986  // Get pointer to node
987  Node* nod_pt=temp_mesh_pt()->boundary_node_pt(ibound,inod);
988 
989  //If we are on the top boundary, set the temperature
990  //to -0.5 (cooled)
991  if(ibound==2) {nod_pt->set_value(0,-0.5);}
992 
993  //If we are on the bottom boundary, set the temperature
994  //to 0.5 (heated)
995  if(ibound==0) {nod_pt->set_value(0,0.5);}
996  }
997  }
998 
999 
1000 
1001  // Loop over all the boundaries on the AD mesh
1002  num_bound=conc_mesh_pt()->nboundary();
1003  for(unsigned ibound=0;ibound<num_bound;ibound++)
1004  {
1005  // Loop over the nodes on boundary
1006  unsigned num_nod=conc_mesh_pt()->nboundary_node(ibound);
1007  for(unsigned inod=0;inod<num_nod;inod++)
1008  {
1009  // Get pointer to node
1010  Node* nod_pt=conc_mesh_pt()->boundary_node_pt(ibound,inod);
1011 
1012  //If we are on the top boundary, set the concentration to be low
1013  //to -0.5 (cooled)
1014  if(ibound==2)
1015  {
1016  nod_pt->set_value(0,-0.5);
1017 
1018  //Add small concentration imperfection if desired
1019  double epsilon = 0.01;
1020 
1021  //Read out the x position
1022  double x = nod_pt->x(0);
1023 
1024  //Set a sinusoidal perturbation in the concentration
1025  double value = sin(2.0*MathematicalConstants::Pi*x/1.5)*
1026  epsilon*time*exp(-10.0*time);
1027  nod_pt->set_value(0, -0.5 + value);
1028  }
1029 
1030  //If we are on the bottom boundary, set the concentration to be high
1031  //to 0.5 (heated)
1032  if(ibound==0) {nod_pt->set_value(0,0.5);}
1033  }
1034  }
1035 
1036 
1037 } // end_of_set_boundary_conditions
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
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
double Pi
Definition: two_d_biharmonic.cc:235
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 exp(const bfloat16 &a)
Definition: BFloat16.h:615
squared absolute value
Definition: GlobalFunctions.h:87
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::SarahBL::epsilon, Eigen::bfloat16_impl::exp(), j, oomph::MathematicalConstants::Pi, oomph::Data::set_value(), sin(), Eigen::value, plotDoE::x, and oomph::Node::x().

◆ temp_mesh_pt()

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<AD_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::temp_mesh_pt ( )
inline

Access function to the Advection-Diffusion mesh.

704  {
705  return dynamic_cast<RectangularQuadMesh<AD_ELEMENT>*>(Temp_mesh_pt);
706  }

Member Data Documentation

◆ Conc_mesh_pt

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<AD_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::Conc_mesh_pt
protected

◆ Doc_info

template<class NST_ELEMENT , class AD_ELEMENT >
DocInfo DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::Doc_info
private

DocInfo object.

◆ Nst_mesh_pt

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<NST_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::Nst_mesh_pt
protected

◆ Temp_mesh_pt

template<class NST_ELEMENT , class AD_ELEMENT >
RectangularQuadMesh<AD_ELEMENT>* DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::Temp_mesh_pt
protected

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