libMesh
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
libMesh::NonlinearSolver< T > Class Template Referenceabstract

This base class can be inherited from to provide interfaces to nonlinear solvers from different packages like PETSc and Trilinos. More...

#include <nonlinear_solver.h>

Inheritance diagram for libMesh::NonlinearSolver< T >:
[legend]

Public Types

typedef NonlinearImplicitSystem sys_type
 The type of system. More...
 

Public Member Functions

 NonlinearSolver (sys_type &s)
 Constructor. More...
 
virtual ~NonlinearSolver ()
 Destructor. More...
 
bool initialized () const
 
virtual void clear ()
 Release all memory and clear data structures. More...
 
virtual void init (const char *name=nullptr)=0
 Initialize data structures if not done so already. More...
 
virtual std::pair< unsigned int, Realsolve (SparseMatrix< T > &, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)=0
 Solves the nonlinear system. More...
 
virtual void print_converged_reason ()
 Prints a useful message about why the latest nonlinear solve con(di)verged. More...
 
virtual int get_total_linear_iterations ()=0
 Get the total number of linear iterations done in the last solve. More...
 
virtual unsigned get_current_nonlinear_iteration_number () const =0
 
const sys_typesystem () const
 
sys_typesystem ()
 
void attach_preconditioner (Preconditioner< T > *preconditioner)
 Attaches a Preconditioner object to be used during the linear solves. More...
 
void set_solver_configuration (SolverConfiguration &solver_configuration)
 Set the solver configuration object. More...
 
virtual bool reuse_preconditioner () const
 Get the reuse_preconditioner flag. More...
 
virtual void set_reuse_preconditioner (bool reuse)
 Set the reuse preconditioner flag. More...
 
virtual unsigned int reuse_preconditioner_max_linear_its () const
 Get the reuse_preconditioner_max_linear_its parameter. More...
 
virtual void set_reuse_preconditioner_max_linear_its (unsigned int i)
 Set the reuse_preconditioner_max_linear_its parameter. More...
 
virtual void force_new_preconditioner ()
 Immediately force a new preconditioner. More...
 
virtual void set_exact_constraint_enforcement (bool enable)
 Enable (or disable; it is true by default) exact enforcement of constraints at the solver level, correcting any constrained DoF coefficients in current_local_solution as well as applying nonlinear residual and Jacobian terms based on constraint equations. More...
 
bool exact_constraint_enforcement ()
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< NonlinearSolver< T > > build (sys_type &s, const SolverPackage solver_package=libMesh::default_solver_package())
 Builds a NonlinearSolver using the nonlinear solver package specified by solver_package. More...
 
static std::string get_info ()
 Gets a string containing the reference information. More...
 
static void print_info (std::ostream &out_stream=libMesh::out)
 Prints the reference information, by default to libMesh::out. More...
 
static unsigned int n_objects ()
 Prints the number of outstanding (created, but not yet destroyed) objects. More...
 
static void enable_print_counter_info ()
 Methods to enable/disable the reference counter output from print_info() More...
 
static void disable_print_counter_info ()
 

Public Attributes

void(* residual )(const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S)
 Function that computes the residual R(X) of the nonlinear system at the input iterate X. More...
 
NonlinearImplicitSystem::ComputeResidualresidual_object
 Object that computes the residual R(X) of the nonlinear system at the input iterate X. More...
 
NonlinearImplicitSystem::ComputeResidualfd_residual_object
 Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose of forming a finite-differenced Jacobian. More...
 
NonlinearImplicitSystem::ComputeResidualmffd_residual_object
 Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose of forming Jacobian-vector products via finite differencing. More...
 
void(* jacobian )(const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S)
 Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X. More...
 
NonlinearImplicitSystem::ComputeJacobianjacobian_object
 Object that computes the Jacobian J(X) of the nonlinear system at the input iterate X. More...
 
void(* matvec )(const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S)
 Function that computes either the residual \( R(X) \) or the Jacobian \( J(X) \) of the nonlinear system at the input iterate \( X \). More...
 
NonlinearImplicitSystem::ComputeResidualandJacobianresidual_and_jacobian_object
 Object that computes either the residual \( R(X) \) or the Jacobian \( J(X) \) of the nonlinear system at the input iterate \( X \). More...
 
void(* bounds )(NumericVector< Number > &XL, NumericVector< Number > &XU, sys_type &S)
 Function that computes the lower and upper bounds XL and XU on the solution of the nonlinear system. More...
 
NonlinearImplicitSystem::ComputeBoundsbounds_object
 Object that computes the bounds vectors \( XL \) and \( XU \). More...
 
void(* nullspace )(std::vector< NumericVector< Number > *> &sp, sys_type &S)
 Function that computes a basis for the Jacobian's nullspace – the kernel or the "zero energy modes" – that can be used in solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP). More...
 
NonlinearImplicitSystem::ComputeVectorSubspacenullspace_object
 A callable object that computes a basis for the Jacobian's nullspace – the kernel or the "zero energy modes" – that can be used in solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP). More...
 
void(* transpose_nullspace )(std::vector< NumericVector< Number > *> &sp, sys_type &S)
 Function that computes a basis for the transpose Jacobian's nullspace – when solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP), it is used to remove contributions outside of R(jac) More...
 
NonlinearImplicitSystem::ComputeVectorSubspacetranspose_nullspace_object
 A callable object that computes a basis for the transpose Jacobian's nullspace – when solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP), it is used to remove contributions outside of R(jac) More...
 
void(* nearnullspace )(std::vector< NumericVector< Number > *> &sp, sys_type &S)
 Function that computes a basis for the Jacobian's near nullspace – the set of "low energy modes" – that can be used for AMG coarsening, if the solver supports it (e.g., ML, PETSc's GAMG). More...
 
NonlinearImplicitSystem::ComputeVectorSubspacenearnullspace_object
 A callable object that computes a basis for the Jacobian's near nullspace – the set of "low energy modes" – that can be used for AMG coarsening, if the solver supports it (e.g., ML, PETSc's GAMG). More...
 
