Rename LinearSolverTerminationType enums.

This increases clarity, drops redundant enums and makes things
cleaner all around.

Change-Id: I761f195ddf17ea6bd8e4e55bf5a72863660c4c3b
diff --git a/internal/ceres/conjugate_gradients_solver.cc b/internal/ceres/conjugate_gradients_solver.cc
index fc1aef4..524cb8a 100644
--- a/internal/ceres/conjugate_gradients_solver.cc
+++ b/internal/ceres/conjugate_gradients_solver.cc
@@ -74,7 +74,7 @@
   CHECK_EQ(A->num_rows(), A->num_cols());
 
   LinearSolver::Summary summary;
-  summary.termination_type = MAX_ITERATIONS;
+  summary.termination_type = LINEAR_SOLVER_NO_CONVERGENCE;
   summary.message = "Maximum number of iterations reached.";
   summary.num_iterations = 0;
 
@@ -85,7 +85,7 @@
   const double norm_b = bref.norm();
   if (norm_b == 0.0) {
     xref.setZero();
-    summary.termination_type = TOLERANCE;
+    summary.termination_type = LINEAR_SOLVER_SUCCESS;
     summary.message = "Convergence. |b| = 0.";
     return summary;
   }
@@ -102,7 +102,7 @@
   r = bref - tmp;
   double norm_r = r.norm();
   if (norm_r <= tol_r) {
-    summary.termination_type = TOLERANCE;
+    summary.termination_type = LINEAR_SOLVER_SUCCESS;
     summary.message =
         StringPrintf("Convergence. |r| = %e <= %e.", norm_r, tol_r);
     return summary;
@@ -127,7 +127,7 @@
     double last_rho = rho;
     rho = r.dot(z);
     if (IsZeroOrInfinity(rho)) {
-      summary.termination_type = FAILURE;
+      summary.termination_type = LINEAR_SOLVER_FAILURE;
       summary.message = StringPrintf("Numerical failure. rho = r'z = %e.", rho);
       break;
     };
@@ -137,7 +137,7 @@
     } else {
       double beta = rho / last_rho;
       if (IsZeroOrInfinity(beta)) {
-        summary.termination_type = FAILURE;
+        summary.termination_type = LINEAR_SOLVER_FAILURE;
         summary.message = StringPrintf(
             "Numerical failure. beta = rho_n / rho_{n-1} = %e.", beta);
         break;
@@ -150,14 +150,14 @@
     A->RightMultiply(p.data(), q.data());
     const double pq = p.dot(q);
     if ((pq <= 0) || IsInfinite(pq))  {
-      summary.termination_type = FAILURE;
+      summary.termination_type = LINEAR_SOLVER_FAILURE;
       summary.message = StringPrintf("Numerical failure. p'q = %e.", pq);
       break;
     }
 
     const double alpha = rho / pq;
     if (IsInfinite(alpha)) {
-      summary.termination_type = FAILURE;
+      summary.termination_type = LINEAR_SOLVER_FAILURE;
       summary.message =
           StringPrintf("Numerical failure. alpha = rho / pq = %e", alpha);
       break;
@@ -208,7 +208,7 @@
     //
     const double zeta = summary.num_iterations * (Q1 - Q0) / Q1;
     if (zeta < per_solve_options.q_tolerance) {
-      summary.termination_type = TOLERANCE;
+      summary.termination_type = LINEAR_SOLVER_SUCCESS;
       summary.message =
           StringPrintf("Convergence: zeta = %e < %e",
                        zeta,
@@ -220,7 +220,7 @@
     // Residual based termination.
     norm_r = r. norm();
     if (norm_r <= tol_r) {
-      summary.termination_type = TOLERANCE;
+      summary.termination_type = LINEAR_SOLVER_SUCCESS;
       summary.message =
           StringPrintf("Convergence. |r| = %e <= %e.", norm_r, tol_r);
       break;
diff --git a/internal/ceres/covariance_impl.cc b/internal/ceres/covariance_impl.cc
index 1bca3c3..91f0393 100644
--- a/internal/ceres/covariance_impl.cc
+++ b/internal/ceres/covariance_impl.cc
@@ -454,7 +454,7 @@
   LinearSolverTerminationType termination_type =
       ss.Cholesky(&cholmod_jacobian_view, factor, &status);
   event_logger.AddEvent("Numeric Factorization");
-  if (termination_type != TOLERANCE) {
+  if (termination_type != LINEAR_SOLVER_SUCCESS) {
     LOG(ERROR) << "Covariance estimation failed. "
                << "CHOLMOD numeric cholesky factorization returned with: "
                << status;
diff --git a/internal/ceres/dense_normal_cholesky_solver.cc b/internal/ceres/dense_normal_cholesky_solver.cc
index b719c94..f44d6da 100644
--- a/internal/ceres/dense_normal_cholesky_solver.cc
+++ b/internal/ceres/dense_normal_cholesky_solver.cc
@@ -95,15 +95,15 @@
 
   LinearSolver::Summary summary;
   summary.num_iterations = 1;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   Eigen::LLT<Matrix, Eigen::Upper> llt =
       lhs.selfadjointView<Eigen::Upper>().llt();
 
   if (llt.info() != Eigen::Success) {
-    summary.termination_type = FAILURE;
+    summary.termination_type = LINEAR_SOLVER_FAILURE;
     summary.message = "Eigen LLT decomposition failed.";
   } else {
-    summary.termination_type = TOLERANCE;
+    summary.termination_type = LINEAR_SOLVER_SUCCESS;
     summary.message = "Success.";
   }
 
diff --git a/internal/ceres/dense_qr_solver.cc b/internal/ceres/dense_qr_solver.cc
index f8927ae..4388357 100644
--- a/internal/ceres/dense_qr_solver.cc
+++ b/internal/ceres/dense_qr_solver.cc
@@ -111,7 +111,7 @@
                                                          rhs_.data(),
                                                          &summary.message);
   event_logger.AddEvent("Solve");
-  if (summary.termination_type == TOLERANCE) {
+  if (summary.termination_type == LINEAR_SOLVER_SUCCESS) {
     VectorRef(x, num_cols) = rhs_.head(num_cols);
   }
 
@@ -159,7 +159,7 @@
   // is good enough or not.
   LinearSolver::Summary summary;
   summary.num_iterations = 1;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.message = "Success.";
 
   event_logger.AddEvent("TearDown");
diff --git a/internal/ceres/dogleg_strategy.cc b/internal/ceres/dogleg_strategy.cc
index d2cab86..f29376d 100644
--- a/internal/ceres/dogleg_strategy.cc
+++ b/internal/ceres/dogleg_strategy.cc
@@ -99,7 +99,7 @@
     }
     TrustRegionStrategy::Summary summary;
     summary.num_iterations = 0;
-    summary.termination_type = TOLERANCE;
+    summary.termination_type = LINEAR_SOLVER_SUCCESS;
     return summary;
   }
 
@@ -135,11 +135,11 @@
   summary.num_iterations = linear_solver_summary.num_iterations;
   summary.termination_type = linear_solver_summary.termination_type;
 
-  if (linear_solver_summary.termination_type == FATAL_ERROR) {
+  if (linear_solver_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
     return summary;
   }
 
-  if (linear_solver_summary.termination_type != FAILURE) {
+  if (linear_solver_summary.termination_type != LINEAR_SOLVER_FAILURE) {
     switch (dogleg_type_) {
       // Interpolate the Cauchy point and the Gauss-Newton step.
       case TRADITIONAL_DOGLEG:
@@ -150,7 +150,7 @@
       // Cauchy point and the (Gauss-)Newton step.
       case SUBSPACE_DOGLEG:
         if (!ComputeSubspaceModel(jacobian)) {
-          summary.termination_type = FAILURE;
+          summary.termination_type = LINEAR_SOLVER_FAILURE;
           break;
         }
         ComputeSubspaceDoglegStep(step);
@@ -517,7 +517,7 @@
     const double* residuals) {
   const int n = jacobian->num_cols();
   LinearSolver::Summary linear_solver_summary;
-  linear_solver_summary.termination_type = FAILURE;
+  linear_solver_summary.termination_type = LINEAR_SOLVER_FAILURE;
 
   // The Jacobian matrix is often quite poorly conditioned. Thus it is
   // necessary to add a diagonal matrix at the bottom to prevent the
@@ -530,7 +530,7 @@
   // If the solve fails, the multiplier to the diagonal is increased
   // up to max_mu_ by a factor of mu_increase_factor_ every time. If
   // the linear solver is still not successful, the strategy returns
-  // with FAILURE.
+  // with LINEAR_SOLVER_FAILURE.
   //
   // Next time when a new Gauss-Newton step is requested, the
   // multiplier starts out from the last successful solve.
@@ -583,21 +583,21 @@
       }
     }
 
-    if (linear_solver_summary.termination_type == FATAL_ERROR) {
+    if (linear_solver_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
       return linear_solver_summary;
     }
 
-    if (linear_solver_summary.termination_type == FAILURE ||
+    if (linear_solver_summary.termination_type == LINEAR_SOLVER_FAILURE ||
         !IsArrayValid(n, gauss_newton_step_.data())) {
       mu_ *= mu_increase_factor_;
       VLOG(2) << "Increasing mu " << mu_;
-      linear_solver_summary.termination_type = FAILURE;
+      linear_solver_summary.termination_type = LINEAR_SOLVER_FAILURE;
       continue;
     }
     break;
   }
 
-  if (linear_solver_summary.termination_type != FAILURE) {
+  if (linear_solver_summary.termination_type != LINEAR_SOLVER_FAILURE) {
     // The scaled Gauss-Newton step is D * GN:
     //
     //     - (D^-1 J^T J D^-1)^-1 (D^-1 g)
diff --git a/internal/ceres/dogleg_strategy_test.cc b/internal/ceres/dogleg_strategy_test.cc
index ace635f..795719d 100644
--- a/internal/ceres/dogleg_strategy_test.cc
+++ b/internal/ceres/dogleg_strategy_test.cc
@@ -144,7 +144,7 @@
                                                               residual_.data(),
                                                               x_.data());
 
-  EXPECT_NE(summary.termination_type, FAILURE);
+  EXPECT_NE(summary.termination_type, LINEAR_SOLVER_FAILURE);
   EXPECT_LE(x_.norm(), options_.initial_radius * (1.0 + 4.0 * kEpsilon));
 }
 
@@ -164,7 +164,7 @@
                                                               residual_.data(),
                                                               x_.data());
 
-  EXPECT_NE(summary.termination_type, FAILURE);
+  EXPECT_NE(summary.termination_type, LINEAR_SOLVER_FAILURE);
   EXPECT_LE(x_.norm(), options_.initial_radius * (1.0 + 4.0 * kEpsilon));
 }
 
@@ -184,7 +184,7 @@
                                                               residual_.data(),
                                                               x_.data());
 
-  EXPECT_NE(summary.termination_type, FAILURE);
+  EXPECT_NE(summary.termination_type, LINEAR_SOLVER_FAILURE);
   EXPECT_NEAR(x_(0), 1.0, kToleranceLoose);
   EXPECT_NEAR(x_(1), 1.0, kToleranceLoose);
   EXPECT_NEAR(x_(2), 1.0, kToleranceLoose);
@@ -246,7 +246,7 @@
                                                               residual_.data(),
                                                               x_.data());
 
-  EXPECT_NE(summary.termination_type, FAILURE);
+  EXPECT_NE(summary.termination_type, LINEAR_SOLVER_FAILURE);
   EXPECT_NEAR(x_(0), 0.0, kToleranceLoose);
   EXPECT_NEAR(x_(1), 0.0, kToleranceLoose);
   EXPECT_NEAR(x_(2), options_.initial_radius, kToleranceLoose);
@@ -274,7 +274,7 @@
                                                               residual_.data(),
                                                               x_.data());
 
-  EXPECT_NE(summary.termination_type, FAILURE);
+  EXPECT_NE(summary.termination_type, LINEAR_SOLVER_FAILURE);
   EXPECT_NEAR(x_(0), 0.0, kToleranceLoose);
   EXPECT_NEAR(x_(1), 0.0, kToleranceLoose);
   EXPECT_NEAR(x_(2), 1.0, kToleranceLoose);
diff --git a/internal/ceres/gmock_gtest_all.cc b/internal/ceres/gmock_gtest_all.cc
index 1a63a8c..e2aabf2 100644
--- a/internal/ceres/gmock_gtest_all.cc
+++ b/internal/ceres/gmock_gtest_all.cc
@@ -4364,7 +4364,7 @@
   int num_disabled = unit_test.reportable_disabled_test_count();
   if (num_disabled && !GTEST_FLAG(also_run_disabled_tests)) {
     if (!num_failures) {
-      printf("\n");  // Add a spacer if no FAILURE banner is displayed.
+      printf("\n");  // Add a spacer if no LINEAR_SOLVER_FAILURE banner is displayed.
     }
     ColoredPrintf(COLOR_YELLOW,
                   "  YOU HAVE %d DISABLED %s\n\n",
diff --git a/internal/ceres/iterative_schur_complement_solver.cc b/internal/ceres/iterative_schur_complement_solver.cc
index c4c77af..6de410b 100644
--- a/internal/ceres/iterative_schur_complement_solver.cc
+++ b/internal/ceres/iterative_schur_complement_solver.cc
@@ -88,7 +88,7 @@
     VLOG(2) << "No parameter blocks left in the schur complement.";
     LinearSolver::Summary cg_summary;
     cg_summary.num_iterations = 0;
-    cg_summary.termination_type = TOLERANCE;
+    cg_summary.termination_type = LINEAR_SOLVER_SUCCESS;
     schur_complement_->BackSubstitute(NULL, x);
     return cg_summary;
   }
@@ -157,7 +157,7 @@
 
   LinearSolver::Summary cg_summary;
   cg_summary.num_iterations = 0;
-  cg_summary.termination_type = FAILURE;
+  cg_summary.termination_type = LINEAR_SOLVER_FAILURE;
 
   // TODO(sameeragarwal): Refactor preconditioners to return a more
   // sane message.
@@ -167,8 +167,8 @@
                                  schur_complement_->rhs().data(),
                                  cg_per_solve_options,
                                  reduced_linear_system_solution_.data());
-    if (cg_summary.termination_type != FAILURE &&
-        cg_summary.termination_type != FATAL_ERROR) {
+    if (cg_summary.termination_type != LINEAR_SOLVER_FAILURE &&
+        cg_summary.termination_type != LINEAR_SOLVER_FATAL_ERROR) {
       schur_complement_->BackSubstitute(
           reduced_linear_system_solution_.data(), x);
     }
diff --git a/internal/ceres/lapack.cc b/internal/ceres/lapack.cc
index 0061433..c4f9302 100644
--- a/internal/ceres/lapack.cc
+++ b/internal/ceres/lapack.cc
@@ -73,7 +73,7 @@
     string* status) {
 #ifdef CERES_NO_LAPACK
   LOG(FATAL) << "Ceres was built without a BLAS library.";
-  return FATAL_ERROR;
+  return LINEAR_SOLVER_FATAL_ERROR;
 #else
   char uplo = 'L';
   int n = num_rows;
@@ -87,7 +87,7 @@
                << "Please report it."
                << "LAPACK::dpotrf fatal error."
                << "Argument: " << -info << " is invalid.";
-    return FATAL_ERROR;
+    return LINEAR_SOLVER_FATAL_ERROR;
   }
 
   if (info > 0) {
@@ -95,7 +95,7 @@
         StringPrintf(
             "LAPACK::dpotrf numerical failure. "
              "The leading minor of order %d  is not positive definite.", info);
-    return FAILURE;
+    return LINEAR_SOLVER_FAILURE;
   }
 
   dpotrs_(&uplo, &n, &nrhs, lhs, &n, rhs_and_solution, &n, &info);
@@ -104,11 +104,11 @@
                << "Please report it."
                << "LAPACK::dpotrs fatal error."
                << "Argument: " << -info << " is invalid.";
-    return FATAL_ERROR;
+    return LINEAR_SOLVER_FATAL_ERROR;
   }
 
   *status = "Success";
-  return TOLERANCE;
+  return LINEAR_SOLVER_SUCCESS;
 #endif
 };
 
@@ -154,7 +154,7 @@
     string* status) {
 #ifdef CERES_NO_LAPACK
   LOG(FATAL) << "Ceres was built without a LAPACK library.";
-  return FATAL_ERROR;
+  return LINEAR_SOLVER_FATAL_ERROR;
 #else
   char trans = 'N';
   int m = num_rows;
@@ -185,7 +185,7 @@
   }
 
   *status = "Success.";
-  return TOLERANCE;
+  return LINEAR_SOLVER_SUCCESS;
 #endif
 }
 
