>> 1. Replace the contents of PetscLinearSolver::print_converged_reason() by
>>
>> KSPConvergedReason reason;
>> KSPGetConvergedReason(_ksp, &reason);
>> libMesh::out << "Linear solver convergence/divergence reason: " <<
>> KSPConvergedReasons[reason] << std::endl;
>>
>> or similar.
>
> This sounds like a good idea, if you want to test it and send us the
> patch.
>
>> 2. Remove the this->clear() call from PetscNonlinearSolver::solve().
>> This prevents the _snes context pointer from being destroyed. In
>> PetscLinearSolver, the _ksp context pointer isn't destroyed either.
>
> I'm up for this unless anyone objects.
>
>> 3. Introduce PetscNonlinearSolver::print_converged_reason()
>> analogously to the linear solver, using the maintained _snes pointer.
>
> Sounds good.
Since there have been no further objections, here's a patch that does
all of the above. I have tested both linear (with example 9 for
instance) and nonlinear print_converged_reason()s, they're working
smoothly in the scope of my testing. Also, no negative side-effects
have occurred to me from keeping the _snes context alive.
You may perhaps commit the new
PetscLinearSolver::print_converged_reason before the rest, if you
prefer.
Roman
Index: include/solvers/petsc_nonlinear_solver.h
===================================================================
--- include/solvers/petsc_nonlinear_solver.h (revision 4219)
+++ include/solvers/petsc_nonlinear_solver.h (working copy)
@@ -95,7 +95,11 @@
const double, // Stopping tolerance
const unsigned int); // N. Iterations
-
+ /**
+ * Prints a useful message about why the latest nonlinear solve
+ * con(di)verged.
+ */
+ virtual void print_converged_reason();
private:
Index: include/solvers/nonlinear_solver.h
===================================================================
--- include/solvers/nonlinear_solver.h (revision 4219)
+++ include/solvers/nonlinear_solver.h (working copy)
@@ -105,6 +105,12 @@
const unsigned int) = 0; // N. Iterations
/**
+ * Prints a useful message about why the latest nonlinear solve
+ * con(di)verged.
+ */
+ virtual void print_converged_reason() { libmesh_not_implemented(); }
+
+ /**
* Function that computes the residual \p R(X) of the nonlinear system
* at the input iterate \p X.
*/
Index: src/solvers/petsc_linear_solver.C
===================================================================
--- src/solvers/petsc_linear_solver.C (revision 4219)
+++ src/solvers/petsc_linear_solver.C (working copy)
@@ -1437,62 +1437,9 @@
template <typename T>
void PetscLinearSolver<T>::print_converged_reason()
{
-#if PETSC_VERSION_LESS_THAN(2,3,1)
- libMesh::out << "This method is currently not supported "
- << "(but may work!) for Petsc 2.3.0 and earlier." << std::endl;
-#else
KSPConvergedReason reason;
KSPGetConvergedReason(_ksp, &reason);
-
- // KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
- // KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
- // KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration)
- // KSP_CONVERGED_STEP_LENGTH
- // KSP_DIVERGED_ITS (required more than its to reach convergence)
- // KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
- // KSP_DIVERGED_NAN (residual norm became Not-a-number likely do to 0/0)
- // KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
-
- switch (reason)
- {
- case KSP_CONVERGED_RTOL:
- {
- libMesh::out << "Linear solver converged, relative tolerance reached." << std::endl;
- break;
- }
- case KSP_CONVERGED_ATOL:
- {
- libMesh::out << "Linear solver converged, absolute tolerance reached." << std::endl;
- break;
- }
-
- // Divergence
- case KSP_DIVERGED_ITS:
- {
- libMesh::out << "Linear solver diverged, max no. of iterations reached." << std::endl;
- break;
- }
- case KSP_DIVERGED_DTOL:
- {
- libMesh::out << "Linear solver diverged, residual norm increase by dtol (default 1.e5)." << std::endl;
- break;
- }
- case KSP_DIVERGED_NAN:
- {
- libMesh::out << "Linear solver diverged, residual norm is NaN." << std::endl;
- break;
- }
- case KSP_DIVERGED_BREAKDOWN:
- {
- libMesh::out << "Linear solver diverged, generic breakdown in the method." << std::endl;
- break;
- }
- default:
- {
- libMesh::out << "Unknown/unsupported con(di)vergence reason: " << reason << std::endl;
- }
- }
-#endif
+ libMesh::out << "Linear solver convergence/divergence reason: " << KSPConvergedReasons[reason] << std::endl;
}
Index: src/solvers/petsc_nonlinear_solver.C
===================================================================
--- src/solvers/petsc_nonlinear_solver.C (revision 4219)
+++ src/solvers/petsc_nonlinear_solver.C (working copy)
@@ -343,8 +343,6 @@
//Based on Petsc 2.3.3 documentation all diverged reasons are negative
this->converged = reason >= 0;
- this->clear();
-
STOP_LOG("solve()", "PetscNonlinearSolver");
// return the # of its. and the final residual norm.
@@ -353,7 +351,17 @@
+template <typename T>
+void PetscNonlinearSolver<T>::print_converged_reason()
+{
+ SNESConvergedReason reason;
+ SNESGetConvergedReason(_snes, &reason);
+ libMesh::out << "Nonlinear solver convergence/divergence reason: " << SNESConvergedReasons[reason] << std::endl;
+}
+
+
+
//------------------------------------------------------------------
// Explicit instantiations
template class PetscNonlinearSolver<Number>;
------------------------------------------------------------------------------
Free Software Download: Index, Search & Analyze Logs and other IT data in
Real-Time with Splunk. Collect, index and harness all the fast moving IT data
generated by your applications, servers and devices whether physical, virtual
or in the cloud. Deliver compliance at lower cost and gain new business
insights. http://p.sf.net/sfu/splunk-dev2dev
_______________________________________________
Libmesh-devel mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-devel