void(* user_presolve )(sys_type &S)
 Customizable function pointer which users can attach to the solver. More...
 
void(* postcheck )(const NumericVector< Number > &old_soln, NumericVector< Number > &search_direction, NumericVector< Number > &new_soln, bool &changed_search_direction, bool &changed_new_soln, sys_type &S)
 Function that performs a "check" on the Newton search direction and solution after each nonlinear step. More...
 
NonlinearImplicitSystem::ComputePostCheckpostcheck_object
 A callable object that is executed after each nonlinear iteration. More...
 
NonlinearImplicitSystem::ComputePreCheckprecheck_object
 
unsigned int max_nonlinear_iterations
 Maximum number of non-linear iterations. More...
 
unsigned int max_function_evaluations
 Maximum number of function evaluations. More...
 
double absolute_residual_tolerance
 The NonlinearSolver should exit after the residual is reduced to either less than absolute_residual_tolerance or less than relative_residual_tolerance times the initial residual. More...
 
double relative_residual_tolerance
 
double divergence_tolerance
 The NonlinearSolver should exit if the residual becomes greater than the initial residual times the divergence_tolerance. More...
 
double absolute_step_tolerance
 The NonlinearSolver should exit after the full nonlinear step norm is reduced to either less than absolute_step_tolerance or less than relative_step_tolerance times the largest nonlinear solution which has been seen so far. More...
 
double relative_step_tolerance
 
unsigned int max_linear_iterations
 Each linear solver step should exit after max_linear_iterations is exceeded. More...
 
double initial_linear_tolerance
 Any required linear solves will at first be done with this tolerance; the NonlinearSolver may tighten the tolerance for later solves. More...
 
double minimum_linear_tolerance
 The tolerance for linear solves is kept above this minimum. More...
 
bool converged
 After a call to solve this will reflect whether or not the nonlinear solve was successful. More...
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 Data structure to log the information. More...
 

Protected Member Functions

void increment_constructor_count (const std::string &name) noexcept
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name) noexcept
 Increments the destruction counter. More...
 

Protected Attributes

bool _reuse_preconditioner
 Whether we should reuse the linear preconditioner. More...
 
bool _exact_constraint_enforcement
 Whether we should enforce exact constraints globally during a solve. More...
 
unsigned int _reuse_preconditioner_max_linear_its
 Number of linear iterations to retain the preconditioner. More...
 
sys_type_system
 A reference to the system we are solving. More...
 
bool _is_initialized
 Flag indicating if the data structures have been initialized. More...
 
Preconditioner< T > * _preconditioner
 Holds the Preconditioner object to be used for the linear solves. More...
 
SolverConfiguration_solver_configuration
 Optionally store a SolverOptions object that can be used to set parameters like solver type, tolerances and iteration limits. More...
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

static Counts _counts
 Actually holds the data. More...
 
static Threads::atomic< unsigned int_n_objects
 The number of objects. More...
 
static Threads::spin_mutex _mutex
 Mutual exclusion object to enable thread-safe reference counting. More...
 
static bool _enable_print_counter = true
 Flag to control whether reference count information is printed when print_info is called. More...
 

Detailed Description

template<typename T>
class libMesh::NonlinearSolver< T >

This base class can be inherited from to provide interfaces to nonlinear solvers from different packages like PETSc and Trilinos.

Author
Benjamin Kirk
Date
2005

Definition at line 52 of file nonlinear_solver.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information.

The log is identified by the class name.

Definition at line 119 of file reference_counter.h.

◆ sys_type

template<typename T>
typedef NonlinearImplicitSystem libMesh::NonlinearSolver< T >::sys_type

The type of system.

Definition at line 59 of file nonlinear_solver.h.

Constructor & Destructor Documentation

◆ NonlinearSolver()

template<typename T >
libMesh::NonlinearSolver< T >::NonlinearSolver ( sys_type s)
inlineexplicit

Constructor.

Initializes Solver data structures

Definition at line 453 of file nonlinear_solver.h.

