Detect large Jacobians and return failure instead of crashing.

Detect when the number of non-zeros overflows when constructing
BlockSparseMatrix and CompressedRowSparseMatrix and return
with an error message instead of crashing.

Change-Id: I45e102f7c0519eef441ce0586b7adf96e4a954a9
diff --git a/internal/ceres/block_jacobian_writer.cc b/internal/ceres/block_jacobian_writer.cc
index b009e96..adeaf41 100644
--- a/internal/ceres/block_jacobian_writer.cc
+++ b/internal/ceres/block_jacobian_writer.cc
@@ -63,7 +63,7 @@
 //
 // TODO(keir): Consider if we should use a boolean for each parameter block
 // instead of num_eliminate_blocks.
-void BuildJacobianLayout(const Program& program,
+bool BuildJacobianLayout(const Program& program,
                          int num_eliminate_blocks,
                          std::vector<int*>* jacobian_layout,
                          std::vector<int>* jacobian_layout_storage) {
@@ -73,8 +73,8 @@
   // Iterate over all the active residual blocks and determine how many E blocks
   // are there. This will determine where the F blocks start in the jacobian
   // matrix. Also compute the number of jacobian blocks.
-  int f_block_pos = 0;
-  int num_jacobian_blocks = 0;
+  unsigned int f_block_pos = 0;
+  unsigned int num_jacobian_blocks = 0;
   for (auto* residual_block : residual_blocks) {
     const int num_residuals = residual_block->NumResiduals();
     const int num_parameter_blocks = residual_block->NumParameterBlocks();
@@ -90,6 +90,11 @@
         }
       }
     }
+    if (num_jacobian_blocks > std::numeric_limits<int>::max()) {
+      LOG(ERROR) << "Overlow error. Too many blocks in the jacobian matrix : "
+                 << num_jacobian_blocks;
+      return false;
+    }
   }
 
   // We now know that the E blocks are laid out starting at zero, and the F
@@ -145,12 +150,18 @@
         jacobian_pos[k] = e_block_pos;
         e_block_pos += jacobian_block_size;
       } else {
-        jacobian_pos[k] = f_block_pos;
+        jacobian_pos[k] = static_cast<int>(f_block_pos);
         f_block_pos += jacobian_block_size;
+        if (f_block_pos > std::numeric_limits<int>::max()) {
+          LOG(ERROR)
+              << "Overlow error. Too many entries in the Jacobian matrix.";
+          return false;
+        }
       }
     }
     jacobian_pos += active_parameter_blocks.size();
   }
+  return true;
 }
 
 }  // namespace
@@ -161,10 +172,10 @@
   CHECK_GE(options.num_eliminate_blocks, 0)
       << "num_eliminate_blocks must be greater than 0.";
 
-  BuildJacobianLayout(*program,
-                      options.num_eliminate_blocks,
-                      &jacobian_layout_,
-                      &jacobian_layout_storage_);
+  jacobian_layout_is_valid_ = BuildJacobianLayout(*program,
+                                                  options.num_eliminate_blocks,
+                                                  &jacobian_layout_,
+                                                  &jacobian_layout_storage_);
 }
 
 // Create evaluate prepareres that point directly into the final jacobian. This
@@ -183,6 +194,12 @@
 }
 
 std::unique_ptr<SparseMatrix> BlockJacobianWriter::CreateJacobian() const {
+  if (!jacobian_layout_is_valid_) {
+    LOG(ERROR) << "Unable to create Jacobian matrix. Too many entries in the "
+                  "Jacobian matrix.";
+    return nullptr;
+  }
+
   auto* bs = new CompressedRowBlockStructure;
 
   const std::vector<ParameterBlock*>& parameter_blocks =
diff --git a/internal/ceres/block_jacobian_writer.h b/internal/ceres/block_jacobian_writer.h
index 30b5f02..627013c 100644
--- a/internal/ceres/block_jacobian_writer.h
+++ b/internal/ceres/block_jacobian_writer.h
@@ -133,6 +133,12 @@
 
   // The pointers in jacobian_layout_ point directly into this vector.
   std::vector<int> jacobian_layout_storage_;
+
+  // The constructor computes the layout of the Jacobian, and this bool keeps
+  // track of whether the computation of the layout completed successfully or
+  // not, if it is false, then jacobian_layout and jacobian_layout_storage are
+  // both in an invalid state.
+  bool jacobian_layout_is_valid_ = false;
 };
 
 }  // namespace ceres::internal
diff --git a/internal/ceres/compressed_row_jacobian_writer.cc b/internal/ceres/compressed_row_jacobian_writer.cc
index 14c5251..5533fcf 100644
--- a/internal/ceres/compressed_row_jacobian_writer.cc
+++ b/internal/ceres/compressed_row_jacobian_writer.cc
@@ -94,7 +94,10 @@
   const int total_num_effective_parameters = program_->NumEffectiveParameters();
 
   // Count the number of jacobian nonzeros.
