Use C++11's inline member initialization syntax

Migrate all Option and Summary structs to use
inline member initialization syntax.

This reduces the amount of code, and collocates the
default values with the documentation for the corresponding
member variable.

Change-Id: I8e6b9ee3b31464699d678667f6166ace5fc137c9
diff --git a/include/ceres/covariance.h b/include/ceres/covariance.h
index a8061b5..0b9f096 100644
--- a/include/ceres/covariance.h
+++ b/include/ceres/covariance.h
@@ -200,26 +200,17 @@
 class CERES_EXPORT Covariance {
  public:
   struct CERES_EXPORT Options {
-    Options() {
-      algorithm_type = SPARSE_QR;
-
-      // Eigen's QR factorization is always available.
-      sparse_linear_algebra_library_type = EIGEN_SPARSE;
-#if !defined(CERES_NO_SUITESPARSE)
-      sparse_linear_algebra_library_type = SUITE_SPARSE;
-#endif
-
-      min_reciprocal_condition_number = 1e-14;
-      null_space_rank = 0;
-      num_threads = 1;
-      apply_loss_function = true;
-    }
-
     // Sparse linear algebra library to use when a sparse matrix
     // factorization is being used to compute the covariance matrix.
     //
     // Currently this only applies to SPARSE_QR.
-    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
+    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
+#if !defined(CERES_NO_SUITESPARSE)
+        SUITE_SPARSE;
+#else
+        // Eigen's QR factorization is always available.
+        EIGEN_SPARSE;
+#endif
 
     // Ceres supports two different algorithms for covariance
     // estimation, which represent different tradeoffs in speed,
@@ -251,7 +242,7 @@
     // Eigen's Sparse QR factorization algorithm will be used or
     // SuiteSparse's high performance SuiteSparseQR algorithm will be
     // used.
-    CovarianceAlgorithmType algorithm_type;
+    CovarianceAlgorithmType algorithm_type = SPARSE_QR;
 
     // If the Jacobian matrix is near singular, then inverting J'J
     // will result in unreliable results, e.g, if
@@ -284,7 +275,7 @@
     //   sparse QR factorization algorithm. It is a fairly reliable
     //   indication of rank deficiency.
     //
-    double min_reciprocal_condition_number;
+    double min_reciprocal_condition_number = 1e-14;
 
     // When using DENSE_SVD, the user has more control in dealing with
     // singular and near singular covariance matrices.
@@ -319,9 +310,9 @@
     //
     // This option has no effect on the SUITE_SPARSE_QR and
     // EIGEN_SPARSE_QR algorithms.
-    int null_space_rank;
+    int null_space_rank = 0;
 
-    int num_threads;
+    int num_threads = 1;
 
     // Even though the residual blocks in the problem may contain loss
     // functions, setting apply_loss_function to false will turn off
@@ -329,7 +320,7 @@
     // function and in turn its effect on the covariance.
     //
     // TODO(sameergaarwal): Expand this based on Jim's experiments.
-    bool apply_loss_function;
+    bool apply_loss_function = true;
   };
 
   explicit Covariance(const Options& options);
diff --git a/include/ceres/gradient_problem_solver.h b/include/ceres/gradient_problem_solver.h
index f9de6fb..5bbc8ce 100644
--- a/include/ceres/gradient_problem_solver.h
+++ b/include/ceres/gradient_problem_solver.h
@@ -54,41 +54,16 @@
   //
   // The constants are defined inside types.h
   struct CERES_EXPORT Options {
-    // Default constructor that sets up a generic sparse problem.
-    Options() {
-      line_search_direction_type = LBFGS;
-      line_search_type = WOLFE;
-      nonlinear_conjugate_gradient_type = FLETCHER_REEVES;
-      max_lbfgs_rank = 20;
-      use_approximate_eigenvalue_bfgs_scaling = false;
-      line_search_interpolation_type = CUBIC;
-      min_line_search_step_size = 1e-9;
-      line_search_sufficient_function_decrease = 1e-4;
-      max_line_search_step_contraction = 1e-3;
-      min_line_search_step_contraction = 0.6;
-      max_num_line_search_step_size_iterations = 20;
-      max_num_line_search_direction_restarts = 5;
-      line_search_sufficient_curvature_decrease = 0.9;
-      max_line_search_step_expansion = 10.0;
-      max_num_iterations = 50;
-      max_solver_time_in_seconds = 1e9;
-      function_tolerance = 1e-6;
-      gradient_tolerance = 1e-10;
-      parameter_tolerance = 1e-8;
-      logging_type = PER_MINIMIZER_ITERATION;
-      minimizer_progress_to_stdout = false;
-      update_state_every_iteration = false;
-    }
-
     // Returns true if the options struct has a valid
     // configuration. Returns false otherwise, and fills in *error
     // with a message describing the problem.
     bool IsValid(std::string* error) const;
 
     // Minimizer options ----------------------------------------
-    LineSearchDirectionType line_search_direction_type;
-    LineSearchType line_search_type;
-    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
+    LineSearchDirectionType line_search_direction_type = LBFGS;
+    LineSearchType line_search_type = WOLFE;
+    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type =
+        FLETCHER_REEVES;
 
     // The LBFGS hessian approximation is a low rank approximation to
     // the inverse of the Hessian matrix. The rank of the
@@ -114,7 +89,7 @@
     //
     // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with
     // Limited Storage". Mathematics of Computation 35 (151): 773–782.
-    int max_lbfgs_rank;
+    int max_lbfgs_rank = 20;
 
     // As part of the (L)BFGS update step (BFGS) / right-multiply step (L-BFGS),
     // the initial inverse Hessian approximation is taken to be the Identity.
@@ -136,18 +111,18 @@
     // Oren S.S., Self-scaling variable metric (SSVM) algorithms
     // Part II: Implementation and experiments, Management Science,
     // 20(5), 863-874, 1974.
-    bool use_approximate_eigenvalue_bfgs_scaling;
+    bool use_approximate_eigenvalue_bfgs_scaling = false;
 
     // Degree of the polynomial used to approximate the objective
     // function. Valid values are BISECTION, QUADRATIC and CUBIC.
     //
     // BISECTION corresponds to pure backtracking search with no
     // interpolation.
-    LineSearchInterpolationType line_search_interpolation_type;
+    LineSearchInterpolationType line_search_interpolation_type = CUBIC;
 
     // If during the line search, the step_size falls below this
     // value, it is truncated to zero.
-    double min_line_search_step_size;
+    double min_line_search_step_size = 1e-9;
 
     // Line search parameters.
 
@@ -161,7 +136,7 @@
     //
     //   f(step_size) <= f(0) + sufficient_decrease * f'(0) * step_size
     //
-    double line_search_sufficient_function_decrease;
+    double line_search_sufficient_function_decrease = 1e-4;
 
     // In each iteration of the line search,
     //
@@ -171,7 +146,7 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double max_line_search_step_contraction;
+    double max_line_search_step_contraction = 1e-3;
 
     // In each iteration of the line search,
     //
@@ -181,19 +156,19 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double min_line_search_step_contraction;
+    double min_line_search_step_contraction = 0.6;
 
     // Maximum number of trial step size iterations during each line search,
     // if a step size satisfying the search conditions cannot be found within
     // this number of trials, the line search will terminate.
-    int max_num_line_search_step_size_iterations;
+    int max_num_line_search_step_size_iterations = 20;
 
     // Maximum number of restarts of the line search direction algorithm before
     // terminating the optimization. Restarts of the line search direction
     // algorithm occur when the current algorithm fails to produce a new descent
     // direction. This typically indicates a numerical failure, or a breakdown
     // in the validity of the approximations used.
-    int max_num_line_search_direction_restarts;
+    int max_num_line_search_direction_restarts = 5;
 
     // The strong Wolfe conditions consist of the Armijo sufficient
     // decrease condition, and an additional requirement that the
@@ -206,7 +181,7 @@
     //
     // Where f() is the line search objective and f'() is the derivative
     // of f w.r.t step_size (d f / d step_size).
-    double line_search_sufficient_curvature_decrease;
+    double line_search_sufficient_curvature_decrease = 0.9;
 
     // During the bracketing phase of the Wolfe search, the step size is
     // increased until either a point satisfying the Wolfe conditions is
@@ -217,49 +192,49 @@
     //   new_step_size <= max_step_expansion * step_size.
     //
     // By definition for expansion, max_step_expansion > 1.0.
-    double max_line_search_step_expansion;
+    double max_line_search_step_expansion = 10.0;
 
     // Maximum number of iterations for the minimizer to run for.
-    int max_num_iterations;
+    int max_num_iterations = 50;
 
     // Maximum time for which the minimizer should run for.
-    double max_solver_time_in_seconds;
+    double max_solver_time_in_seconds = 1e9;
 
     // Minimizer terminates when
     //
     //   (new_cost - old_cost) < function_tolerance * old_cost;
     //
-    double function_tolerance;
+    double function_tolerance = 1e-6;
 
     // Minimizer terminates when
     //
     //   max_i |x - Project(Plus(x, -g(x))| < gradient_tolerance
     //
     // This value should typically be 1e-4 * function_tolerance.
-    double gradient_tolerance;
+    double gradient_tolerance = 1e-10;
 
     // Minimizer terminates when
     //
     //   |step|_2 <= parameter_tolerance * ( |x|_2 +  parameter_tolerance)
     //
-    double parameter_tolerance;
+    double parameter_tolerance = 1e-8;
 
     // Logging options ---------------------------------------------------------
 
-    LoggingType logging_type;
+    LoggingType logging_type = PER_MINIMIZER_ITERATION;
 
     // By default the Minimizer progress is logged to VLOG(1), which
     // is sent to STDERR depending on the vlog level. If this flag is
     // set to true, and logging_type is not SILENT, the logging output
     // is sent to STDOUT.
-    bool minimizer_progress_to_stdout;
+    bool minimizer_progress_to_stdout = false;
 
     // If true, the user's parameter blocks are updated at the end of
     // every Minimizer iteration, otherwise they are updated when the
     // Minimizer terminates. This is useful if, for example, the user
     // wishes to visualize the state of the optimization every
     // iteration.
-    bool update_state_every_iteration;
+    bool update_state_every_iteration = false;
 
     // Callbacks that are executed at the end of each iteration of the
     // Minimizer. An iteration may terminate midway, either due to
@@ -280,8 +255,6 @@
   };
 
   struct CERES_EXPORT Summary {
-    Summary();
-
     // A brief one line description of the state of the solver after
     // termination.
     std::string BriefReport() const;
@@ -293,65 +266,65 @@
     bool IsSolutionUsable() const;
 
     // Minimizer summary -------------------------------------------------
-    TerminationType termination_type;
+    TerminationType termination_type = FAILURE;
 
     // Reason why the solver terminated.
-    std::string message;
+    std::string message = "ceres::GradientProblemSolve was not called.";
 
     // Cost of the problem (value of the objective function) before
     // the optimization.
-    double initial_cost;
+    double initial_cost = -1.0;
 
     // Cost of the problem (value of the objective function) after the
     // optimization.
-    double final_cost;
+    double final_cost = -1.0;
 
     // IterationSummary for each minimizer iteration in order.
     std::vector<IterationSummary> iterations;
 
     // Number of times the cost (and not the gradient) was evaluated.
-    int num_cost_evaluations;
+    int num_cost_evaluations = -1;
 
     // Number of times the gradient (and the cost) were evaluated.
-    int num_gradient_evaluations;
+    int num_gradient_evaluations = -1;
 
     // Sum total of all time spent inside Ceres when Solve is called.
-    double total_time_in_seconds;
+    double total_time_in_seconds = -1.0;
 
     // Time (in seconds) spent evaluating the cost.
-    double cost_evaluation_time_in_seconds;
+    double cost_evaluation_time_in_seconds = -1.0;
 
     // Time (in seconds) spent evaluating the gradient.
-    double gradient_evaluation_time_in_seconds;
+    double gradient_evaluation_time_in_seconds = -1.0;
 
     // Time (in seconds) spent minimizing the interpolating polynomial
     // to compute the next candidate step size as part of a line search.
-    double line_search_polynomial_minimization_time_in_seconds;
+    double line_search_polynomial_minimization_time_in_seconds = -1.0;
 
     // Number of parameters in the probem.
-    int num_parameters;
+    int num_parameters = -1;
 
     // Dimension of the tangent space of the problem.
-    int num_local_parameters;
+    int num_local_parameters = -1;
 
     // Type of line search direction used.
-    LineSearchDirectionType line_search_direction_type;
+    LineSearchDirectionType line_search_direction_type = LBFGS;
 
     // Type of the line search algorithm used.
-    LineSearchType line_search_type;
+    LineSearchType line_search_type = WOLFE;
 
     //  When performing line search, the degree of the polynomial used
     //  to approximate the objective function.
-    LineSearchInterpolationType line_search_interpolation_type;
+    LineSearchInterpolationType line_search_interpolation_type = CUBIC;
 
     // If the line search direction is NONLINEAR_CONJUGATE_GRADIENT,
     // then this indicates the particular variant of non-linear
     // conjugate gradient used.
-    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
+    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type = FLETCHER_REEVES;
 
     // If the type of the line search direction is LBFGS, then this
     // indicates the rank of the Hessian approximation.
-    int max_lbfgs_rank;
+    int max_lbfgs_rank = -1;
   };
 
   // Once a least squares problem has been built, this function takes
diff --git a/include/ceres/numeric_diff_options.h b/include/ceres/numeric_diff_options.h
index 119c8a8..3356be9 100644
--- a/include/ceres/numeric_diff_options.h
+++ b/include/ceres/numeric_diff_options.h
@@ -37,18 +37,10 @@
 // Options pertaining to numeric differentiation (e.g., convergence criteria,
 // step sizes).
 struct CERES_EXPORT NumericDiffOptions {
-  NumericDiffOptions() {
-    relative_step_size = 1e-6;
-    ridders_relative_initial_step_size = 1e-2;
-    max_num_ridders_extrapolations = 10;
-    ridders_epsilon = 1e-12;
-    ridders_step_shrink_factor = 2.0;
-  }
-
   // Numeric differentiation step size (multiplied by parameter block's
   // order of magnitude). If parameters are close to zero, the step size
   // is set to sqrt(machine_epsilon).
-  double relative_step_size;
+  double relative_step_size = 1e-6;
 
   // Initial step size for Ridders adaptive numeric differentiation (multiplied
   // by parameter block's order of magnitude).
@@ -59,19 +51,19 @@
   // Note: For Ridders' method to converge, the step size should be initialized
   // to a value that is large enough to produce a significant change in the
   // function. As the derivative is estimated, the step size decreases.
-  double ridders_relative_initial_step_size;
+  double ridders_relative_initial_step_size = 1e-2;
 
   // Maximal number of adaptive extrapolations (sampling) in Ridders' method.
-  int max_num_ridders_extrapolations;
+  int max_num_ridders_extrapolations = 10;
 
   // Convergence criterion on extrapolation error for Ridders adaptive
   // differentiation. The available error estimation methods are defined in
   // NumericDiffErrorType and set in the "ridders_error_method" field.
-  double ridders_epsilon;
+  double ridders_epsilon = 1e-12;
 
   // The factor in which to shrink the step size with each extrapolation in
   // Ridders' method.
-  double ridders_step_shrink_factor;
+  double ridders_step_shrink_factor = 2.0;
 };
 
 }  // namespace ceres
