oomph::FpPressureAdvectionDiffusionProblem< ELEMENT > Class Template Reference

#include <navier_stokes_preconditioners.h>

+ Inheritance diagram for oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >:

Public Member Functions

 FpPressureAdvectionDiffusionProblem (Mesh *navier_stokes_mesh_pt, Problem *orig_problem_pt)
 Constructor: Pass Navier-Stokes mesh and pointer to orig problem. More...
 
void get_pressure_advection_diffusion_jacobian (CRDoubleMatrix &jacobian)
 Get the pressure advection diffusion Jacobian. More...
 
void reset_pin_status ()
 Reset pin status of all values. More...
 
void pin_all_non_pressure_dofs (const bool &set_pressure_bc_for_validation=false)
 
void validate (DocInfo &doc_info)
 
void doc_solution (DocInfo &doc_info)
 
- 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

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

Detailed Description

template<class ELEMENT>
class oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >

//////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// Auxiliary Problem that can be used to assemble the pressure advection diffusion matrix needed by the FpPreconditoner

Constructor & Destructor Documentation

◆ FpPressureAdvectionDiffusionProblem()

template<class ELEMENT >
oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::FpPressureAdvectionDiffusionProblem ( Mesh navier_stokes_mesh_pt,
Problem orig_problem_pt 
)
inline

Constructor: Pass Navier-Stokes mesh and pointer to orig problem.

160  {
161  // Pass across mesh -- boundary conditions have already
162  // been applied and the equations have been enumerated.
163  mesh_pt() = navier_stokes_mesh_pt;
164 
165  // Store pointer to orig problem so we can re-assign the
166  // orig eqn numbers when we're done
167  Orig_problem_pt = orig_problem_pt;
168 
169  // Get the spatial dimension
170  Ndim = mesh_pt()->finite_element_pt(0)->dim();
171 
172  // Set new assembly handler
173  this->assembly_handler_pt() = new FpPreconditionerAssemblyHandler(Ndim);
174  }
unsigned dim() const
Definition: elements.h:2611
unsigned Ndim
The spatial dimension.
Definition: navier_stokes_preconditioners.h:529
Problem * Orig_problem_pt
Definition: navier_stokes_preconditioners.h:533
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
AssemblyHandler *& assembly_handler_pt()
Return a pointer to the assembly handler object.
Definition: problem.h:1570
Mesh *& mesh_pt()
Return a pointer to the global mesh.
Definition: problem.h:1280