-  int num_jacobian_nonzeros = 0;
+  //
+  // We used an unsigned int here, so that we can compare it INT_MAX without
+  // triggering overflow behaviour.
+  unsigned int num_jacobian_nonzeros = total_num_effective_parameters;
   for (auto* residual_block : residual_blocks) {
     const int num_residuals = residual_block->NumResiduals();
     const int num_parameter_blocks = residual_block->NumParameterBlocks();
@@ -102,6 +105,12 @@
       auto parameter_block = residual_block->parameter_blocks()[j];
       if (!parameter_block->IsConstant()) {
         num_jacobian_nonzeros += num_residuals * parameter_block->TangentSize();
+        if (num_jacobian_nonzeros > std::numeric_limits<int>::max()) {
+          LOG(ERROR) << "Unable to create Jacobian matrix: Too many entries in "
+                        "the Jacobian matrix. num_jacobian_nonzeros = "
+                     << num_jacobian_nonzeros;
+          return nullptr;
+        }
       }
     }
   }
@@ -113,7 +122,7 @@
   auto jacobian = std::make_unique<CompressedRowSparseMatrix>(
       total_num_residuals,
       total_num_effective_parameters,
-      num_jacobian_nonzeros + total_num_effective_parameters);
+      static_cast<int>(num_jacobian_nonzeros));
 
   // At this stage, the CompressedRowSparseMatrix is an invalid state. But
   // this seems to be the only way to construct it without doing a memory
diff --git a/internal/ceres/evaluator_test.cc b/internal/ceres/evaluator_test.cc
index 34ec78c..773cbcd 100644
--- a/internal/ceres/evaluator_test.cc
+++ b/internal/ceres/evaluator_test.cc
@@ -676,5 +676,51 @@
   }
 }
 
+class HugeCostFunction : public SizedCostFunction<46341, 46345> {
+  bool Evaluate(double const* const* parameters,
+                double* residuals,
+                double** jacobians) const override {
+    return true;
+  }
+};
+
+TEST(Evaluator, LargeProblemDoesNotCauseCrashBlockJacobianWriter) {
+  ProblemImpl problem;
+  std::vector<double> x(46345);
+
+  problem.AddResidualBlock(new HugeCostFunction, nullptr, x.data());
+  Evaluator::Options options;
+  options.linear_solver_type = SPARSE_NORMAL_CHOLESKY;
+  options.context = problem.context();
+  options.num_eliminate_blocks = 0;
+  options.dynamic_sparsity = false;
+  std::string error;
+  auto program = problem.mutable_program();
+  program->SetParameterOffsetsAndIndex();
+  auto evaluator = Evaluator::Create(options, program, &error);
+  auto jacobian = evaluator->CreateJacobian();
+  EXPECT_EQ(jacobian, nullptr);
+}
+
+TEST(Evaluator, LargeProblemDoesNotCauseCrashCompressedRowJacobianWriter) {
+  ProblemImpl problem;
+  std::vector<double> x(46345);
+
+  problem.AddResidualBlock(new HugeCostFunction, nullptr, x.data());
+  Evaluator::Options options;
+  // CGNR on CUDA_SPARSE is the only combination that triggers a
+  // CompressedRowJacobianWriter.
+  options.linear_solver_type = CGNR;
+  options.sparse_linear_algebra_library_type = CUDA_SPARSE;
+  options.context = problem.context();
+  options.num_eliminate_blocks = 0;
+  std::string error;
+  auto program = problem.mutable_program();
+  program->SetParameterOffsetsAndIndex();
+  auto evaluator = Evaluator::Create(options, program, &error);
+  auto jacobian = evaluator->CreateJacobian();
+  EXPECT_EQ(jacobian, nullptr);
+}
+
 }  // namespace internal
 }  // namespace ceres
diff --git a/internal/ceres/trust_region_preprocessor.cc b/internal/ceres/trust_region_preprocessor.cc
index 519f851..2b036a5 100644
--- a/internal/ceres/trust_region_preprocessor.cc
+++ b/internal/ceres/trust_region_preprocessor.cc
@@ -336,13 +336,19 @@
 }
 
 // Configure and create a TrustRegionMinimizer object.
-void SetupMinimizerOptions(PreprocessedProblem* pp) {
+bool SetupMinimizerOptions(PreprocessedProblem* pp) {
   const Solver::Options& options = pp->options;
 
   SetupCommonMinimizerOptions(pp);
   pp->minimizer_options.is_constrained =
       pp->reduced_program->IsBoundsConstrained();
   pp->minimizer_options.jacobian = pp->evaluator->CreateJacobian();
+  if (pp->minimizer_options.jacobian == nullptr) {
+    pp->error =
+        "Unable to create Jacobian matrix. Likely because it is too large.";
+    return false;
+  }
+
   pp->minimizer_options.inner_iteration_minimizer =
       pp->inner_iteration_minimizer;
 
@@ -360,6 +366,7 @@
   pp->minimizer_options.trust_region_strategy =
       TrustRegionStrategy::Create(strategy_options);
   CHECK(pp->minimizer_options.trust_region_strategy != nullptr);
+  return true;
 }
 
 }  // namespace
@@ -395,8 +402,7 @@
     return false;
   }
 
-  SetupMinimizerOptions(pp);
-  return true;
+  return SetupMinimizerOptions(pp);
 }
 
 }  // namespace ceres::internal