diff --git a/include/ceres/problem.h b/include/ceres/problem.h
index 344685d..b63395a 100644
--- a/include/ceres/problem.h
+++ b/include/ceres/problem.h
@@ -121,14 +121,6 @@
 class CERES_EXPORT Problem {
  public:
   struct CERES_EXPORT Options {
-    Options()
-        : cost_function_ownership(TAKE_OWNERSHIP),
-          loss_function_ownership(TAKE_OWNERSHIP),
-          local_parameterization_ownership(TAKE_OWNERSHIP),
-          enable_fast_removal(false),
-          disable_all_safety_checks(false),
-          context(NULL) {}
-
     // These flags control whether the Problem object owns the cost
     // functions, loss functions, and parameterizations passed into
     // the Problem. If set to TAKE_OWNERSHIP, then the problem object
@@ -136,9 +128,9 @@
     // destruction. The destructor is careful to delete the pointers
     // only once, since sharing cost/loss/parameterizations is
     // allowed.
-    Ownership cost_function_ownership;
-    Ownership loss_function_ownership;
-    Ownership local_parameterization_ownership;
+    Ownership cost_function_ownership = TAKE_OWNERSHIP;
+    Ownership loss_function_ownership = TAKE_OWNERSHIP;
+    Ownership local_parameterization_ownership = TAKE_OWNERSHIP;
 
     // If true, trades memory for faster RemoveResidualBlock() and
     // RemoveParameterBlock() operations.
@@ -154,7 +146,7 @@
     // The increase in memory usage is twofold: an additonal hash set per
     // parameter block containing all the residuals that depend on the parameter
     // block; and a hash set in the problem containing all residuals.
-    bool enable_fast_removal;
+    bool enable_fast_removal = false;
 
     // By default, Ceres performs a variety of safety checks when constructing
     // the problem. There is a small but measurable performance penalty to
@@ -165,14 +157,14 @@
     //
     // WARNING: Do not set this to true, unless you are absolutely sure of what
     // you are doing.
-    bool disable_all_safety_checks;
+    bool disable_all_safety_checks = false;
 
     // A Ceres global context to use for solving this problem. This may help to
     // reduce computation time as Ceres can reuse expensive objects to create.
     // The context object can be NULL, in which case Ceres may create one.
     //
     // Ceres does NOT take ownership of the pointer.
-    Context* context;
+    Context* context = nullptr;
   };
 
   // The default constructor is equivalent to the
@@ -401,11 +393,6 @@
 
   // Options struct to control Problem::Evaluate.
   struct EvaluateOptions {
-    EvaluateOptions()
-        : apply_loss_function(true),
-          num_threads(1) {
-    }
-
     // The set of parameter blocks for which evaluation should be
     // performed. This vector determines the order that parameter
     // blocks occur in the gradient vector and in the columns of the
@@ -438,9 +425,9 @@
     // function. This is of use for example if the user wishes to
     // analyse the solution quality by studying the distribution of
     // residuals before and after the solve.
-    bool apply_loss_function;
+    bool apply_loss_function = true;
 
-    int num_threads;
+    int num_threads = 1;
   };
 
   // Evaluate Problem. Any of the output pointers can be NULL. Which
diff --git a/include/ceres/solver.h b/include/ceres/solver.h
index 09157cc..592d329 100644
--- a/include/ceres/solver.h
+++ b/include/ceres/solver.h
@@ -59,88 +59,6 @@
   //
   // The constants are defined inside types.h
   struct CERES_EXPORT Options {
-    // Default constructor that sets up a generic sparse problem.
-    Options() {
-      minimizer_type = TRUST_REGION;
-      line_search_direction_type = LBFGS;
-      line_search_type = WOLFE;
-      nonlinear_conjugate_gradient_type = FLETCHER_REEVES;
-      max_lbfgs_rank = 20;
-      use_approximate_eigenvalue_bfgs_scaling = false;
-      line_search_interpolation_type = CUBIC;
-      min_line_search_step_size = 1e-9;
-      line_search_sufficient_function_decrease = 1e-4;
-      max_line_search_step_contraction = 1e-3;
-      min_line_search_step_contraction = 0.6;
-      max_num_line_search_step_size_iterations = 20;
-      max_num_line_search_direction_restarts = 5;
-      line_search_sufficient_curvature_decrease = 0.9;
-      max_line_search_step_expansion = 10.0;
-      trust_region_strategy_type = LEVENBERG_MARQUARDT;
-      dogleg_type = TRADITIONAL_DOGLEG;
-      use_nonmonotonic_steps = false;
-      max_consecutive_nonmonotonic_steps = 5;
-      max_num_iterations = 50;
-      max_solver_time_in_seconds = 1e9;
-      num_threads = 1;
-      initial_trust_region_radius = 1e4;
-      max_trust_region_radius = 1e16;
-      min_trust_region_radius = 1e-32;
-      min_relative_decrease = 1e-3;
-      min_lm_diagonal = 1e-6;
-      max_lm_diagonal = 1e32;
-      max_num_consecutive_invalid_steps = 5;
-      function_tolerance = 1e-6;
-      gradient_tolerance = 1e-10;
-      parameter_tolerance = 1e-8;
-
-#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE) && !defined(CERES_ENABLE_LGPL_CODE)  // NOLINT
-      linear_solver_type = DENSE_QR;
-#else
-      linear_solver_type = SPARSE_NORMAL_CHOLESKY;
-#endif
-
-      preconditioner_type = JACOBI;
-      visibility_clustering_type = CANONICAL_VIEWS;
-      dense_linear_algebra_library_type = EIGEN;
-
-      // Choose a default sparse linear algebra library in the order:
-      //
-      //   SUITE_SPARSE > CX_SPARSE > EIGEN_SPARSE > NO_SPARSE
-      sparse_linear_algebra_library_type = NO_SPARSE;
-#if !defined(CERES_NO_SUITESPARSE)
-      sparse_linear_algebra_library_type = SUITE_SPARSE;
-#else
-  #if !defined(CERES_NO_CXSPARSE)
-      sparse_linear_algebra_library_type = CX_SPARSE;
-  #else
-    #if defined(CERES_USE_EIGEN_SPARSE)
-      sparse_linear_algebra_library_type = EIGEN_SPARSE;
-    #endif
-  #endif
-#endif
-
-      num_linear_solver_threads = -1;
-      use_explicit_schur_complement = false;
-      use_postordering = false;
-      dynamic_sparsity = false;
-      min_linear_solver_iterations = 0;
-      max_linear_solver_iterations = 500;
-      eta = 1e-1;
-      jacobi_scaling = true;
-      use_inner_iterations = false;
-      inner_iteration_tolerance = 1e-3;
-      logging_type = PER_MINIMIZER_ITERATION;
-      minimizer_progress_to_stdout = false;
-      trust_region_problem_dump_directory = "/tmp";
-      trust_region_problem_dump_format_type = TEXTFILE;
-      check_gradients = false;
-      gradient_check_relative_precision = 1e-8;
-      gradient_check_numeric_derivative_relative_step_size = 1e-6;
-      update_state_every_iteration = false;
-      evaluation_callback = NULL;
-    }
-
     // Returns true if the options struct has a valid
     // configuration. Returns false otherwise, and fills in *error
     // with a message describing the problem.
@@ -171,11 +89,12 @@
     // trust region methods first choose a step size (the size of the
     // trust region) and then a step direction while line search methods
     // first choose a step direction and then a step size.
-    MinimizerType minimizer_type;
+    MinimizerType minimizer_type = TRUST_REGION;
 
-    LineSearchDirectionType line_search_direction_type;
-    LineSearchType line_search_type;
-    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
+    LineSearchDirectionType line_search_direction_type = LBFGS;
+    LineSearchType line_search_type = WOLFE;
+    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type =
+        FLETCHER_REEVES;
 
     // The LBFGS hessian approximation is a low rank approximation to
     // the inverse of the Hessian matrix. The rank of the
@@ -201,7 +120,7 @@
     //
     // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with
     // Limited Storage". Mathematics of Computation 35 (151): 773–782.
-    int max_lbfgs_rank;
+    int max_lbfgs_rank = 20;
 
     // As part of the (L)BFGS update step (BFGS) / right-multiply step (L-BFGS),
     // the initial inverse Hessian approximation is taken to be the Identity.
@@ -223,18 +142,18 @@
     // Oren S.S., Self-scaling variable metric (SSVM) algorithms
     // Part II: Implementation and experiments, Management Science,
     // 20(5), 863-874, 1974.
-    bool use_approximate_eigenvalue_bfgs_scaling;
+    bool use_approximate_eigenvalue_bfgs_scaling = false;
 
     // Degree of the polynomial used to approximate the objective
     // function. Valid values are BISECTION, QUADRATIC and CUBIC.
     //
     // BISECTION corresponds to pure backtracking search with no
     // interpolation.
-    LineSearchInterpolationType line_search_interpolation_type;
+    LineSearchInterpolationType line_search_interpolation_type = CUBIC;
 
     // If during the line search, the step_size falls below this
     // value, it is truncated to zero.
-    double min_line_search_step_size;
+    double min_line_search_step_size = 1e-9;
 
     // Line search parameters.
 
@@ -248,7 +167,7 @@
     //
     //   f(step_size) <= f(0) + sufficient_decrease * f'(0) * step_size
     //
-    double line_search_sufficient_function_decrease;
+    double line_search_sufficient_function_decrease = 1e-4;
 
     // In each iteration of the line search,
     //
@@ -258,7 +177,7 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double max_line_search_step_contraction;
+    double max_line_search_step_contraction = 1e-3;
 
     // In each iteration of the line search,
     //
@@ -268,19 +187,19 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double min_line_search_step_contraction;
+    double min_line_search_step_contraction = 0.6;
 
     // Maximum number of trial step size iterations during each line search,
     // if a step size satisfying the search conditions cannot be found within
     // this number of trials, the line search will terminate.
-    int max_num_line_search_step_size_iterations;
+    int max_num_line_search_step_size_iterations = 20;
 
     // Maximum number of restarts of the line search direction algorithm before
     // terminating the optimization. Restarts of the line search direction
     // algorithm occur when the current algorithm fails to produce a new descent
     // direction. This typically indicates a numerical failure, or a breakdown
     // in the validity of the approximations used.
-    int max_num_line_search_direction_restarts;
+    int max_num_line_search_direction_restarts = 5;
 
     // The strong Wolfe conditions consist of the Armijo sufficient
     // decrease condition, and an additional requirement that the
@@ -293,7 +212,7 @@
     //
     // Where f() is the line search objective and f'() is the derivative
     // of f w.r.t step_size (d f / d step_size).
-    double line_search_sufficient_curvature_decrease;
+    double line_search_sufficient_curvature_decrease = 0.9;
 
     // During the bracketing phase of the Wolfe search, the step size is
     // increased until either a point satisfying the Wolfe conditions is
@@ -304,12 +223,12 @@
     //   new_step_size <= max_step_expansion * step_size.
     //
     // By definition for expansion, max_step_expansion > 1.0.
-    double max_line_search_step_expansion;
+    double max_line_search_step_expansion = 10.0;
 
-    TrustRegionStrategyType trust_region_strategy_type;
+    TrustRegionStrategyType trust_region_strategy_type = LEVENBERG_MARQUARDT;
 
     // Type of dogleg strategy to use.
-    DoglegType dogleg_type;
+    DoglegType dogleg_type = TRADITIONAL_DOGLEG;
 
     // The classical trust region methods are descent methods, in that
     // they only accept a point if it strictly reduces the value of
@@ -336,30 +255,30 @@
     // than the minimum value encountered over the course of the
     // optimization, the final parameters returned to the user are the
     // ones corresponding to the minimum cost over all iterations.
-    bool use_nonmonotonic_steps;
-    int max_consecutive_nonmonotonic_steps;
+    bool use_nonmonotonic_steps = false;
+    int max_consecutive_nonmonotonic_steps = 5;
 
     // Maximum number of iterations for the minimizer to run for.
-    int max_num_iterations;
+    int max_num_iterations = 50;
 
     // Maximum time for which the minimizer should run for.
-    double max_solver_time_in_seconds;
+    double max_solver_time_in_seconds = 1e9;
 
     // Number of threads used by Ceres for evaluating the cost and
     // jacobians.
-    int num_threads;
+    int num_threads = 1;
 
     // Trust region minimizer settings.
-    double initial_trust_region_radius;
-    double max_trust_region_radius;
+    double initial_trust_region_radius = 1e4;
+    double max_trust_region_radius = 1e16;
 
     // Minimizer terminates when the trust region radius becomes
     // smaller than this value.
-    double min_trust_region_radius;
+    double min_trust_region_radius = 1e-32;
 
     // Lower bound for the relative decrease before a step is
     // accepted.
-    double min_relative_decrease;
+    double min_relative_decrease = 1e-3;
 
     // For the Levenberg-Marquadt algorithm, the scaled diagonal of
     // the normal equations J'J is used to control the size of the
@@ -368,46 +287,52 @@
     // fail. max_lm_diagonal and min_lm_diagonal, clamp the values of
     // diag(J'J) from above and below. In the normal course of
     // operation, the user should not have to modify these parameters.
-    double min_lm_diagonal;
-    double max_lm_diagonal;
+    double min_lm_diagonal = 1e-6;
+    double max_lm_diagonal = 1e32;
 
     // Sometimes due to numerical conditioning problems or linear
     // solver flakiness, the trust region strategy may return a
     // numerically invalid step that can be fixed by reducing the
     // trust region size. So the TrustRegionMinimizer allows for a few
     // successive invalid steps before it declares NUMERICAL_FAILURE.
-    int max_num_consecutive_invalid_steps;
+    int max_num_consecutive_invalid_steps = 5;
 
     // Minimizer terminates when
     //
     //   (new_cost - old_cost) < function_tolerance * old_cost;
     //
-    double function_tolerance;
+    double function_tolerance = 1e-6;
 
     // Minimizer terminates when
     //
     //   max_i |x - Project(Plus(x, -g(x))| < gradient_tolerance
     //
     // This value should typically be 1e-4 * function_tolerance.
-    double gradient_tolerance;
+    double gradient_tolerance = 1e-10;
 
     // Minimizer terminates when
     //
     //   |step|_2 <= parameter_tolerance * ( |x|_2 +  parameter_tolerance)
     //
-    double parameter_tolerance;
+    double parameter_tolerance = 1e-8;
 
     // Linear least squares solver options -------------------------------------
 
-    LinearSolverType linear_solver_type;
+    LinearSolverType linear_solver_type =
+#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE) && \
+    !defined(CERES_USE_EIGEN_SPARSE)  // NOLINT
+        DENSE_QR;
+#else
+        SPARSE_NORMAL_CHOLESKY;
+#endif
 
     // Type of preconditioner to use with the iterative linear solvers.
-    PreconditionerType preconditioner_type;
+    PreconditionerType preconditioner_type = JACOBI;
 
     // Type of clustering algorithm to use for visibility based
     // preconditioning. This option is used only when the
     // preconditioner_type is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
-    VisibilityClusteringType visibility_clustering_type;
+    VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
 
     // Ceres supports using multiple dense linear algebra libraries
     // for dense matrix factorizations. Currently EIGEN and LAPACK are
@@ -420,20 +345,33 @@
     // is a fine choice but for large problems, an optimized LAPACK +
     // BLAS implementation can make a substantial difference in
     // performance.
-    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
+    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type = EIGEN;
 
     // Ceres supports using multiple sparse linear algebra libraries
     // for sparse matrix ordering and factorizations. Currently,
     // SUITE_SPARSE and CX_SPARSE are the valid choices, depending on
     // whether they are linked into Ceres at build time.
-    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
+    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
+#if !defined(CERES_NO_SUITESPARSE)
+        SUITE_SPARSE;
+#else
+  #if !defined(CERES_NO_CXSPARSE)
+        CX_SPARSE;
+  #else
+    #if defined(CERES_USE_EIGEN_SPARSE)
+        EIGEN_SPARSE;
+    #else
+        NO_SPARSE;
+    #endif
+  #endif
+#endif
 
     // NOTE: This field is deprecated, and is ignored by
     // Ceres. Solver::Options::num_threads controls threading for all
     // of Ceres Solver.
     //
     // This setting is scheduled to be removed in 1.15.0.
-    int num_linear_solver_threads;
+    int num_linear_solver_threads = -1;
 
     // The order in which variables are eliminated in a linear solver
     // can have a significant of impact on the efficiency and accuracy
@@ -524,7 +462,7 @@
     //
     // NOTE: This option can only be used with the SCHUR_JACOBI
     // preconditioner.
-    bool use_explicit_schur_complement;
+    bool use_explicit_schur_complement = false;
 
     // Sparse Cholesky factorization algorithms use a fill-reducing
     // ordering to permute the columns of the Jacobian matrix. There
@@ -545,7 +483,7 @@
     // reordering algorithm which has slightly better runtime
     // performance at the expense of an extra copy of the Jacobian
     // matrix. Setting use_postordering to true enables this tradeoff.
-    bool use_postordering;
+    bool use_postordering = false;
 
     // Some non-linear least squares problems are symbolically dense but
     // numerically sparse. i.e. at any given state only a small number
@@ -560,7 +498,7 @@
     // likely lead to worse performance.
 
     // This settings affects the SPARSE_NORMAL_CHOLESKY solver.
-    bool dynamic_sparsity;
+    bool dynamic_sparsity = false;
 
     // Some non-linear least squares problems have additional
     // structure in the way the parameter blocks interact that it is
@@ -620,7 +558,7 @@
     // displays better convergence behaviour per iteration. Setting
     // Solver::Options::num_threads to the maximum number possible is
     // highly recommended.
-    bool use_inner_iterations;
+    bool use_inner_iterations = false;
 
     // If inner_iterations is true, then the user has two choices.
     //
@@ -643,17 +581,17 @@
     // inner iterations drops below inner_iteration_tolerance, the use
     // of inner iterations in subsequent trust region minimizer
     // iterations is disabled.
-    double inner_iteration_tolerance;
+    double inner_iteration_tolerance = 1e-3;
 
     // Minimum number of iterations for which the linear solver should
     // run, even if the convergence criterion is satisfied.
-    int min_linear_solver_iterations;
+    int min_linear_solver_iterations = 0;
 
     // Maximum number of iterations for which the linear solver should
     // run. If the solver does not converge in less than
     // max_linear_solver_iterations, then it returns MAX_ITERATIONS,
     // as its termination type.
-    int max_linear_solver_iterations;
+    int max_linear_solver_iterations = 500;
 
     // Forcing sequence parameter. The truncated Newton solver uses
     // this number to control the relative accuracy with which the
@@ -663,21 +601,21 @@
     // it to terminate the iterations when
     //
     //  (Q_i - Q_{i-1})/Q_i < eta/i
-    double eta;
+    double eta = 1e-1;
 
     // Normalize the jacobian using Jacobi scaling before calling
     // the linear least squares solver.
-    bool jacobi_scaling;
+    bool jacobi_scaling = true;
 
     // Logging options ---------------------------------------------------------
 
-    LoggingType logging_type;
+    LoggingType logging_type = PER_MINIMIZER_ITERATION;
 
     // By default the Minimizer progress is logged to VLOG(1), which
     // is sent to STDERR depending on the vlog level. If this flag is
     // set to true, and logging_type is not SILENT, the logging output
     // is sent to STDOUT.
-    bool minimizer_progress_to_stdout;
+    bool minimizer_progress_to_stdout = false;
 
     // List of iterations at which the minimizer should dump the trust
     // region problem. Useful for testing and benchmarking. If empty
@@ -688,8 +626,8 @@
     // non-empty if trust_region_minimizer_iterations_to_dump is
     // non-empty and trust_region_problem_dump_format_type is not
     // CONSOLE.
-    std::string trust_region_problem_dump_directory;
-    DumpFormatType trust_region_problem_dump_format_type;
+    std::string trust_region_problem_dump_directory = "/tmp";
+    DumpFormatType trust_region_problem_dump_format_type = TEXTFILE;
 
     // Finite differences options ----------------------------------------------
 
@@ -699,12 +637,12 @@
     // etc), then also computing it using finite differences. The
     // results are compared, and if they differ substantially, details
     // are printed to the log.
-    bool check_gradients;
+    bool check_gradients = false;
 
     // Relative precision to check for in the gradient checker. If the
     // relative difference between an element in a jacobian exceeds
     // this number, then the jacobian for that cost term is dumped.
-    double gradient_check_relative_precision;
+    double gradient_check_relative_precision = 1e-8;
 
     // WARNING: This option only applies to the to the numeric
     // differentiation used for checking the user provided derivatives
@@ -738,7 +676,7 @@
     // theory a good choice is sqrt(eps) * x, which for doubles means
     // about 1e-8 * x. However, I have found this number too
     // optimistic. This number should be exposed for users to change.
-    double gradient_check_numeric_derivative_relative_step_size;
+    double gradient_check_numeric_derivative_relative_step_size = 1e-6;
 
     // If true, the user's parameter blocks are updated at the end of
     // every Minimizer iteration, otherwise they are updated when the
@@ -760,7 +698,7 @@
     // BUT the solver will ensure that before the user provided
     // IterationCallbacks are called, the user visible state will be
     // updated to the current best point found by the solver.
-    bool update_state_every_iteration;
+    bool update_state_every_iteration = false;
 
     // Callbacks that are executed at the end of each iteration of the
     // Minimizer. An iteration may terminate midway, either due to
@@ -790,12 +728,10 @@
     // the documentation for that option for more details.
     //
     // The solver does NOT take ownership of the pointer.
-    EvaluationCallback* evaluation_callback;
+    EvaluationCallback* evaluation_callback = nullptr;
   };
 
   struct CERES_EXPORT Summary {
-    Summary();
-
     // A brief one line description of the state of the solver after
     // termination.
     std::string BriefReport() const;
@@ -807,25 +743,25 @@
     bool IsSolutionUsable() const;
 
     // Minimizer summary -------------------------------------------------
-    MinimizerType minimizer_type;
+    MinimizerType minimizer_type = TRUST_REGION;
 
-    TerminationType termination_type;
+    TerminationType termination_type = FAILURE;
 
     // Reason why the solver terminated.
-    std::string message;
+    std::string message = "ceres::Solve was not called.";
 
     // Cost of the problem (value of the objective function) before
     // the optimization.
-    double initial_cost;
+    double initial_cost = -1.0;
 
     // Cost of the problem (value of the objective function) after the
     // optimization.
-    double final_cost;
+    double final_cost = -1.0;
 
     // The part of the total cost that comes from residual blocks that
     // were held fixed by the preprocessor because all the parameter
     // blocks that they depend on were fixed.
-    double fixed_cost;
+    double fixed_cost = -1.0;
 
     // IterationSummary for each minimizer iteration in order.
     std::vector<IterationSummary> iterations;
@@ -834,22 +770,22 @@
     // accepted. Unless use_non_monotonic_steps is true this is also
     // the number of steps in which the objective function value/cost
     // went down.
-    int num_successful_steps;
+    int num_successful_steps = -1.0;
 
     // Number of minimizer iterations in which the step was rejected
     // either because it did not reduce the cost enough or the step
     // was not numerically valid.
-    int num_unsuccessful_steps;
+    int num_unsuccessful_steps = -1.0;
 
     // Number of times inner iterations were performed.
-    int num_inner_iteration_steps;
+    int num_inner_iteration_steps = -1.0;
 
     // Total number of iterations inside the line search algorithm
     // across all invocations. We call these iterations "steps" to
     // distinguish them from the outer iterations of the line search
     // and trust region minimizer algorithms which call the line
     // search algorithm as a subroutine.
-    int num_line_search_steps;
+    int num_line_search_steps = -1.0;
 
     // All times reported below are wall times.
 
@@ -857,42 +793,42 @@
     // occurs, Ceres performs a number of preprocessing steps. These
     // include error checks, memory allocations, and reorderings. This
     // time is accounted for as preprocessing time.
-    double preprocessor_time_in_seconds;
+    double preprocessor_time_in_seconds = -1.0;
 
     // Time spent in the TrustRegionMinimizer.
-    double minimizer_time_in_seconds;
+    double minimizer_time_in_seconds = -1.0;
 
     // After the Minimizer is finished, some time is spent in
     // re-evaluating residuals etc. This time is accounted for in the
     // postprocessor time.
-    double postprocessor_time_in_seconds;
+    double postprocessor_time_in_seconds = -1.0;
 
     // Some total of all time spent inside Ceres when Solve is called.
-    double total_time_in_seconds;
+    double total_time_in_seconds = -1.0;
 
     // Time (in seconds) spent in the linear solver computing the
     // trust region step.
-    double linear_solver_time_in_seconds;
+    double linear_solver_time_in_seconds = -1.0;
 
     // Number of times the Newton step was computed by solving a
     // linear system. This does not include linear solves used by
     // inner iterations.
-    int num_linear_solves;
+    int num_linear_solves = -1;
 
     // Time (in seconds) spent evaluating the residual vector.
-    double residual_evaluation_time_in_seconds;
+    double residual_evaluation_time_in_seconds = 1.0;
 
     // Number of residual only evaluations.
-    int num_residual_evaluations;
+    int num_residual_evaluations = -1;
 
     // Time (in seconds) spent evaluating the jacobian matrix.
-    double jacobian_evaluation_time_in_seconds;
+    double jacobian_evaluation_time_in_seconds = -1.0;
 
     // Number of Jacobian (and residual) evaluations.
-    int num_jacobian_evaluations;
+    int num_jacobian_evaluations = -1;
 
     // Time (in seconds) spent doing inner iterations.
-    double inner_iteration_time_in_seconds;
+    double inner_iteration_time_in_seconds = -1.0;
 
     // Cumulative timing information for line searches performed as part of the
     // solve.  Note that in addition to the case when the Line Search minimizer
@@ -901,69 +837,69 @@
 
     // Time (in seconds) spent evaluating the univariate cost function as part
     // of a line search.
-    double line_search_cost_evaluation_time_in_seconds;
+    double line_search_cost_evaluation_time_in_seconds = -1.0;
 
     // Time (in seconds) spent evaluating the gradient of the univariate cost
     // function as part of a line search.
-    double line_search_gradient_evaluation_time_in_seconds;
+    double line_search_gradient_evaluation_time_in_seconds = -1.0;
 
     // Time (in seconds) spent minimizing the interpolating polynomial
     // to compute the next candidate step size as part of a line search.
-    double line_search_polynomial_minimization_time_in_seconds;
+    double line_search_polynomial_minimization_time_in_seconds = -1.0;
 
     // Total time (in seconds) spent performing line searches.
-    double line_search_total_time_in_seconds;
+    double line_search_total_time_in_seconds = -1.0;
 
     // Number of parameter blocks in the problem.
-    int num_parameter_blocks;
+    int num_parameter_blocks = -1;
 
     // Number of parameters in the probem.
-    int num_parameters;
+    int num_parameters = -1;
 
     // Dimension of the tangent space of the problem (or the number of
     // columns in the Jacobian for the problem). This is different
     // from num_parameters if a parameter block is associated with a
     // LocalParameterization
-    int num_effective_parameters;
+    int num_effective_parameters = -1;
 
     // Number of residual blocks in the problem.
-    int num_residual_blocks;
+    int num_residual_blocks = -1;
 
     // Number of residuals in the problem.
-    int num_residuals;
+    int num_residuals = -1;
 
     // Number of parameter blocks in the problem after the inactive
     // and constant parameter blocks have been removed. A parameter
     // block is inactive if no residual block refers to it.
-    int num_parameter_blocks_reduced;
+    int num_parameter_blocks_reduced = -1;
 
     // Number of parameters in the reduced problem.
-    int num_parameters_reduced;
+    int num_parameters_reduced = -1;
 
     // Dimension of the tangent space of the reduced problem (or the
     // number of columns in the Jacobian for the reduced
     // problem). This is different from num_parameters_reduced if a
     // parameter block in the reduced problem is associated with a
     // LocalParameterization.
-    int num_effective_parameters_reduced;
+    int num_effective_parameters_reduced = -1;
 
     // Number of residual blocks in the reduced problem.
-    int num_residual_blocks_reduced;
+    int num_residual_blocks_reduced = -1;
 
     //  Number of residuals in the reduced problem.
-    int num_residuals_reduced;
+    int num_residuals_reduced = -1;
 
     // Is the reduced problem bounds constrained.
-    bool is_constrained;
+    bool is_constrained = false;
 
     //  Number of threads specified by the user for Jacobian and
     //  residual evaluation.
-    int num_threads_given;
+    int num_threads_given = -1;
 
     // Number of threads actually used by the solver for Jacobian and
     // residual evaluation. This number is not equal to
     // num_threads_given if OpenMP is not available.
-    int num_threads_used;
+    int num_threads_used = -1;
 
     // NOTE: This field is deprecated,
     // Solver::Summary::num_threads_given should be used instead.
@@ -974,7 +910,7 @@
     //
     // Number of threads specified by the user for solving the trust
     // region problem.
-    int num_linear_solver_threads_given;
+    int num_linear_solver_threads_given = -1;
 
     // NOTE: This field is deprecated,
     // Solver::Summary::num_threads_used should be used instead.
@@ -986,18 +922,29 @@
     // Number of threads actually used by the solver for solving the
     // trust region problem. This number is not equal to
     // num_threads_given if OpenMP is not available.
-    int num_linear_solver_threads_used;
+    int num_linear_solver_threads_used = -1;
 
     // Type of the linear solver requested by the user.
-    LinearSolverType linear_solver_type_given;
-
+    LinearSolverType linear_solver_type_given =
+#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE) && \
+    !defined(CERES_USE_EIGEN_SPARSE)  // NOLINT
+        DENSE_QR;
+#else
+        SPARSE_NORMAL_CHOLESKY;
+#endif
     // Type of the linear solver actually used. This may be different
     // from linear_solver_type_given if Ceres determines that the
     // problem structure is not compatible with the linear solver
     // requested or if the linear solver requested by the user is not
     // available, e.g. The user requested SPARSE_NORMAL_CHOLESKY but
     // no sparse linear algebra library was available.
-    LinearSolverType linear_solver_type_used;
+    LinearSolverType linear_solver_type_used =
+#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE) && \
+    !defined(CERES_USE_EIGEN_SPARSE)  // NOLINT
+        DENSE_QR;
+#else
+        SPARSE_NORMAL_CHOLESKY;
+#endif
 
     // Size of the elimination groups given by the user as hints to
     // the linear solver.
@@ -1027,13 +974,13 @@
 
     // True if the user asked for inner iterations to be used as part
     // of the optimization.
-    bool inner_iterations_given;
+    bool inner_iterations_given = false;
 
     // True if the user asked for inner iterations to be used as part
     // of the optimization and the problem structure was such that
     // they were actually performed. e.g., in a problem with just one
     // parameter block, inner iterations are not performed.
-    bool inner_iterations_used;
+    bool inner_iterations_used = false;
 
     // Size of the parameter groups given by the user for performing
     // inner iterations.
@@ -1048,51 +995,53 @@
     std::vector<int> inner_iteration_ordering_used;
 
     // Type of the preconditioner requested by the user.
-    PreconditionerType preconditioner_type_given;
+    PreconditionerType preconditioner_type_given = IDENTITY;
 
     // Type of the preconditioner actually used. This may be different
     // from linear_solver_type_given if Ceres determines that the
     // problem structure is not compatible with the linear solver
     // requested or if the linear solver requested by the user is not
     // available.
-    PreconditionerType preconditioner_type_used;
+    PreconditionerType preconditioner_type_used = IDENTITY;
 
     // Type of clustering algorithm used for visibility based
     // preconditioning. Only meaningful when the preconditioner_type
     // is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
-    VisibilityClusteringType visibility_clustering_type;
+    VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
 
     //  Type of trust region strategy.
-    TrustRegionStrategyType trust_region_strategy_type;
+    TrustRegionStrategyType trust_region_strategy_type = LEVENBERG_MARQUARDT;
 
     //  Type of dogleg strategy used for solving the trust region
     //  problem.
-    DoglegType dogleg_type;
+    DoglegType dogleg_type = TRADITIONAL_DOGLEG;
 
     //  Type of the dense linear algebra library used.
-    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
+    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type = EIGEN;
 
     // Type of the sparse linear algebra library used.
-    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
+    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
+        NO_SPARSE;
 
     // Type of line search direction used.
-    LineSearchDirectionType line_search_direction_type;
+    LineSearchDirectionType line_search_direction_type = LBFGS;
 
     // Type of the line search algorithm used.
-    LineSearchType line_search_type;
+    LineSearchType line_search_type = WOLFE;
 
     //  When performing line search, the degree of the polynomial used
     //  to approximate the objective function.
-    LineSearchInterpolationType line_search_interpolation_type;
+    LineSearchInterpolationType line_search_interpolation_type = CUBIC;
 
     // If the line search direction is NONLINEAR_CONJUGATE_GRADIENT,
     // then this indicates the particular variant of non-linear
     // conjugate gradient used.
-    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
+    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type =
+        FLETCHER_REEVES;
 
     // If the type of the line search direction is LBFGS, then this
     // indicates the rank of the Hessian approximation.
-    int max_lbfgs_rank;
+    int max_lbfgs_rank = -1;
   };
 
   // Once a least squares problem has been built, this function takes
diff --git a/include/ceres/tiny_solver.h b/include/ceres/tiny_solver.h
index 3cf1c32..d447665 100644
--- a/include/ceres/tiny_solver.h
+++ b/include/ceres/tiny_solver.h
@@ -148,32 +148,20 @@
   };
 
   struct Options {
-    Options()
-        : gradient_tolerance(1e-10),
-          parameter_tolerance(1e-8),
-          cost_threshold(std::numeric_limits<Scalar>::epsilon()),
-          initial_trust_region_radius(1e4),
-          max_num_iterations(50) {}
-    Scalar gradient_tolerance;   // eps > max(J'*f(x))
-    Scalar parameter_tolerance;  // eps > ||dx|| / ||x||
-    Scalar cost_threshold;       // eps > ||f(x)||
-    Scalar initial_trust_region_radius;
-    int max_num_iterations;
+    Scalar gradient_tolerance = 1e-10;  // eps > max(J'*f(x))
+    Scalar parameter_tolerance = 1e-8;  // eps > ||dx|| / ||x||
+    Scalar cost_threshold =             // eps > ||f(x)||
+        std::numeric_limits<Scalar>::epsilon();
+    Scalar initial_trust_region_radius = 1e4;
+    int max_num_iterations = 50;
   };
 
   struct Summary {
-    Summary()
-        : initial_cost(-1),
-          final_cost(-1),
-          gradient_max_norm(-1),
-          iterations(0),
-          status(HIT_MAX_ITERATIONS) {}
-
-    Scalar initial_cost;       // 1/2 ||f(x)||^2
-    Scalar final_cost;         // 1/2 ||f(x)||^2
-    Scalar gradient_max_norm;  // max(J'f(x))
-    int iterations;
-    Status status;
+    Scalar initial_cost = -1;       // 1/2 ||f(x)||^2
+    Scalar final_cost = -1;         // 1/2 ||f(x)||^2
+    Scalar gradient_max_norm = -1;  // max(J'f(x))
+    int iterations = -1;
+    Status status = HIT_MAX_ITERATIONS;
   };
 
   bool Update(const Function& function, const Parameters &x) {
@@ -211,6 +199,7 @@
     assert(x_and_min);
     Parameters& x = *x_and_min;
     summary = Summary();
+    summary.iterations = 0;
 
     // TODO(sameeragarwal): Deal with failure here.
     Update(function, x);
diff --git a/internal/ceres/block_sparse_matrix.h b/internal/ceres/block_sparse_matrix.h
index b93a6fa..aa74cfa 100644
--- a/internal/ceres/block_sparse_matrix.h
+++ b/internal/ceres/block_sparse_matrix.h
@@ -96,27 +96,17 @@
       const std::vector<Block>& column_blocks);
 
   struct RandomMatrixOptions {
-    RandomMatrixOptions()
-        : num_row_blocks(0),
-          min_row_block_size(0),
-          max_row_block_size(0),
-          num_col_blocks(0),
-          min_col_block_size(0),
-          max_col_block_size(0),
-          block_density(0.0) {
-    }
-
-    int num_row_blocks;
-    int min_row_block_size;
-    int max_row_block_size;
-    int num_col_blocks;
-    int min_col_block_size;
-    int max_col_block_size;
+    int num_row_blocks = 0;
+    int min_row_block_size = 0;
+    int max_row_block_size = 0;
+    int num_col_blocks = 0;
+    int min_col_block_size = 0;
+    int max_col_block_size = 0;
 
     // 0 < block_density <= 1 is the probability of a block being
     // present in the matrix. A given random matrix will not have
     // precisely this density.
-    double block_density;
+    double block_density = 0.0;
 
     // If col_blocks is non-empty, then the generated random matrix
     // has this block structure and the column related options in this
diff --git a/internal/ceres/canonical_views_clustering.h b/internal/ceres/canonical_views_clustering.h
index 6514827..e7d2472 100644
--- a/internal/ceres/canonical_views_clustering.h
+++ b/internal/ceres/canonical_views_clustering.h
@@ -101,27 +101,21 @@
     std::unordered_map<int, int>* membership);
 
 struct CanonicalViewsClusteringOptions {
-  CanonicalViewsClusteringOptions()
-      : min_views(3),
-        size_penalty_weight(5.75),
-        similarity_penalty_weight(100.0),
-        view_score_weight(0.0) {
-  }
   // The minimum number of canonical views to compute.
-  int min_views;
+  int min_views = 3;
 
   // Penalty weight for the number of canonical views.  A higher
   // number will result in fewer canonical views.
-  double size_penalty_weight;
+  double size_penalty_weight = 5.75;
 
   // Penalty weight for the diversity (orthogonality) of the
   // canonical views.  A higher number will encourage less similar
   // canonical views.
-  double similarity_penalty_weight;
+  double similarity_penalty_weight = 100;
 
   // Weight for per-view scores.  Lower weight places less
   // confidence in the view scores.
-  double view_score_weight;
+  double view_score_weight = 0.0;
 };
 
 }  // namespace internal