diff --git a/internal/ceres/levenberg_marquardt_strategy.cc b/internal/ceres/levenberg_marquardt_strategy.cc
index 2227da6..ce3b69a 100644
--- a/internal/ceres/levenberg_marquardt_strategy.cc
+++ b/internal/ceres/levenberg_marquardt_strategy.cc
@@ -106,12 +106,12 @@
   LinearSolver::Summary linear_solver_summary =
       linear_solver_->Solve(jacobian, residuals, solve_options, step);
 
-  if (linear_solver_summary.termination_type == FATAL_ERROR) {
+  if (linear_solver_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
     LOG(WARNING) << "Linear solver fatal error.";
-  } else if (linear_solver_summary.termination_type == FAILURE ||
+  } else if (linear_solver_summary.termination_type == LINEAR_SOLVER_FAILURE ||
              !IsArrayValid(num_parameters, step)) {
     LOG(WARNING) << "Linear solver failure. Failed to compute a finite step.";
-    linear_solver_summary.termination_type = FAILURE;
+    linear_solver_summary.termination_type = LINEAR_SOLVER_FAILURE;
   } else {
     VectorRef(step, num_parameters) *= -1.0;
   }
diff --git a/internal/ceres/levenberg_marquardt_strategy_test.cc b/internal/ceres/levenberg_marquardt_strategy_test.cc
index 86302b7..ac7ddbc 100644
--- a/internal/ceres/levenberg_marquardt_strategy_test.cc
+++ b/internal/ceres/levenberg_marquardt_strategy_test.cc
@@ -150,7 +150,7 @@
 
     TrustRegionStrategy::Summary summary =
         lms.ComputeStep(pso, &dsm, &residual, x);
-    EXPECT_EQ(summary.termination_type, FAILURE);
+    EXPECT_EQ(summary.termination_type, LINEAR_SOLVER_FAILURE);
   }
 }
 