References oomph::Problem::assembly_handler_pt(), oomph::FiniteElement::dim(), oomph::Mesh::finite_element_pt(), oomph::Problem::mesh_pt(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Ndim, and oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Orig_problem_pt.

Member Function Documentation

◆ doc_solution()

template<class ELEMENT >
void oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::doc_solution ( DocInfo doc_info)
inline

Doc solution (only needed during development – kept alive in case validation is required during code maintenance)

479  {
480  std::ofstream some_file;
481  std::ostringstream filename;
482 
483  // Number of plot points
484  unsigned npts;
485  npts = 5;
486 
487  // Check value of FE solution in first element
488  Vector<double> s(Ndim, 0.0);
489  Vector<double> x(Ndim);
490  double p = dynamic_cast<ELEMENT*>(mesh_pt()->finite_element_pt(0))
491  ->interpolated_p_nst(s);
492  dynamic_cast<ELEMENT*>(mesh_pt()->finite_element_pt(0))
493  ->interpolated_x(s, x);
494 
495  // Get offset-free exact solution
496  double u_exact = 0;
498 
499  // Adjust offset
500  unsigned nnode = mesh_pt()->nnode();
501  for (unsigned j = 0; j < nnode; j++)
502  {
503  if (mesh_pt()->node_pt(j)->nvalue() == 3)
504  {
505  *(mesh_pt()->node_pt(j)->value_pt(2)) -= p - u_exact;
506  }
507  }
508 
509  // Output solution
510  filename << doc_info.directory() << "/fp_soln" << doc_info.number()
511  << ".dat";
512  some_file.open(filename.str().c_str());
513  some_file.precision(20);
514  mesh_pt()->output(some_file, npts);
515  some_file.close();
516 
517  filename.str("");
518  filename << doc_info.directory() << "/fp_exact_soln" << doc_info.number()
519  << ".dat";
520  some_file.open(filename.str().c_str());
521  some_file.precision(20);
522  mesh_pt()->output_fct(
524  some_file.close();
525  }
float * p
Definition: Tutorial_Map_using.cpp:9
double * value_pt(const unsigned &i) const
Definition: nodes.h:324
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
Output a given Vector function at f(n_plot) points in each element.
Definition: mesh.cc:2199
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
void output(std::ostream &outfile)
Output for all elements.
Definition: mesh.cc:2027
Node *& node_pt(const unsigned long &n)
Return pointer to global node n.
Definition: mesh.h:436
RealScalar s
Definition: level1_cplx_impl.h:130
string filename
Definition: MergeRestartFiles.py:39
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
Definition: navier_stokes_preconditioners.cc:63
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::DocInfo::directory(), MergeRestartFiles::filename, oomph::Mesh::finite_element_pt(), oomph::PressureAdvectionDiffusionValidation::get_exact_u(), j, oomph::Problem::mesh_pt(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Ndim, oomph::Mesh::nnode(), oomph::Mesh::node_pt(), oomph::DocInfo::number(), oomph::Mesh::output(), oomph::Mesh::output_fct(), p, s, oomph::Data::value_pt(), and plotDoE::x.

Referenced by oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::validate().

◆ get_pressure_advection_diffusion_jacobian()

template<class ELEMENT >
void oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::get_pressure_advection_diffusion_jacobian ( CRDoubleMatrix jacobian)
inline

Get the pressure advection diffusion Jacobian.

179  {
180  // Pin all non-pressure dofs
182 
183  // Re-assign the equation numbers for all elements in Navier-Stokes mesh
184  unsigned n_dof = assign_eqn_numbers();
185  oomph_info << "Eqn numbers after pinning veloc dofs: " << n_dof
186  << std::endl;
187 
188  // Get "Jacobian" of the modified system
189  DoubleVector dummy_residuals;
190  this->get_jacobian(dummy_residuals, jacobian);
191 
192  // Reset pin status
194 
195  // Reassign equation numbering of original problem
196  oomph_info << "Eqn numbers in orig problem after re-assignment: "
197  << Orig_problem_pt->assign_eqn_numbers() << std::endl;
198  }
std::vector< double > DoubleVector
loads clump configuration
Definition: ClumpInput.h:26
void pin_all_non_pressure_dofs(const bool &set_pressure_bc_for_validation=false)
Definition: navier_stokes_preconditioners.h:257
void reset_pin_status()
Reset pin status of all values.
Definition: navier_stokes_preconditioners.h:202
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Definition: problem.cc:1989
virtual void get_jacobian(DoubleVector &residuals, DenseDoubleMatrix &jacobian)
Definition: problem.cc:3890
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References oomph::Problem::assign_eqn_numbers(), oomph::Problem::get_jacobian(), oomph::oomph_info, oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Orig_problem_pt, oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), and oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::reset_pin_status().

◆ pin_all_non_pressure_dofs()

template<class ELEMENT >
void oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs ( const bool set_pressure_bc_for_validation = false)
inline

Pin all non-pressure dofs and (if boolean flag is set to true also all pressure dofs along domain boundaries – only used for validation)

259  {
260  // Backup pin status and then pin all non-pressure degrees of freedom
261  unsigned nelem = mesh_pt()->nelement();
262  for (unsigned e = 0; e < nelem; e++)
263  {
264  // Get actual element (needed because different elements
265  // store nodal pressure in different places)
266  ELEMENT* el_pt =
267  dynamic_cast<ELEMENT*>(mesh_pt()->finite_element_pt(e));
268 
269 
270 #ifdef PARANOID
271  if (el_pt == 0)
272  {
273  std::ostringstream error_message;
274  error_message << "Navier Stokes mesh must contain only Navier Stokes"
275  << "bulk elements\n";
276  throw OomphLibError(error_message.str(),
279  }
280 #endif
281 
282  // Check if element has internal pressure representation
283  // usually discontinuous -- preconditioner doesn't work for that case!
284  if (el_pt->p_nodal_index_nst() < 0)
285  {
286  std::ostringstream error_message;
287  error_message << "Cannot use Fp preconditioner with discontinuous\n"
288  << "pressures.\n";
289  throw OomphLibError(error_message.str(),
292  }
293 
294  // Loop over internal data and pin the values (having established that
295  // pressure dofs aren't amongst those)
296  unsigned nint = el_pt->ninternal_data();
297  for (unsigned j = 0; j < nint; j++)
298  {
299  Data* data_pt = el_pt->internal_data_pt(j);
300  if (Eqn_number_backup[data_pt].size() == 0)
301  {
302  unsigned nvalue = data_pt->nvalue();
303  Eqn_number_backup[data_pt].resize(nvalue);
304  for (unsigned i = 0; i < nvalue; i++)
305  {
306  // Backup
307  Eqn_number_backup[data_pt][i] = data_pt->eqn_number(i);
308 
309  // Pin everything
310  data_pt->pin(i);
311  }
312  }
313  }
314 
315  // Now deal with nodal values
316  unsigned nnod = el_pt->nnode();
317  for (unsigned j = 0; j < nnod; j++)
318  {
319  Node* nod_pt = el_pt->node_pt(j);
320  if (Eqn_number_backup[nod_pt].size() == 0)
321  {
322  unsigned nvalue = nod_pt->nvalue();
323  Eqn_number_backup[nod_pt].resize(nvalue);
324  for (unsigned i = 0; i < nvalue; i++)
325  {
326  // Pin everything apart from the nodal pressure
327  // value
328  if (int(i) != el_pt->p_nodal_index_nst())
329  {
330  // Backup and pin
331  Eqn_number_backup[nod_pt][i] = nod_pt->eqn_number(i);
332  nod_pt->pin(i);
333  }
334  // Else it's a pressure value
335  else
336  {
337  // Exclude non-nodal pressure based elements
338  if (el_pt->p_nodal_index_nst() >= 0)
339  {
340  // Backup
341  Eqn_number_backup[nod_pt][i] = nod_pt->eqn_number(i);
342  }
343  }
344  }
345  }
346 
347  // Set wind
348  if (set_pressure_bc_for_validation)
349  {
350  Vector<double> x(2);
351  x[0] = nod_pt->x(0);
352  x[1] = nod_pt->x(1);
353  Vector<double> u(2);
355  nod_pt->set_value(el_pt->u_index_nst(0), u[0]);
356  nod_pt->set_value(el_pt->u_index_nst(1), u[1]);
357  }
358  }
359  }
360  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
std::map< Data *, std::vector< int > > Eqn_number_backup
Map to store original equation numbers.
Definition: navier_stokes_preconditioners.h:536
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
void wind_function(const Vector< double > &x, Vector< double > &wind)
Wind.
Definition: navier_stokes_preconditioners.cc:48
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References e(), oomph::Data::eqn_number(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Eqn_number_backup, oomph::Mesh::finite_element_pt(), i, j, oomph::Problem::mesh_pt(), oomph::Mesh::nelement(), oomph::Data::nvalue(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Data::pin(), oomph::Data::set_value(), size, oomph::PressureAdvectionDiffusionValidation::wind_function(), plotDoE::x, and oomph::Node::x().

Referenced by oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::get_pressure_advection_diffusion_jacobian(), and oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::validate().

◆ reset_pin_status()

template<class ELEMENT >
void oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::reset_pin_status ( )
inline

Reset pin status of all values.

203  {
204  // Reset pin status for nodes
205  unsigned nnod = mesh_pt()->nnode();
206  for (unsigned j = 0; j < nnod; j++)
207  {
208  Node* nod_pt = mesh_pt()->node_pt(j);
209  unsigned nval = nod_pt->nvalue();
210  for (unsigned i = 0; i < nval; i++)
211  {
212  nod_pt->eqn_number(i) = Eqn_number_backup[nod_pt][i];
213  }
214  }
215 
216  // ... and internal data
217  unsigned nelem = mesh_pt()->nelement();
218  for (unsigned e = 0; e < nelem; e++)
219  {
220  // Get actual element
221  ELEMENT* el_pt =
222  dynamic_cast<ELEMENT*>(mesh_pt()->finite_element_pt(e));
223 
224 
225 #ifdef PARANOID
226  if (el_pt == 0)
227  {
228  std::ostringstream error_message;
229  error_message << "Navier Stokes mesh must contain only Navier Stokes"
230  << "bulk elements\n";
231  throw OomphLibError(error_message.str(),
234  }
235 #endif
236 
237  unsigned nint = el_pt->ninternal_data();
238  for (unsigned j = 0; j < nint; j++)
239  {
240  Data* data_pt = el_pt->internal_data_pt(j);
241  unsigned nvalue = data_pt->nvalue();
242  for (unsigned i = 0; i < nvalue; i++)
243  {
244  data_pt->eqn_number(i) = Eqn_number_backup[data_pt][i];
245  }
246  }
247  }
248 
249  // Free up storage
250  Eqn_number_backup.clear();
251  }
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483

References e(), oomph::Data::eqn_number(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Eqn_number_backup, oomph::Mesh::finite_element_pt(), i, j, oomph::Problem::mesh_pt(), oomph::Mesh::nelement(), oomph::Mesh::nnode(), oomph::Mesh::node_pt(), oomph::Data::nvalue(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::get_pressure_advection_diffusion_jacobian(), and oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::validate().

◆ validate()

template<class ELEMENT >
void oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::validate ( DocInfo doc_info)
inline

Validate pressure advection diffusion problem and doc exact and computed solution in directory specified in DocInfo object

366  {
367  oomph_info << "\n\n==============================================\n\n";
368  oomph_info << "Doing validation for pressure adv diff problem\n";
369 
370  // Choose exact solution
372 
373  // Pin all non-pressure dofs and pressure dofs along boundaries for
374  // validation
375  bool set_pressure_bc_for_validation = true;
376  pin_all_non_pressure_dofs(set_pressure_bc_for_validation);
377 
378 
379  // Set Peclet number
381  dynamic_cast<NavierStokesEquations<2>*>(mesh_pt()->element_pt(0))->re();
382 
383  // Loop over all elements and set source function pointer
384  unsigned nel = mesh_pt()->nelement();
385  for (unsigned e = 0; e < nel; e++)
386  {
387  dynamic_cast<NavierStokesEquations<2>*>(mesh_pt()->element_pt(e))
388  ->source_fct_for_pressure_adv_diff() =
390  }
391 
392  // Re-assign the equation numbers for all elements in Navier-Stokes mesh
393  oomph_info << "Eqn numbers after pinning veloc dofs: "
394  << assign_eqn_numbers() << std::endl;
395 
396  // Attach Robin BC elements
397  unsigned nbound = mesh_pt()->nboundary();
399  {
400  // Loop over all boundaries of Navier Stokes mesh
401  for (unsigned b = 0; b < nbound; b++)
402  {
403  // How many bulk elements are adjacent to boundary b?
404  unsigned n_element = mesh_pt()->nboundary_element(b);
405 
406  // Loop over the bulk elements adjacent to boundary b?
407  for (unsigned e = 0; e < n_element; e++)
408  {
409  TemplateFreeNavierStokesEquationsBase* bulk_elem_pt =
410  dynamic_cast<TemplateFreeNavierStokesEquationsBase*>(
412 
413  // What is the index of the face of the bulk element e on bondary b
414  int face_index = mesh_pt()->face_index_at_boundary(b, e);
415 
416  // Build face element
417  bulk_elem_pt->build_fp_press_adv_diff_robin_bc_element(face_index);
418 
419  } // end of loop over bulk elements adjacent to boundary b
420  }
421  }
422 
423 
424  // Loop over all elements and set source function pointer
425  std::ofstream outfile;
426  outfile.open("robin_elements.dat");
427  for (unsigned e = 0; e < nel; e++)
428  {
429  dynamic_cast<NavierStokesEquations<2>*>(mesh_pt()->element_pt(e))
430  ->output_pressure_advection_diffusion_robin_elements(outfile);
431  }
432  outfile.close();
433 
434  // Solve it
435  newton_solve();
436 
437  // And output the solution...
438  doc_solution(doc_info);
439 
440  // Kill Robin BC elements
442  {
443  // Loop over all boundaries of Navier Stokes mesh
444  for (unsigned b = 0; b < nbound; b++)
445  {
446  // How many bulk elements are adjacent to boundary b?
447  unsigned n_element = mesh_pt()->nboundary_element(b);
448 
449  // Loop over the bulk elements adjacent to boundary b?
450  for (unsigned e = 0; e < n_element; e++)
451  {
452  TemplateFreeNavierStokesEquationsBase* bulk_elem_pt =
453  dynamic_cast<TemplateFreeNavierStokesEquationsBase*>(
455 
456  // Kill associated Robin elements
457  bulk_elem_pt->delete_pressure_advection_diffusion_robin_elements();
458 
459  } // end of loop over bulk elements adjacent to boundary b
460  }
461  }
462 
463  // Reset pin status
465 
466  // Reassign equation numbering of original problem
467  oomph_info << "Eqn numbers in orig problem after re-assignment: "
468  << Orig_problem_pt->assign_eqn_numbers() << std::endl;
469 
470 
471  oomph_info << "Done validation for pressure adv diff problem\n";
472  oomph_info << "\n\n==============================================\n\n";
473  }
Scalar * b
Definition: benchVecAdd.cpp:17
void doc_solution(DocInfo &doc_info)
Definition: navier_stokes_preconditioners.h:478
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
unsigned nboundary() const
Return number of boundaries.
Definition: mesh.h:827
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
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
void newton_solve()
Use Newton method to solve the problem.
Definition: problem.cc:8783
double source_function(const Vector< double > &x_vect)
Source function required to make the solution above an exact solution.
Definition: navier_stokes_preconditioners.cc:93
unsigned Flag
Flag for solution.
Definition: navier_stokes_preconditioners.cc:42
double Peclet
Peclet number – overwrite with actual Reynolds number.
Definition: navier_stokes_preconditioners.cc:45

References oomph::Problem::assign_eqn_numbers(), b, oomph::Mesh::boundary_element_pt(), oomph::TemplateFreeNavierStokesEquationsBase::build_fp_press_adv_diff_robin_bc_element(), oomph::TemplateFreeNavierStokesEquationsBase::delete_pressure_advection_diffusion_robin_elements(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::doc_solution(), e(), oomph::Mesh::element_pt(), oomph::Mesh::face_index_at_boundary(), oomph::PressureAdvectionDiffusionValidation::Flag, oomph::Problem::mesh_pt(), oomph::Mesh::nboundary(), oomph::Mesh::nboundary_element(), oomph::Mesh::nelement(), oomph::Problem::newton_solve(), oomph::oomph_info, oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Orig_problem_pt, oomph::PressureAdvectionDiffusionValidation::Peclet, oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::reset_pin_status(), and oomph::PressureAdvectionDiffusionValidation::source_function().

Referenced by oomph::NavierStokesSchurComplementPreconditioner::validate().

Member Data Documentation

◆ Eqn_number_backup

template<class ELEMENT >
std::map<Data*, std::vector<int> > oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::Eqn_number_backup
private

◆ Ndim

◆ Orig_problem_pt


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