453  :
454  ParallelObject (s),
455  residual (nullptr),
456  residual_object (nullptr),
457  fd_residual_object (nullptr),
458  mffd_residual_object (nullptr),
459  jacobian (nullptr),
460  jacobian_object (nullptr),
461  matvec (nullptr),
463  bounds (nullptr),
464  bounds_object (nullptr),
465  nullspace (nullptr),
466  nullspace_object (nullptr),
467  transpose_nullspace (nullptr),
468  transpose_nullspace_object (nullptr),
469  nearnullspace (nullptr),
470  nearnullspace_object (nullptr),
471  user_presolve (nullptr),
472  postcheck (nullptr),
473  postcheck_object (nullptr),
474  precheck_object (nullptr),
485  converged(false),
486  _reuse_preconditioner(false),
489  _system(s),
490  _is_initialized (false),
491  _preconditioner (nullptr),
492  _solver_configuration(nullptr)
493 {
494 }
ParallelObject(const Parallel::Communicator &comm_in)
Constructor.
unsigned int _reuse_preconditioner_max_linear_its
Number of linear iterations to retain the preconditioner.
void(* residual)(const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S)
Function that computes the residual R(X) of the nonlinear system at the input iterate X...
NonlinearImplicitSystem::ComputeResidualandJacobian * residual_and_jacobian_object
Object that computes either the residual or the Jacobian of the nonlinear system at the input itera...
sys_type & _system
A reference to the system we are solving.
bool converged
After a call to solve this will reflect whether or not the nonlinear solve was successful.
unsigned int max_function_evaluations
Maximum number of function evaluations.
Preconditioner< T > * _preconditioner
Holds the Preconditioner object to be used for the linear solves.
void(* user_presolve)(sys_type &S)
Customizable function pointer which users can attach to the solver.
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
SolverConfiguration * _solver_configuration
Optionally store a SolverOptions object that can be used to set parameters like solver type...
NonlinearImplicitSystem::ComputeResidual * residual_object
Object that computes the residual R(X) of the nonlinear system at the input iterate X...
void(* jacobian)(const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S)
Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X...
bool _is_initialized
Flag indicating if the data structures have been initialized.
NonlinearImplicitSystem::ComputeJacobian * jacobian_object
Object that computes the Jacobian J(X) of the nonlinear system at the input iterate X...
double minimum_linear_tolerance
The tolerance for linear solves is kept above this minimum.
NonlinearImplicitSystem::ComputeBounds * bounds_object
Object that computes the bounds vectors and .
bool _exact_constraint_enforcement
Whether we should enforce exact constraints globally during a solve.
void(* transpose_nullspace)(std::vector< NumericVector< Number > *> &sp, sys_type &S)
Function that computes a basis for the transpose Jacobian&#39;s nullspace – when solving a degenerate pr...
void(* nullspace)(std::vector< NumericVector< Number > *> &sp, sys_type &S)
Function that computes a basis for the Jacobian&#39;s nullspace – the kernel or the "zero energy modes" ...
void(* nearnullspace)(std::vector< NumericVector< Number > *> &sp, sys_type &S)
Function that computes a basis for the Jacobian&#39;s near nullspace – the set of "low energy modes" – ...
double divergence_tolerance
The NonlinearSolver should exit if the residual becomes greater than the initial residual times the d...
NonlinearImplicitSystem::ComputePostCheck * postcheck_object
A callable object that is executed after each nonlinear iteration.
double absolute_step_tolerance
The NonlinearSolver should exit after the full nonlinear step norm is reduced to either less than abs...
bool _reuse_preconditioner
Whether we should reuse the linear preconditioner.
void(* bounds)(NumericVector< Number > &XL, NumericVector< Number > &XU, sys_type &S)
Function that computes the lower and upper bounds XL and XU on the solution of the nonlinear system...
NonlinearImplicitSystem::ComputeResidual * mffd_residual_object
Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose...
double absolute_residual_tolerance
The NonlinearSolver should exit after the residual is reduced to either less than absolute_residual_t...
void(* postcheck)(const NumericVector< Number > &old_soln, NumericVector< Number > &search_direction, NumericVector< Number > &new_soln, bool &changed_search_direction, bool &changed_new_soln, sys_type &S)
Function that performs a "check" on the Newton search direction and solution after each nonlinear ste...
NonlinearImplicitSystem::ComputeVectorSubspace * nearnullspace_object
A callable object that computes a basis for the Jacobian&#39;s near nullspace – the set of "low energy m...
NonlinearImplicitSystem::ComputePreCheck * precheck_object
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the NonlinearSolver may tighten...
unsigned int max_nonlinear_iterations
Maximum number of non-linear iterations.
NonlinearImplicitSystem::ComputeVectorSubspace * transpose_nullspace_object
A callable object that computes a basis for the transpose Jacobian&#39;s nullspace – when solving a dege...
void(* matvec)(const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S)
Function that computes either the residual or the Jacobian of the nonlinear system at the input ite...
NonlinearImplicitSystem::ComputeResidual * fd_residual_object
Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose...
NonlinearImplicitSystem::ComputeVectorSubspace * nullspace_object
A callable object that computes a basis for the Jacobian&#39;s nullspace – the kernel or the "zero energ...

◆ ~NonlinearSolver()

template<typename T >
libMesh::NonlinearSolver< T >::~NonlinearSolver ( )
inlinevirtual

Destructor.

Definition at line 500 of file nonlinear_solver.h.

501 {
502  this->NonlinearSolver::clear ();
503 }
virtual void clear()
Release all memory and clear data structures.

Member Function Documentation

◆ attach_preconditioner()

template<typename T>
void libMesh::NonlinearSolver< T >::attach_preconditioner ( Preconditioner< T > *  preconditioner)

Attaches a Preconditioner object to be used during the linear solves.

Definition at line 74 of file nonlinear_solver.C.

75 {
76  libmesh_error_msg_if(this->_is_initialized, "Preconditioner must be attached before the solver is initialized!");
77 
78  _preconditioner = preconditioner;
79 }
Preconditioner< T > * _preconditioner
Holds the Preconditioner object to be used for the linear solves.
bool _is_initialized
Flag indicating if the data structures have been initialized.

◆ build()

