denoising.cc 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2015 Google Inc. All rights reserved.
  3. // 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
  36. // by substituting another function for QuadraticCostFunction.
  37. //
  38. // [1] S. Roth and M.J. Black. "Fields of Experts." International Journal of
  39. // Computer Vision, 82(2):205--229, 2009.
  40. #include <algorithm>
  41. #include <cmath>
  42. #include <iostream>
  43. #include <random>
  44. #include <sstream>
  45. #include <string>
  46. #include <vector>
  47. #include "ceres/ceres.h"
  48. #include "fields_of_experts.h"
  49. #include "gflags/gflags.h"
  50. #include "glog/logging.h"
  51. #include "pgm_image.h"
  52. DEFINE_string(input, "", "File to which the output image should be written");
  53. DEFINE_string(foe_file, "", "FoE file to use");
  54. DEFINE_string(output, "", "File to which the output image should be written");
  55. DEFINE_double(sigma, 20.0, "Standard deviation of noise");
  56. DEFINE_string(trust_region_strategy,
  57. "levenberg_marquardt",
  58. "Options are: levenberg_marquardt, dogleg.");
  59. DEFINE_string(dogleg,
  60. "traditional_dogleg",
  61. "Options are: traditional_dogleg,"
  62. "subspace_dogleg.");
  63. DEFINE_string(linear_solver,
  64. "sparse_normal_cholesky",
  65. "Options are: "
  66. "sparse_normal_cholesky and cgnr.");
  67. DEFINE_string(preconditioner,
  68. "jacobi",
  69. "Options are: "
  70. "identity, jacobi, subset");
  71. DEFINE_string(sparse_linear_algebra_library,
  72. "suite_sparse",
  73. "Options are: suite_sparse, cx_sparse and eigen_sparse");
  74. DEFINE_double(eta,
  75. 1e-2,
  76. "Default value for eta. Eta determines the "
  77. "accuracy of each linear solve of the truncated newton step. "
  78. "Changing this parameter can affect solve performance.");
  79. DEFINE_int32(num_threads, 1, "Number of threads.");
  80. DEFINE_int32(num_iterations, 10, "Number of iterations.");
  81. DEFINE_bool(nonmonotonic_steps,
  82. false,
  83. "Trust region algorithm can use"
  84. " nonmonotic steps.");
  85. DEFINE_bool(inner_iterations,
  86. false,
  87. "Use inner iterations to non-linearly "
  88. "refine each successful trust region step.");
  89. DEFINE_bool(mixed_precision_solves, false, "Use mixed precision solves.");
  90. DEFINE_int32(max_num_refinement_iterations,
  91. 0,
  92. "Iterative refinement iterations");
  93. DEFINE_bool(line_search,
  94. false,
  95. "Use a line search instead of trust region "
  96. "algorithm.");
  97. DEFINE_double(subset_fraction,
  98. 0.2,
  99. "The fraction of residual blocks to use for the"
  100. " subset preconditioner.");
  101. namespace ceres::examples {
  102. namespace {
  103. // This cost function is used to build the data term.
  104. //
  105. // f_i(x) = a * (x_i - b)^2
  106. //
  107. class QuadraticCostFunction : public ceres::SizedCostFunction<1, 1> {
  108. public:
  109. QuadraticCostFunction(double a, double b) : sqrta_(std::sqrt(a)), b_(b) {}
  110. bool Evaluate(double const* const* parameters,
  111. double* residuals,
  112. double** jacobians) const override {
  113. const double x = parameters[0][0];
  114. residuals[0] = sqrta_ * (x - b_);
  115. if (jacobians != nullptr && jacobians[0] != nullptr) {
  116. jacobians[0][0] = sqrta_;
  117. }
  118. return true;
  119. }
  120. private:
  121. double sqrta_, b_;
  122. };
  123. // Creates a Fields of Experts MAP inference problem.
  124. void CreateProblem(const FieldsOfExperts& foe,
  125. const PGMImage<double>& image,
  126. Problem* problem,
  127. PGMImage<double>* solution) {
  128. // Create the data term
  129. CHECK_GT(CERES_GET_FLAG(FLAGS_sigma), 0.0);
  130. const double coefficient =
  131. 1 / (2.0 * CERES_GET_FLAG(FLAGS_sigma) * CERES_GET_FLAG(FLAGS_sigma));
  132. for (int index = 0; index < image.NumPixels(); ++index) {
  133. ceres::CostFunction* cost_function = new QuadraticCostFunction(
  134. coefficient, image.PixelFromLinearIndex(index));
  135. problem->AddResidualBlock(
  136. cost_function, nullptr, solution->MutablePixelFromLinearIndex(index));
  137. }
  138. // Create Ceres cost and loss functions for regularization. One is needed for
  139. // each filter.
  140. std::vector<ceres::LossFunction*> loss_function(foe.NumFilters());
  141. std::vector<ceres::CostFunction*> cost_function(foe.NumFilters());
  142. for (int alpha_index = 0; alpha_index < foe.NumFilters(); ++alpha_index) {
  143. loss_function[alpha_index] = foe.NewLossFunction(alpha_index);
  144. cost_function[alpha_index] = foe.NewCostFunction(alpha_index);
  145. }
  146. // Add FoE regularization for each patch in the image.
  147. for (int x = 0; x < image.width() - (foe.Size() - 1); ++x) {
  148. for (int y = 0; y < image.height() - (foe.Size() - 1); ++y) {
  149. // Build a vector with the pixel indices of this patch.
  150. std::vector<double*> pixels;
  151. const std::vector<int>& x_delta_indices = foe.GetXDeltaIndices();
  152. const std::vector<int>& y_delta_indices = foe.GetYDeltaIndices();
  153. for (int i = 0; i < foe.NumVariables(); ++i) {
  154. double* pixel = solution->MutablePixel(x + x_delta_indices[i],
  155. y + y_delta_indices[i]);
  156. pixels.push_back(pixel);
  157. }
  158. // For this patch with coordinates (x, y), we will add foe.NumFilters()
  159. // terms to the objective function.
  160. for (int alpha_index = 0; alpha_index < foe.NumFilters(); ++alpha_index) {
  161. problem->AddResidualBlock(
  162. cost_function[alpha_index], loss_function[alpha_index], pixels);
  163. }
  164. }
  165. }
  166. }
  167. void SetLinearSolver(Solver::Options* options) {
  168. CHECK(StringToLinearSolverType(CERES_GET_FLAG(FLAGS_linear_solver),
  169. &options->linear_solver_type));
  170. CHECK(StringToPreconditionerType(CERES_GET_FLAG(FLAGS_preconditioner),
  171. &options->preconditioner_type));
  172. CHECK(StringToSparseLinearAlgebraLibraryType(
  173. CERES_GET_FLAG(FLAGS_sparse_linear_algebra_library),
  174. &options->sparse_linear_algebra_library_type));
  175. options->use_mixed_precision_solves =
  176. CERES_GET_FLAG(FLAGS_mixed_precision_solves);
  177. options->max_num_refinement_iterations =
  178. CERES_GET_FLAG(FLAGS_max_num_refinement_iterations);
  179. }
  180. void SetMinimizerOptions(Solver::Options* options) {
  181. options->max_num_iterations = CERES_GET_FLAG(FLAGS_num_iterations);
  182. options->minimizer_progress_to_stdout = true;
  183. options->num_threads = CERES_GET_FLAG(FLAGS_num_threads);
  184. options->eta = CERES_GET_FLAG(FLAGS_eta);
  185. options->use_nonmonotonic_steps = CERES_GET_FLAG(FLAGS_nonmonotonic_steps);
  186. if (CERES_GET_FLAG(FLAGS_line_search)) {
  187. options->minimizer_type = ceres::LINE_SEARCH;
  188. }
  189. CHECK(StringToTrustRegionStrategyType(
  190. CERES_GET_FLAG(FLAGS_trust_region_strategy),
  191. &options->trust_region_strategy_type));
  192. CHECK(
  193. StringToDoglegType(CERES_GET_FLAG(FLAGS_dogleg), &options->dogleg_type));
  194. options->use_inner_iterations = CERES_GET_FLAG(FLAGS_inner_iterations);
  195. }
  196. // Solves the FoE problem using Ceres and post-processes it to make sure the
  197. // solution stays within [0, 255].
  198. void SolveProblem(Problem* problem, PGMImage<double>* solution) {
  199. // These parameters may be experimented with. For example, ceres::DOGLEG tends
  200. // to be faster for 2x2 filters, but gives solutions with slightly higher
  201. // objective function value.
  202. ceres::Solver::Options options;
  203. SetMinimizerOptions(&options);
  204. SetLinearSolver(&options);
  205. options.function_tolerance = 1e-3; // Enough for denoising.
  206. if (options.linear_solver_type == ceres::CGNR &&
  207. options.preconditioner_type == ceres::SUBSET) {
  208. std::vector<ResidualBlockId> residual_blocks;
  209. problem->GetResidualBlocks(&residual_blocks);
  210. // To use the SUBSET preconditioner we need to provide a list of
  211. // residual blocks (rows of the Jacobian). The denoising problem
  212. // has fairly general sparsity, and there is no apriori reason to
  213. // select one residual block over another, so we will randomly
  214. // subsample the residual blocks with probability subset_fraction.
  215. std::default_random_engine engine;
  216. std::uniform_real_distribution<> distribution(0, 1); // rage 0 - 1
  217. for (auto residual_block : residual_blocks) {
  218. if (distribution(engine) <= CERES_GET_FLAG(FLAGS_subset_fraction)) {
  219. options.residual_blocks_for_subset_preconditioner.insert(
  220. residual_block);
  221. }
  222. }
  223. }
  224. ceres::Solver::Summary summary;
  225. ceres::Solve(options, problem, &summary);
  226. std::cout << summary.FullReport() << "\n";
  227. // Make the solution stay in [0, 255].
  228. for (int x = 0; x < solution->width(); ++x) {
  229. for (int y = 0; y < solution->height(); ++y) {
  230. *solution->MutablePixel(x, y) =
  231. std::min(255.0, std::max(0.0, solution->Pixel(x, y)));
  232. }
  233. }
  234. }
  235. } // namespace
  236. } // namespace ceres::examples
  237. int main(int argc, char** argv) {
  238. using namespace ceres::examples;
  239. GFLAGS_NAMESPACE::ParseCommandLineFlags(&argc, &argv, true);
  240. google::InitGoogleLogging(argv[0]);
  241. if (CERES_GET_FLAG(FLAGS_input).empty()) {
  242. std::cerr << "Please provide an image file name using -input.\n";
  243. return 1;
  244. }
  245. if (CERES_GET_FLAG(FLAGS_foe_file).empty()) {
  246. std::cerr << "Please provide a Fields of Experts file name using -foe_file."
  247. "\n";
  248. return 1;
  249. }
  250. // Load the Fields of Experts filters from file.
  251. FieldsOfExperts foe;
  252. if (!foe.LoadFromFile(CERES_GET_FLAG(FLAGS_foe_file))) {
  253. std::cerr << "Loading \"" << CERES_GET_FLAG(FLAGS_foe_file)
  254. << "\" failed.\n";
  255. return 2;
  256. }
  257. // Read the images
  258. PGMImage<double> image(CERES_GET_FLAG(FLAGS_input));
  259. if (image.width() == 0) {
  260. std::cerr << "Reading \"" << CERES_GET_FLAG(FLAGS_input) << "\" failed.\n";
  261. return 3;
  262. }
  263. PGMImage<double> solution(image.width(), image.height());
  264. solution.Set(0.0);
  265. ceres::Problem problem;
  266. CreateProblem(foe, image, &problem, &solution);
  267. SolveProblem(&problem, &solution);
  268. if (!CERES_GET_FLAG(FLAGS_output).empty()) {
  269. CHECK(solution.WriteToFile(CERES_GET_FLAG(FLAGS_output)))
  270. << "Writing \"" << CERES_GET_FLAG(FLAGS_output) << "\" failed.";
  271. }
  272. return 0;
  273. }