diff --git a/internal/ceres/linear_solver.h b/internal/ceres/linear_solver.h
index fa151e0..cb26356 100644
--- a/internal/ceres/linear_solver.h
+++ b/internal/ceres/linear_solver.h
@@ -51,26 +51,19 @@
 namespace internal {
 
 enum LinearSolverTerminationType {
-  // Termination criterion was met. For factorization based solvers
-  // the tolerance is assumed to be zero. Any user provided values are
-  // ignored.
-  TOLERANCE,
+  // Termination criterion was met.
+  LINEAR_SOLVER_SUCCESS,
 
   // Solver ran for max_num_iterations and terminated before the
-  // termination tolerance could be satified.
-  MAX_ITERATIONS,
+  // termination tolerance could be satisfied.
+  LINEAR_SOLVER_NO_CONVERGENCE,
 
-  // Solver is stuck and further iterations will not result in any
-  // measurable progress.
-  STAGNATION,
-
-  // Solver failed. Solver was terminated due to numerical errors. The
-  // exact cause of failure depends on the particular solver being
-  // used.
-  FAILURE,
+  // Solver was terminated due to numerical problems, generally due to
+  // the linear system being poorly conditioned.
+  LINEAR_SOLVER_FAILURE,
 
   // Solver failed with a fatal error that cannot be recovered from.
-  FATAL_ERROR
+  LINEAR_SOLVER_FATAL_ERROR
 };
 
 
@@ -267,7 +260,7 @@
     Summary()
         : residual_norm(0.0),
           num_iterations(-1),
-          termination_type(FAILURE) {
+          termination_type(LINEAR_SOLVER_FAILURE) {
     }
 
     double residual_norm;
diff --git a/internal/ceres/schur_complement_solver.cc b/internal/ceres/schur_complement_solver.cc
index 9bf9d9d..2f407fd 100644
--- a/internal/ceres/schur_complement_solver.cc
+++ b/internal/ceres/schur_complement_solver.cc
@@ -83,7 +83,7 @@
       SolveReducedLinearSystem(reduced_solution);
   event_logger.AddEvent("ReducedSolve");
 
-  if (summary.termination_type == TOLERANCE) {
+  if (summary.termination_type == LINEAR_SOLVER_SUCCESS) {
     eliminator_->BackSubstitute(A, b, per_solve_options.D, reduced_solution, x);
     event_logger.AddEvent("BackSubstitute");
   }
@@ -116,7 +116,7 @@
 DenseSchurComplementSolver::SolveReducedLinearSystem(double* solution) {
   LinearSolver::Summary summary;
   summary.num_iterations = 0;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.message = "Success.";
 
   const BlockRandomAccessDenseMatrix* m =
@@ -137,7 +137,7 @@
         .selfadjointView<Eigen::Upper>()
         .llt();
     if (llt.info() != Eigen::Success) {
-      summary.termination_type = FAILURE;
+      summary.termination_type = LINEAR_SOLVER_FAILURE;
       summary.message = "Eigen LLT decomposition failed.";
       return summary;
     }
@@ -275,7 +275,7 @@
     double* solution) {
   LinearSolver::Summary summary;
   summary.num_iterations = 0;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.message = "Success.";
 
   TripletSparseMatrix* tsm =
@@ -325,13 +325,13 @@
 
   if (factor_ == NULL) {
     ss_.Free(cholmod_lhs);
-    summary.termination_type = FATAL_ERROR;
+    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
     return summary;
   }
 
   summary.termination_type =
       ss_.Cholesky(cholmod_lhs, factor_, &summary.message);
-  if (summary.termination_type != TOLERANCE) {
+  if (summary.termination_type != LINEAR_SOLVER_SUCCESS) {
     return summary;
   }
 
@@ -344,7 +344,7 @@
   ss_.Free(cholmod_rhs);
 
   if (cholmod_solution == NULL) {
-    summary.termination_type = FAILURE;
+    summary.termination_type = LINEAR_SOLVER_FAILURE;
     return summary;
   }
 
@@ -371,7 +371,7 @@
     double* solution) {
   LinearSolver::Summary summary;
   summary.num_iterations = 0;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.message = "Success.";
 
   // Extract the TripletSparseMatrix that is used for actually storing S.
@@ -395,11 +395,11 @@
   }
 
   if (cxsparse_factor_ == NULL) {
-    summary.termination_type = FATAL_ERROR;
+    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
     summary.message =
         "CXSparse failure. Unable to find symbolic factorization.";
   } else if (!cxsparse_.SolveCholesky(lhs, cxsparse_factor_, solution)) {
-    summary.termination_type = FAILURE;
+    summary.termination_type = LINEAR_SOLVER_FAILURE;
     summary.message = "CXSparse::SolveCholesky failed.";
   }
 
diff --git a/internal/ceres/sparse_normal_cholesky_solver.cc b/internal/ceres/sparse_normal_cholesky_solver.cc
index ceeb654..1ead8f7 100644
--- a/internal/ceres/sparse_normal_cholesky_solver.cc
+++ b/internal/ceres/sparse_normal_cholesky_solver.cc
@@ -102,7 +102,7 @@
 
   LinearSolver::Summary summary;
   summary.num_iterations = 1;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.message = "Success.";
 
   const int num_cols = A->num_cols();
@@ -150,13 +150,13 @@
   event_logger.AddEvent("Analysis");
 
   if (cxsparse_factor_ == NULL) {
-    summary.termination_type = FATAL_ERROR;
+    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
     summary.message =
         "CXSparse failure. Unable to find symbolic factorization.";
   } else if (cxsparse_.SolveCholesky(AtA, cxsparse_factor_, Atb.data())) {
     VectorRef(x, Atb.rows()) = Atb;
   } else {
-    summary.termination_type = FAILURE;
+    summary.termination_type = LINEAR_SOLVER_FAILURE;
   }
   event_logger.AddEvent("Solve");
 
@@ -185,7 +185,7 @@
     double * x) {
   EventLogger event_logger("SparseNormalCholeskySolver::SuiteSparse::Solve");
   LinearSolver::Summary summary;
-  summary.termination_type = TOLERANCE;
+  summary.termination_type = LINEAR_SOLVER_SUCCESS;
   summary.num_iterations = 1;
   summary.message = "Success.";
 
@@ -220,12 +220,12 @@
     if (per_solve_options.D != NULL) {
       A->DeleteRows(num_cols);
     }
-    summary.termination_type = FATAL_ERROR;
+    summary.termination_type = LINEAR_SOLVER_FATAL_ERROR;
     return summary;
   }
 
   summary.termination_type = ss_.Cholesky(&lhs, factor_, &summary.message);
-  if (summary.termination_type != TOLERANCE) {
+  if (summary.termination_type != LINEAR_SOLVER_SUCCESS) {
     if (per_solve_options.D != NULL) {
       A->DeleteRows(num_cols);
     }
@@ -245,7 +245,7 @@
     memcpy(x, sol->x, num_cols * sizeof(*x));
     ss_.Free(sol);
   } else {
-    summary.termination_type = FAILURE;
+    summary.termination_type = LINEAR_SOLVER_FAILURE;
   }
 
   event_logger.AddEvent("Teardown");
diff --git a/internal/ceres/suitesparse.cc b/internal/ceres/suitesparse.cc
index 69a07ce..06cc0a8 100644
--- a/internal/ceres/suitesparse.cc
+++ b/internal/ceres/suitesparse.cc
@@ -269,41 +269,41 @@
   switch (cc_.status) {
     case CHOLMOD_NOT_INSTALLED:
       *status = "CHOLMOD failure: Method not installed.";
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
     case CHOLMOD_OUT_OF_MEMORY:
       *status = "CHOLMOD failure: Out of memory.";
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
     case CHOLMOD_TOO_LARGE:
       *status = "CHOLMOD failure: Integer overflow occured.";
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
     case CHOLMOD_INVALID:
       *status = "CHOLMOD failure: Invalid input.";
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
     case CHOLMOD_NOT_POSDEF:
       *status = "CHOLMOD warning: Matrix not positive definite.";
-      return FAILURE;
+      return LINEAR_SOLVER_FAILURE;
     case CHOLMOD_DSMALL:
       *status = "CHOLMOD warning: D for LDL' or diag(L) or "
                 "LL' has tiny absolute value.";
-      return FAILURE;
+      return LINEAR_SOLVER_FAILURE;
     case CHOLMOD_OK:
       if (cholmod_status != 0) {
-        return TOLERANCE;
+        return LINEAR_SOLVER_SUCCESS;
       }
 
       *status = "CHOLMOD failure: cholmod_factorize returned false "
                 "but cholmod_common::status is CHOLMOD_OK."
                 "Please report this to ceres-solver@googlegroups.com.";
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
     default:
       *status =
           StringPrintf("Unknown cholmod return code: %d. "
                        "Please report this to ceres-solver@googlegroups.com.",
                        cc_.status);
-      return FATAL_ERROR;
+      return LINEAR_SOLVER_FATAL_ERROR;
   }
 
-  return FATAL_ERROR;
+  return LINEAR_SOLVER_FATAL_ERROR;
 }
 
 cholmod_dense* SuiteSparse::Solve(cholmod_factor* L,
diff --git a/internal/ceres/symmetric_linear_solver_test.cc b/internal/ceres/symmetric_linear_solver_test.cc
index f33adb4..ac5a774 100644
--- a/internal/ceres/symmetric_linear_solver_test.cc
+++ b/internal/ceres/symmetric_linear_solver_test.cc
@@ -71,7 +71,7 @@
   LinearSolver::Summary summary =
       solver.Solve(A.get(), b.data(), per_solve_options, x.data());
 
-  EXPECT_EQ(summary.termination_type, TOLERANCE);
+  EXPECT_EQ(summary.termination_type, LINEAR_SOLVER_SUCCESS);
   ASSERT_EQ(summary.num_iterations, 1);
 
   ASSERT_DOUBLE_EQ(1, x(0));
@@ -128,7 +128,7 @@
   LinearSolver::Summary summary =
       solver.Solve(A.get(), b.data(), per_solve_options, x.data());
 
-  EXPECT_EQ(summary.termination_type, TOLERANCE);
+  EXPECT_EQ(summary.termination_type, LINEAR_SOLVER_SUCCESS);
 
   ASSERT_DOUBLE_EQ(0, x(0));
   ASSERT_DOUBLE_EQ(1, x(1));
diff --git a/internal/ceres/trust_region_minimizer.cc b/internal/ceres/trust_region_minimizer.cc
index 9b630e1..6094408 100644
--- a/internal/ceres/trust_region_minimizer.cc
+++ b/internal/ceres/trust_region_minimizer.cc
@@ -227,6 +227,15 @@
                               residuals.data(),
                               trust_region_step.data());
 
+    if (strategy_summary.termination_type == LINEAR_SOLVER_FATAL_ERROR) {
+      summary->error =
+          "Terminating. Linear solver failed due to unrecoverable "
+          "non-numeric causes. Please see the error log for clues. ";
+      summary->termination_type = NUMERICAL_FAILURE;
+      LOG_IF(WARNING, is_not_silent) << summary->error;
+      return;
+    }
+
     iteration_summary = IterationSummary();
     iteration_summary.iteration = summary->iterations.back().iteration + 1;
     iteration_summary.step_solver_time_in_seconds =
@@ -246,7 +255,7 @@
     }
 
     double model_cost_change = 0.0;
-    if (strategy_summary.termination_type != FAILURE) {
+    if (strategy_summary.termination_type != LINEAR_SOLVER_FAILURE) {
       // new_model_cost
       //  = 1/2 [f + J * step]^2
       //  = 1/2 [ f'f + 2f'J * step + step' * J' * J * step ]
diff --git a/internal/ceres/trust_region_strategy.h b/internal/ceres/trust_region_strategy.h
index 3f078d5..998514f 100644
--- a/internal/ceres/trust_region_strategy.h
+++ b/internal/ceres/trust_region_strategy.h
@@ -106,7 +106,7 @@
     Summary()
         : residual_norm(0.0),
           num_iterations(-1),
-          termination_type(FAILURE) {
+          termination_type(LINEAR_SOLVER_FAILURE) {
     }
 
     // If the trust region problem is,
diff --git a/internal/ceres/unsymmetric_linear_solver_test.cc b/internal/ceres/unsymmetric_linear_solver_test.cc
index af9dffe..9784b46 100644
--- a/internal/ceres/unsymmetric_linear_solver_test.cc
+++ b/internal/ceres/unsymmetric_linear_solver_test.cc
@@ -99,13 +99,13 @@
                       per_solve_options,
                       x_regularized.data());
 
-    EXPECT_EQ(unregularized_solve_summary.termination_type, TOLERANCE);
+    EXPECT_EQ(unregularized_solve_summary.termination_type, LINEAR_SOLVER_SUCCESS);
 
     for (int i = 0; i < A_->num_cols(); ++i) {
       EXPECT_NEAR(sol_unregularized_[i], x_unregularized[i], 1e-8);
     }
 
-    EXPECT_EQ(regularized_solve_summary.termination_type, TOLERANCE);
+    EXPECT_EQ(regularized_solve_summary.termination_type, LINEAR_SOLVER_SUCCESS);
     for (int i = 0; i < A_->num_cols(); ++i) {
       EXPECT_NEAR(sol_regularized_[i], x_regularized[i], 1e-8);
     }
diff --git a/internal/ceres/visibility_based_preconditioner.cc b/internal/ceres/visibility_based_preconditioner.cc
index 13fe27b..a3bebed 100644
--- a/internal/ceres/visibility_based_preconditioner.cc
+++ b/internal/ceres/visibility_based_preconditioner.cc
@@ -370,7 +370,7 @@
   // scaling is not needed, which is quite often in our experience.
   LinearSolverTerminationType status = Factorize();
 
-  if (status == FATAL_ERROR) {
+  if (status == LINEAR_SOLVER_FATAL_ERROR) {
     return false;
   }
 
@@ -379,7 +379,7 @@
   // belong to the edges of the degree-2 forest. In the CLUSTER_JACOBI
   // case, the preconditioner is guaranteed to be positive
   // semidefinite.
-  if (status == FAILURE && options_.type == CLUSTER_TRIDIAGONAL) {
+  if (status == LINEAR_SOLVER_FAILURE && options_.type == CLUSTER_TRIDIAGONAL) {
     VLOG(1) << "Unscaled factorization failed. Retrying with off-diagonal "
             << "scaling";
     ScaleOffDiagonalCells();
@@ -387,7 +387,7 @@
   }
 
   VLOG(2) << "Compute time: " << time(NULL) - start_time;
-  return (status == TOLERANCE);
+  return (status == LINEAR_SOLVER_SUCCESS);
 }
 
 // Consider the preconditioner matrix as meta-block matrix, whose
@@ -444,7 +444,9 @@
   }
 
   const LinearSolverTerminationType termination_type =
-      (factor_ != NULL) ? ss_.Cholesky(lhs, factor_, &status) : FATAL_ERROR;
+      (factor_ != NULL)
+      ? ss_.Cholesky(lhs, factor_, &status)
+      : LINEAR_SOLVER_FATAL_ERROR;
 
   ss_.Free(lhs);
   return termination_type;