template<typename T >
std::unique_ptr< NonlinearSolver< T > > libMesh::NonlinearSolver< T >::build ( sys_type s,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
static

Builds a NonlinearSolver using the nonlinear solver package specified by solver_package.

Definition at line 40 of file nonlinear_solver.C.

41 {
42  // Build the appropriate solver
43  switch (solver_package)
44  {
45 
46 #ifdef LIBMESH_HAVE_PETSC
47  case PETSC_SOLVERS:
48  return std::make_unique<PetscNonlinearSolver<T>>(s);
49 #endif // LIBMESH_HAVE_PETSC
50 
51 #if defined(LIBMESH_TRILINOS_HAVE_NOX) && defined(LIBMESH_TRILINOS_HAVE_EPETRA)
52  case TRILINOS_SOLVERS:
53  return std::make_unique<NoxNonlinearSolver<T>>(s);
54 #endif
55 
56  default:
57  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
58  }
59 }

◆ clear()

template<typename T>
virtual void libMesh::NonlinearSolver< T >::clear ( )
inlinevirtual

Release all memory and clear data structures.

Reimplemented in libMesh::PetscNonlinearSolver< T >, libMesh::PetscNonlinearSolver< Number >, libMesh::NoxNonlinearSolver< T >, and libMesh::NoxNonlinearSolver< Number >.

Definition at line 88 of file nonlinear_solver.h.

88 {}

◆ comm()

const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 97 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult_add(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::DofMap::add_constraints_to_send_list(), add_cube_convex_hull_to_mesh(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::TransientRBConstruction::add_IC_to_RB_space(), libMesh::EigenSystem::add_matrices(), libMesh::System::add_matrix(), libMesh::RBConstruction::add_scaled_matrix_and_vector(), libMesh::System::add_variable(), libMesh::System::add_variables(), libMesh::System::add_vector(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::TransientRBConstruction::allocate_data_structures(), libMesh::RBConstruction::allocate_data_structures(), libMesh::TransientRBConstruction::assemble_affine_expansion(), libMesh::FEMSystem::assemble_qoi(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::MeshCommunication::assign_global_indices(), libMesh::Partitioner::assign_partitioning(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_data(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::RBConstruction::compute_Fq_representor_innerprods(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::RBConstruction::compute_output_dual_innerprods(), libMesh::RBConstruction::compute_residual_dual_norm_slow(), libMesh::RBSCMConstruction::compute_SCM_bounds_on_training_set(), libMesh::DofMap::computed_sparsity_already(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ContinuationSystem::ContinuationSystem(), libMesh::MeshBase::copy_constraint_rows(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::DofMap::create_dof_constraints(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::PetscMatrix< libMesh::Number >::create_submatrix_nosort(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::RBEIMEvaluation::distribute_bfs(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::DTKSolutionTransfer::DTKSolutionTransfer(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_interiors(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_nodes(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_sides(), libMesh::TransientRBConstruction::enrich_RB_space(), libMesh::EpetraVector< T >::EpetraVector(), AssembleOptimization::equality_constraints(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_error_tolerance(), libMesh::MeshRefinement::flag_elements_by_mean_stddev(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::DofMap::gather_constraints(), libMesh::MeshfreeInterpolation::gather_remote_data(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::RBEIMEvaluation::get_eim_basis_function_node_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_side_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_value(), libMesh::MeshBase::get_info(), libMesh::System::get_info(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::RBEIMConstruction::get_max_abs_value(), libMesh::RBEIMConstruction::get_node_max_abs_value(), libMesh::RBEIMEvaluation::get_parametrized_function_node_value(), libMesh::RBEIMEvaluation::get_parametrized_function_side_value(), libMesh::RBEIMEvaluation::get_parametrized_function_value(), libMesh::RBEIMConstruction::get_random_point(), AssembleOptimization::inequality_constraints(), AssembleOptimization::inequality_constraints_jacobian(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::RBEIMConstruction::initialize_parametrized_functions_in_training_set(), libMesh::RBEIMConstruction::inner_product(), integrate_function(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_equal_connectivity(), libMesh::MeshTools::libmesh_assert_equal_points(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_linesearch_shellfunc(), libMesh::libmesh_petsc_preconditioner_apply(), libMesh::libmesh_petsc_recalculate_monitor(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_interface(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_precheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::LinearImplicitSystem::LinearImplicitSystem(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_bcids_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::TransientRBConstruction::mass_matrix_scaled_matvec(), libMesh::FEMSystem::mesh_position_set(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), LinearElasticityWithContact::move_mesh(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::DofMap::n_constrained_dofs(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), MixedOrderTest::n_neighbor_links(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SparsityPattern::Build::n_nonzeros(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_distribute_bfs(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::RBEIMConstruction::node_inner_product(), libMesh::MeshBase::operator==(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::MeshTools::paranoid_n_levels(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Partitioner::processor_pairs_to_interface_nodes(), libMesh::InterMeshProjection::project_system_vectors(), FEMParameters::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::RBEIMEvaluation::read_in_interior_basis_functions(), libMesh::RBEIMEvaluation::read_in_node_basis_functions(), libMesh::RBEIMEvaluation::read_in_side_basis_functions(), libMesh::RBEvaluation::read_in_vectors_from_multiple_files(), libMesh::System::read_legacy_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), LinearElasticityWithContact::residual_and_jacobian(), OverlappingAlgebraicGhostingTest::run_ghosting_test(), OverlappingCouplingGhostingTest::run_sparsity_pattern_test(), scale_mesh_and_plot(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::send_and_insert_dof_values(), libMesh::TransientRBConstruction::set_error_temporal_data(), libMesh::Partitioner::set_interface_node_processor_ids_BFS(), libMesh::Partitioner::set_interface_node_processor_ids_linear(), libMesh::Partitioner::set_interface_node_processor_ids_petscpartitioner(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::RBEIMEvaluation::side_distribute_bfs(), libMesh::RBEIMEvaluation::side_gather_bfs(), libMesh::RBEIMConstruction::side_inner_product(), libMesh::Partitioner::single_partition(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::RBEIMConstruction::store_eim_solutions_for_training_set(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), ConstraintOperatorTest::test1DCoarseningNewNodes(), ConstraintOperatorTest::test1DCoarseningOperator(), libMesh::MeshRefinement::test_level_one(), MeshfunctionDFEM::test_mesh_function_dfem(), MeshfunctionDFEM::test_mesh_function_dfem_grad(), MeshFunctionTest::test_p_level(), libMesh::MeshRefinement::test_unflagged(), DofMapTest::testBadElemFECombo(), SystemsTest::testBlockRestrictedVarNDofs(), BoundaryInfoTest::testBoundaryOnChildrenErrors(), ConstraintOperatorTest::testCoreform(), MeshInputTest::testExodusIGASidesets(), MeshTriangulationTest::testFoundCenters(), PointLocatorTest::testLocator(), BoundaryInfoTest::testMesh(), PointLocatorTest::testPlanar(), MeshTriangulationTest::testPoly2TriRefinementBase(), SystemsTest::testProjectCubeWithMeshFunction(), BoundaryInfoTest::testRenumber(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), MeshTriangulationTest::testTriangulatorInterp(), MeshTriangulationTest::testTriangulatorMeshedHoles(), MeshTriangulationTest::testTriangulatorRoundHole(), libMesh::MeshTools::total_weight(), libMesh::RBConstruction::train_reduced_basis_with_POD(), libMesh::MeshFunctionSolutionTransfer::transfer(), libMesh::MeshfreeSolutionTransfer::transfer(), libMesh::Poly2TriTriangulator::triangulate(), libMesh::TransientRBConstruction::truth_assembly(), libMesh::RBConstruction::truth_assembly(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::TransientRBConstruction::update_RB_initial_condition_all_N(), libMesh::TransientRBConstruction::update_RB_system_matrices(), libMesh::RBConstruction::update_RB_system_matrices(), libMesh::TransientRBConstruction::update_residual_terms(), libMesh::RBConstruction::update_residual_terms(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::VTKIO::write_nodal_data(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), libMesh::RBEvaluation::write_out_vectors(), libMesh::TransientRBConstruction::write_riesz_representors_to_files(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::RBDataSerialization::RBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::TransientRBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::RBEIMEvaluationSerialization::write_to_file(), and libMesh::RBDataSerialization::RBSCMEvaluationSerialization::write_to_file().

98  { return _communicator; }
const Parallel::Communicator & _communicator

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = false;
103  return;
104 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 94 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

95 {
96  _enable_print_counter = true;
97  return;
98 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ exact_constraint_enforcement()

template<typename T>
bool libMesh::NonlinearSolver< T >::exact_constraint_enforcement ( )
inline

Definition at line 403 of file nonlinear_solver.h.

404  {
406  }
bool _exact_constraint_enforcement
Whether we should enforce exact constraints globally during a solve.

◆ force_new_preconditioner()

template<typename T>
virtual void libMesh::NonlinearSolver< T >::force_new_preconditioner ( )
inlinevirtual

Immediately force a new preconditioner.

Reimplemented in libMesh::PetscNonlinearSolver< T >, and libMesh::PetscNonlinearSolver< Number >.

Definition at line 382 of file nonlinear_solver.h.

382 {};

◆ get_current_nonlinear_iteration_number()

template<typename T>
virtual unsigned libMesh::NonlinearSolver< T >::get_current_nonlinear_iteration_number ( ) const
pure virtual
Returns
The current nonlinear iteration number if called during the solve(), for example by the user-specified residual or Jacobian function.

Must be overridden in derived classes.

Implemented in libMesh::PetscNonlinearSolver< T >, libMesh::PetscNonlinearSolver< Number >, libMesh::NoxNonlinearSolver< T >, and libMesh::NoxNonlinearSolver< Number >.

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & [name, cd] : _counts)
59  oss << "| " << name << " reference count information:\n"
60  << "| Creations: " << cd.first << '\n'
61  << "| Destructions: " << cd.second << '\n';
62 
63  oss << " ---------------------------------------------------------------------------- \n";
64 
65  return oss.str();
66 
67 #else
68 
69  return "";
70 
71 #endif
72 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
static Counts _counts
Actually holds the data.

◆ get_total_linear_iterations()

template<typename T>
virtual int libMesh::NonlinearSolver< T >::get_total_linear_iterations ( )
pure virtual

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectednoexceptinherited

Increments the construction counter.

Should be called in the constructor of any derived class that will be reference counted.

Definition at line 183 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

184 {
185  libmesh_try
186  {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189  p.first++;
190  }
191  libmesh_catch (...)
192  {
193  auto stream = libMesh::err.get();
194  stream->exceptions(stream->goodbit); // stream must not throw
195  libMesh::err << "Encountered unrecoverable error while calling "
196  << "ReferenceCounter::increment_constructor_count() "
197  << "for a(n) " << name << " object." << std::endl;
198  std::terminate();
199  }
200 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectednoexceptinherited

Increments the destruction counter.

Should be called in the destructor of any derived class that will be reference counted.

Definition at line 207 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

208 {
209  libmesh_try
210  {
211  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
212  std::pair<unsigned int, unsigned int> & p = _counts[name];
213  p.second++;
214  }
215  libmesh_catch (...)
216  {
217  auto stream = libMesh::err.get();
218  stream->exceptions(stream->goodbit); // stream must not throw
219  libMesh::err << "Encountered unrecoverable error while calling "
220  << "ReferenceCounter::increment_destructor_count() "
221  << "for a(n) " << name << " object." << std::endl;
222  std::terminate();
223  }
224 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ init()

template<typename T>
virtual void libMesh::NonlinearSolver< T >::init ( const char *  name = nullptr)
pure virtual

Initialize data structures if not done so already.

May assign a name to the solver in some implementations

Implemented in libMesh::PetscNonlinearSolver< T >, libMesh::PetscNonlinearSolver< Number >, libMesh::NoxNonlinearSolver< T >, and libMesh::NoxNonlinearSolver< Number >.

◆ initialized()

template<typename T>
bool libMesh::NonlinearSolver< T >::initialized ( ) const
inline
Returns
true if the data structures are initialized, false otherwise.

Definition at line 83 of file nonlinear_solver.h.

83 { return _is_initialized; }
bool _is_initialized
Flag indicating if the data structures have been initialized.

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 85 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

Referenced by libMesh::LibMeshInit::~LibMeshInit().

86  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects
The number of objects.

◆ n_processors()

processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 103 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, libMesh::libmesh_assert(), and TIMPI::Communicator::size().

Referenced by libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DofMap::add_constraints_to_send_list(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::System::add_vector(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::DofMap::prepare_send_list(), libMesh::DofMap::print_dof_constraints(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), OverlappingFunctorTest::run_partitioner_test(), libMesh::DofMap::scatter_constraints(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), WriteVecAndScalar::setupTests(), libMesh::RBEIMEvaluation::side_gather_bfs(), DistributedMeshTest::testRemoteElemError(), CheckpointIOTest::testSplitter(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

104  {
105  processor_id_type returnval =
106  cast_int<processor_id_type>(_communicator.size());
107  libmesh_assert(returnval); // We never have an empty comm
108  return returnval;
109  }
const Parallel::Communicator & _communicator
processor_id_type size() const
uint8_t processor_id_type
libmesh_assert(ctx)

◆ print_converged_reason()

template<typename T>
virtual void libMesh::NonlinearSolver< T >::print_converged_reason ( )
inlinevirtual

Prints a useful message about why the latest nonlinear solve con(di)verged.

Reimplemented in libMesh::PetscNonlinearSolver< T >, and libMesh::PetscNonlinearSolver< Number >.

Definition at line 109 of file nonlinear_solver.h.

109 { libmesh_not_implemented(); }

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out_stream = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 81 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

Referenced by libMesh::LibMeshInit::~LibMeshInit().

82 {
84  out_stream << ReferenceCounter::get_info();
85 }
static std::string get_info()
Gets a string containing the reference information.
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ processor_id()

processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 114 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and TIMPI::Communicator::rank().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Partitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::find_dofs_to_send(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::DofMap::get_local_constraints(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::LaplaceMeshSmoother::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), HeatSystem::init_data(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::TransientRBEvaluation::legacy_write_offline_data_to_files(), libMesh::RBSCMEvaluation::legacy_write_offline_data_to_files(), libMesh::RBEvaluation::legacy_write_offline_data_to_files(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), AugmentSparsityOnInterface::mesh_reinit(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::DynaIO::read_mesh(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::RBEIMEvaluation::side_gather_bfs(), ExodusTest< elem_type >::test_read_gold(), ExodusTest< elem_type >::test_write(), MeshInputTest::testAbaqusRead(), MeshInputTest::testCopyElementSolutionImpl(), MeshInputTest::testCopyElementVectorImpl(), MeshInputTest::testCopyNodalSolutionImpl(), DefaultCouplingTest::testCoupling(), PointNeighborCouplingTest::testCoupling(), MeshInputTest::testDynaFileMappings(), MeshInputTest::testDynaNoSplines(), MeshInputTest::testDynaReadElem(), MeshInputTest::testDynaReadPatch(), MeshInputTest::testExodusFileMappings(), MeshInputTest::testExodusIGASidesets(), MeshInputTest::testExodusWriteElementDataFromDiscontinuousNodalData(), MeshInputTest::testLowOrderEdgeBlocks(), SystemsTest::testProjectMatrix1D(), SystemsTest::testProjectMatrix2D(), SystemsTest::testProjectMatrix3D(), BoundaryInfoTest::testShellFaceConstraints(), MeshInputTest::testSingleElementImpl(), WriteVecAndScalar::testSolution(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::DTKAdapter::update_variable_values(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_element_values_element_major(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO_Helper::write_elemset_data(), libMesh::ExodusII_IO_Helper::write_elemsets(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::ExodusII_IO_Helper::write_nodeset_data(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), write_output_solvedata(), libMesh::System::write_parallel_data(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::ExodusII_IO_Helper::write_sideset_data(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().

115  { return cast_int<processor_id_type>(_communicator.rank()); }
processor_id_type rank() const
const Parallel::Communicator & _communicator

◆ reuse_preconditioner()

template<typename T >
bool libMesh::NonlinearSolver< T >::reuse_preconditioner ( ) const
virtual

Get the reuse_preconditioner flag.

Reimplemented in libMesh::PetscNonlinearSolver< T >, and libMesh::PetscNonlinearSolver< Number >.

Definition at line 88 of file nonlinear_solver.C.

89 {
90  libmesh_not_implemented();
91 }

◆ reuse_preconditioner_max_linear_its()

template<typename T >
unsigned int libMesh::NonlinearSolver< T >::reuse_preconditioner_max_linear_its ( ) const
virtual

Get the reuse_preconditioner_max_linear_its parameter.

Reimplemented in libMesh::PetscNonlinearSolver< T >, and libMesh::PetscNonlinearSolver< Number >.

Definition at line 100 of file nonlinear_solver.C.

101 {
102  libmesh_not_implemented();
103 }

◆ set_exact_constraint_enforcement()

template<typename T>
virtual void libMesh::NonlinearSolver< T >::set_exact_constraint_enforcement ( bool  enable)
inlinevirtual

Enable (or disable; it is true by default) exact enforcement of constraints at the solver level, correcting any constrained DoF coefficients in current_local_solution as well as applying nonlinear residual and Jacobian terms based on constraint equations.

This is probably only safe to disable if user code is setting nonlinear residual and Jacobian terms based on constraint equations at an element-by-element level, by combining the asymmetric_constraint_rows option with the residual_constrain_element_vector processing option in DofMap.

Definition at line 398 of file nonlinear_solver.h.

399  {
401  }
bool _exact_constraint_enforcement
Whether we should enforce exact constraints globally during a solve.

◆ set_reuse_preconditioner()

template<typename T >
void libMesh::NonlinearSolver< T >::set_reuse_preconditioner ( bool  reuse)
virtual

Set the reuse preconditioner flag.

Definition at line 94 of file nonlinear_solver.C.

95 {
96  _reuse_preconditioner = reuse;
97 }
bool _reuse_preconditioner
Whether we should reuse the linear preconditioner.

◆ set_reuse_preconditioner_max_linear_its()

template<typename T >
void libMesh::NonlinearSolver< T >::set_reuse_preconditioner_max_linear_its ( unsigned int  i)
virtual

Set the reuse_preconditioner_max_linear_its parameter.

Definition at line 106 of file nonlinear_solver.C.

107 {
109 }
unsigned int _reuse_preconditioner_max_linear_its
Number of linear iterations to retain the preconditioner.

◆ set_solver_configuration()

template<typename T >
void libMesh::NonlinearSolver< T >::set_solver_configuration ( SolverConfiguration solver_configuration)

Set the solver configuration object.

Definition at line 82 of file nonlinear_solver.C.

83 {
84  _solver_configuration = &solver_configuration;
85 }
SolverConfiguration * _solver_configuration
Optionally store a SolverOptions object that can be used to set parameters like solver type...

◆ solve()

template<typename T>
virtual std::pair<unsigned int, Real> libMesh::NonlinearSolver< T >::solve ( SparseMatrix< T > &  ,
NumericVector< T > &  ,
NumericVector< T > &  ,
const double  ,
const unsigned  int 
)
pure virtual

◆ system() [1/2]

template<typename T>
const sys_type& libMesh::NonlinearSolver< T >::system ( ) const
inline

◆ system() [2/2]

template<typename T>
sys_type& libMesh::NonlinearSolver< T >::system ( )
inline
Returns
A writable reference to the system we are solving.

Definition at line 278 of file nonlinear_solver.h.

278 { return _system; }
sys_type & _system
A reference to the system we are solving.

Member Data Documentation

◆ _communicator

const Parallel::Communicator& libMesh::ParallelObject::_communicator
protectedinherited

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

Actually holds the data.

Definition at line 124 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::get_info().

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 143 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _exact_constraint_enforcement

template<typename T>
bool libMesh::NonlinearSolver< T >::_exact_constraint_enforcement
protected

◆ _is_initialized

template<typename T>
bool libMesh::NonlinearSolver< T >::_is_initialized
protected

Flag indicating if the data structures have been initialized.

Definition at line 433 of file nonlinear_solver.h.

Referenced by libMesh::NonlinearSolver< Number >::initialized().

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

◆ _n_objects

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects.

Print the reference count information when the number returns to 0.

Definition at line 132 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _preconditioner

template<typename T>
Preconditioner<T>* libMesh::NonlinearSolver< T >::_preconditioner
protected

Holds the Preconditioner object to be used for the linear solves.

Definition at line 438 of file nonlinear_solver.h.

◆ _reuse_preconditioner

template<typename T>
bool libMesh::NonlinearSolver< T >::_reuse_preconditioner
protected

Whether we should reuse the linear preconditioner.

Definition at line 412 of file nonlinear_solver.h.

◆ _reuse_preconditioner_max_linear_its

template<typename T>
unsigned int libMesh::NonlinearSolver< T >::_reuse_preconditioner_max_linear_its
protected

Number of linear iterations to retain the preconditioner.

Definition at line 423 of file nonlinear_solver.h.

◆ _solver_configuration

template<typename T>
SolverConfiguration* libMesh::NonlinearSolver< T >::_solver_configuration
protected

Optionally store a SolverOptions object that can be used to set parameters like solver type, tolerances and iteration limits.

Definition at line 444 of file nonlinear_solver.h.

◆ _system

template<typename T>
sys_type& libMesh::NonlinearSolver< T >::_system
protected

A reference to the system we are solving.

Definition at line 428 of file nonlinear_solver.h.

Referenced by libMesh::NonlinearSolver< Number >::system().

◆ absolute_residual_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::absolute_residual_tolerance

The NonlinearSolver should exit after the residual is reduced to either less than absolute_residual_tolerance or less than relative_residual_tolerance times the initial residual.

Users should increase any of these tolerances that they want to use for a stopping condition.

Definition at line 305 of file nonlinear_solver.h.

◆ absolute_step_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::absolute_step_tolerance

The NonlinearSolver should exit after the full nonlinear step norm is reduced to either less than absolute_step_tolerance or less than relative_step_tolerance times the largest nonlinear solution which has been seen so far.

Users should increase any of these tolerances that they want to use for a stopping condition.

Note
Not all NonlinearSolvers support relative_step_tolerance!

Definition at line 328 of file nonlinear_solver.h.

◆ bounds

template<typename T>
void(* libMesh::NonlinearSolver< T >::bounds) (NumericVector< Number > &XL, NumericVector< Number > &XU, sys_type &S)

Function that computes the lower and upper bounds XL and XU on the solution of the nonlinear system.

Definition at line 190 of file nonlinear_solver.h.

◆ bounds_object

template<typename T>
NonlinearImplicitSystem::ComputeBounds* libMesh::NonlinearSolver< T >::bounds_object

Object that computes the bounds vectors \( XL \) and \( XU \).

Definition at line 196 of file nonlinear_solver.h.

◆ converged

template<typename T>
bool libMesh::NonlinearSolver< T >::converged

After a call to solve this will reflect whether or not the nonlinear solve was successful.

Definition at line 352 of file nonlinear_solver.h.

◆ divergence_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::divergence_tolerance

The NonlinearSolver should exit if the residual becomes greater than the initial residual times the divergence_tolerance.

Users should adjust this tolerances to prevent divergence of the NonlinearSolver.

Definition at line 315 of file nonlinear_solver.h.

◆ fd_residual_object

template<typename T>
NonlinearImplicitSystem::ComputeResidual* libMesh::NonlinearSolver< T >::fd_residual_object

Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose of forming a finite-differenced Jacobian.

Definition at line 143 of file nonlinear_solver.h.

Referenced by libMesh::libmesh_petsc_snes_fd_residual().

◆ initial_linear_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::initial_linear_tolerance

Any required linear solves will at first be done with this tolerance; the NonlinearSolver may tighten the tolerance for later solves.

Definition at line 341 of file nonlinear_solver.h.

◆ jacobian

template<typename T>
void(* libMesh::NonlinearSolver< T >::jacobian) (const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S)

Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X.

Definition at line 156 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), and libMesh::libmesh_petsc_snes_jacobian().

◆ jacobian_object

template<typename T>
NonlinearImplicitSystem::ComputeJacobian* libMesh::NonlinearSolver< T >::jacobian_object

Object that computes the Jacobian J(X) of the nonlinear system at the input iterate X.

Definition at line 164 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), and libMesh::libmesh_petsc_snes_jacobian().

◆ matvec

template<typename T>
void(* libMesh::NonlinearSolver< T >::matvec) (const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S)

Function that computes either the residual \( R(X) \) or the Jacobian \( J(X) \) of the nonlinear system at the input iterate \( X \).

Note
Either R or J could be nullptr.

Definition at line 173 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::libmesh_petsc_snes_jacobian(), and libMesh::libmesh_petsc_snes_residual().

◆ max_function_evaluations

template<typename T>
unsigned int libMesh::NonlinearSolver< T >::max_function_evaluations

Maximum number of function evaluations.

Definition at line 293 of file nonlinear_solver.h.

◆ max_linear_iterations

template<typename T>
unsigned int libMesh::NonlinearSolver< T >::max_linear_iterations

Each linear solver step should exit after max_linear_iterations is exceeded.

Definition at line 335 of file nonlinear_solver.h.

◆ max_nonlinear_iterations

template<typename T>
unsigned int libMesh::NonlinearSolver< T >::max_nonlinear_iterations

Maximum number of non-linear iterations.

Definition at line 288 of file nonlinear_solver.h.

◆ mffd_residual_object

template<typename T>
NonlinearImplicitSystem::ComputeResidual* libMesh::NonlinearSolver< T >::mffd_residual_object

Object that computes the residual R(X) of the nonlinear system at the input iterate X for the purpose of forming Jacobian-vector products via finite differencing.

Definition at line 150 of file nonlinear_solver.h.

Referenced by libMesh::libmesh_petsc_snes_mffd_residual().

◆ minimum_linear_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::minimum_linear_tolerance

The tolerance for linear solves is kept above this minimum.

Definition at line 346 of file nonlinear_solver.h.

◆ nearnullspace

template<typename T>
void(* libMesh::NonlinearSolver< T >::nearnullspace) (std::vector< NumericVector< Number > * > &sp, sys_type &S)

Function that computes a basis for the Jacobian's near nullspace – the set of "low energy modes" – that can be used for AMG coarsening, if the solver supports it (e.g., ML, PETSc's GAMG).

Definition at line 233 of file nonlinear_solver.h.

◆ nearnullspace_object

template<typename T>
NonlinearImplicitSystem::ComputeVectorSubspace* libMesh::NonlinearSolver< T >::nearnullspace_object

A callable object that computes a basis for the Jacobian's near nullspace – the set of "low energy modes" – that can be used for AMG coarsening, if the solver supports it (e.g., ML, PETSc's GAMG).

Definition at line 240 of file nonlinear_solver.h.

◆ nullspace

template<typename T>
void(* libMesh::NonlinearSolver< T >::nullspace) (std::vector< NumericVector< Number > * > &sp, sys_type &S)

Function that computes a basis for the Jacobian's nullspace – the kernel or the "zero energy modes" – that can be used in solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP).

Definition at line 204 of file nonlinear_solver.h.

◆ nullspace_object

template<typename T>
NonlinearImplicitSystem::ComputeVectorSubspace* libMesh::NonlinearSolver< T >::nullspace_object

A callable object that computes a basis for the Jacobian's nullspace – the kernel or the "zero energy modes" – that can be used in solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP).

Definition at line 212 of file nonlinear_solver.h.

◆ postcheck

template<typename T>
void(* libMesh::NonlinearSolver< T >::postcheck) (const NumericVector< Number > &old_soln, NumericVector< Number > &search_direction, NumericVector< Number > &new_soln, bool &changed_search_direction, bool &changed_new_soln, sys_type &S)

Function that performs a "check" on the Newton search direction and solution after each nonlinear step.

See documentation for the NonlinearImplicitSystem::ComputePostCheck object for more information about the calling sequence.

Definition at line 254 of file nonlinear_solver.h.

Referenced by libMesh::libmesh_petsc_snes_postcheck().

◆ postcheck_object

template<typename T>
NonlinearImplicitSystem::ComputePostCheck* libMesh::NonlinearSolver< T >::postcheck_object

A callable object that is executed after each nonlinear iteration.

Allows the user to modify both the search direction and the solution vector in an application-specific way.

Definition at line 266 of file nonlinear_solver.h.

Referenced by libMesh::libmesh_petsc_snes_postcheck().

◆ precheck_object

template<typename T>
NonlinearImplicitSystem::ComputePreCheck* libMesh::NonlinearSolver< T >::precheck_object

Definition at line 268 of file nonlinear_solver.h.

Referenced by libMesh::libmesh_petsc_snes_precheck().

◆ relative_residual_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::relative_residual_tolerance

Definition at line 306 of file nonlinear_solver.h.

◆ relative_step_tolerance

template<typename T>
double libMesh::NonlinearSolver< T >::relative_step_tolerance

Definition at line 329 of file nonlinear_solver.h.

◆ residual

template<typename T>
void(* libMesh::NonlinearSolver< T >::residual) (const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S)

Function that computes the residual R(X) of the nonlinear system at the input iterate X.

Definition at line 129 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeF(), and libMesh::libmesh_petsc_snes_residual().

◆ residual_and_jacobian_object

template<typename T>
NonlinearImplicitSystem::ComputeResidualandJacobian* libMesh::NonlinearSolver< T >::residual_and_jacobian_object

Object that computes either the residual \( R(X) \) or the Jacobian \( J(X) \) of the nonlinear system at the input iterate \( X \).

Note
Either R or J could be nullptr.

Definition at line 185 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::libmesh_petsc_snes_jacobian(), and libMesh::libmesh_petsc_snes_residual().

◆ residual_object

template<typename T>
NonlinearImplicitSystem::ComputeResidual* libMesh::NonlinearSolver< T >::residual_object

Object that computes the residual R(X) of the nonlinear system at the input iterate X.

Definition at line 137 of file nonlinear_solver.h.

Referenced by libMesh::Problem_Interface::computeF(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_mffd_residual(), and libMesh::libmesh_petsc_snes_residual().

◆ transpose_nullspace

template<typename T>
void(* libMesh::NonlinearSolver< T >::transpose_nullspace) (std::vector< NumericVector< Number > * > &sp, sys_type &S)

Function that computes a basis for the transpose Jacobian's nullspace – when solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP), it is used to remove contributions outside of R(jac)

Definition at line 219 of file nonlinear_solver.h.

◆ transpose_nullspace_object

template<typename T>
NonlinearImplicitSystem::ComputeVectorSubspace* libMesh::NonlinearSolver< T >::transpose_nullspace_object

A callable object that computes a basis for the transpose Jacobian's nullspace – when solving a degenerate problem iteratively, if the solver supports it (e.g., PETSc's KSP), it is used to remove contributions outside of R(jac)

Definition at line 226 of file nonlinear_solver.h.

◆ user_presolve

template<typename T>
void(* libMesh::NonlinearSolver< T >::user_presolve) (sys_type &S)

Customizable function pointer which users can attach to the solver.

Gets called prior to every call to solve().

Definition at line 246 of file nonlinear_solver.h.


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