diff --git a/internal/ceres/compressed_row_sparse_matrix.h b/internal/ceres/compressed_row_sparse_matrix.h
index 54b3053..3e262e2 100644
--- a/internal/ceres/compressed_row_sparse_matrix.h
+++ b/internal/ceres/compressed_row_sparse_matrix.h
@@ -173,15 +173,6 @@
   // zero or not. If the answer is no, then we generate entries for the
   // block which are distributed normally.
   struct RandomMatrixOptions {
-    RandomMatrixOptions()
-        : storage_type(UNSYMMETRIC),
-          min_row_block_size(0),
-          max_row_block_size(0),
-          num_col_blocks(0),
-          min_col_block_size(0),
-          max_col_block_size(0),
-          block_density(0.0) {}
-
     // Type of matrix to create.
     //
     // If storage_type is UPPER_TRIANGULAR (LOWER_TRIANGULAR), then
@@ -189,19 +180,19 @@
     // (lower triangular) part. In this case, num_col_blocks,
     // min_col_block_size and max_col_block_size will be ignored and
     // assumed to be equal to the corresponding row settings.
-    StorageType storage_type;
+    StorageType storage_type = UNSYMMETRIC;
 
-    int num_row_blocks;
-    int min_row_block_size;
-    int max_row_block_size;
-    int num_col_blocks;
-    int min_col_block_size;
-    int max_col_block_size;
+    int num_row_blocks = 0;
+    int min_row_block_size = 0;
+    int max_row_block_size = 0;
+    int num_col_blocks = 0;
+    int min_col_block_size = 0;
+    int max_col_block_size = 0;
 
     // 0 < block_density <= 1 is the probability of a block being
     // present in the matrix. A given random matrix will not have
     // precisely this density.
-    double block_density;
+    double block_density = 0.0;
   };
 
   // Create a random CompressedRowSparseMatrix whose entries are
diff --git a/internal/ceres/evaluator.h b/internal/ceres/evaluator.h
index 532f437..b820958 100644
--- a/internal/ceres/evaluator.h
+++ b/internal/ceres/evaluator.h
@@ -60,20 +60,12 @@
   virtual ~Evaluator();
 
   struct Options {
-    Options()
-        : num_threads(1),
-          num_eliminate_blocks(-1),
-          linear_solver_type(DENSE_QR),
-          dynamic_sparsity(false),
-          context(NULL),
-          evaluation_callback(NULL) {}
-
-    int num_threads;
-    int num_eliminate_blocks;
-    LinearSolverType linear_solver_type;
-    bool dynamic_sparsity;
-    ContextImpl* context;
-    EvaluationCallback* evaluation_callback;
+    int num_threads = 1;
+    int num_eliminate_blocks = -1;
+    LinearSolverType linear_solver_type = DENSE_QR;
+    bool dynamic_sparsity = false;
+    ContextImpl* context = nullptr;
+    EvaluationCallback* evaluation_callback = nullptr;
   };
 
   static Evaluator* Create(const Options& options,
@@ -100,17 +92,12 @@
 
   // Options struct to control Evaluator::Evaluate;
   struct EvaluateOptions {
-    EvaluateOptions()
-        : apply_loss_function(true),
-          new_evaluation_point(true) {
-    }
-
     // If false, the loss function correction is not applied to the
     // residual blocks.
-    bool apply_loss_function;
+    bool apply_loss_function = true;
 
     // If false, this evaluation point is the same as the last one.
-    bool new_evaluation_point;
+    bool new_evaluation_point = true;
   };
 
   // Evaluate the cost function for the given state. Returns the cost,
diff --git a/internal/ceres/gradient_problem_solver.cc b/internal/ceres/gradient_problem_solver.cc
index 5ef36ad..47d8ade 100644
--- a/internal/ceres/gradient_problem_solver.cc
+++ b/internal/ceres/gradient_problem_solver.cc
@@ -194,26 +194,6 @@
   summary->total_time_in_seconds = WallTimeInSeconds() - start_time;
 }
 
-// Invalid values for most fields, to ensure that we are not
-// accidentally reporting default values.
-GradientProblemSolver::Summary::Summary()
-    : termination_type(FAILURE),
-      message("ceres::GradientProblemSolve was not called."),
-      initial_cost(-1.0),
-      final_cost(-1.0),
-      total_time_in_seconds(-1.0),
-      cost_evaluation_time_in_seconds(-1.0),
-      gradient_evaluation_time_in_seconds(-1.0),
-      line_search_polynomial_minimization_time_in_seconds(-1.0),
-      num_parameters(-1),
-      num_local_parameters(-1),
-      line_search_direction_type(LBFGS),
-      line_search_type(ARMIJO),
-      line_search_interpolation_type(BISECTION),
-      nonlinear_conjugate_gradient_type(FLETCHER_REEVES),
-      max_lbfgs_rank(-1) {
-}
-
 bool GradientProblemSolver::Summary::IsSolutionUsable() const {
   return internal::IsSolutionUsable(*this);
 }
diff --git a/internal/ceres/line_search.h b/internal/ceres/line_search.h
index 6e39a65..1c849a0 100644
--- a/internal/ceres/line_search.h
+++ b/internal/ceres/line_search.h
@@ -61,21 +61,9 @@
   struct Summary;
 
   struct Options {
-    Options()
-        : interpolation_type(CUBIC),
-          sufficient_decrease(1e-4),
-          max_step_contraction(1e-3),
-          min_step_contraction(0.9),
-          min_step_size(1e-9),
-          max_num_iterations(20),
-          sufficient_curvature_decrease(0.9),
-          max_step_expansion(10.0),
-          is_silent(false),
-          function(NULL) {}
-
     // Degree of the polynomial used to approximate the objective
     // function.
-    LineSearchInterpolationType interpolation_type;
+    LineSearchInterpolationType interpolation_type = CUBIC;
 
     // Armijo and Wolfe line search parameters.
 
@@ -88,7 +76,7 @@
     // s.t.
     //
     //  f(step_size) <= f(0) + sufficient_decrease * f'(0) * step_size
-    double sufficient_decrease;
+    double sufficient_decrease = 1e-4;
 
     // In each iteration of the Armijo / Wolfe line search,
     //
@@ -98,7 +86,7 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double max_step_contraction;
+    double max_step_contraction = 1e-3;
 
     // In each iteration of the Armijo / Wolfe line search,
     //
@@ -107,16 +95,16 @@
     //
     //  0 < max_step_contraction < min_step_contraction < 1
     //
-    double min_step_contraction;
+    double min_step_contraction = 0.9;
 
     // If during the line search, the step_size falls below this
     // value, it is truncated to zero.
-    double min_step_size;
+    double min_step_size = 1e-9;
 
     // Maximum number of trial step size iterations during each line search,
     // if a step size satisfying the search conditions cannot be found within
     // this number of trials, the line search will terminate.
-    int max_num_iterations;
+    int max_num_iterations = 20;
 
     // Wolfe-specific line search parameters.
 
@@ -131,7 +119,7 @@
     //
     // Where f() is the line search objective and f'() is the derivative
     // of f w.r.t step_size (d f / d step_size).
-    double sufficient_curvature_decrease;
+    double sufficient_curvature_decrease = 0.9;
 
     // During the bracketing phase of the Wolfe search, the step size is
     // increased until either a point satisfying the Wolfe conditions is
@@ -142,42 +130,32 @@
     //   new_step_size <= max_step_expansion * step_size.
     //
     // By definition for expansion, max_step_expansion > 1.0.
-    double max_step_expansion;
+    double max_step_expansion = 10;
 
-    bool is_silent;
+    bool is_silent = false;
 
     // The one dimensional function that the line search algorithm
     // minimizes.
-    LineSearchFunction* function;
+    LineSearchFunction* function = nullptr;
   };
 
   // Result of the line search.
   struct Summary {
-    Summary()
-        : success(false),
-          num_function_evaluations(0),
-          num_gradient_evaluations(0),
-          num_iterations(0),
-          cost_evaluation_time_in_seconds(-1.0),
-          gradient_evaluation_time_in_seconds(-1.0),
-          polynomial_minimization_time_in_seconds(-1.0),
-          total_time_in_seconds(-1.0) {}
-
-    bool success;
+    bool success = false;
     FunctionSample optimal_point;
-    int num_function_evaluations;
-    int num_gradient_evaluations;
-    int num_iterations;
+    int num_function_evaluations = 0;
+    int num_gradient_evaluations = 0;
+    int num_iterations = 0;
     // Cumulative time spent evaluating the value of the cost function across
     // all iterations.
-    double cost_evaluation_time_in_seconds;
+    double cost_evaluation_time_in_seconds = 0.0;
     // Cumulative time spent evaluating the gradient of the cost function across
     // all iterations.
-    double gradient_evaluation_time_in_seconds;
+    double gradient_evaluation_time_in_seconds = 0.0;
     // Cumulative time spent minimizing the interpolating polynomial to compute
     // the next candidate step size across all iterations.
-    double polynomial_minimization_time_in_seconds;
-    double total_time_in_seconds;
+    double polynomial_minimization_time_in_seconds = 0.0;
+    double total_time_in_seconds = 0.0;
     std::string error;
   };
 
