////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// More...
#include <segregated_fsi_solver.h>
Public Types | |
enum | convergence_criteria { Assess_convergence_based_on_absolute_solid_change , Assess_convergence_based_on_relative_solid_change , Assess_convergence_based_on_max_global_residual } |
Enumerated flags for convergence criteria. More... | |
enum | solve_type { Full_solve , Fluid_solve , Solid_solve } |
Enumerated flags to indicate which solve is taking place. More... | |
Public Types inherited from oomph::Problem | |
enum | Distributed_problem_matrix_distribution { Default_matrix_distribution , Problem_matrix_distribution , Uniform_matrix_distribution } |
enum for distribution of distributed jacobians. 1 - Automatic - the Problem distribution is employed, unless any processor has number of rows equal to 110% of N/P, in which case uniform distribution is employed. 2 - Problem - the jacobian on processor p only contains rows that correspond to equations that are on this processor. (minimises communication) 3 - Uniform - each processor holds as close to N/P matrix rows as possible. (very well load balanced) More... | |
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 Member Functions | |
SegregatableFSIProblem () | |
Constructor. Set default values for solver parameters: More... | |
virtual | ~SegregatableFSIProblem () |
Empty virtual destructor. More... | |
virtual void | identify_fluid_and_solid_dofs (Vector< Data * > &fluid_data_pt, Vector< Data * > &solid_data_pt, Mesh *&fluid_mesh_pt, Mesh *&solid_mesh_pt)=0 |
Identify the fluid and solid Data. This is a pure virtual function that MUST be implemented for every specific problem that is to be solved by the segregated solver. The two mesh pointers identify meshes that contain elements and nodes used during the fluid or solid solves respectively. Elements that feature during both phases of the segretated solution must be included in both. These pointers may be set to NULL. In this case, all elements in the global mesh (set up during the monolithic discretisation of the problem) contribute to the global Jacobian matrix and the residual vector, even if their contributions only contain zero entries. This can be costly, though the code will still compute the correct results. More... | |
void | setup_segregated_solver (const bool &full_setup_of_fluid_and_solid_dofs=true) |
Setup the segregated solver: Backup the pinned status of the fluid and solid dofs and allocate the internal storage based on the input provided by identify_fluid_and_solid_dofs(...) In addition, reset storage associated with convergence acceleration techniques. If the problem and degrees of freedom has not changed between calls to the solver then it is wasteful to call identify_fluid_and_solid_dofs(...) again and again. If the optional boolean flag is set to false then the storage for convergence acceleration techniques is reset, but the fluid and solid dofs are not altered. More... | |
PicardConvergenceData | segregated_solve () |
Segregated solver. Peform a segregated step from the present state of the system. Returns PicardConvergenceData object that contains the vital stats of the iteration. More... | |
PicardConvergenceData | steady_segregated_solve () |
Steady version of segregated solver. Makes all timesteppers steady before solving. Returns PicardConvergenceData object that contains the vital stats of the iteration. More... | |
PicardConvergenceData | unsteady_segregated_solve (const double &dt) |
Unsteady segregated solver, advance time by dt and solve by the segregated solver. The time values are always shifted by this function. Returns PicardConvergenceData object that contains the vital stats of the iteration. More... | |
PicardConvergenceData | unsteady_segregated_solve (const double &dt, const bool &shift_values) |
Unsteady segregated solver. Advance time by dt and solve the system by a segregated method. The boolean flag is used to control whether the time values should be shifted. If it is true the current data values will be shifted (stored as previous timesteps) before the solution step. Returns PicardConvergenceData object that contains the vital stats of the iteration. More... | |
void | assess_convergence_based_on_max_global_residual (const double &tol) |
Assess convergence based on max. residual of coupled system of eqns. The argument specifies the convergence tolerance. More... | |
void | assess_convergence_based_on_max_global_residual () |
Assess convergence based on max. residuals of coupled system of eqns. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used. More... | |
void | assess_convergence_based_on_absolute_solid_change (const double &tol) |
Assess convergence based on max. absolute change of solid dofs. The argument specifies the convergence tolerance. More... | |
void | assess_convergence_based_on_absolute_solid_change () |
Assess convergence based on max. absolute change of solid dofs. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used. More... | |
void | assess_convergence_based_on_relative_solid_change (const double &tol) |
Assess convergence based on max. relative change of solid dofs. The argument specifies the convergence tolerance. More... | |
void | assess_convergence_based_on_relative_solid_change () |
Assess convergence based on max. relative change of solid dofs. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used. More... | |
void | enable_pointwise_aitken (const unsigned &pointwise_aitken_start) |
Use pointwise Aitken extrapolation. The argument is used to specify the Picard iteration after which pointwise Aitken extrapolation is to be used for the first time. More... | |
void | enable_pointwise_aitken () |
Use pointwise Aitken extrapolation. This interface has no argument and the current value of Pointwise_aitken_start will be used. The default is zero, extrapolation starts immediately. More... | |
void | disable_pointwise_aitken () |
Disable the use of pointwise Aitken extrapolation. More... | |
void | enable_under_relaxation (const double &omega=1.0) |
Use under-relaxation and (optionally) specify under-relaxation parameter. Default: omega=1.0 (i.e. no actual under-relaxation; Other extreme: omega=0.0 (freeze wall shape). Under-relaxation parameter can also be computed dynamically by setting use_irons_and_tuck_extrapolation() More... | |
void | enable_irons_and_tuck_extrapolation () |
Use Irons and Tuck extrapolation for solid dofs. More... | |
void | disable_irons_and_tuck_extrapolation () |
Do not use Irons and Tuck extrapolation for solid dofs. More... | |
void | get_solid_change (double &rms_change, double &max_change, double &rms_norm) |
Get rms of change in the solid dofs; the max. change of the solid dofs and the rms norm of the solid dofs themselves. Change is computed relative to the reference values stored when store_solid_dofs() was last called. More... | |
void | store_solid_dofs () |
Store the current solid values as reference values for future convergence check. Also add another entry to pointwise Aitken history if required. More... | |
void | reset_timer () |
Reset timer. More... | |
void | restart_timer () |
(Re-)start timer (e.g. after completing non-essential parts of the code such as documentation of the iteration's progress) More... | |
void | halt_timer () |
Halt timer (e.g. before performing non-essential parts of the code such as documentation of the iteration's progress) More... | |
double | t_spent_on_actual_solve () |
Total elapsed time since start of solve. More... | |
Public Member Functions inherited from oomph::Problem | |
virtual void | debug_hook_fct (const unsigned &i) |
Hook for debugging. Can be overloaded in driver code; argument allows identification of where we're coming from. More... | |
void | set_analytic_dparameter (double *const ¶meter_pt) |
Function to turn on analytic calculation of the parameter derivatives in continuation and bifurcation detection problems. More... | |
void | unset_analytic_dparameter (double *const ¶meter_pt) |
Function to turn off analytic calculation of the parameter derivatives in continuation and bifurcation detection problems. More... | |
bool | is_dparameter_calculated_analytically (double *const ¶meter_pt) |
Function to determine whether the parameter derivatives are calculated analytically. More... | |
void | set_analytic_hessian_products () |
Function to turn on analytic calculation of the parameter derivatives in continuation and bifurcation detection problems. More... | |
void | unset_analytic_hessian_products () |
Function to turn off analytic calculation of the parameter derivatives in continuation and bifurcation detection problems. More... | |
bool | are_hessian_products_calculated_analytically () |
Function to determine whether the hessian products are calculated analytically. More... | |
void | set_pinned_values_to_zero () |
Set all pinned values to zero. Used to set boundary conditions to be homogeneous in the copy of the problem used in adaptive bifurcation tracking (ALH: TEMPORARY HACK, WILL BE FIXED) More... | |
bool | distributed () const |
If we have MPI return the "problem has been distributed" flag, otherwise it can't be distributed so return false. More... | |
Distributed_problem_matrix_distribution & | distributed_problem_matrix_distribution () |
accesss function to the distributed matrix distribution method 1 - Automatic - the Problem distribution is employed, unless any processor has number of rows equal to 110% of N/P, in which case uniform distribution is employed. 2 - Problem - the jacobian on processor p only contains rows that correspond to equations that are on this processor. (minimises communication) 3 - Uniform - each processor holds as close to N/P matrix rows as possible. (very well load balanced) More... | |
void | enable_doc_imbalance_in_parallel_assembly () |
Enable doc of load imbalance in parallel assembly of distributed problem. More... | |
void | disable_doc_imbalance_in_parallel_assembly () |
Disable doc of load imbalance in parallel assembly of distributed problem. More... | |
Vector< double > | elemental_assembly_time () |
Return vector of most-recent elemental assembly times (used for load balancing). Zero sized if no Jacobian has been computed since last re-assignment of equation numbers. More... | |
void | clear_elemental_assembly_time () |
Clear storage of most-recent elemental assembly times (used for load balancing). Next load balancing operation will be based on the number of elements associated with a tree root. More... | |
void | enable_problem_distributed () |
Enable problem distributed. More... | |
void | disable_problem_distributed () |
Disable problem distributed. More... | |
void | set_default_first_and_last_element_for_assembly () |
Set default first and last elements for parallel assembly of non-distributed problem. More... | |
void | set_first_and_last_element_for_assembly (Vector< unsigned > &first_el_for_assembly, Vector< unsigned > &last_el_for_assembly) |
Manually set first and last elements for parallel assembly of non-distributed problem. More... | |
void | get_first_and_last_element_for_assembly (Vector< unsigned > &first_el_for_assembly, Vector< unsigned > &last_el_for_assembly) const |
Get first and last elements for parallel assembly of non-distributed problem. More... | |
virtual void | actions_before_adapt () |
Actions that are to be performed before a mesh adaptation. These might include removing any additional elements, such as traction boundary elements before the adaptation. More... | |
virtual void | actions_after_adapt () |
Actions that are to be performed after a mesh adaptation. More... | |
OomphCommunicator * | communicator_pt () |
access function to the oomph-lib communicator More... | |
const OomphCommunicator * | communicator_pt () const |
access function to the oomph-lib communicator, const version More... | |
Problem () | |
Constructor: Allocate space for one time stepper and set all pointers to NULL and set defaults for all parameters. More... | |
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) |
Return a pointer to the i-th submesh. If there are no submeshes, the 0-th submesh is the global mesh itself. More... | |
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) |
Add a submesh to the problem and return its number, i, by which it can be accessed via mesh_pt(i). More... | |
void | flush_sub_meshes () |
Flush the problem's collection of sub-meshes. Must be followed by call to rebuild_global_mesh(). More... | |
void | build_global_mesh () |
Build the global mesh by combining the all the submeshes. Note: The nodes boundary information refers to the boundary numbers within the submesh! More... | |
void | rebuild_global_mesh () |
If one of the submeshes has changed (e.g. by mesh adaptation) we need to update the global mesh. Note: The nodes boundary information refers to the boundary numbers within the submesh! More... | |
virtual void | build_mesh () |
Function to build the Problem's base mesh; this must be supplied by the user if they wish to use the load_balance() functionality, which is only available to problems that have already been distributed. If the problem has multiple meshes, each mesh must be built, added as as a submesh, and a call to build_global_mesh() must be made in this function. On return from this function all meshes must have been refined to the same level that they were in the when Problem::distribute() was first called. More... | |
void | load_balance () |
Balance the load of a (possibly non-uniformly refined) problem that has already been distributed, by re-distributing elements over processors. More... | |
void | load_balance (const bool &report_stats) |
Balance the load of a (possibly non-uniformly refined) problem that has already been distributed, by re-distributing elements over processors. Produce explicit stats of load balancing process if boolean, report_stats, is set to true. More... | |
void | load_balance (DocInfo &doc_info, const bool &report_stats) |
Balance the load of a (possibly non-uniformly refined) problem that has already been distributed, by re-distributing elements over processors. Produce explicit stats of load balancing process if boolean, report_stats, is set to true. More... | |
void | load_balance (DocInfo &doc_info, const bool &report_stats, const Vector< unsigned > &input_target_domain_for_local_non_halo_element) |
Balance the load of a (possibly non-uniformly refined) problem that has already been distributed, by re-distributing elements over processors. Produce explicit stats of load balancing process if boolean, report_stats, is set to true and doc various bits of data (mainly for debugging) in directory specified by DocInfo object. If final input vector is non-zero-sized it provides an imposed partitioning. More... | |
void | set_default_partition_in_load_balance () |
Set the use of the default partition in the load balance. More... | |
void | unset_default_partition_in_load_balance () |
Do not use the default partition in the load balance. More... | |
void | refine_distributed_base_mesh (Vector< Vector< Vector< unsigned >>> &to_be_refined_on_each_root, const unsigned &max_level_overall) |
Load balance helper routine: refine each new base (sub)mesh based upon the elements to be refined within each tree at each root on the current processor. More... | |
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 () |
Return a pointer to the linear solver object used for explicit time stepping. More... | |
LinearSolver * | mass_matrix_solver_for_explicit_timestepper_pt () const |
Return a pointer to the linear solver object used for explicit time stepping (const version). More... | |
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... | |
Time * | time_pt () const |
Return a pointer to the global time object (const version). More... | |
double & | time () |
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 () |
Access function for the pointer to the first (presumably only) timestepper. More... | |
const TimeStepper * | time_stepper_pt () const |
Access function for the pointer to the first (presumably only) timestepper. More... | |
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) |
Set all problem data to have the same timestepper (timestepper_pt) Return the new number of dofs in the problem. More... | |
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... | |
double & | minimum_dt () |
Access function to min timestep in adaptive timestepping. More... | |
double & | maximum_dt () |
Access function to max timestep in adaptive timestepping. More... | |
unsigned & | max_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... | |
double & | max_residuals () |
Access function to max residuals in Newton iterations before giving up. More... | |
bool & | time_adaptive_newton_crash_on_solve_fail () |
Access function for Time_adaptive_newton_crash_on_solve_fail. More... | |
double & | newton_solver_tolerance () |
Access function to tolererance of the Newton solver, i.e. the maximum value of the residuals that will be accepted. More... | |
void | add_time_stepper_pt (TimeStepper *const &time_stepper_pt) |
Add a timestepper to the problem. The function will automatically create or resize the Time object so that it contains the appropriate number of levels of storage. More... | |
void | set_explicit_time_stepper_pt (ExplicitTimeStepper *const &explicit_time_stepper_pt) |
Set the explicit timestepper for the problem. The function will automatically create or resize the Time object so that it contains the appropriate number of levels of storage. More... | |
void | initialise_dt (const double &dt) |
Set all timesteps to the same value, dt, and assign weights for all timesteppers in the problem. More... | |
void | initialise_dt (const Vector< double > &dt) |
Set the value of the timesteps to be equal to the values passed in a vector and assign weights for all timesteppers in the problem. More... | |
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) |
Add Data to the Problem's global data – the Problem will perform equation numbering etc. for such Data. More... | |
void | flush_global_data () |
Flush the Problem's global data – resizes container to zero. Data objects are not deleted! More... | |
void | create_new_linear_algebra_distribution (LinearAlgebraDistribution *&dist_pt) |
Get new linear algebra distribution (you're in charge of deleting it!) More... | |
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 () |
Insist that local dof pointers are set up in each element when equation numbering takes place. More... | |
void | disable_store_local_dof_pt_in_elements () |
Insist that local dof pointers are NOT set up in each element when equation numbering takes place (the default) More... | |
unsigned long | assign_eqn_numbers (const bool &assign_local_eqn_numbers=true) |
Assign all equation numbers for problem: Deals with global data (= data that isn't attached to any elements) and then does the equation numbering for the elements. Virtual so it can be overloaded in MPI problems. Bool argument can be set to false to ignore assigning local equation numbers (found to be necessary in the parallel implementation of locate_zeta between multiple meshes). More... | |
void | describe_dofs (std::ostream &out= *(oomph_info.stream_pt())) const |
Function to describe the dofs in terms of the global equation number, i.e. what type of value (nodal value of a Node; value in a Data object; value of internal Data in an element; etc) is the unknown with a certain global equation number. Output stream defaults to oomph_info. More... | |
void | enable_discontinuous_formulation () |
Indicate that the problem involves discontinuous elements This allows for a more efficiently assembly and inversion of the mass matrix. More... | |
void | disable_discontinuous_formulation () |
Disable the use of a discontinuous-element formulation. Note that the methods will all still work even if the elements are discontinuous, we will just be assembling a larger system matrix than necessary. More... | |
void | get_dofs (DoubleVector &dofs) const |
Return the vector of dofs, i.e. a vector containing the current values of all unknowns. More... | |
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) |
Set history values of dofs from the type of vector stored in problem::Dof_pt. More... | |
void | add_to_dofs (const double &lambda, const DoubleVector &increment_dofs) |
Add lambda x incremenet_dofs[l] to the l-th dof. More... | |
double * | global_dof_pt (const unsigned &i) |
Return a pointer to the dof, indexed by global equation number which may be haloed or stored locally. If it is haloed, a local copy must have been requested on this processor in the Halo_scheme_pt. More... | |
double & | dof (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... | |
double * | dof_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) |
Return the residual vector multiplied by the inverse mass matrix Virtual so that it can be overloaded for mpi problems. More... | |
virtual void | get_dvaluesdt (DoubleVector &f) |
Get the time derivative of all values (using get_inverse_mass_matrix_times_residuals(..) with all time steppers set to steady) e.g. for use in explicit time steps. The approach used is slighty hacky, beware if you have a residual which is non-linear or implicit in the derivative or if you have overloaded get_jacobian(...). More... | |
virtual void | get_residuals (DoubleVector &residuals) |
Return the fully-assembled residuals Vector for the problem: Virtual so it can be overloaded in for mpi problems. More... | |
virtual void | get_jacobian (DoubleVector &residuals, DenseDoubleMatrix &jacobian) |
Return the fully-assembled Jacobian and residuals for the problem Interface for the case when the Jacobian matrix is dense. This is virtual so, if we feel like it (e.g. for testing iterative solvers with specific test matrices, we can bypass the default assembly procedure for the Jacobian and the residual vector. More... | |
virtual void | get_jacobian (DoubleVector &residuals, CRDoubleMatrix &jacobian) |
Return the fully-assembled Jacobian and residuals for the problem. Interface for the case when the Jacobian is in row-compressed storage format. This is virtual so, if we feel like it (e.g. for testing iterative solvers with specific test matrices), we can bypass the default assembly procedure for the Jacobian and the residual vector. More... | |
virtual void | get_jacobian (DoubleVector &residuals, CCDoubleMatrix &jacobian) |
Return the fully-assembled Jacobian and residuals for the problem. Interface for the case when the Jacobian is in column-compressed storage format. This is virtual so, if we feel like it (e.g. for testing iterative solvers with specific test matrices), we can bypass the default assembly procedure for the Jacobian and the residual vector. More... | |
virtual void | get_jacobian (DoubleVector &residuals, SumOfMatrices &jacobian) |
Dummy virtual function that must be overloaded by the problem to specify which matrices should be summed to give the final Jacobian. More... | |
void | get_fd_jacobian (DoubleVector &residuals, DenseMatrix< double > &jacobian) |
Return the fully-assembled Jacobian and residuals, generated by finite differences. More... | |
void | get_derivative_wrt_global_parameter (double *const ¶meter_pt, DoubleVector &result) |
Get the derivative of the entire residuals vector wrt a global parameter, used in continuation problems. More... | |
void | get_hessian_vector_products (DoubleVectorWithHaloEntries const &Y, Vector< DoubleVectorWithHaloEntries > const &C, Vector< DoubleVectorWithHaloEntries > &product) |
Return the product of the global hessian (derivative of Jacobian matrix with respect to all variables) with an eigenvector, Y, and any number of other specified vectors C (d(J_{ij})/d u_{k}) Y_{j} C_{k}. This function is used in assembling and solving the augmented systems associated with bifurcation tracking. The default implementation is to use finite differences at the global level. More... | |
void | solve_eigenproblem_legacy (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &steady=true) |
Get derivative of an element in the problem wrt a global parameter, used in continuation problems. More... | |
void | solve_eigenproblem_legacy (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true) |
Solve an eigenproblem as assembled by the Problem's constituent elements. Calculate (at least) n_eval eigenvalues. The boolean flag (default true) specifies whether the steady jacobian should be assembled. If the flag is false then the weighted mass-matrix terms from the timestepper will be included in the jacobian — this is almost certainly never wanted. Legacy version. More... | |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &alpha, Vector< double > &beta, Vector< DoubleVector > &eigenvector_real, Vector< DoubleVector > &eigenvector_imag, const bool &steady=true) |
Solve an eigenproblem as assembled by the Problem's constituent elements. Calculate (at least) n_eval eigenvalues and return the corresponding eigenvectors. The boolean flag (default true) specifies whether the steady jacobian should be assembled. If the flag is false then the weighted mass-matrix terms from the timestepper will be included in the jacobian — this is almost certainly never wanted. The eigenvalues and eigenvectors are, in general, complex. Eigenvalues may be infinite and are therefore returned as where is complex while is real. The actual eigenvalues may then be computed by doing the division, checking for zero betas to avoid NaNs. There's a convenience wrapper to this function that simply computes these eigenvalues regardless. That version may die in NaN checking is enabled (via the fenv.h header and the associated feenable function). More... | |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector_real, Vector< DoubleVector > &eigenvector_imag, const bool &steady=true) |
Solve an eigenproblem as assembled by the Problem's constituent elements. Calculate (at least) n_eval eigenvalues and return the corresponding eigenvectors. The boolean flag (default true) specifies whether the steady jacobian should be assembled. If the flag is false then the weighted mass-matrix terms from the timestepper will be included in the jacobian — this is almost certainly never wanted. Note that the eigenvalues and eigenvectors are, in general, complex and the eigenvalues may be infinite. In this case it's safer to use the other version of this function which returns the eigenvalues in terms of a fractional representation. More... | |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &make_timesteppers_steady=true) |
Solve an eigenproblem as assembled by the Problem's constituent elements but only return the eigenvalues, not the eigenvectors. At least n_eval eigenvalues are computed. The boolean flag (default true) is used to specify whether the weighted mass-matrix terms from the timestepping scheme should be included in the jacobian — this is almost certainly never wanted. Note that the eigenvalues may be infinite. In this case it's safer to use the other version of this function which returns the eigenvalues in terms of a fractional representation. More... | |
void | solve_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &alpha, Vector< double > &beta, const bool &steady=true) |
Solve an eigenproblem as assembled by the Problem's constituent elements but only return the eigenvalues, not the eigenvectors. At least n_eval eigenvalues are computed. The boolean flag (default true) is used to specify whether the weighted mass-matrix terms from the timestepping scheme should be included in the jacobian — this is almost certainly never wanted. Note that the eigenvalues may be infinite and are therefore returned as where is complex while is real. The actual eigenvalues may then be computed by doing the division, checking for zero betas to avoid NaNs. More... | |
void | solve_adjoint_eigenproblem_legacy (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &make_timesteppers_steady=true) |
Solve an adjoint eigenvalue problem using the same procedure as solve_eigenproblem. See the documentation on that function for more details. Note: this is a legacy version of this function that stores re & imag parts of eigenvectors in some solver-specific collection of real vectors. More... | |
void | solve_adjoint_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector_real, Vector< DoubleVector > &eigenvector_imag, const bool &steady=true) |
Solve an adjoint eigenvalue problem using the same procedure as solve_eigenproblem. See the documentation on that function for more details. More... | |
void | solve_adjoint_eigenproblem (const unsigned &n_eval, Vector< std::complex< double >> &eigenvalue, const bool &steady=true) |
Solve an adjoint eigenvalue problem using the same procedure as solve_eigenproblem but only return the eigenvalues, not the eigenvectors. At least n_eval eigenvalues are computed. See the documentation on that function for more details. More... | |
virtual void | get_eigenproblem_matrices (CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0) |
Get the matrices required by a eigensolver. If the shift parameter is non-zero the second matrix will be shifted. More... | |
void | assign_eigenvector_to_dofs (DoubleVector &eigenvector) |
Assign the eigenvector passed to the function to the dofs in the problem so that it can be output by the usual routines. More... | |
void | add_eigenvector_to_dofs (const double &epsilon, const DoubleVector &eigenvector) |
Add the eigenvector passed to the function scaled by the constat epsilon to the dofs in the problem so that it can be output by the usual routines. More... | |
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 () |
Enable recycling of Jacobian in Newton iteration (if the linear solver allows it). Useful for linear problems with constant Jacobians or nonlinear problems where you're willing to risk the trade-off between faster solve times and degraded Newton convergence rate. More... | |
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... | |
bool & | use_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) |
Adaptive Newton solve: up to max_adapt adaptations of all refineable submeshes are performed to achieve the the error targets specified in the refineable submeshes. More... | |
void | steady_newton_solve (unsigned const &max_adapt=0) |
Solve a steady problem using adaptive Newton's method, but in the context of an overall unstady problem, perhaps to determine an initial condition. This is achieved by setting the weights in the timesteppers to be zero which has the effect of rendering them steady timesteppers. The optional argument max_adapt specifies the max. number of adaptations of all refineable submeshes are performed to achieve the the error targets specified in the refineable submeshes. More... | |
void | copy (Problem *orig_problem_pt) |
Copy Data values, nodal positions etc from specified problem. Note: This is not a copy constructor. We assume that the current and the "original" problem have both been created by calling the same problem constructor so that all Data objects, time steppers etc. in the two problems are completely independent. This function copies the nodal, internal and global values, and the time parameters from the original problem into "this" one. This functionality is required, e.g. for multigrid computations. More... | |
virtual Problem * | make_copy () |
Make and return a pointer to the copy of the problem. A virtual function that must be filled in by the user is they wish to perform adaptive refinement in bifurcation tracking or in multigrid problems. ALH: WILL NOT BE NECESSARY IN BIFURCATION TRACKING IN LONG RUN... More... | |
virtual void | read (std::ifstream &restart_file, bool &unsteady_restart) |
Read refinement pattern of all refineable meshes and refine them accordingly, then read all Data and nodal position info from file for restart. Return flag to indicate if the restart was from steady or unsteady solution. More... | |
virtual void | read (std::ifstream &restart_file) |
Read refinement pattern of all refineable meshes and refine them accordingly, then read all Data and nodal position info from file for restart. More... | |
virtual void | dump (std::ofstream &dump_file) const |
Dump refinement pattern of all refineable meshes and all generic Problem data to file for restart. More... | |
void | dump (const std::string &dump_file_name) const |
Dump refinement pattern of all refineable meshes and all generic Problem data to file for restart. More... | |
void | get_all_halo_data (std::map< unsigned, double * > &map_of_halo_data) |
Get pointers to all possible halo data indexed by global equation number in a map. More... | |
long | synchronise_eqn_numbers (const bool &assign_local_eqn_numbers=true) |
Classify any non-classified nodes into halo/haloed and synchronise equation numbers. Return the total number of degrees of freedom in the overall problem. More... | |
void | synchronise_dofs (const bool &do_halos, const bool &do_external_halos) |
Synchronise the degrees of freedom by overwriting the haloed values with their non-halo counterparts held on other processors. Bools control if we deal with data associated with external halo/ed elements/nodes or the "normal" halo/ed ones. More... | |
void | synchronise_all_dofs () |
Perform all required synchronisation in solvers. More... | |
void | check_halo_schemes (DocInfo &doc_info) |
Check the halo/haloed node/element schemes. More... | |
void | check_halo_schemes () |
Check the halo/haloed node/element schemes. More... | |
Vector< unsigned > | distribute (const Vector< unsigned > &element_partition, DocInfo &doc_info, const bool &report_stats=false) |
Distribute the problem and doc, using the specified partition; returns a vector which details the partitioning. More... | |
Vector< unsigned > | distribute (DocInfo &doc_info, const bool &report_stats=false) |
Distribute the problem; returns a vector which details the partitioning. More... | |
Vector< unsigned > | distribute (const Vector< unsigned > &element_partition, const bool &report_stats=false) |
Distribute the problem using the specified partition; returns a vector which details the partitioning. More... | |
Vector< unsigned > | distribute (const bool &report_stats=false) |
Distribute the problem; returns a vector which details the partitioning. More... | |
virtual void | partition_global_mesh (Mesh *&global_mesh_pt, DocInfo &doc_info, Vector< unsigned > &element_domain, const bool &report_stats=false) |
Partition the global mesh, return vector specifying the processor number for each element. Virtual so that it can be overloaded by any user; the default is to use METIS to perform the partitioning (with a bit of cleaning up afterwards to sort out "special cases"). More... | |
void | prune_halo_elements_and_nodes (DocInfo &doc_info, const bool &report_stats) |
(Irreversibly) prune halo(ed) elements and nodes, usually after another round of refinement, to get rid of excessively wide halo layers. Note that the current mesh will be now regarded as the base mesh and no unrefinement relative to it will be possible once this function has been called. More... | |
void | prune_halo_elements_and_nodes (const bool &report_stats=false) |
(Irreversibly) prune halo(ed) elements and nodes, usually after another round of refinement, to get rid of excessively wide halo layers. Note that the current mesh will be now regarded as the base mesh and no unrefinement relative to it will be possible once this function has been called. More... | |
bool | problem_has_been_distributed () |
Access to Problem_has_been_distributed flag. More... | |
void | delete_all_external_storage () |
Wrapper function to delete external storage for each submesh of the problem. More... | |
virtual void | symmetrise_eigenfunction_for_adaptive_pitchfork_tracking () |
Virtual function that is used to symmetrise the problem so that the current solution exactly satisfies any symmetries within the system. Used when adpativly solving pitchfork detection problems when small asymmetries in the coarse solution can be magnified leading to very inaccurate answers on the fine mesh. This is always problem-specific and must be filled in by the user The default issues a warning. More... | |
double * | bifurcation_parameter_pt () const |
Return pointer to the parameter that is used in the bifurcation detection. If we are not tracking a bifurcation then an error will be thrown by the AssemblyHandler. More... | |
void | get_bifurcation_eigenfunction (Vector< DoubleVector > &eigenfunction) |
Return the eigenfunction calculated as part of a bifurcation tracking process. If we are not tracking a bifurcation then an error will be thrown by the AssemblyHandler. More... | |
void | activate_fold_tracking (double *const ¶meter_pt, const bool &block_solve=true) |
Turn on fold tracking using the augmented system specified in the FoldHandler class. After a call to this function subsequent calls of the standard solution methods will converge to a fold (limit) point at a particular value of the variable addressed by parameter_pt. The system may not converge if the initial guess is sufficiently poor or, alternatively, if finite differencing is used to calculate the jacobian matrix in the elements. If the boolean flag block_solver is true (the default) then a block factorisation is used to solve the augmented system which is both faster and uses less memory. More... | |
void | activate_bifurcation_tracking (double *const ¶meter_pt, const DoubleVector &eigenvector, const bool &block_solve=true) |
Activate generic bifurcation tracking for a single (real) eigenvalue where the initial guess for the eigenvector can be specified. More... | |
void | activate_bifurcation_tracking (double *const ¶meter_pt, const DoubleVector &eigenvector, const DoubleVector &normalisation, const bool &block_solve=true) |
Activate generic bifurcation tracking for a single (real) eigenvalue where the initial guess for the eigenvector can be specified and the normalisation condition can also be specified. More... | |
void | activate_pitchfork_tracking (double *const ¶meter_pt, const DoubleVector &symmetry_vector, const bool &block_solve=true) |
Turn on pitchfork tracking using the augmented system specified in the PitchForkHandler class. After a call to this function subsequent calls of the standard solution methods will converge to a pitchfork bifurcation at a particular value of the variable addressed by parameter_pt. The symmetry that is to be broken must be specified by supplying a symmetry_vector(ndof). The easiest way to determine such a vector is to solve the associated eigenproblem and pass in the eigenvector. This is not always necessary however, if the symmetry is easy to construct. The system may not converge if the initial guess is sufficiently poor or, alternatively, if finite differencing is used to calculate the jacobian matrix in the elements. If the boolean flag block_solver is true (the default) then a block factorisation is used to solve the augmented system which is both faster and requires less memory. More... | |
void | activate_hopf_tracking (double *const ¶meter_pt, const bool &block_solve=true) |
Turn on Hopf bifurcation tracking using the augmented system specified in the HopfHandler class. After a call to this function subsequent calls of the standard solution methods will converge to a Hopf bifuraction at a particular value of the variable addressed by parameter_pt. The system may not converge if the initial guess is sufficiently poor or, alternatively, if finite differencing is used to calculate the jacobian matrix in the elements. More... | |
void | activate_hopf_tracking (double *const ¶meter_pt, const double &omega, const DoubleVector &null_real, const DoubleVector &null_imag, const bool &block_solve=true) |
Turn on Hopf bifurcation tracking using the augmented system specified in the HopfHandler class. After a call to this function subsequent calls of the standard solution methods will converge to a Hopf bifuraction at a particular value of the variable addressed by parameter_pt. The system may not converge if the initial guess is sufficiently poor or, alternatively, if finite differencing is used to calculate the jacobian matrix in the elements. This interface allows specification of an inital guess for the frequency and real and imaginary parts of the null vector, such as might be obtained from an eigensolve. More... | |
void | deactivate_bifurcation_tracking () |
Deactivate all bifuraction tracking, by reseting the assembly handler to the default. More... | |
void | reset_assembly_handler_to_default () |
Reset the system to the standard non-augemented state. More... | |
double | arc_length_step_solve (double *const ¶meter_pt, const double &ds, const unsigned &max_adapt=0) |
Solve a steady problem using arc-length continuation, when the parameter that becomes a variable corresponding to the arc-length constraint equation is an external double: parameter_pt is a pointer to that double, ds is the desired arc_length and max_adapt is the maximum number of spatial adaptations (default zero, no adaptation). More... | |
double | arc_length_step_solve (Data *const &data_pt, const unsigned &data_index, const double &ds, const unsigned &max_adapt=0) |
Solve a steady problem using arc-length continuation, when the variable corresponding to the arc-length constraint equation is already stored in data used in the problem: data_pt is a pointer to the appropriate data object, data_index is the index of the value that will be traded for the constriant, ds is the desired arc_length and max_adapt is the maximum number of spatial adaptations (default zero, no adaptation). Note that the value must be pinned in order for this formulation to work. More... | |
void | reset_arc_length_parameters () |
Reset the "internal" arc-length continuation parameters, so as to allow continuation in another parameter. N.B. The parameters that are reset are the "minimum" that are required, others should perhaps be reset, depending upon the application. More... | |
int & | sign_of_jacobian () |
Access function for the sign of the global jacobian matrix. This will be set by the linear solver, if possible (direct solver). If not alternative methods must be used to detect bifurcations (solving the associated eigenproblem). More... | |
void | explicit_timestep (const double &dt, const bool &shift_values=true) |
Take an explicit timestep of size dt and optionally shift any stored values of the time history. More... | |
void | unsteady_newton_solve (const double &dt) |
Advance time by dt and solve by Newton's method. This version always shifts time values. More... | |
void | unsteady_newton_solve (const double &dt, const bool &shift_values) |
Advance time by dt and solve the system, using Newton's method. The boolean flag is used to control whether the time values should be shifted. If it is true the current data values will be shifted (copied to the locations where there are stored as previous timesteps) before solution. More... | |
void | unsteady_newton_solve (const double &dt, const unsigned &max_adapt, const bool &first, const bool &shift=true) |
Unsteady adaptive Newton solve: up to max_adapt adaptations of all refineable submeshes are performed to achieve the the error targets specified in the refineable submeshes. If first==true, the initial conditions are re-assigned after the mesh adaptations. Shifting of time can be suppressed by overwriting the default value of shift (true). [Shifting must be done if first_timestep==true because we're constantly re-assigning the initial conditions; if first_timestep==true and shift==false shifting is performed anyway and a warning is issued. More... | |
double | doubly_adaptive_unsteady_newton_solve (const double &dt, const double &epsilon, const unsigned &max_adapt, const bool &first, const bool &shift=true) |
Unsteady "doubly" adaptive Newton solve: Does temporal adaptation first, i.e. we try to do a timestep with an increment of dt, and adjusting dt until the solution on the given mesh satisfies the temporal error measure with tolerance epsilon. Following this, we do up to max_adapt spatial adaptions (without re-examining the temporal error). If first==true, the initial conditions are re-assigned after the mesh adaptations. Shifting of time can be suppressed by overwriting the default value of shift (true). [Shifting must be done if first_timestep==true because we're constantly re-assigning the initial conditions; if first_timestep==true and shift==false shifting is performed anyway and a warning is issued. More... | |
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) |
Unsteady "doubly" adaptive Newton solve: Does temporal adaptation first, i.e. we try to do a timestep with an increment of dt, and adjusting dt until the solution on the given mesh satisfies the temporal error measure with tolerance epsilon. Following this, we do up to max_adapt spatial adaptions (without re-examining the temporal error). If first==true, the initial conditions are re-assigned after the mesh adaptations. Shifting of time can be suppressed by overwriting the default value of shift (true). [Shifting must be done if first_timestep==true because we're constantly re-assigning the initial conditions; if first_timestep==true and shift==false shifting is performed anyway and a warning is issued. Pseudo-Boolean flag suppress_resolve_after_spatial_adapt [0: false; 1: true] does what it says.]. More... | |
double | adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon) |
Attempt to advance timestep by dt_desired. If the solution fails the timestep will be halved until convergence is achieved, or the timestep falls below NewtonSolverParameters::Minimum_time_step. The error control parameter epsilon represents the (approximate) desired magnitude of the global error at each timestep. The routine returns a double that is the suggested next timestep and should be passed as dt_desired the next time the routine is called. This version always shifts the time values. More... | |
double | adaptive_unsteady_newton_solve (const double &dt_desired, const double &epsilon, const bool &shift_values) |
Attempt to advance timestep by dt_desired. If the solution fails the timestep will be halved until convergence is achieved, or the timestep falls below NewtonSolverParameters::Minimum_time_step. The error control parameter epsilon represents the (approximate) desired magnitude of the global error at each timestep. The routine returns a double that is the suggested next timestep and should be passed as dt_desired the next time the routine is called. Once again the boolean flag, shift_values, is used to control whether the time values are shifted. More... | |
void | assign_initial_values_impulsive () |
Initialise data and nodal positions to simulate impulsive start from initial configuration/solution. More... | |
void | assign_initial_values_impulsive (const double &dt) |
Initialise data and nodal positions to simulate an impulsive start and also set the initial and previous values of dt. More... | |
void | calculate_predictions () |
Calculate predictions. More... | |
void | enable_mass_matrix_reuse () |
Enable recycling of the mass matrix in explicit timestepping schemes. Useful for timestepping on fixed meshes when you want to avoid the linear solve phase. More... | |
void | disable_mass_matrix_reuse () |
Turn off recyling of the mass matrix in explicit timestepping schemes. More... | |
bool | mass_matrix_reuse_is_enabled () |
Return whether the mass matrix is being reused. More... | |
void | refine_uniformly (const Vector< unsigned > &nrefine_for_mesh) |
Refine refineable sub-meshes, each as many times as specified in the vector and rebuild problem. More... | |
void | refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
Refine refineable sub-meshes, each as many times as specified in the vector and rebuild problem; doc refinement process. More... | |
void | refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh) |
Refine refineable sub-meshes, each as many times as specified in the vector and rebuild problem. Prune after refinements. More... | |
void | refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
Refine refineable sub-meshes, each as many times as specified in the vector and rebuild problem; doc refinement process. More... | |
void | refine_uniformly (DocInfo &doc_info) |
Refine (all) refineable (sub)mesh(es) uniformly and rebuild problem; doc refinement process. More... | |
void | refine_uniformly_and_prune (DocInfo &doc_info) |
Refine (all) refineable (sub)mesh(es) uniformly and rebuild problem; doc refinement process. More... | |
void | refine_uniformly () |
Refine (all) refineable (sub)mesh(es) uniformly and rebuild problem. More... | |
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) |
p-refine p-refineable sub-meshes, each as many times as specified in the vector and rebuild problem More... | |
void | p_refine_uniformly (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
p-refine p-refineable sub-meshes, each as many times as specified in the vector and rebuild problem; doc refinement process More... | |
void | p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh) |
p-refine p-refineable sub-meshes, each as many times as specified in the vector and rebuild problem. Prune after refinements More... | |
void | p_refine_uniformly_and_prune (const Vector< unsigned > &nrefine_for_mesh, DocInfo &doc_info) |
p-refine p-refineable sub-meshes, each as many times as specified in the vector and rebuild problem; doc refinement process More... | |
void | p_refine_uniformly (DocInfo &doc_info) |
p-refine (all) p-refineable (sub)mesh(es) uniformly and rebuild problem; doc refinement process. More... | |
void | p_refine_uniformly_and_prune (DocInfo &doc_info) |
p-refine (all) p-refineable (sub)mesh(es) uniformly and rebuild problem; doc refinement process. More... | |
void | p_refine_uniformly () |
p-refine (all) p-refineable (sub)mesh(es) uniformly and rebuild problem More... | |
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) |
Refine (one and only!) mesh by splitting the elements identified by their numbers relative to the problems' only mesh, then rebuild the problem. More... | |
void | refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined_pt) |
Refine (one and only!) mesh by splitting the elements identified by their pointers, then rebuild the problem. More... | |
void | refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined) |
Refine specified submesh by splitting the elements identified by their numbers relative to the submesh, then rebuild the problem. More... | |
void | refine_selected_elements (const unsigned &i_mesh, const Vector< RefineableElement * > &elements_to_be_refined_pt) |
Refine specified submesh by splitting the elements identified by their pointers, then rebuild the problem. More... | |
void | refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined) |
Refine all submeshes by splitting the elements identified by their numbers relative to each submesh in a Vector of Vectors, then rebuild the problem. More... | |
void | refine_selected_elements (const Vector< Vector< RefineableElement * >> &elements_to_be_refined_pt) |
Refine all submeshes by splitting the elements identified by their pointers within each submesh in a Vector of Vectors, then rebuild the problem. More... | |
void | p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined) |
p-refine (one and only!) mesh by refining the elements identified by their numbers relative to the problems' only mesh, then rebuild the problem. More... | |
void | p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt) |
p-refine (one and only!) mesh by refining the elements identified by their pointers, then rebuild the problem. More... | |
void | p_refine_selected_elements (const unsigned &i_mesh, const Vector< unsigned > &elements_to_be_refined) |
p-refine specified submesh by refining the elements identified by their numbers relative to the submesh, then rebuild the problem. More... | |
void | p_refine_selected_elements (const unsigned &i_mesh, const Vector< PRefineableElement * > &elements_to_be_refined_pt) |
p-refine specified submesh by refining the elements identified by their pointers, then rebuild the problem. More... | |
void | p_refine_selected_elements (const Vector< Vector< unsigned >> &elements_to_be_refined) |
p-refine all submeshes by refining the elements identified by their numbers relative to each submesh in a Vector of Vectors, then rebuild the problem. More... | |
void | p_refine_selected_elements (const Vector< Vector< PRefineableElement * >> &elements_to_be_refined_pt) |
p-refine all submeshes by refining the elements identified by their pointers within each submesh in a Vector of Vectors, then rebuild the problem. More... | |
unsigned | unrefine_uniformly () |
Refine (all) refineable (sub)mesh(es) uniformly and rebuild problem. Return 0 for success, 1 for failure (if unrefinement has reached the coarsest permitted level) More... | |
unsigned | unrefine_uniformly (const unsigned &i_mesh) |
Do uniform refinement for submesh i_mesh without documentation. Return 0 for success, 1 for failure (if unrefinement has reached the coarsest permitted level) More... | |
void | p_unrefine_uniformly (DocInfo &doc_info) |
p-unrefine (all) p-refineable (sub)mesh(es) uniformly and rebuild problem. More... | |
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) |
Adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on their own error estimates and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. Return # of refined/unrefined elements. On return from this function, Problem can immediately be solved again. More... | |
void | adapt () |
Adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on their own error estimates and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. On return from this function, Problem can immediately be solved again. [Argument-free wrapper]. More... | |
void | p_adapt (unsigned &n_refined, unsigned &n_unrefined) |
p-adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on their own error estimates and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. Return # of refined/unrefined elements. On return from this function, Problem can immediately be solved again. More... | |
void | p_adapt () |
p-adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on their own error estimates and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. On return from this function, Problem can immediately be solved again. [Argument-free wrapper] More... | |
void | adapt_based_on_error_estimates (unsigned &n_refined, unsigned &n_unrefined, Vector< Vector< double >> &elemental_error) |
Adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on the error estimates in elemental_error and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. Return # of refined/unrefined elements. On return from this function, Problem can immediately be solved again. More... | |
void | adapt_based_on_error_estimates (Vector< Vector< double >> &elemental_error) |
Adapt problem: Perform mesh adaptation for (all) refineable (sub)mesh(es), based on the error estimates in elemental_error and the target errors specified in the mesh(es). Following mesh adaptation, update global mesh, and re-assign equation numbers. Return # of refined/unrefined elements. On return from this function, Problem can immediately be solved again. [Wrapper without n_refined and n_unrefined arguments]. More... | |
void | get_all_error_estimates (Vector< Vector< double >> &elemental_error) |
Return the error estimates computed by (all) refineable (sub)mesh(es) in the elemental_error structure, which consists of a vector of vectors of elemental errors, one vector for each (sub)mesh. More... | |
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 () |
Enable the output of information when in the newton solver (Default) More... | |
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 () |
Empty virtual function to do anything needed before a stage of an explicit time step (Runge-Kutta steps contain multiple stages per time step, most others only contain one). More... | |
virtual void | actions_after_explicit_stage () |
Empty virtual function that should be overloaded to update any dependent data or boundary conditions that should be advanced after each stage of an explicit time step (Runge-Kutta steps contain multiple stages per time step, most others only contain one). More... | |
Protected Member Functions | |
virtual void | actions_before_segregated_solve () |
This function is called once at the start of each segregated solve. More... | |
virtual void | actions_after_segregated_solve () |
This function is called once at the end of each segregated solve. More... | |
virtual void | actions_before_segregated_convergence_check () |
This function is to be filled with actions that take place before the check for convergence of the entire segregated solve. More... | |
void | rebuild_monolithic_mesh () |
Rebuild global mesh for monolithic discretisation. More... | |
void | restore_fluid_dofs () |
Restore pinned status of fluid dofs. More... | |
void | pin_solid_dofs () |
Pin solid dofs. More... | |
void | restore_solid_dofs () |
Restore pinned status of solid dofs. More... | |
void | pointwise_aitken_extrapolate () |
Do pointwise Aitken extrapolation for solid. More... | |
Protected Member Functions inherited from oomph::Problem | |
unsigned | setup_element_count_per_dof () |
Function that populates the Element_counter_per_dof vector with the number of elements that contribute to each dof. For example, with linear elements in 1D each dof contains contributions from two elements apart from those on the boundary. Returns the total number of elements in the problem. More... | |
void | setup_dof_halo_scheme () |
Function that is used to setup the halo scheme. More... | |
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) |
Protected helper function that is used to assemble the Jacobian matrix in the case when the storage is row or column compressed. The boolean Flag indicates if we want compressed row format (true) or compressed column. More... | |
virtual void | actions_before_newton_solve () |
Any actions that are to be performed before a complete Newton solve (e.g. adjust boundary conditions). CAREFUL: This step should (and if the FD-based LinearSolver FD_LU is used, must) only update values that are pinned! More... | |
virtual void | actions_after_newton_solve () |
Any actions that are to be performed after a complete Newton solve, e.g. post processing. CAREFUL: This step should (and if the FD-based LinearSolver FD_LU is used, must) only update values that are pinned! More... | |
virtual void | actions_before_newton_convergence_check () |
Any actions that are to be performed before the residual is checked in the Newton method, e.g. update any boundary conditions that depend upon variables of the problem; update any ‘dependent’ variables; or perform an update of the nodal positions in SpineMeshes etc. CAREFUL: This step should (and if the FD-based LinearSolver FD_LU is used, must) only update values that are pinned! More... | |
virtual void | actions_before_newton_step () |
Any actions that are to be performed before each individual Newton step. Most likely to be used for diagnostic purposes to doc the solution during a non-converging iteration, say. More... | |
virtual void | actions_after_newton_step () |
Any actions that are to be performed after each individual Newton step. Most likely to be used for diagnostic purposes to doc the solution during a non-converging iteration, say. More... | |
virtual void | actions_before_implicit_timestep () |
Actions that should be performed before each implicit time step. This is needed when one wants to solve a steady problem before timestepping and needs to distinguish between the two cases. More... | |
virtual void | actions_after_implicit_timestep () |
Actions that should be performed after each implicit time step. This is needed when one wants to solve a steady problem before timestepping and needs to distinguish between the two cases. More... | |
virtual void | actions_after_implicit_timestep_and_error_estimation () |
Actions that should be performed after each implicit time step. This is needed if your actions_after_implicit_timestep() modify the solution in a way that affects the error estimate. More... | |
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 () |
Actions that are to be performed before reading in restart data for problems involving unstructured bulk meshes. If the problem contains such meshes we need to strip out any face elements that are attached to them because restart of unstructured meshes re-creates their elements and nodes from scratch, leading to dangling pointers from the face elements to the old elements and nodes. This function is virtual and (practically) empty but toggles a flag to indicate that it has been called. This is used to issue a warning, prompting the user to consider overloading it if the problem is found to contain unstructured bulk meshes during restarts. More... | |
virtual void | actions_after_read_unstructured_meshes () |
Actions that are to be performed before reading in restart data for problems involving unstructured bulk meshes. Typically used to re-attach FaceElements, say, that were stripped out in actions_before_read_unstructured_meshes(). This function is virtual and (practically) empty but toggles a flag to indicate that it has been called. This is used to issue a warning, prompting the user to consider overloading it if the problem is found to contain unstructured bulk meshes during restarts. More... | |
virtual void | actions_before_distribute () |
Actions to be performed before a (mesh) distribution. More... | |
virtual void | actions_after_distribute () |
Actions to be performed after a (mesh) distribution. More... | |
virtual void | actions_after_change_in_global_parameter (double *const ¶meter_pt) |
Actions that are to be performed when the global parameter addressed by parameter_pt has been changed in the function get_derivative_wrt_global_parameter() The default is to call actions_before_newton_solve(), actions_before_newton_convergence_check() and actions_after_newton_solve(). This could be amazingly inefficient in certain problems and should be overloaded in such cases. An example would be when a remesh is required in general, but the global parameter does not affect the mesh directly. More... | |
virtual void | actions_after_change_in_bifurcation_parameter () |
Actions that are to be performed after a change in the parameter that is being varied as part of the solution of a bifurcation detection problem. The default is to call actions_before_newton_solve(), actions_before_newton_convergence_check() and actions_after_newton_solve(). This could be amazingly inefficient in certain problems and should be overloaded in such cases. An example would be when a remesh is required in general, but the global parameter does not affect the mesh directly. More... | |
virtual void | actions_after_parameter_increase (double *const ¶meter_pt) |
Empty virtual function; provides hook to perform actions after the increase in the arclength parameter (during continuation) More... | |
double & | dof_derivative (const unsigned &i) |
Access function to the derivative of the i-th (local) dof with respect to the arc length, used in continuation. More... | |
double & | dof_current (const unsigned &i) |
Access function to the current value of the i-th (local) dof at the start of a continuation step. More... | |
virtual void | set_initial_condition () |
Set initial condition (incl previous timesteps). We need to establish this interface because I.C. needs to be reset when problem is adapted during the first timestep. More... | |
virtual double | global_temporal_error_norm () |
Function to calculate a global error norm, used in adaptive timestepping to control the change in timestep. Individual errors for each data object can be obtained via the data timestepper's temporal_error_in_value or temporal_error_in_position functions and should be combined to construct a global norm. For example, in fluids problems a suitable norm is usually the weighted sum of the errors in the velocities; for moving mesh problems is it usually better to use the weighted sum of the errors in position. More... | |
unsigned | newton_solve_continuation (double *const ¶meter_pt) |
Perform a basic arc-length continuation step using Newton's method. Returns number of Newton steps taken. More... | |
unsigned | newton_solve_continuation (double *const ¶meter_pt, DoubleVector &z) |
This function performs a basic continuation step using the Newton method. The number of Newton steps taken is returned, to be used in any external step-size control routines. The governing parameter of the problem is passed as a pointer to the routine, as is a vector in which to store the derivatives wrt the parameter, if required. More... | |
void | calculate_continuation_derivatives (double *const ¶meter_pt) |
A function to calculate the derivatives wrt the arc-length. This version of the function actually does a linear solve so that the derivatives are calculated "exactly" rather than using the values at the Newton step just before convergence. This is necessary in spatially adaptive problems, in which the number of degrees of freedom changes and so the appropriate derivatives must be calculated for the new variables. This function is called if no Newton steps were taken in the continuation routine ... i.e. the initial residuals were sufficiently small. More... | |
void | calculate_continuation_derivatives (const DoubleVector &z) |
A function to calculate the derivatives with respect to the arc-length required for continuation. The arguments is the solution of the linear system, Jz = dR/dparameter, that gives du/dparameter and the direction output from the newton_solve_continuation function. The derivatives are stored in the ContinuationParameters namespace. More... | |
void | calculate_continuation_derivatives_fd (double *const ¶meter_pt) |
A function to calculate the derivatives with respect to the arc-length required for continuation by finite differences, using the previous values of the solution. The derivatives are stored in the ContinuationParameters namespace. More... | |
bool | does_pointer_correspond_to_problem_data (double *const ¶meter_pt) |
Return a boolean flag to indicate whether the pointer parameter_pt refers to values stored in a Data object that is contained within the problem. More... | |
void | set_consistent_pinned_values_for_continuation () |
Private helper function that is used to set the appropriate pinned values for continuation. More... | |
Protected Attributes | |
int | Pointwise_aitken_counter |
Number of Aitken histories available (int because after extrapolation it's re-initialised to -1 to force the computation of three new genuine iterates). More... | |
bool | Use_pointwise_aitken |
Use pointwise Aitken extrapolation? More... | |
unsigned | Pointwise_aitken_start |
Start pointwise Aitken extrpolation after specified number of Picard iterations. More... | |
int | Solve_type |
Solve that is taking place (enumerated flag) More... | |
double | Convergence_tolerance |
Convergence tolerance for Picard iteration. More... | |
unsigned | Max_picard |
Max. number of Picard iterations. More... | |
bool | Doc_max_global_residual |
Doc maximum global residual during iteration? (default: false) More... | |
Vector< Data * > | Fluid_data_pt |
Vector storing the Data objects associated with the fluid problem: Tyically the nodal and internal data of the elements in the fluid bulk mesh. More... | |
Vector< std::vector< bool > > | Fluid_value_is_pinned |
Vector of vectors that store the pinned status of fluid Data values. More... | |
Vector< Data * > | Solid_data_pt |
Vector storing the Data objects associated with the solid problem: Typically the positional data of solid nodes and any quantities associated with displacement control, say. More... | |
Vector< std::vector< bool > > | Solid_value_is_pinned |
Vector of vectors that store the pinned status of solid Data values. More... | |
Vector< double > | Previous_solid_value |
Vector storing the previous solid values – used for convergence check. More... | |
Mesh * | Fluid_mesh_pt |
Mesh containing only fluid elements – the elements in this Mesh will be excluded from the assembly process when the solid problem is solved. More... | |
Mesh * | Solid_mesh_pt |
Mesh containing only solid elements – the elements in this mesh will be excluded from the assembly process when the fluid problem is solved. More... | |
Vector< Mesh * > | Orig_sub_mesh_pt |
Backup for the pointers to the submeshes in the original problem. More... | |
Vector< double > | Del_irons_and_tuck |
Vector of changes in Irons and Tuck under-relaxation. More... | |
double | R_irons_and_tuck |
Irons and Tuck relaxation factor. More... | |
Vector< Vector< double > > | Pointwise_aitken_solid_value |
Vector of Vectors containing up to three previous iterates for the solid dofs; used for pointwise Aitken extrapolation. More... | |
bool | Recheck_convergence_after_pointwise_aitken |
Have we just done a pointwise Aitken step. More... | |
Protected Attributes inherited from oomph::Problem | |
Vector< Problem * > | Copy_of_problem_pt |
Vector of pointers to copies of the problem used in adaptive bifurcation tracking problems (ALH: TEMPORARY HACK, WILL BE FIXED) More... | |
std::map< double *, bool > | Calculate_dparameter_analytic |
Map used to determine whether the derivatives with respect to a parameter should be finite differenced. The default is that finite differences should be used. More... | |
bool | Calculate_hessian_products_analytic |
Map used to determine whether the hessian products should be computed using finite differences. The default is that finite differences will be used. More... | |
LinearAlgebraDistribution * | Dof_distribution_pt |
The distribution of the DOFs in this problem. This object is created in the Problem constructor and setup when assign_eqn_numbers(...) is called. If this problem is distributed then this distribution will match the distribution of the equation numbers. If this problem is not distributed then this distribution will be uniform over all processors. More... | |
Vector< double * > | Dof_pt |
Vector of pointers to dofs. More... | |
DoubleVectorWithHaloEntries | Element_count_per_dof |
Counter that records how many elements contribute to each dof. Used to determine the (discrete) arc-length automatically. It really should be an integer, but is a double so that the distribution information can be used for distributed problems. More... | |
Vector< double > | Elemental_assembly_time |
Storage for assembly times (used for load balancing) More... | |
DoubleVectorHaloScheme * | Halo_scheme_pt |
Pointer to the halo scheme for any global vectors that have the Dof_distribution. More... | |
Vector< double * > | Halo_dof_pt |
Storage for the halo degrees of freedom (only required) when accessing via the global equation number in distributed problems. More... | |
double | Relaxation_factor |
Relaxation fator for Newton method (only a fractional Newton correction is applied if this is less than 1; default 1). More... | |
double | Newton_solver_tolerance |
The Tolerance below which the Newton Method is deemed to have converged. More... | |
unsigned | Max_newton_iterations |
Maximum number of Newton iterations. More... | |
unsigned | Nnewton_iter_taken |
Actual number of Newton iterations taken during the most recent iteration. More... | |
Vector< double > | Max_res |
Maximum residuals at start and after each newton iteration. More... | |
double | Max_residuals |
Maximum desired residual: if the maximum residual exceeds this value, the program will exit. More... | |
bool | Time_adaptive_newton_crash_on_solve_fail |
Bool to specify what to do if a Newton solve fails within a time adaptive solve. Default (false) is to half the step and try again. If true then crash instead. More... | |
bool | Jacobian_reuse_is_enabled |
Is re-use of Jacobian in Newton iteration enabled? Default: false. More... | |
bool | Jacobian_has_been_computed |
Has a Jacobian been computed (and can therefore be re-used if required)? Default: false. More... | |
bool | Problem_is_nonlinear |
Boolean flag indicating if we're dealing with a linear or nonlinear Problem – if set to false the Newton solver will not check the residual before or after the linear solve. Set to true by default; can be over-written in specific Problem class. More... | |
bool | Pause_at_end_of_sparse_assembly |
Protected boolean flag to halt program execution during sparse assemble process to assess peak memory usage. Initialised to false (obviously!) More... | |
bool | Doc_time_in_distribute |
Protected boolean flag to provide comprehensive timimings during problem distribution. Initialised to false. More... | |
unsigned | Sparse_assembly_method |
Flag to determine which sparse assembly method to use By default we use assembly by vectors of pairs. More... | |
unsigned | Sparse_assemble_with_arrays_initial_allocation |
the number of elements to initially allocate for a matrix row within the sparse_assembly_with_two_arrays(...) method (if a matrix of this size has not been assembled already). If a matrix of this size has been assembled then the number of elements in each row in that matrix is used as the initial allocation More... | |
unsigned | Sparse_assemble_with_arrays_allocation_increment |
the number of elements to add to a matrix row when the initial allocation is exceeded within the sparse_assemble_with_two_arrays(...) method. More... | |
Vector< Vector< unsigned > > | Sparse_assemble_with_arrays_previous_allocation |
the number of elements in each row of a compressed matrix in the previous matrix assembly. More... | |
double | Numerical_zero_for_sparse_assembly |
A tolerance used to determine whether the entry in a sparse matrix is zero. If it is then storage need not be allocated. More... | |
double | FD_step_used_in_get_hessian_vector_products |
bool | Mass_matrix_reuse_is_enabled |
Is re-use of the mass matrix in explicit timestepping enabled Default:false. More... | |
bool | Mass_matrix_has_been_computed |
Has the mass matrix been computed (and can therefore be reused) Default: false. More... | |
bool | Discontinuous_element_formulation |
Is the problem a discontinuous one, i.e. can the elemental contributions be treated independently. Default: false. More... | |
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 |
Maximum possible increase of dt between time-steps in adaptive schemes. More... | |
double | DTSF_min_decrease |
Minimum allowed decrease of dt between time-steps in adaptive schemes. Lower scaling values will reject the time-step (and retry with a smaller dt). More... | |
double | Minimum_dt_but_still_proceed |
If Minimum_dt_but_still_proceed positive, then dt will not be reduced below this value during adaptive timestepping and the computation will continue with this value, accepting the larger errors that this will incur). Note: This option is disabled by default as this value is initialised to -1.0. More... | |
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 |
Value of the scaling parameter required so that the parameter occupies the desired proportion of the arc length. NOTE: If you wish to change this then make sure to set the value of Scale_arc_length to false to ensure the value of this isn't overwritten during the arc-length process. Instead of changing this variable, it's better to actually update the Desired_proportion_of_arc_length value. More... | |
int | Sign_of_jacobian |
Storage for the sign of the global Jacobian. More... | |
double | Continuation_direction |
The direction of the change in parameter that will ensure that a branch is followed in one direction only. More... | |
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 |
Storage for the offset for the continuation derivatives from the original dofs (default is 1, but this will be differnet when continuing bifurcations and periodic orbits) More... | |
unsigned | Dof_current_offset |
Storage for the offset for the current values of dofs from the original dofs (default is 2, but this will be differnet when continuing bifurcations and periodic orbits) More... | |
Vector< double > | Dof_derivative |
Storage for the derivative of the problem variables wrt arc-length. More... | |
Vector< double > | Dof_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 |
The desired number of Newton Steps to reach convergence at each step along the arc. More... | |
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 |
Boolean to specify which scheme to use to calculate the continuation derivatievs. More... | |
bool | Problem_has_been_distributed |
Has the problem been distributed amongst multiple processors? More... | |
bool | Bypass_increase_in_dof_check_during_pruning |
Boolean to bypass check of increase in dofs during pruning. More... | |
OomphCommunicator * | Communicator_pt |
The communicator for this problem. More... | |
bool | Always_take_one_newton_step |
Boolean to indicate whether a Newton step should be taken even if the initial residuals are below the required tolerance. More... | |
double | Timestep_reduction_factor_after_nonconvergence |
What it says: If temporally adaptive Newton solver fails to to converge, reduce timestep by this factor and try again; defaults to 1/2; can be over-written by user in derived problem. More... | |
bool | Keep_temporal_error_below_tolerance |
Boolean to decide if a timestep is to be rejected if the error estimate post-solve (computed by global_temporal_error_norm()) exceeds the tolerance required in the call to adaptive_unsteady_newton_solve(...). Defaults to true. More... | |
Private Member Functions | |
void | extrapolate_solid_data () |
Extrapolate solid data and update fluid mesh during unsteady run. More... | |
void | under_relax_solid () |
Under-relax the most recently computed solid variables, either by classical relaxation or by Irons & Tuck. More... | |
void | use_only_fluid_elements () |
Only include fluid elements in the Problem's mesh. This is called before the segregated fluid solve. The fluid elements are identified by the user via the fluid_mesh_pt argument in the pure virtual function identify_fluid_and_solid_dofs(...). If a NULL pointer is returned by this function (i.e. if the user hasn't bothered to identify the fluids elements in a submesh, then no stripping of non-fluid elements is performed. The result of the computation will be correct but it won't be very efficient. More... | |
void | use_only_solid_elements () |
Only include solid elements in the Problem's mesh. This is called before the segregated solid solve. The solid elements are identified by the user via the solid_mesh_pt argument in the pure virtual function identify_fluid_and_solid_dofs(...). If a NULL pointer is returned by this function (i.e. if the user hasn't bothered to identify the solid elements in a submesh, then no stripping of non-solid elements is performed. The result of the computation will be correct but it won't be very efficient. More... | |
void | pin_fluid_dofs () |
Pin fluid dofs. More... | |
Private Attributes | |
double | Omega_relax |
Under-relaxation parameter. (1.0: no under-relaxation; 0.0: Freeze wall shape) More... | |
bool | Use_irons_and_tuck_extrapolation |
Boolean flag to indicate use of Irons and Tuck's extrapolation for solid values. More... | |
int | Convergence_criterion |
Convergence criterion (enumerated flag) More... | |
clock_t | T_ref |
Reference time for segregated solve. Can be re-initialised whenever total elapsed time has been stored (before entering non-essential doc sections of the code) More... | |
double | T_spent_on_actual_solve |
Total elapsed time since start of solve, can be accumulated by adding bits of time spent in relevant parts of code (bypassing sections that only document the progress) More... | |
bool | Timer_has_been_halted |
boolean flag to indicate if timer has been halted More... | |
Additional Inherited Members | |
Public Attributes inherited from oomph::Problem | |
double | Max_permitted_error_for_halo_check |
Threshold for error throwing in Problem::check_halo_schemes() More... | |
bool | Shut_up_in_newton_solve |
Boolean to indicate if all output is suppressed in Problem::newton_solve(). Defaults to false. More... | |
Static Public Attributes inherited from oomph::Problem | |
static bool | Suppress_warning_about_actions_before_read_unstructured_meshes |
Flag to allow suppression of warning messages re reading in unstructured meshes during restart. More... | |
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... | |
Static Protected Attributes inherited from oomph::Problem | |
static ContinuationStorageScheme | Continuation_time_stepper |
Storage for the single static continuation timestorage object. More... | |
////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////
Base class for problems that can be solved by segregated FSI solver
Definition at line 167 of file segregated_fsi_solver.h.
Enumerated flags for convergence criteria.
Enumerator | |
---|---|
Assess_convergence_based_on_absolute_solid_change | |
Assess_convergence_based_on_relative_solid_change | |
Assess_convergence_based_on_max_global_residual |
Definition at line 409 of file segregated_fsi_solver.h.
Enumerated flags to indicate which solve is taking place.
Enumerator | |
---|---|
Full_solve | |
Fluid_solve | |
Solid_solve |
Definition at line 418 of file segregated_fsi_solver.h.
|
inline |
Constructor. Set default values for solver parameters:
boolean flag to indicate if timer has been halted
Definition at line 191 of file segregated_fsi_solver.h.
References Assess_convergence_based_on_max_global_residual, Convergence_criterion, Convergence_tolerance, Doc_max_global_residual, Fluid_mesh_pt, Full_solve, Max_picard, oomph::Problem::Newton_solver_tolerance, Omega_relax, Pointwise_aitken_start, Recheck_convergence_after_pointwise_aitken, Solid_mesh_pt, Solve_type, T_ref, T_spent_on_actual_solve, Timer_has_been_halted, Use_irons_and_tuck_extrapolation, and Use_pointwise_aitken.
|
inlinevirtual |
Empty virtual destructor.
Definition at line 242 of file segregated_fsi_solver.h.
|
inlineprotectedvirtual |
This function is called once at the end of each segregated solve.
Definition at line 176 of file segregated_fsi_solver.h.
Referenced by segregated_solve().
|
inlineprotectedvirtual |
This function is to be filled with actions that take place before the check for convergence of the entire segregated solve.
Definition at line 180 of file segregated_fsi_solver.h.
Referenced by segregated_solve().
|
inlineprotectedvirtual |
This function is called once at the start of each segregated solve.
Definition at line 172 of file segregated_fsi_solver.h.
Referenced by segregated_solve().
|
inline |
Assess convergence based on max. absolute change of solid dofs. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used.
Definition at line 338 of file segregated_fsi_solver.h.
References oomph::Problem::Newton_solver_tolerance.
|
inline |
Assess convergence based on max. absolute change of solid dofs. The argument specifies the convergence tolerance.
Definition at line 328 of file segregated_fsi_solver.h.
References Assess_convergence_based_on_absolute_solid_change, Convergence_criterion, and Convergence_tolerance.
|
inline |
Assess convergence based on max. residuals of coupled system of eqns. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used.
Definition at line 320 of file segregated_fsi_solver.h.
References oomph::Problem::Newton_solver_tolerance.
|
inline |
Assess convergence based on max. residual of coupled system of eqns. The argument specifies the convergence tolerance.
Definition at line 310 of file segregated_fsi_solver.h.
References Assess_convergence_based_on_max_global_residual, Convergence_criterion, and Convergence_tolerance.
|
inline |
Assess convergence based on max. relative change of solid dofs. This interface has no argument and the default convergence tolerance for the Newton solver, Problem::Newton_solver_tolerance is used.
Definition at line 356 of file segregated_fsi_solver.h.
References oomph::Problem::Newton_solver_tolerance.
|
inline |
Assess convergence based on max. relative change of solid dofs. The argument specifies the convergence tolerance.
Definition at line 346 of file segregated_fsi_solver.h.
References Assess_convergence_based_on_relative_solid_change, Convergence_criterion, and Convergence_tolerance.
|
inline |
Do not use Irons and Tuck extrapolation for solid dofs.
Definition at line 403 of file segregated_fsi_solver.h.
References Use_irons_and_tuck_extrapolation.
|
inline |
Disable the use of pointwise Aitken extrapolation.
Definition at line 381 of file segregated_fsi_solver.h.
References Use_pointwise_aitken.
|
inline |
Use Irons and Tuck extrapolation for solid dofs.
Definition at line 397 of file segregated_fsi_solver.h.
References Use_irons_and_tuck_extrapolation.
|
inline |
Use pointwise Aitken extrapolation. This interface has no argument and the current value of Pointwise_aitken_start will be used. The default is zero, extrapolation starts immediately.
Definition at line 375 of file segregated_fsi_solver.h.
References Use_pointwise_aitken.
|
inline |
Use pointwise Aitken extrapolation. The argument is used to specify the Picard iteration after which pointwise Aitken extrapolation is to be used for the first time.
Definition at line 366 of file segregated_fsi_solver.h.
References Pointwise_aitken_start, and Use_pointwise_aitken.
|
inline |
Use under-relaxation and (optionally) specify under-relaxation parameter. Default: omega=1.0 (i.e. no actual under-relaxation; Other extreme: omega=0.0 (freeze wall shape). Under-relaxation parameter can also be computed dynamically by setting use_irons_and_tuck_extrapolation()
Definition at line 391 of file segregated_fsi_solver.h.
References Omega_relax.
|
private |
Extrapolate solid data and update fluid mesh during unsteady run.
Extrapolate solid data and update fluid mesh. Extrapolation is based on previous history values and is therefore only called in time-dependent runs.
Definition at line 39 of file segregated_fsi_solver.cc.
References Fluid_mesh_pt, i, oomph::Problem::mesh_pt(), oomph::Mesh::node_update(), and Solid_data_pt.
Referenced by unsteady_segregated_solve().
void oomph::SegregatableFSIProblem::get_solid_change | ( | double & | rms_change, |
double & | max_change, | ||
double & | rms_norm | ||
) |
Get rms of change in the solid dofs; the max. change of the solid dofs and the rms norm of the solid dofs themselves. Change is computed relative to the reference values stored when store_solid_dofs() was last called.
Definition at line 303 of file segregated_fsi_solver.cc.
References i, Previous_solid_value, and Solid_data_pt.
Referenced by segregated_solve().
|
inline |
Halt timer (e.g. before performing non-essential parts of the code such as documentation of the iteration's progress)
Definition at line 460 of file segregated_fsi_solver.h.
References T_ref, T_spent_on_actual_solve, and Timer_has_been_halted.
Referenced by t_spent_on_actual_solve().
|
pure virtual |
Identify the fluid and solid Data. This is a pure virtual function that MUST be implemented for every specific problem that is to be solved by the segregated solver. The two mesh pointers identify meshes that contain elements and nodes used during the fluid or solid solves respectively. Elements that feature during both phases of the segretated solution must be included in both. These pointers may be set to NULL. In this case, all elements in the global mesh (set up during the monolithic discretisation of the problem) contribute to the global Jacobian matrix and the residual vector, even if their contributions only contain zero entries. This can be costly, though the code will still compute the correct results.
Referenced by setup_segregated_solver().
|
private |
Pin fluid dofs.
Pin all fluid dofs.
Definition at line 153 of file segregated_fsi_solver.cc.
References Fluid_data_pt, and i.
Referenced by segregated_solve().
|
protected |
Pin solid dofs.
Pin all solid dofs.
Definition at line 206 of file segregated_fsi_solver.cc.
References i, and Solid_data_pt.
Referenced by segregated_solve().
|
protected |
Do pointwise Aitken extrapolation for solid.
Pointwise Aitken extrapolation for solid variables.
Definition at line 490 of file segregated_fsi_solver.cc.
References i, Pointwise_aitken_counter, Pointwise_aitken_solid_value, and Solid_data_pt.
Referenced by segregated_solve().
|
protected |
Rebuild global mesh for monolithic discretisation.
Rebuild global mesh for monolithic problem.
Definition at line 87 of file segregated_fsi_solver.cc.
References oomph::Problem::add_sub_mesh(), oomph::Problem::flush_sub_meshes(), i, Orig_sub_mesh_pt, and oomph::Problem::rebuild_global_mesh().
Referenced by segregated_solve().
|
inline |
Reset timer.
Definition at line 439 of file segregated_fsi_solver.h.
References T_ref, T_spent_on_actual_solve, and Timer_has_been_halted.
Referenced by segregated_solve().
|
inline |
(Re-)start timer (e.g. after completing non-essential parts of the code such as documentation of the iteration's progress)
Definition at line 450 of file segregated_fsi_solver.h.
References T_ref, and Timer_has_been_halted.
Referenced by t_spent_on_actual_solve().
|
protected |
Restore pinned status of fluid dofs.
Restore the pinned status of all fluid dofs.
Definition at line 176 of file segregated_fsi_solver.cc.
References Fluid_data_pt, Fluid_value_is_pinned, and i.
Referenced by segregated_solve().
|
protected |
Restore pinned status of solid dofs.
Restore the pinned status of all solid dofs.
Definition at line 229 of file segregated_fsi_solver.cc.
References i, Solid_data_pt, and Solid_value_is_pinned.
Referenced by segregated_solve().
PicardConvergenceData oomph::SegregatableFSIProblem::segregated_solve | ( | ) |
Segregated solver. Peform a segregated step from the present state of the system. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Segregated fixed point solver with optional pointwise Aitken acceleration on selected solid dofs. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Definition at line 538 of file segregated_fsi_solver.cc.
References actions_after_segregated_solve(), actions_before_segregated_convergence_check(), actions_before_segregated_solve(), Assess_convergence_based_on_absolute_solid_change, Assess_convergence_based_on_max_global_residual, Assess_convergence_based_on_relative_solid_change, oomph::Problem::assign_eqn_numbers(), Convergence_criterion, Convergence_tolerance, oomph::PicardConvergenceData::cpu_for_global_residual(), oomph::PicardConvergenceData::cpu_total(), Doc_max_global_residual, oomph::PicardConvergenceData::essential_cpu_total(), Fluid_solve, Full_solve, oomph::Problem::get_residuals(), get_solid_change(), oomph::DoubleVector::max(), Max_picard, oomph::Problem::newton_solve(), oomph::PicardConvergenceData::niter(), oomph::oomph_info, pin_fluid_dofs(), pin_solid_dofs(), Pointwise_aitken_counter, pointwise_aitken_extrapolate(), Pointwise_aitken_start, rebuild_monolithic_mesh(), Recheck_convergence_after_pointwise_aitken, reset_timer(), restore_fluid_dofs(), restore_solid_dofs(), oomph::PicardConvergenceData::set_solver_converged(), Solid_solve, Solve_type, store_solid_dofs(), t_spent_on_actual_solve(), oomph::PicardConvergenceData::tol_achieved(), under_relax_solid(), use_only_fluid_elements(), use_only_solid_elements(), and Use_pointwise_aitken.
Referenced by steady_segregated_solve(), and unsteady_segregated_solve().
void oomph::SegregatableFSIProblem::setup_segregated_solver | ( | const bool & | full_setup_of_fluid_and_solid_dofs = true | ) |
Setup the segregated solver: Backup the pinned status of the fluid and solid dofs and allocate the internal storage based on the input provided by identify_fluid_and_solid_dofs(...) In addition, reset storage associated with convergence acceleration techniques. If the problem and degrees of freedom has not changed between calls to the solver then it is wasteful to call identify_fluid_and_solid_dofs(...) again and again. If the optional boolean flag is set to false then the storage for convergence acceleration techniques is reset, but the fluid and solid dofs are not altered.
Setup segregated solver, using the information provided by the user in his/her implementation of the pure virtual function identify_fluid_and_solid_dofs(...).
Definition at line 1111 of file segregated_fsi_solver.cc.
References Del_irons_and_tuck, Fluid_data_pt, Fluid_mesh_pt, Fluid_value_is_pinned, i, identify_fluid_and_solid_dofs(), oomph::Problem::mesh_pt(), oomph::Problem::nsub_mesh(), Omega_relax, oomph::oomph_info, Orig_sub_mesh_pt, Pointwise_aitken_counter, Pointwise_aitken_solid_value, Previous_solid_value, R_irons_and_tuck, Solid_data_pt, Solid_mesh_pt, Solid_value_is_pinned, Use_irons_and_tuck_extrapolation, and Use_pointwise_aitken.
PicardConvergenceData oomph::SegregatableFSIProblem::steady_segregated_solve | ( | ) |
Steady version of segregated solver. Makes all timesteppers steady before solving. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Segregated fixed point solver with optional pointwise Aitken acceleration on selected solid dofs. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Definition at line 964 of file segregated_fsi_solver.cc.
References oomph::Problem::assign_initial_values_impulsive(), i, oomph::TimeStepper::is_steady(), oomph::TimeStepper::make_steady(), oomph::Problem::ntime_stepper(), oomph::oomph_info, oomph::SegregatedSolverError::Ran_out_of_iterations, segregated_solve(), oomph::Problem::time_stepper_pt(), and oomph::TimeStepper::undo_make_steady().
void oomph::SegregatableFSIProblem::store_solid_dofs | ( | ) |
Store the current solid values as reference values for future convergence check. Also add another entry to pointwise Aitken history if required.
Store the current solid values as reference values for future convergence check. Also add another entry to pointwise Aitken history.
Definition at line 260 of file segregated_fsi_solver.cc.
References i, Pointwise_aitken_counter, Pointwise_aitken_solid_value, Previous_solid_value, Solid_data_pt, and Use_pointwise_aitken.
Referenced by segregated_solve().
|
inline |
Total elapsed time since start of solve.
Definition at line 471 of file segregated_fsi_solver.h.
References halt_timer(), restart_timer(), T_spent_on_actual_solve, and oomph::Problem::time().
Referenced by segregated_solve().
|
private |
Under-relax the most recently computed solid variables, either by classical relaxation or by Irons & Tuck.
Under-relax solid, either by classical relaxation or by Irons & Tuck.
Definition at line 354 of file segregated_fsi_solver.cc.
References Del_irons_and_tuck, i, Omega_relax, Previous_solid_value, R_irons_and_tuck, Solid_data_pt, and Use_irons_and_tuck_extrapolation.
Referenced by segregated_solve().
PicardConvergenceData oomph::SegregatableFSIProblem::unsteady_segregated_solve | ( | const double & | dt | ) |
Unsteady segregated solver, advance time by dt and solve by the segregated solver. The time values are always shifted by this function. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Segregated fixed point solver with optional pointwise Aitken acceleration on selected solid dofs. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Definition at line 1033 of file segregated_fsi_solver.cc.
PicardConvergenceData oomph::SegregatableFSIProblem::unsteady_segregated_solve | ( | const double & | dt, |
const bool & | shift_values | ||
) |
Unsteady segregated solver. Advance time by dt and solve the system by a segregated method. The boolean flag is used to control whether the time values should be shifted. If it is true the current data values will be shifted (stored as previous timesteps) before the solution step. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Segregated fixed point solver with optional pointwise Aitken acceleration on selected solid dofs. Returns PicardConvergenceData object that contains the vital stats of the iteration.
Definition at line 1045 of file segregated_fsi_solver.cc.
References oomph::Problem::actions_after_implicit_timestep(), oomph::Problem::actions_before_implicit_timestep(), oomph::Time::dt(), extrapolate_solid_data(), i, oomph::Problem::ntime_stepper(), oomph::oomph_info, oomph::SegregatedSolverError::Ran_out_of_iterations, segregated_solve(), oomph::TimeStepper::set_weights(), oomph::Problem::shift_time_values(), oomph::Time::time(), oomph::Problem::time_pt(), and oomph::Problem::time_stepper_pt().
|
private |
Only include fluid elements in the Problem's mesh. This is called before the segregated fluid solve. The fluid elements are identified by the user via the fluid_mesh_pt argument in the pure virtual function identify_fluid_and_solid_dofs(...). If a NULL pointer is returned by this function (i.e. if the user hasn't bothered to identify the fluids elements in a submesh, then no stripping of non-fluid elements is performed. The result of the computation will be correct but it won't be very efficient.
Definition at line 138 of file segregated_fsi_solver.cc.
References oomph::Problem::add_sub_mesh(), Fluid_mesh_pt, oomph::Problem::flush_sub_meshes(), and oomph::Problem::rebuild_global_mesh().
Referenced by segregated_solve().
|
private |
Only include solid elements in the Problem's mesh. This is called before the segregated solid solve. The solid elements are identified by the user via the solid_mesh_pt argument in the pure virtual function identify_fluid_and_solid_dofs(...). If a NULL pointer is returned by this function (i.e. if the user hasn't bothered to identify the solid elements in a submesh, then no stripping of non-solid elements is performed. The result of the computation will be correct but it won't be very efficient.
Only include solid elements in the Problem's mesh. This is called before the segregated solid solve. The solid elements are identified by the user via the solid_mesh_pt argument in the pure virtual function identify_fluid_and_solid_dofs(...). If a NULL pointer is returned by this function (i.e. if the user hasn't bothered to identify the solids elements in a submesh, then no stripping of non-solid elements is performed. The result of the computation will be correct but it won't be very efficient.
Definition at line 115 of file segregated_fsi_solver.cc.
References oomph::Problem::add_sub_mesh(), oomph::Problem::flush_sub_meshes(), oomph::Problem::rebuild_global_mesh(), and Solid_mesh_pt.
Referenced by segregated_solve().
|
private |
Convergence criterion (enumerated flag)
Definition at line 611 of file segregated_fsi_solver.h.
Referenced by assess_convergence_based_on_absolute_solid_change(), assess_convergence_based_on_max_global_residual(), assess_convergence_based_on_relative_solid_change(), SegregatableFSIProblem(), and segregated_solve().
|
protected |
Convergence tolerance for Picard iteration.
Definition at line 500 of file segregated_fsi_solver.h.
Referenced by assess_convergence_based_on_absolute_solid_change(), assess_convergence_based_on_max_global_residual(), assess_convergence_based_on_relative_solid_change(), SegregatableFSIProblem(), and segregated_solve().
|
protected |
Vector of changes in Irons and Tuck under-relaxation.
Definition at line 557 of file segregated_fsi_solver.h.
Referenced by setup_segregated_solver(), and under_relax_solid().
|
protected |
Doc maximum global residual during iteration? (default: false)
Definition at line 506 of file segregated_fsi_solver.h.
Referenced by SegregatableFSIProblem(), and segregated_solve().
Vector storing the Data objects associated with the fluid problem: Tyically the nodal and internal data of the elements in the fluid bulk mesh.
Definition at line 524 of file segregated_fsi_solver.h.
Referenced by pin_fluid_dofs(), restore_fluid_dofs(), and setup_segregated_solver().
|
protected |
Mesh containing only fluid elements – the elements in this Mesh will be excluded from the assembly process when the solid problem is solved.
Definition at line 546 of file segregated_fsi_solver.h.
Referenced by extrapolate_solid_data(), SegregatableFSIProblem(), setup_segregated_solver(), and use_only_fluid_elements().
|
protected |
Vector of vectors that store the pinned status of fluid Data values.
Definition at line 528 of file segregated_fsi_solver.h.
Referenced by restore_fluid_dofs(), and setup_segregated_solver().
|
protected |
Max. number of Picard iterations.
Definition at line 503 of file segregated_fsi_solver.h.
Referenced by SegregatableFSIProblem(), and segregated_solve().
|
private |
Under-relaxation parameter. (1.0: no under-relaxation; 0.0: Freeze wall shape)
Definition at line 604 of file segregated_fsi_solver.h.
Referenced by enable_under_relaxation(), SegregatableFSIProblem(), setup_segregated_solver(), and under_relax_solid().
Backup for the pointers to the submeshes in the original problem.
Definition at line 554 of file segregated_fsi_solver.h.
Referenced by rebuild_monolithic_mesh(), and setup_segregated_solver().
|
protected |
Number of Aitken histories available (int because after extrapolation it's re-initialised to -1 to force the computation of three new genuine iterates).
Definition at line 487 of file segregated_fsi_solver.h.
Referenced by pointwise_aitken_extrapolate(), segregated_solve(), setup_segregated_solver(), and store_solid_dofs().
Vector of Vectors containing up to three previous iterates for the solid dofs; used for pointwise Aitken extrapolation.
Definition at line 564 of file segregated_fsi_solver.h.
Referenced by pointwise_aitken_extrapolate(), setup_segregated_solver(), and store_solid_dofs().
|
protected |
Start pointwise Aitken extrpolation after specified number of Picard iterations.
Definition at line 494 of file segregated_fsi_solver.h.
Referenced by enable_pointwise_aitken(), SegregatableFSIProblem(), and segregated_solve().
|
protected |
Vector storing the previous solid values – used for convergence check.
Definition at line 541 of file segregated_fsi_solver.h.
Referenced by get_solid_change(), setup_segregated_solver(), store_solid_dofs(), and under_relax_solid().
|
protected |
Irons and Tuck relaxation factor.
Definition at line 560 of file segregated_fsi_solver.h.
Referenced by setup_segregated_solver(), and under_relax_solid().
|
protected |
Have we just done a pointwise Aitken step.
Definition at line 567 of file segregated_fsi_solver.h.
Referenced by SegregatableFSIProblem(), and segregated_solve().
Vector storing the Data objects associated with the solid problem: Typically the positional data of solid nodes and any quantities associated with displacement control, say.
Definition at line 533 of file segregated_fsi_solver.h.
Referenced by extrapolate_solid_data(), get_solid_change(), pin_solid_dofs(), pointwise_aitken_extrapolate(), restore_solid_dofs(), setup_segregated_solver(), store_solid_dofs(), and under_relax_solid().
|
protected |
Mesh containing only solid elements – the elements in this mesh will be excluded from the assembly process when the fluid problem is solved.
Definition at line 551 of file segregated_fsi_solver.h.
Referenced by SegregatableFSIProblem(), setup_segregated_solver(), and use_only_solid_elements().
|
protected |
Vector of vectors that store the pinned status of solid Data values.
Definition at line 537 of file segregated_fsi_solver.h.
Referenced by restore_solid_dofs(), and setup_segregated_solver().
|
protected |
Solve that is taking place (enumerated flag)
Definition at line 497 of file segregated_fsi_solver.h.
Referenced by SegregatableFSIProblem(), and segregated_solve().
|
private |
Reference time for segregated solve. Can be re-initialised whenever total elapsed time has been stored (before entering non-essential doc sections of the code)
Definition at line 616 of file segregated_fsi_solver.h.
Referenced by halt_timer(), reset_timer(), restart_timer(), and SegregatableFSIProblem().
|
private |
Total elapsed time since start of solve, can be accumulated by adding bits of time spent in relevant parts of code (bypassing sections that only document the progress)
Definition at line 621 of file segregated_fsi_solver.h.
Referenced by halt_timer(), reset_timer(), SegregatableFSIProblem(), and t_spent_on_actual_solve().
|
private |
boolean flag to indicate if timer has been halted
Definition at line 624 of file segregated_fsi_solver.h.
Referenced by halt_timer(), reset_timer(), restart_timer(), and SegregatableFSIProblem().
|
private |
Boolean flag to indicate use of Irons and Tuck's extrapolation for solid values.
Definition at line 608 of file segregated_fsi_solver.h.
Referenced by disable_irons_and_tuck_extrapolation(), enable_irons_and_tuck_extrapolation(), SegregatableFSIProblem(), setup_segregated_solver(), and under_relax_solid().
|
protected |
Use pointwise Aitken extrapolation?
Definition at line 490 of file segregated_fsi_solver.h.
Referenced by disable_pointwise_aitken(), enable_pointwise_aitken(), SegregatableFSIProblem(), segregated_solve(), setup_segregated_solver(), and store_solid_dofs().