blob: dc13d19646eb93013cb3d236b2639a5ffcb5369a [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// http://ceres-solver.org/
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7//
8// * Redistributions of source code must retain the above copyright notice,
9// this list of conditions and the following disclaimer.
10// * Redistributions in binary form must reproduce the above copyright notice,
11// this list of conditions and the following disclaimer in the documentation
12// and/or other materials provided with the distribution.
13// * Neither the name of Google Inc. nor the names of its contributors may be
14// used to endorse or promote products derived from this software without
15// specific prior written permission.
16//
17// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27// POSSIBILITY OF SUCH DAMAGE.
28//
29// Author: strandmark@google.com (Petter Strandmark)
30//
31// Denoising using Fields of Experts and the Ceres minimizer.
32//
33// Note that for good denoising results the weighting between the data term
34// and the Fields of Experts term needs to be adjusted. This is discussed
35// in [1]. This program assumes Gaussian noise. The noise model can be changed
Austin Schuh3de38b02024-06-25 18:25:10 -070036// by substituting another function for QuadraticCostFunction.
Austin Schuh70cc9552019-01-21 19:46:48 -080037//
38// [1] S. Roth and M.J. Black. "Fields of Experts." International Journal of
39// Computer Vision, 82(2):205--229, 2009.
40
41#include <algorithm>
42#include <cmath>
43#include <iostream>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080044#include <random>
Austin Schuh70cc9552019-01-21 19:46:48 -080045#include <sstream>
46#include <string>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080047#include <vector>
Austin Schuh70cc9552019-01-21 19:46:48 -080048
49#include "ceres/ceres.h"
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080050#include "fields_of_experts.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080051#include "gflags/gflags.h"
52#include "glog/logging.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080053#include "pgm_image.h"
54
55DEFINE_string(input, "", "File to which the output image should be written");
56DEFINE_string(foe_file, "", "FoE file to use");
57DEFINE_string(output, "", "File to which the output image should be written");
58DEFINE_double(sigma, 20.0, "Standard deviation of noise");
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080059DEFINE_string(trust_region_strategy,
60 "levenberg_marquardt",
61 "Options are: levenberg_marquardt, dogleg.");
62DEFINE_string(dogleg,
63 "traditional_dogleg",
64 "Options are: traditional_dogleg,"
65 "subspace_dogleg.");
66DEFINE_string(linear_solver,
67 "sparse_normal_cholesky",
68 "Options are: "
69 "sparse_normal_cholesky and cgnr.");
70DEFINE_string(preconditioner,
71 "jacobi",
72 "Options are: "
73 "identity, jacobi, subset");
74DEFINE_string(sparse_linear_algebra_library,
75 "suite_sparse",
76 "Options are: suite_sparse, cx_sparse and eigen_sparse");
77DEFINE_double(eta,
78 1e-2,
79 "Default value for eta. Eta determines the "
80 "accuracy of each linear solve of the truncated newton step. "
81 "Changing this parameter can affect solve performance.");
82DEFINE_int32(num_threads, 1, "Number of threads.");
83DEFINE_int32(num_iterations, 10, "Number of iterations.");
84DEFINE_bool(nonmonotonic_steps,
85 false,
86 "Trust region algorithm can use"
87 " nonmonotic steps.");
88DEFINE_bool(inner_iterations,
89 false,
90 "Use inner iterations to non-linearly "
91 "refine each successful trust region step.");
92DEFINE_bool(mixed_precision_solves, false, "Use mixed precision solves.");
93DEFINE_int32(max_num_refinement_iterations,
94 0,
95 "Iterative refinement iterations");
96DEFINE_bool(line_search,
97 false,
98 "Use a line search instead of trust region "
Austin Schuh70cc9552019-01-21 19:46:48 -080099 "algorithm.");
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800100DEFINE_double(subset_fraction,
101 0.2,
102 "The fraction of residual blocks to use for the"
103 " subset preconditioner.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800104
Austin Schuh3de38b02024-06-25 18:25:10 -0700105namespace ceres::examples {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800106namespace {
Austin Schuh70cc9552019-01-21 19:46:48 -0800107
108// This cost function is used to build the data term.
109//
110// f_i(x) = a * (x_i - b)^2
111//
112class QuadraticCostFunction : public ceres::SizedCostFunction<1, 1> {
113 public:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800114 QuadraticCostFunction(double a, double b) : sqrta_(std::sqrt(a)), b_(b) {}
Austin Schuh3de38b02024-06-25 18:25:10 -0700115 bool Evaluate(double const* const* parameters,
116 double* residuals,
117 double** jacobians) const override {
Austin Schuh70cc9552019-01-21 19:46:48 -0800118 const double x = parameters[0][0];
119 residuals[0] = sqrta_ * (x - b_);
Austin Schuh3de38b02024-06-25 18:25:10 -0700120 if (jacobians != nullptr && jacobians[0] != nullptr) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800121 jacobians[0][0] = sqrta_;
122 }
123 return true;
124 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800125
Austin Schuh70cc9552019-01-21 19:46:48 -0800126 private:
127 double sqrta_, b_;
128};
129
130// Creates a Fields of Experts MAP inference problem.
131void CreateProblem(const FieldsOfExperts& foe,
132 const PGMImage<double>& image,
133 Problem* problem,
134 PGMImage<double>* solution) {
135 // Create the data term
Austin Schuh3de38b02024-06-25 18:25:10 -0700136 CHECK_GT(CERES_GET_FLAG(FLAGS_sigma), 0.0);
137 const double coefficient =
138 1 / (2.0 * CERES_GET_FLAG(FLAGS_sigma) * CERES_GET_FLAG(FLAGS_sigma));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800139 for (int index = 0; index < image.NumPixels(); ++index) {
140 ceres::CostFunction* cost_function = new QuadraticCostFunction(
141 coefficient, image.PixelFromLinearIndex(index));
142 problem->AddResidualBlock(
Austin Schuh3de38b02024-06-25 18:25:10 -0700143 cost_function, nullptr, solution->MutablePixelFromLinearIndex(index));
Austin Schuh70cc9552019-01-21 19:46:48 -0800144 }
145
146 // Create Ceres cost and loss functions for regularization. One is needed for
147 // each filter.
148 std::vector<ceres::LossFunction*> loss_function(foe.NumFilters());
149 std::vector<ceres::CostFunction*> cost_function(foe.NumFilters());
150 for (int alpha_index = 0; alpha_index < foe.NumFilters(); ++alpha_index) {
151 loss_function[alpha_index] = foe.NewLossFunction(alpha_index);
152 cost_function[alpha_index] = foe.NewCostFunction(alpha_index);
153 }
154
155 // Add FoE regularization for each patch in the image.
156 for (int x = 0; x < image.width() - (foe.Size() - 1); ++x) {
157 for (int y = 0; y < image.height() - (foe.Size() - 1); ++y) {
158 // Build a vector with the pixel indices of this patch.
159 std::vector<double*> pixels;
160 const std::vector<int>& x_delta_indices = foe.GetXDeltaIndices();
161 const std::vector<int>& y_delta_indices = foe.GetYDeltaIndices();
162 for (int i = 0; i < foe.NumVariables(); ++i) {
163 double* pixel = solution->MutablePixel(x + x_delta_indices[i],
164 y + y_delta_indices[i]);
165 pixels.push_back(pixel);
166 }
167 // For this patch with coordinates (x, y), we will add foe.NumFilters()
168 // terms to the objective function.
169 for (int alpha_index = 0; alpha_index < foe.NumFilters(); ++alpha_index) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800170 problem->AddResidualBlock(
171 cost_function[alpha_index], loss_function[alpha_index], pixels);
Austin Schuh70cc9552019-01-21 19:46:48 -0800172 }
173 }
174 }
175}
176
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800177void SetLinearSolver(Solver::Options* options) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700178 CHECK(StringToLinearSolverType(CERES_GET_FLAG(FLAGS_linear_solver),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800179 &options->linear_solver_type));
Austin Schuh3de38b02024-06-25 18:25:10 -0700180 CHECK(StringToPreconditionerType(CERES_GET_FLAG(FLAGS_preconditioner),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800181 &options->preconditioner_type));
182 CHECK(StringToSparseLinearAlgebraLibraryType(
Austin Schuh3de38b02024-06-25 18:25:10 -0700183 CERES_GET_FLAG(FLAGS_sparse_linear_algebra_library),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800184 &options->sparse_linear_algebra_library_type));
Austin Schuh3de38b02024-06-25 18:25:10 -0700185 options->use_mixed_precision_solves =
186 CERES_GET_FLAG(FLAGS_mixed_precision_solves);
187 options->max_num_refinement_iterations =
188 CERES_GET_FLAG(FLAGS_max_num_refinement_iterations);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800189}
190
191void SetMinimizerOptions(Solver::Options* options) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700192 options->max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800193 options->minimizer_progress_to_stdout = true;
Austin Schuh3de38b02024-06-25 18:25:10 -0700194 options->num_threads = CERES_GET_FLAG(FLAGS_num_threads);
195 options->eta = CERES_GET_FLAG(FLAGS_eta);
196 options->use_nonmonotonic_steps = CERES_GET_FLAG(FLAGS_nonmonotonic_steps);
197 if (CERES_GET_FLAG(FLAGS_line_search)) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800198 options->minimizer_type = ceres::LINE_SEARCH;
199 }
200
Austin Schuh3de38b02024-06-25 18:25:10 -0700201 CHECK(StringToTrustRegionStrategyType(
202 CERES_GET_FLAG(FLAGS_trust_region_strategy),
203 &options->trust_region_strategy_type));
204 CHECK(
205 StringToDoglegType(CERES_GET_FLAG(FLAGS_dogleg), &options->dogleg_type));
206 options->use_inner_iterations = CERES_GET_FLAG(FLAGS_inner_iterations);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800207}
208
Austin Schuh70cc9552019-01-21 19:46:48 -0800209// Solves the FoE problem using Ceres and post-processes it to make sure the
210// solution stays within [0, 255].
211void SolveProblem(Problem* problem, PGMImage<double>* solution) {
212 // These parameters may be experimented with. For example, ceres::DOGLEG tends
213 // to be faster for 2x2 filters, but gives solutions with slightly higher
214 // objective function value.
215 ceres::Solver::Options options;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800216 SetMinimizerOptions(&options);
217 SetLinearSolver(&options);
Austin Schuh70cc9552019-01-21 19:46:48 -0800218 options.function_tolerance = 1e-3; // Enough for denoising.
219
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800220 if (options.linear_solver_type == ceres::CGNR &&
221 options.preconditioner_type == ceres::SUBSET) {
222 std::vector<ResidualBlockId> residual_blocks;
223 problem->GetResidualBlocks(&residual_blocks);
224
225 // To use the SUBSET preconditioner we need to provide a list of
226 // residual blocks (rows of the Jacobian). The denoising problem
227 // has fairly general sparsity, and there is no apriori reason to
228 // select one residual block over another, so we will randomly
229 // subsample the residual blocks with probability subset_fraction.
230 std::default_random_engine engine;
231 std::uniform_real_distribution<> distribution(0, 1); // rage 0 - 1
232 for (auto residual_block : residual_blocks) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700233 if (distribution(engine) <= CERES_GET_FLAG(FLAGS_subset_fraction)) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800234 options.residual_blocks_for_subset_preconditioner.insert(
235 residual_block);
236 }
237 }
238 }
239
Austin Schuh70cc9552019-01-21 19:46:48 -0800240 ceres::Solver::Summary summary;
241 ceres::Solve(options, problem, &summary);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800242 std::cout << summary.FullReport() << "\n";
Austin Schuh70cc9552019-01-21 19:46:48 -0800243
244 // Make the solution stay in [0, 255].
245 for (int x = 0; x < solution->width(); ++x) {
246 for (int y = 0; y < solution->height(); ++y) {
247 *solution->MutablePixel(x, y) =
248 std::min(255.0, std::max(0.0, solution->Pixel(x, y)));
249 }
250 }
251}
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800252
253} // namespace
Austin Schuh3de38b02024-06-25 18:25:10 -0700254} // namespace ceres::examples
Austin Schuh70cc9552019-01-21 19:46:48 -0800255
256int main(int argc, char** argv) {
257 using namespace ceres::examples;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800258 GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
Austin Schuh70cc9552019-01-21 19:46:48 -0800259 google::InitGoogleLogging(argv[0]);
260
Austin Schuh3de38b02024-06-25 18:25:10 -0700261 if (CERES_GET_FLAG(FLAGS_input).empty()) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800262 std::cerr << "Please provide an image file name using -input.\n";
Austin Schuh70cc9552019-01-21 19:46:48 -0800263 return 1;
264 }
265
Austin Schuh3de38b02024-06-25 18:25:10 -0700266 if (CERES_GET_FLAG(FLAGS_foe_file).empty()) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800267 std::cerr << "Please provide a Fields of Experts file name using -foe_file."
268 "\n";
Austin Schuh70cc9552019-01-21 19:46:48 -0800269 return 1;
270 }
271
272 // Load the Fields of Experts filters from file.
273 FieldsOfExperts foe;
Austin Schuh3de38b02024-06-25 18:25:10 -0700274 if (!foe.LoadFromFile(CERES_GET_FLAG(FLAGS_foe_file))) {
275 std::cerr << "Loading \"" << CERES_GET_FLAG(FLAGS_foe_file)
276 << "\" failed.\n";
Austin Schuh70cc9552019-01-21 19:46:48 -0800277 return 2;
278 }
279
280 // Read the images
Austin Schuh3de38b02024-06-25 18:25:10 -0700281 PGMImage<double> image(CERES_GET_FLAG(FLAGS_input));
Austin Schuh70cc9552019-01-21 19:46:48 -0800282 if (image.width() == 0) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700283 std::cerr << "Reading \"" << CERES_GET_FLAG(FLAGS_input) << "\" failed.\n";
Austin Schuh70cc9552019-01-21 19:46:48 -0800284 return 3;
285 }
286 PGMImage<double> solution(image.width(), image.height());
287 solution.Set(0.0);
288
289 ceres::Problem problem;
290 CreateProblem(foe, image, &problem, &solution);
291
292 SolveProblem(&problem, &solution);
293
Austin Schuh3de38b02024-06-25 18:25:10 -0700294 if (!CERES_GET_FLAG(FLAGS_output).empty()) {
295 CHECK(solution.WriteToFile(CERES_GET_FLAG(FLAGS_output)))
296 << "Writing \"" << CERES_GET_FLAG(FLAGS_output) << "\" failed.";
Austin Schuh70cc9552019-01-21 19:46:48 -0800297 }
298
299 return 0;
300}