diff --git a/internal/ceres/linear_solver.h b/internal/ceres/linear_solver.h
index 3f58cfb..a7c61bc 100644
--- a/internal/ceres/linear_solver.h
+++ b/internal/ceres/linear_solver.h
@@ -102,43 +102,25 @@
 class LinearSolver {
  public:
   struct Options {
-    Options()
-        : type(SPARSE_NORMAL_CHOLESKY),
-          preconditioner_type(JACOBI),
-          visibility_clustering_type(CANONICAL_VIEWS),
-          dense_linear_algebra_library_type(EIGEN),
-          sparse_linear_algebra_library_type(SUITE_SPARSE),
-          use_postordering(false),
-          dynamic_sparsity(false),
-          use_explicit_schur_complement(false),
-          min_num_iterations(1),
-          max_num_iterations(1),
-          num_threads(1),
-          residual_reset_period(10),
-          row_block_size(Eigen::Dynamic),
-          e_block_size(Eigen::Dynamic),
-          f_block_size(Eigen::Dynamic),
-          context(NULL) {
-    }
-
-    LinearSolverType type;
-    PreconditionerType preconditioner_type;
-    VisibilityClusteringType visibility_clustering_type;
-    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
-    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
+    LinearSolverType type = SPARSE_NORMAL_CHOLESKY;
+    PreconditionerType preconditioner_type = JACOBI;
+    VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
+    DenseLinearAlgebraLibraryType dense_linear_algebra_library_type = EIGEN;
+    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
+        SUITE_SPARSE;
 
     // See solver.h for information about these flags.
-    bool use_postordering;
-    bool dynamic_sparsity;
-    bool use_explicit_schur_complement;
+    bool use_postordering = false;
+    bool dynamic_sparsity = false;
+    bool use_explicit_schur_complement = false;
 
     // Number of internal iterations that the solver uses. This
     // parameter only makes sense for iterative solvers like CG.
-    int min_num_iterations;
-    int max_num_iterations;
+    int min_num_iterations = 1;
+    int max_num_iterations = 1;
 
     // If possible, how many threads can the solver use.
-    int num_threads;
+    int num_threads = 1;
 
     // Hints about the order in which the parameter blocks should be
     // eliminated by the linear solver.
@@ -163,7 +145,7 @@
     // inaccurate over time. Thus for non-zero values of this
     // parameter, the solver can be told to recalculate the value of
     // the residual using a |b - Ax| evaluation.
-    int residual_reset_period;
+    int residual_reset_period = 10;
 
     // If the block sizes in a BlockSparseMatrix are fixed, then in
     // some cases the Schur complement based solvers can detect and
@@ -174,22 +156,15 @@
     //
     // Please see schur_complement_solver.h and schur_eliminator.h for
     // more details.
-    int row_block_size;
-    int e_block_size;
-    int f_block_size;
+    int row_block_size = Eigen::Dynamic;
+    int e_block_size = Eigen::Dynamic;
+    int f_block_size = Eigen::Dynamic;
 
-    ContextImpl* context;
+    ContextImpl* context = nullptr;
   };
 
   // Options for the Solve method.
   struct PerSolveOptions {
-    PerSolveOptions()
-        : D(NULL),
-          preconditioner(NULL),
-          r_tolerance(0.0),
-          q_tolerance(0.0) {
-    }
-
     // This option only makes sense for unsymmetric linear solvers
     // that can solve rectangular linear systems.
     //
@@ -213,7 +188,7 @@
     // does not have full column rank, the results returned by the
     // solver cannot be relied on. D, if it is not null is an array of
     // size n.  b is an array of size m and x is an array of size n.
-    double * D;
+    double* D = nullptr;
 
     // This option only makes sense for iterative solvers.
     //
@@ -235,7 +210,7 @@
     //
     // A null preconditioner is equivalent to an identity matrix being
     // used a preconditioner.
-    LinearOperator* preconditioner;
+    LinearOperator* preconditioner = nullptr;
 
 
     // The following tolerance related options only makes sense for
@@ -247,7 +222,7 @@
     //
     // This is the most commonly used termination criterion for
     // iterative solvers.
-    double r_tolerance;
+    double r_tolerance = 0.0;
 
     // For PSD matrices A, let
     //
@@ -271,22 +246,16 @@
     //      Journal of Computational and Applied Mathematics,
     //      124(1-2), 45-59, 2000.
     //
-    double q_tolerance;
+    double q_tolerance = 0.0;
   };
 
   // Summary of a call to the Solve method. We should move away from
   // the true/false method for determining solver success. We should
   // let the summary object do the talking.
   struct Summary {
-    Summary()
-        : residual_norm(0.0),
-          num_iterations(-1),
-          termination_type(LINEAR_SOLVER_FAILURE) {
-    }
-
-    double residual_norm;
-    int num_iterations;
-    LinearSolverTerminationType termination_type;
+    double residual_norm = -1.0;
+    int num_iterations = -1;
+    LinearSolverTerminationType termination_type = LINEAR_SOLVER_FAILURE;
     std::string message;
   };
 
diff --git a/internal/ceres/preconditioner.h b/internal/ceres/preconditioner.h
index 37b221e..47fcdaf 100644
--- a/internal/ceres/preconditioner.h
+++ b/internal/ceres/preconditioner.h
@@ -48,22 +48,9 @@
 class Preconditioner : public LinearOperator {
  public:
   struct Options {
-    Options()
-        : type(JACOBI),
-          visibility_clustering_type(CANONICAL_VIEWS),
-          sparse_linear_algebra_library_type(SUITE_SPARSE),
-          subset_preconditioner_start_row_block(-1),
-          use_postordering(false),
-          num_threads(1),
-          row_block_size(Eigen::Dynamic),
-          e_block_size(Eigen::Dynamic),
-          f_block_size(Eigen::Dynamic),
-          context(NULL) {
-    }
-
-    PreconditionerType type;
-    VisibilityClusteringType visibility_clustering_type;
-    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
+    PreconditionerType type = JACOBI;
+    VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
+    SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type = SUITE_SPARSE;
 
     // When using the subset preconditioner, all row blocks starting
     // from this row block are used to construct the preconditioner.
@@ -75,13 +62,13 @@
     //
     // where P has subset_preconditioner_start_row_block row blocks,
     // and the preconditioner is the inverse of the matrix Q'Q.
-    int subset_preconditioner_start_row_block;
+    int subset_preconditioner_start_row_block = -1;
 
     // See solver.h for information about these flags.
-    bool use_postordering;
+    bool use_postordering = false;
 
     // If possible, how many threads the preconditioner can use.
-    int num_threads;
+    int num_threads = 1;
 
     // Hints about the order in which the parameter blocks should be
     // eliminated by the linear solver.
@@ -110,11 +97,11 @@
     //
     // Please see schur_complement_solver.h and schur_eliminator.h for
     // more details.
-    int row_block_size;
-    int e_block_size;
-    int f_block_size;
+    int row_block_size = Eigen::Dynamic;
+    int e_block_size = Eigen::Dynamic;
+    int f_block_size = Eigen::Dynamic;
 
-    ContextImpl* context;
+    ContextImpl* context = nullptr;
   };
 
   // If the optimization problem is such that there are no remaining
diff --git a/internal/ceres/single_linkage_clustering.h b/internal/ceres/single_linkage_clustering.h
index 374125b..ccd6f8e 100644
--- a/internal/ceres/single_linkage_clustering.h
+++ b/internal/ceres/single_linkage_clustering.h
@@ -38,13 +38,9 @@
 namespace internal {
 
 struct SingleLinkageClusteringOptions {
-  SingleLinkageClusteringOptions()
-      : min_similarity(0.99) {
-  }
-
   // Graph edges with edge weight less than min_similarity are ignored
   // during the clustering process.
-  double min_similarity;
+  double min_similarity = 0.99;
 };
 
 // Compute a partitioning of the vertices of the graph using the
diff --git a/internal/ceres/solver.cc b/internal/ceres/solver.cc
index de3329e..f66adce 100644
--- a/internal/ceres/solver.cc
+++ b/internal/ceres/solver.cc
@@ -635,66 +635,6 @@
   solver.Solve(options, problem, summary);
 }
 
-Solver::Summary::Summary()
-    // Invalid values for most fields, to ensure that we are not
-    // accidentally reporting default values.
-    : minimizer_type(TRUST_REGION),
-      termination_type(FAILURE),
-      message("ceres::Solve was not called."),
-      initial_cost(-1.0),
-      final_cost(-1.0),
-      fixed_cost(-1.0),
-      num_successful_steps(-1),
-      num_unsuccessful_steps(-1),
-      num_inner_iteration_steps(-1),
-      num_line_search_steps(-1),
-      preprocessor_time_in_seconds(-1.0),
-      minimizer_time_in_seconds(-1.0),
-      postprocessor_time_in_seconds(-1.0),
-      total_time_in_seconds(-1.0),
-      linear_solver_time_in_seconds(-1.0),
-      num_linear_solves(-1),
-      residual_evaluation_time_in_seconds(-1.0),
-      num_residual_evaluations(-1),
-      jacobian_evaluation_time_in_seconds(-1.0),
-      num_jacobian_evaluations(-1),
-      inner_iteration_time_in_seconds(-1.0),
-      line_search_cost_evaluation_time_in_seconds(-1.0),
-      line_search_gradient_evaluation_time_in_seconds(-1.0),
-      line_search_polynomial_minimization_time_in_seconds(-1.0),
-      line_search_total_time_in_seconds(-1.0),
-      num_parameter_blocks(-1),
-      num_parameters(-1),
-      num_effective_parameters(-1),
-      num_residual_blocks(-1),
-      num_residuals(-1),
-      num_parameter_blocks_reduced(-1),
-      num_parameters_reduced(-1),
-      num_effective_parameters_reduced(-1),
-      num_residual_blocks_reduced(-1),
-      num_residuals_reduced(-1),
-      is_constrained(false),
-      num_threads_given(-1),
-      num_threads_used(-1),
-      num_linear_solver_threads_given(-1),
-      num_linear_solver_threads_used(-1),
-      linear_solver_type_given(SPARSE_NORMAL_CHOLESKY),
-      linear_solver_type_used(SPARSE_NORMAL_CHOLESKY),
-      inner_iterations_given(false),
-      inner_iterations_used(false),
-      preconditioner_type_given(IDENTITY),
-      preconditioner_type_used(IDENTITY),
-      visibility_clustering_type(CANONICAL_VIEWS),
-      trust_region_strategy_type(LEVENBERG_MARQUARDT),
-      dense_linear_algebra_library_type(EIGEN),
-      sparse_linear_algebra_library_type(SUITE_SPARSE),
-      line_search_direction_type(LBFGS),
-      line_search_type(ARMIJO),
-      line_search_interpolation_type(BISECTION),
-      nonlinear_conjugate_gradient_type(FLETCHER_REEVES),
-      max_lbfgs_rank(-1) {
-}
-
 using internal::StringAppendF;
 using internal::StringPrintf;
 
diff --git a/internal/ceres/trust_region_strategy.h b/internal/ceres/trust_region_strategy.h
index 36e8e98..b3b2e5d 100644
--- a/internal/ceres/trust_region_strategy.h
+++ b/internal/ceres/trust_region_strategy.h
@@ -56,43 +56,29 @@
 class TrustRegionStrategy {
  public:
   struct Options {
-    Options()
-        : trust_region_strategy_type(LEVENBERG_MARQUARDT),
-          initial_radius(1e4),
-          max_radius(1e32),
-          min_lm_diagonal(1e-6),
-          max_lm_diagonal(1e32),
-          dogleg_type(TRADITIONAL_DOGLEG) {
-    }
-
-    TrustRegionStrategyType trust_region_strategy_type;
+    TrustRegionStrategyType trust_region_strategy_type = LEVENBERG_MARQUARDT;
     // Linear solver used for actually solving the trust region step.
-    LinearSolver* linear_solver;
-    double initial_radius;
-    double max_radius;
+    LinearSolver* linear_solver = nullptr;
+    double initial_radius = 1e4;
+    double max_radius = 1e32;
 
     // Minimum and maximum values of the diagonal damping matrix used
     // by LevenbergMarquardtStrategy. The DoglegStrategy also uses
     // these bounds to construct a regularizing diagonal to ensure
     // that the Gauss-Newton step computation is of full rank.
-    double min_lm_diagonal;
-    double max_lm_diagonal;
+    double min_lm_diagonal = 1e-6;
+    double max_lm_diagonal = 1e32;
 
     // Further specify which dogleg method to use
-    DoglegType dogleg_type;
+    DoglegType dogleg_type = TRADITIONAL_DOGLEG;
   };
 
   // Per solve options.
   struct PerSolveOptions {
-    PerSolveOptions()
-        : eta(0),
-          dump_format_type(TEXTFILE) {
-    }
-
     // Forcing sequence for inexact solves.
-    double eta;
+    double eta = 1e-1;
 
-    DumpFormatType dump_format_type;
+    DumpFormatType dump_format_type = TEXTFILE;
 
     // If non-empty and dump_format_type is not CONSOLE, the trust
     // regions strategy will write the linear system to file(s) with
@@ -103,12 +89,6 @@
   };
 
   struct Summary {
-    Summary()
-        : residual_norm(0.0),
-          num_iterations(-1),
-          termination_type(LINEAR_SOLVER_FAILURE) {
-    }
-
     // If the trust region problem is,
     //
     //   1/2 x'Ax + b'x + c,
@@ -116,15 +96,15 @@
     // then
     //
     //   residual_norm = |Ax -b|
-    double residual_norm;
+    double residual_norm = -1;
 
     // Number of iterations used by the linear solver. If a linear
     // solver was not called (e.g., DogLegStrategy after an
     // unsuccessful step), then this would be zero.
-    int num_iterations;
+    int num_iterations = -1;
 
     // Status of the linear solver used to solve the Newton system.
-    LinearSolverTerminationType termination_type;
+    LinearSolverTerminationType termination_type = LINEAR_SOLVER_FAILURE;
   };
 
   virtual ~TrustRegionStrategy();