libMesh::PetscDMNonlinearSolver< T > Class Template Reference
#include <petsc_dm_nonlinear_solver.h>

Detailed Description
template<typename T>
class libMesh::PetscDMNonlinearSolver< T >
This class provides an interface to PETSc DM-based iterative solvers that is compatible with the libMesh NonlinearSolver<>
Definition at line 60 of file petsc_dm_nonlinear_solver.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
| typedef NonlinearImplicitSystem libMesh::PetscDMNonlinearSolver< T >::sys_type |
The type of system
Reimplemented from libMesh::PetscNonlinearSolver< T >.
Definition at line 66 of file petsc_dm_nonlinear_solver.h.
Constructor & Destructor Documentation
| libMesh::PetscDMNonlinearSolver< T >::PetscDMNonlinearSolver | ( | sys_type & | system | ) | [inline, explicit] |
Constructor. Initializes PETSc data structures
Definition at line 55 of file petsc_dm_nonlinear_solver.C.
References libMesh::PetscDMRegister().
00055 : 00056 PetscNonlinearSolver<T>(system_in) 00057 { 00058 PetscDMRegister(); 00059 }
| libMesh::PetscDMNonlinearSolver< T >::~PetscDMNonlinearSolver | ( | ) | [inline] |
Destructor.
Definition at line 63 of file petsc_dm_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::clear().
00064 { 00065 this->clear (); 00066 }
Member Function Documentation
| void libMesh::NonlinearSolver< T >::attach_preconditioner | ( | Preconditioner< T > * | preconditioner | ) | [inline, inherited] |
Attaches a Preconditioner object to be used during the linear solves.
Definition at line 91 of file nonlinear_solver.C.
References libMesh::NonlinearSolver< T >::_is_initialized, libMesh::NonlinearSolver< T >::_preconditioner, and libMesh::err.
00092 { 00093 if(this->_is_initialized) 00094 { 00095 libMesh::err << "Preconditioner must be attached before the solver is initialized!"<<std::endl; 00096 libmesh_error(); 00097 } 00098 00099 _preconditioner = preconditioner; 00100 }
| AutoPtr< NonlinearSolver< T > > libMesh::NonlinearSolver< T >::build | ( | sys_type & | s, | |
| const SolverPackage | solver_package = libMesh::default_solver_package() | |||
| ) | [inline, static, inherited] |
Builds a NonlinearSolver using the nonlinear solver package specified by solver_package
Definition at line 38 of file nonlinear_solver.C.
References libMesh::err, libMesh::on_command_line(), libMeshEnums::PETSC_SOLVERS, libMesh::AutoPtr< Tp >::reset(), and libMesh::TRILINOS_SOLVERS.
00039 { 00040 AutoPtr<NonlinearSolver<T> > ap; 00041 00042 // Build the appropriate solver 00043 switch (solver_package) 00044 { 00045 00046 #ifdef LIBMESH_HAVE_PETSC 00047 case PETSC_SOLVERS: 00048 #if PETSC_VERSION_LESS_THAN(3,3,0) 00049 ap.reset(new PetscNonlinearSolver<T>(s)); 00050 break; 00051 #else 00052 if (libMesh::on_command_line ("--use-petsc-dm")){ 00053 ap.reset(new PetscDMNonlinearSolver<T>(s)); 00054 } 00055 else { 00056 ap.reset(new PetscNonlinearSolver<T>(s)); 00057 } 00058 break; 00059 #endif 00060 #endif // LIBMESH_HAVE_PETSC 00061 00062 #ifdef LIBMESH_HAVE_NOX 00063 case TRILINOS_SOLVERS: 00064 ap.reset(new NoxNonlinearSolver<T>(s)); 00065 break; 00066 #endif 00067 00068 default: 00069 libMesh::err << "ERROR: Unrecognized solver package: " 00070 << solver_package 00071 << std::endl; 00072 libmesh_error(); 00073 } 00074 00075 return ap; 00076 }
| void libMesh::PetscNonlinearSolver< T >::clear | ( | ) | [inline, virtual, inherited] |
Release all memory and clear data structures.
Reimplemented from libMesh::NonlinearSolver< T >.
Definition at line 261 of file petsc_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::_current_nonlinear_iteration_number, libMesh::NonlinearSolver< T >::_is_initialized, libMesh::PetscNonlinearSolver< T >::_snes, libMesh::COMM_WORLD, and libMesh::NonlinearSolver< T >::initialized().
Referenced by libMesh::PetscNonlinearSolver< T >::solve(), libMesh::PetscDMNonlinearSolver< T >::solve(), libMesh::PetscDMNonlinearSolver< T >::~PetscDMNonlinearSolver(), and libMesh::PetscNonlinearSolver< T >::~PetscNonlinearSolver().
00262 { 00263 if (this->initialized()) 00264 { 00265 this->_is_initialized = false; 00266 00267 int ierr=0; 00268 00269 ierr = LibMeshSNESDestroy(&_snes); 00270 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00271 00272 // Reset the nonlinear iteration counter. This information is only relevant 00273 // *during* the solve(). After the solve is completed it should return to 00274 // the default value of 0. 00275 _current_nonlinear_iteration_number = 0; 00276 } 00277 }
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00107 { 00108 _enable_print_counter = false; 00109 return; 00110 }
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00101 { 00102 _enable_print_counter = true; 00103 return; 00104 }
| SNESConvergedReason libMesh::PetscNonlinearSolver< T >::get_converged_reason | ( | ) | [inline, inherited] |
Returns the currently-available (or most recently obtained, if the SNES object has been destroyed) convergence reason. Refer to PETSc docs for the meaning of different SNESConvergedReasons.
Definition at line 567 of file petsc_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::_reason, libMesh::PetscNonlinearSolver< T >::_snes, libMesh::COMM_WORLD, and libMesh::NonlinearSolver< T >::initialized().
Referenced by libMesh::PetscNonlinearSolver< T >::print_converged_reason().
00568 { 00569 int ierr=0; 00570 00571 if (this->initialized()) 00572 { 00573 ierr = SNESGetConvergedReason(_snes, &_reason); 00574 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00575 } 00576 00577 return _reason; 00578 }
| virtual unsigned libMesh::PetscNonlinearSolver< T >::get_current_nonlinear_iteration_number | ( | ) | const [inline, virtual, inherited] |
If called *during* the solve(), for example by the user-specified residual or Jacobian function, returns the current nonlinear iteration number.
Implements libMesh::NonlinearSolver< T >.
Definition at line 126 of file petsc_nonlinear_solver.h.
References libMesh::PetscNonlinearSolver< T >::_current_nonlinear_iteration_number.
00126 { return _current_nonlinear_iteration_number; }
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
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().
00048 { 00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00050 00051 std::ostringstream oss; 00052 00053 oss << '\n' 00054 << " ---------------------------------------------------------------------------- \n" 00055 << "| Reference count information |\n" 00056 << " ---------------------------------------------------------------------------- \n"; 00057 00058 for (Counts::iterator it = _counts.begin(); 00059 it != _counts.end(); ++it) 00060 { 00061 const std::string name(it->first); 00062 const unsigned int creations = it->second.first; 00063 const unsigned int destructions = it->second.second; 00064 00065 oss << "| " << name << " reference count information:\n" 00066 << "| Creations: " << creations << '\n' 00067 << "| Destructions: " << destructions << '\n'; 00068 } 00069 00070 oss << " ---------------------------------------------------------------------------- \n"; 00071 00072 return oss.str(); 00073 00074 #else 00075 00076 return ""; 00077 00078 #endif 00079 }
| int libMesh::PetscNonlinearSolver< T >::get_total_linear_iterations | ( | ) | [inline, virtual, inherited] |
Get the total number of linear iterations done in the last solve
Implements libMesh::NonlinearSolver< T >.
Definition at line 581 of file petsc_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::_n_linear_iterations.
00582 { 00583 return _n_linear_iterations; 00584 }
| void libMesh::ReferenceCounter::increment_constructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
00164 { 00165 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00166 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00167 00168 p.first++; 00169 }
| void libMesh::ReferenceCounter::increment_destructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
00177 { 00178 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00179 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00180 00181 p.second++; 00182 }
| void libMesh::PetscDMNonlinearSolver< T >::init | ( | ) | [inline, virtual] |
Initialize data structures if not done so already.
Reimplemented from libMesh::PetscNonlinearSolver< T >.
Definition at line 71 of file petsc_dm_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::_snes, libMesh::NonlinearSolver< T >::absolute_residual_tolerance, libMesh::NonlinearSolver< T >::absolute_step_tolerance, libMesh::COMM_WORLD, DMCreateLibMesh(), libMesh::NonlinearSolver< T >::initial_linear_tolerance, libMesh::NonlinearSolver< T >::max_function_evaluations, libMesh::NonlinearSolver< T >::max_linear_iterations, libMesh::NonlinearSolver< T >::max_nonlinear_iterations, libMesh::NonlinearSolver< T >::relative_residual_tolerance, and libMesh::NonlinearSolver< T >::system().
Referenced by libMesh::PetscDMNonlinearSolver< T >::solve().
00072 { 00073 PetscErrorCode ierr; 00074 DM dm; 00075 this->PetscNonlinearSolver<T>::init(); 00076 00077 // Attaching a DM with the function and Jacobian callbacks to SNES. 00078 ierr = DMCreateLibMesh(libMesh::COMM_WORLD, this->system(), &dm); CHKERRABORT(libMesh::COMM_WORLD, ierr); 00079 ierr = DMSetFromOptions(dm); CHKERRABORT(libMesh::COMM_WORLD, ierr); 00080 ierr = DMSetUp(dm); CHKERRABORT(libMesh::COMM_WORLD, ierr); 00081 ierr = SNESSetDM(this->_snes, dm); CHKERRABORT(libMesh::COMM_WORLD, ierr); 00082 // SNES now owns the reference to dm. 00083 ierr = DMDestroy(&dm); CHKERRABORT(libMesh::COMM_WORLD, ierr); 00084 KSP ksp; 00085 ierr = SNESGetKSP (this->_snes, &ksp); CHKERRABORT(libMesh::COMM_WORLD,ierr); 00086 00087 // Set the tolerances for the iterative solver. Use the user-supplied 00088 // tolerance for the relative residual & leave the others at default values 00089 ierr = KSPSetTolerances (ksp, this->initial_linear_tolerance, PETSC_DEFAULT,PETSC_DEFAULT, this->max_linear_iterations); CHKERRABORT(libMesh::COMM_WORLD,ierr); 00090 00091 // Set the tolerances for the non-linear solver. 00092 ierr = SNESSetTolerances(this->_snes, 00093 this->absolute_residual_tolerance, 00094 this->relative_residual_tolerance, 00095 this->absolute_step_tolerance, 00096 this->max_nonlinear_iterations, 00097 this->max_function_evaluations); 00098 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00099 00100 //Pull in command-line options 00101 KSPSetFromOptions(ksp); 00102 SNESSetFromOptions(this->_snes); 00103 }
| bool libMesh::NonlinearSolver< T >::initialized | ( | ) | const [inline, inherited] |
- Returns:
- true if the data structures are initialized, false otherwise.
Definition at line 84 of file nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::clear(), libMesh::PetscNonlinearSolver< T >::get_converged_reason(), libMesh::NoxNonlinearSolver< T >::init(), and libMesh::PetscNonlinearSolver< T >::init().
00084 { return _is_initialized; }
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
00080 { return _n_objects; }
| void libMesh::PetscNonlinearSolver< T >::print_converged_reason | ( | ) | [inline, virtual, inherited] |
Prints a useful message about why the latest nonlinear solve con(di)verged.
Reimplemented from libMesh::NonlinearSolver< T >.
Definition at line 557 of file petsc_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::get_converged_reason(), and libMesh::out.
00558 { 00559 00560 libMesh::out << "Nonlinear solver convergence/divergence reason: " 00561 << SNESConvergedReasons[this->get_converged_reason()] << std::endl; 00562 }
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out |
) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
00089 { 00090 if( _enable_print_counter ) out_stream << ReferenceCounter::get_info(); 00091 }
| void libMesh::PetscNonlinearSolver< T >::set_current_nonlinear_iteration_number | ( | unsigned | num | ) | [inline, inherited] |
This public setter is necessary since the value is computed in the __libmesh_petsc_snes_residual()/jacobian() function and must be stored somehow.
Definition at line 133 of file petsc_nonlinear_solver.h.
References libMesh::PetscNonlinearSolver< T >::_current_nonlinear_iteration_number.
00133 { _current_nonlinear_iteration_number = num; }
| SNES libMesh::PetscNonlinearSolver< T >::snes | ( | ) | [inline, inherited] |
Returns the raw PETSc snes context pointer.
Definition at line 91 of file petsc_nonlinear_solver.h.
References libMesh::PetscNonlinearSolver< T >::_snes, and libMesh::PetscNonlinearSolver< T >::init().
| std::pair< unsigned int, Real > libMesh::PetscDMNonlinearSolver< T >::solve | ( | SparseMatrix< T > & | jac_in, | |
| NumericVector< T > & | x_in, | |||
| NumericVector< T > & | r_in, | |||
| const | double, | |||
| const unsigned | int | |||
| ) | [inline, virtual] |
Call the Petsc solver. It calls the method below, using the same matrix for the system and preconditioner matrices.
Reimplemented from libMesh::PetscNonlinearSolver< T >.
Definition at line 108 of file petsc_dm_nonlinear_solver.C.
References libMesh::PetscNonlinearSolver< T >::_n_linear_iterations, libMesh::NonlinearSolver< T >::_preconditioner, libMesh::PetscNonlinearSolver< T >::_reason, libMesh::PetscNonlinearSolver< T >::_snes, libMesh::PetscNonlinearSolver< T >::clear(), libMesh::COMM_WORLD, libMesh::NonlinearSolver< T >::converged, libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::Real, libMesh::NonlinearSolver< T >::system(), and libMesh::NonlinearSolver< T >::user_presolve.
00113 { 00114 START_LOG("solve()", "PetscNonlinearSolver"); 00115 this->init (); 00116 00117 // Make sure the data passed in are really of Petsc types 00118 libmesh_cast_ptr<PetscMatrix<T>*>(&jac_in); 00119 libmesh_cast_ptr<PetscVector<T>*>(&r_in); 00120 00121 // Extract solution vector 00122 PetscVector<T>* x = libmesh_cast_ptr<PetscVector<T>*>(&x_in); 00123 00124 int ierr=0; 00125 int n_iterations =0; 00126 00127 // Should actually be a PetscReal, but I don't know which version of PETSc first introduced PetscReal 00128 Real final_residual_norm=0.; 00129 00130 if (this->user_presolve) 00131 this->user_presolve(this->system()); 00132 00133 //Set the preconditioning matrix 00134 if (this->_preconditioner) 00135 this->_preconditioner->set_matrix(jac_in); 00136 00137 ierr = SNESSolve (this->_snes, PETSC_NULL, x->vec()); 00138 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00139 00140 ierr = SNESGetIterationNumber(this->_snes,&n_iterations); 00141 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00142 00143 ierr = SNESGetLinearSolveIterations(this->_snes, &this->_n_linear_iterations); 00144 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00145 00146 ierr = SNESGetFunctionNorm(this->_snes,&final_residual_norm); 00147 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00148 00149 // Get and store the reason for convergence 00150 SNESGetConvergedReason(this->_snes, &this->_reason); 00151 00152 //Based on Petsc 2.3.3 documentation all diverged reasons are negative 00153 this->converged = (this->_reason >= 0); 00154 00155 this->clear(); 00156 00157 STOP_LOG("solve()", "PetscNonlinearSolver"); 00158 00159 // return the # of its. and the final residual norm. 00160 return std::make_pair(n_iterations, final_residual_norm); 00161 }
| sys_type& libMesh::NonlinearSolver< T >::system | ( | ) | [inline, inherited] |
- Returns:
- a writeable reference to the system we are solving.
Definition at line 222 of file nonlinear_solver.h.
00222 { return _system; }
| const sys_type& libMesh::NonlinearSolver< T >::system | ( | ) | const [inline, inherited] |
- Returns:
- a constant reference to the system we are solving.
Definition at line 217 of file nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::build_mat_null_space(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
00217 { return _system; }
Member Data Documentation
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
unsigned libMesh::PetscNonlinearSolver< T >::_current_nonlinear_iteration_number [inherited] |
Stores the current nonlinear iteration number
Definition at line 157 of file petsc_nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::clear(), libMesh::PetscNonlinearSolver< T >::get_current_nonlinear_iteration_number(), and libMesh::PetscNonlinearSolver< T >::set_current_nonlinear_iteration_number().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
bool libMesh::NonlinearSolver< T >::_is_initialized [protected, inherited] |
Flag indicating if the data structures have been initialized.
Definition at line 298 of file nonlinear_solver.h.
Referenced by libMesh::NonlinearSolver< T >::attach_preconditioner(), libMesh::PetscNonlinearSolver< T >::clear(), libMesh::PetscNonlinearSolver< T >::init(), and libMesh::NonlinearSolver< Number >::initialized().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
int libMesh::PetscNonlinearSolver< T >::_n_linear_iterations [inherited] |
Stores the total number of linear iterations from the last solve.
Definition at line 152 of file petsc_nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::get_total_linear_iterations(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
Preconditioner<T>* libMesh::NonlinearSolver< T >::_preconditioner [protected, inherited] |
Holds the Preconditioner object to be used for the linear solves.
Definition at line 303 of file nonlinear_solver.h.
Referenced by libMesh::NonlinearSolver< T >::attach_preconditioner(), libMesh::PetscNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
SNESConvergedReason libMesh::PetscNonlinearSolver< T >::_reason [inherited] |
Store the reason for SNES convergence/divergence for use even after the _snes has been cleared. Note that print_converged_reason() will always *try* to get the current reason with SNESGetConvergedReason(), but if the SNES object has already been cleared, it will fall back on this stored value. Note that this value is therefore necessarily *not* cleared by the clear() function.
Definition at line 147 of file petsc_nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::get_converged_reason(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
SNES libMesh::PetscNonlinearSolver< T >::_snes [inherited] |
Nonlinear solver context
Definition at line 138 of file petsc_nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::clear(), libMesh::PetscNonlinearSolver< T >::get_converged_reason(), libMesh::PetscNonlinearSolver< T >::init(), libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::PetscNonlinearSolver< T >::snes(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
sys_type& libMesh::NonlinearSolver< T >::_system [protected, inherited] |
A reference to the system we are solving.
Definition at line 293 of file nonlinear_solver.h.
Referenced by libMesh::NonlinearSolver< Number >::system().
Real libMesh::NonlinearSolver< T >::absolute_residual_tolerance [inherited] |
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 249 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
Real libMesh::NonlinearSolver< T >::absolute_step_tolerance [inherited] |
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 that not all NonlinearSolvers support relative_step_tolerance!
Definition at line 263 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), and libMesh::NoxNonlinearSolver< T >::solve().
void(* libMesh::NonlinearSolver< T >::bounds)(NumericVector< Number > &XL, NumericVector< Number > &XU, sys_type &S) [inherited] |
Function that computes the lower and upper bounds XL and XU on the solution of the nonlinear system.
NonlinearImplicitSystem::ComputeBounds* libMesh::NonlinearSolver< T >::bounds_object [inherited] |
Object that computes the bounds vectors
and
.
Definition at line 179 of file nonlinear_solver.h.
bool libMesh::NonlinearSolver< T >::converged [inherited] |
After a call to solve this will reflect whether or not the nonlinear solve was successful.
Definition at line 287 of file nonlinear_solver.h.
Referenced by libMesh::NoxNonlinearSolver< T >::solve(), libMesh::PetscNonlinearSolver< T >::solve(), and libMesh::PetscDMNonlinearSolver< T >::solve().
Real libMesh::NonlinearSolver< T >::initial_linear_tolerance [inherited] |
Any required linear solves will at first be done with this tolerance; the NonlinearSolver may tighten the tolerance for later solves.
Definition at line 276 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
void(* libMesh::NonlinearSolver< T >::jacobian)(const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S) [inherited] |
Function that computes the Jacobian J(X) of the nonlinear system at the input iterate X.
Referenced by libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
NonlinearImplicitSystem::ComputeJacobian* libMesh::NonlinearSolver< T >::jacobian_object [inherited] |
Object that computes the Jacobian J(X) of the nonlinear system at the input iterate X.
Definition at line 149 of file nonlinear_solver.h.
Referenced by libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
void(* libMesh::NonlinearSolver< T >::matvec)(const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S) [inherited] |
Function that computes either the residual
or the Jacobian
of the nonlinear system at the input iterate
. Note that either R or J could be XSNULL.
Referenced by libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), and libMesh::Problem_Interface::computePreconditioner().
unsigned int libMesh::NonlinearSolver< T >::max_function_evaluations [inherited] |
Maximum number of function evaluations.
Definition at line 237 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), and libMesh::PetscNonlinearSolver< T >::solve().
unsigned int libMesh::NonlinearSolver< T >::max_linear_iterations [inherited] |
Each linear solver step should exit after max_linear_iterations is exceeded.
Definition at line 270 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
unsigned int libMesh::NonlinearSolver< T >::max_nonlinear_iterations [inherited] |
Maximum number of non-linear iterations.
Definition at line 232 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
Real libMesh::NonlinearSolver< T >::minimum_linear_tolerance [inherited] |
The tolerance for linear solves is kept above this minimum
Definition at line 281 of file nonlinear_solver.h.
void(* libMesh::NonlinearSolver< T >::nearnullspace)(std::vector< NumericVector< Number > * > &sp, sys_type &S) [inherited] |
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).
Referenced by libMesh::PetscNonlinearSolver< T >::solve().
NonlinearImplicitSystem::ComputeVectorSubspace* libMesh::NonlinearSolver< T >::nearnullspace_object [inherited] |
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 209 of file nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::solve().
void(* libMesh::NonlinearSolver< T >::nullspace)(std::vector< NumericVector< Number > * > &sp, sys_type &S) [inherited] |
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).
Referenced by libMesh::PetscNonlinearSolver< T >::solve().
NonlinearImplicitSystem::ComputeVectorSubspace* libMesh::NonlinearSolver< T >::nullspace_object [inherited] |
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 195 of file nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::solve().
Real libMesh::NonlinearSolver< T >::relative_residual_tolerance [inherited] |
Definition at line 250 of file nonlinear_solver.h.
Referenced by libMesh::PetscDMNonlinearSolver< T >::init(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
Real libMesh::NonlinearSolver< T >::relative_step_tolerance [inherited] |
Definition at line 264 of file nonlinear_solver.h.
Referenced by libMesh::PetscNonlinearSolver< T >::solve().
void(* libMesh::NonlinearSolver< T >::residual)(const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S) [inherited] |
Function that computes the residual R(X) of the nonlinear system at the input iterate X.
Referenced by libMesh::Problem_Interface::computeF().
NonlinearImplicitSystem::ComputeResidualandJacobian* libMesh::NonlinearSolver< T >::residual_and_jacobian_object [inherited] |
Object that computes either the residual
or the Jacobian
of the nonlinear system at the input iterate
. Note that either R or J could be XSNULL.
Definition at line 168 of file nonlinear_solver.h.
Referenced by libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::NoxNonlinearSolver< T >::solve(), and libMesh::PetscNonlinearSolver< T >::solve().
NonlinearImplicitSystem::ComputeResidual* libMesh::NonlinearSolver< T >::residual_object [inherited] |
Object that computes the residual R(X) of the nonlinear system at the input iterate X.
Definition at line 135 of file nonlinear_solver.h.
Referenced by libMesh::Problem_Interface::computeF().
void(* libMesh::NonlinearSolver< T >::user_presolve)(sys_type &S) [inherited] |
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:33 UTC
Hosted By: