GammaLib  2.0.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GOptimizerLM.cpp
Go to the documentation of this file.
1 /***************************************************************************
2  * GOptimizerLM.cpp - Levenberg Marquardt optimizer *
3  * ----------------------------------------------------------------------- *
4  * copyright (C) 2009-2022 by Juergen Knoedlseder *
5  * ----------------------------------------------------------------------- *
6  * *
7  * This program is free software: you can redistribute it and/or modify *
8  * it under the terms of the GNU General Public License as published by *
9  * the Free Software Foundation, either version 3 of the License, or *
10  * (at your option) any later version. *
11  * *
12  * This program is distributed in the hope that it will be useful, *
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15  * GNU General Public License for more details. *
16  * *
17  * You should have received a copy of the GNU General Public License *
18  * along with this program. If not, see <http://www.gnu.org/licenses/>. *
19  * *
20  ***************************************************************************/
21 /**
22  * @file GOptimizerLM.cpp
23  * @brief Levenberg-Marquardt optimizer class implementation
24  * @author Juergen Knoedlseder
25  */
26 
27 /* __ Includes ___________________________________________________________ */
28 #ifdef HAVE_CONFIG_H
29 #include <config.h>
30 #endif
31 #include "GOptimizerLM.hpp"
32 #include "GTools.hpp"
33 #include "GException.hpp"
34 
35 /* __ Method name definitions ____________________________________________ */
36 
37 /* __ Macros _____________________________________________________________ */
38 
39 /* __ Coding definitions _________________________________________________ */
40 
41 /* __ Debug definitions __________________________________________________ */
42 //#define G_DEBUG_OPT //!< Define to debug optimize() method
43 //#define G_DEBUG_ITER //!< Define to debug iteration() method
44 //#define G_DEBUG_SHOW_GRAD_COVAR //!< Define to show grad and curvature
45 
46 
47 /*==========================================================================
48  = =
49  = Constructors/destructors =
50  = =
51  ==========================================================================*/
52 
53 /***********************************************************************//**
54  * @brief Void constructor
55  ***************************************************************************/
57 {
58  // Initialise private members for clean destruction
59  init_members();
60 
61  // Return
62  return;
63 }
64 
65 
66 /***********************************************************************//**
67  * @brief Constructor with logger
68  *
69  * @param[in] log Logger to use in optimizer.
70  ***************************************************************************/
72 {
73  // Initialise private members for clean destruction
74  init_members();
75 
76  // Set pointer to logger
77  m_logger = log;
78 
79  // Return
80  return;
81 }
82 
83 
84 /***********************************************************************//**
85  * @brief Copy constructor
86  *
87  * @param[in] opt Optimizer from which the instance should be built.
88  ***************************************************************************/
90 {
91  // Initialise private members for clean destruction
92  init_members();
93 
94  // Copy members
95  copy_members(opt);
96 
97  // Return
98  return;
99 }
100 
101 
102 /***********************************************************************//**
103  * @brief Destructor
104  ***************************************************************************/
106 {
107  // Free members
108  free_members();
109 
110  // Return
111  return;
112 }
113 
114 
115 /*==========================================================================
116  = =
117  = Operators =
118  = =
119  ==========================================================================*/
120 
121 /***********************************************************************//**
122  * @brief Assignment operator
123  *
124  * @param[in] opt Optimizer to be assigned.
125  ***************************************************************************/
127 {
128  // Execute only if object is not identical
129  if (this != &opt) {
130 
131  // Copy base class members
132  this->GOptimizer::operator=(opt);
133 
134  // Free members
135  free_members();
136 
137  // Initialise private members for clean destruction
138  init_members();
139 
140  // Copy members
141  copy_members(opt);
142 
143  } // endif: object was not identical
144 
145  // Return
146  return *this;
147 }
148 
149 
150 /*==========================================================================
151  = =
152  = Public methods =
153  = =
154  ==========================================================================*/
155 
156 /***********************************************************************//**
157  * @brief Clear object
158  *
159  * This method properly resets the object to an initial state.
160  ***************************************************************************/
162 {
163  // Free class members (base and derived classes, derived class first)
164  free_members();
165  this->GOptimizer::free_members();
166 
167  // Initialise members
168  this->GOptimizer::init_members();
169  init_members();
170 
171  // Return
172  return;
173 }
174 
175 
176 /***********************************************************************//**
177  * @brief Clone object
178 ***************************************************************************/
180 {
181  return new GOptimizerLM(*this);
182 }
183 
184 
185 /***********************************************************************//**
186  * @brief Optimize function parameters
187  *
188  * @param[in] fct Optimization function.
189  * @param[in] pars Function parameters.
190  ***************************************************************************/
192 {
193  // Initialize optimizer parameters
195  m_value = 0.0;
196  m_delta = 0.0;
198  m_iter = 0;
199  m_num_dec = 0;
200 
201  // Get number of parameters. Continue only if there are free parameters
202  m_npars = pars.size();
203  m_nfree = pars.nfree();
204  if (m_nfree > 0) {
205 
206  // Initialise bookkeeping arrays
207  m_hit_boundary.clear();
208  m_hit_minimum.clear();
209  m_hit_maximum.clear();
210  m_par_freeze.clear();
211  m_par_remove.clear();
212  m_hit_boundary.reserve(m_npars);
213  m_hit_minimum.reserve(m_npars);
214  m_hit_maximum.reserve(m_npars);
215  m_par_freeze.reserve(m_npars);
216  m_par_remove.reserve(m_npars);
217  for (int i = 0; i < m_npars; ++i) {
218  m_hit_boundary.push_back(false);
219  m_hit_minimum.push_back(0);
220  m_hit_maximum.push_back(0);
221  m_par_freeze.push_back(false);
222  m_par_remove.push_back(false);
223  }
224 
225  // Initial function evaluation
226  fct.eval(pars);
227 
228  // If a free parameter has a zero diagonal element in the curvature
229  // matrix then remove this parameter definitely from the fit as it
230  // otherwise will block the fit. The problem appears in the unbinned
231  // fitting where parameter gradients may be zero (due to the truncation
232  // of the PSF), but the Npred gradient is not zero. In principle we
233  // could use the Npred gradient for fitting (I guess), but I still
234  // have to figure out how ... (the diagonal loading was not so
235  // successful as it faked early convergence)
236  for (int ipar = 0; ipar < m_npars; ++ipar) {
237  if (pars[ipar]->is_free()) {
238  if ((*fct.curvature())(ipar,ipar) == 0.0) {
239  if (m_logger != NULL) {
240  *m_logger << " Parameter \"" << pars[ipar]->name();
241  *m_logger << "\" has zero curvature.";
242  *m_logger << " Fix parameter." << std::endl;
243  }
244  m_par_remove[ipar] = true;
245  pars[ipar]->fix();
246  }
247  }
248  }
249 
250  // Save function value
251  m_value = fct.value();
252 
253  // Save initial statistics and lambda values
254  double lambda_old = m_lambda;
255  int lambda_inc = 0;
256 
257  // Optionally write initial iteration into logger
258  if (m_logger != NULL) {
259  (*m_logger)(">Iteration %3d: -logL=%.3f, Lambda=%.1e",
260  0, m_value, m_lambda);
261  }
262  #if defined(G_DEBUG_OPT)
263  std::cout << "Initial iteration: func=" << m_value << ", Lambda="
264  << m_lambda << std::endl;
265  #endif
266  #if defined(G_DEBUG_SHOW_GRAD_COVAR)
267  if (m_logger != NULL) {
268  *m_logger << *fct.gradient() << std::endl;
269  *m_logger << *fct.curvature() << std::endl;
270  }
271  #endif
272 
273  // Initialise iteration flags
274  bool check_for_freeze = true;
275 
276  // Iterative fitting
277  for (m_iter = 1; m_iter <= m_max_iter; ++m_iter) {
278 
279  // Perform one iteration
280  double step = iteration(fct, pars);
281 
282  // Determine maximum (scaled) gradient
283  double grad_max = 0.0;
284  int grad_imax = -1;
285  for (int ipar = 0; ipar < m_npars; ++ipar) {
286  if (pars[ipar]->is_free()) {
287  double grad = pars[ipar]->factor_gradient();
288  if (std::abs(grad) > std::abs(grad_max)) {
289  grad_max = grad;
290  grad_imax = ipar;
291  }
292  if (grad == 0.0) {
293  if (m_logger != NULL) {
294  *m_logger << "Parameter " << ipar;
295  *m_logger << " (" << pars[ipar]->name() << ")";
296  *m_logger << " has a zero gradient." << std::endl;
297  }
298  }
299  }
300  }
301 
302  // Optionally write iteration results into logger
303  if (m_logger != NULL) {
304  std::string stalled = "";
305  std::string status = "";
306  if (m_lambda > lambda_old) {
307  status = " ";
308  stalled = " (stalled)";
309  }
310  else {
311  status = ">";
312  }
313  std::string parname = "";
314  if (grad_imax != -1) {
315  parname = " [" + pars[grad_imax]->name() + ":" +
316  gammalib::str(grad_imax) + "]";
317  }
318  (*m_logger)("%sIteration %3d: -logL=%.3f, Lambda=%.1e,"
319  " delta=%.3f, step=%.1e, max(|grad|)=%f%s%s",
320  status.c_str(), m_iter, m_value, lambda_old,
321  m_delta, step, grad_max,
322  parname.c_str(), stalled.c_str());
323  }
324  #if defined(G_DEBUG_OPT)
325  std::cout << "Iteration " << m_iter << ": func="
326  << m_value << ", Lambda=" << lambda_old
327  << ", delta=" << m_delta << std::endl;
328  #endif
329  #if defined(G_DEBUG_SHOW_GRAD_COVAR)
330  if (m_logger != NULL) {
331  *m_logger << *fct.gradient() << std::endl;
332  *m_logger << *fct.curvature() << std::endl;
333  }
334  #endif
335 
336  // Reset lambda increment if we had success
337  if (m_lambda < lambda_old) {
338  lambda_inc = 0;
339  }
340 
341  // Stop if convergence was reached. Before stopping, check
342  // if some parameters were frozen, and if this was the case,
343  // free them now and continue. We do this only once, i.e.
344  // the next time a parameter is frozen and convergence is
345  // reached we really stop.
346  if ((m_lambda <= lambda_old) &&
347  (m_delta >= 0.0) &&
348  (m_delta < step*m_eps)) {
349 
350  // Check for frozen parameters, and if some exist, free
351  // them and start over
352  if (check_for_freeze) {
353 
354  // Signal that we won't check frozen parameters
355  // anymore
356  check_for_freeze = false;
357 
358  // Free frozen parameters and determine how many
359  // have been frozen
360  int nfrozen = 0;
361  for (int ipar = 0; ipar < m_npars; ++ipar) {
362  if (m_par_freeze[ipar]) {
363  nfrozen++;
364  pars[ipar]->free();
365  if (m_logger != NULL) {
366  *m_logger << " Free parameter \""
367  << pars[ipar]->name()
368  << "\" after convergence was"
369  " reached with frozen"
370  " parameter." << std::endl;
371  }
372  }
373  }
374 
375  // If there were frozen parameters then start over
376  // again (initialise optimizer)
377  if (nfrozen > 0) {
379  lambda_old = m_lambda;
380  lambda_inc = 0;
381  continue;
382  }
383  }
384 
385  // ... otherwise, convergence was reached and we can stop
386  // now
387  break;
388 
389  } // endif: convergence check
390 
391  // Monitor the number of subsequent increases of lambda and
392  // stop if the number of increases exceeds threshold
393  lambda_inc = (m_lambda > lambda_old) ? lambda_inc + 1 : 0;
394  if (lambda_inc >= m_max_stall) {
396  break;
397  }
398 
399  // Bookkeeping of actual result (we always store the last
400  // lambda to detect turn arounds in the lambda tendency; however
401  // we always keep the best function value)
402  lambda_old = m_lambda;
403 
404  } // endfor: iterations
405 
406  // Free now all temporarily frozen parameters so that the resulting
407  // model has the same attributes as the initial model
408  for (int ipar = 0; ipar < m_npars; ++ipar) {
409  if (m_par_freeze[ipar] || m_par_remove[ipar]) {
410  pars[ipar]->free();
411  if (m_logger != NULL) {
412  *m_logger << " Free parameter \""
413  << pars[ipar]->name()
414  << "\" after convergence was"
415  " reached with frozen"
416  " parameter." << std::endl;
417  }
418  }
419  }
420 
421  } // endif: there were free parameters to fit
422 
423  // ... otherwise just evaluate function and store the parameters
424  else {
425 
426  // Evaluate function
427  fct.eval(pars);
428 
429  // Save function value
430  m_value = fct.value();
431 
432  // Optionally write initial iteration into logger
433  if (m_logger != NULL) {
434  (*m_logger)(">Iteration %3d: -logL=%.3f, Lambda=%.1e (no free parameters)",
435  0, m_value, m_lambda);
436  }
437 
438  }
439 
440  // Return
441  return;
442 }
443 
444 
445 /***********************************************************************//**
446  * @brief Compute parameter uncertainties
447  *
448  * @param[in] fct Optimizer function.
449  * @param[in] pars Function parameters.
450  *
451  * Compute parameter uncertainties from the diagonal elements of the
452  * curvature matrix.
453  ***************************************************************************/
455 {
456  // Get number of parameters
457  int npars = pars.size();
458 
459  // Perform final parameter evaluation
460  fct.eval(pars);
461 
462  // Fetch sparse matrix pointer. We have to do this after the eval()
463  // method since eval() will allocate new memory for the curvature
464  // matrix!
465  GMatrixSparse* curvature = fct.curvature();
466 
467  // Save best fitting value
468  m_value = fct.value();
469 
470  // Save curvature matrix
471  GMatrixSparse save_curvature = GMatrixSparse(*curvature);
472 
473  // Signal no diagonal element loading
474  bool diag_loaded = false;
475 
476  // Loop over error computation (maximum 2 turns)
477  for (int i = 0; i < 2; ++i) {
478 
479  // Solve: curvature * X = unit
480  try {
481  GMatrixSparse decomposition = curvature->cholesky_decompose(true);
482  GVector unit(npars);
483  for (int ipar = 0; ipar < npars; ++ipar) {
484  unit[ipar] = 1.0;
485  GVector x = decomposition.cholesky_solver(unit, true);
486  if (x[ipar] >= 0.0) {
487  pars[ipar]->factor_error(sqrt(x[ipar]));
488  }
489  else {
490  pars[ipar]->factor_error(0.0);
492  }
493  unit[ipar] = 0.0;
494  }
495  }
496  catch (GException::runtime_error &e) {
497 
498  // Load diagonal if this has not yet been tried
499  if (!diag_loaded) {
500 
501  // Flag errors as inaccurate
503  if (m_logger != NULL) {
504  *m_logger << "Non-Positive definite curvature matrix encountered."
505  << std::endl;
506  *m_logger << "Load diagonal elements with 1e-10."
507  << " Fit errors may be inaccurate."
508  << std::endl;
509  }
510 
511  // Try now with diagonal loaded matrix
512  *curvature = save_curvature;
513  for (int ipar = 0; ipar < npars; ++ipar) {
514  (*curvature)(ipar,ipar) += 1.0e-10;
515  }
516 
517  // Signal loading
518  diag_loaded = true;
519 
520  // Try again
521  continue;
522 
523  } // endif: diagonal has not yet been loaded
524 
525  // ... otherwise signal an error
526  else {
528  if (m_logger != NULL) {
529  *m_logger << "Non-Positive definite curvature matrix encountered,"
530  << " even after diagonal loading." << std::endl;
531  }
532  break;
533  }
534  }
535  catch (std::exception &e) {
536  throw;
537  }
538 
539  // If no error occured then break now
540  break;
541 
542  } // endfor: looped over error computation
543 
544  // Return
545  return;
546 }
547 
548 
549 /***********************************************************************//**
550  * @brief Print optimizer information
551  *
552  * @param[in] chatter Chattiness.
553  * @return String containing optimizer information.
554  ***************************************************************************/
555 std::string GOptimizerLM::print(const GChatter& chatter) const
556 {
557  // Initialise result string
558  std::string result;
559 
560  // Continue only if chatter is not silent
561  if (chatter != SILENT) {
562 
563  // Append header
564  result.append("=== GOptimizerLM ===");
565 
566  // Append information
567  result.append("\n"+gammalib::parformat("Optimized function value"));
568  result.append(gammalib::str(m_value, 3));
569  result.append("\n"+gammalib::parformat("Absolute precision"));
570  result.append(gammalib::str(m_eps));
571  result.append("\n"+gammalib::parformat("Acceptable value decrease"));
572  result.append(gammalib::str(std::abs(m_accept_dec)));
573 
574  // Append status
575  result.append("\n"+gammalib::parformat("Optimization status"));
576  switch (m_status) {
577  case G_LM_CONVERGED:
578  result.append("converged");
579  break;
580  case G_LM_STALLED:
581  result.append("stalled");
582  break;
583  case G_LM_SINGULAR:
584  result.append("singular curvature matrix encountered");
585  break;
587  result.append("curvature matrix not positive definite");
588  break;
589  case G_LM_BAD_ERRORS:
590  result.append("errors are inaccurate");
591  break;
592  default:
593  result.append("unknown");
594  break;
595  }
596 
597  // Append further information
598  result.append("\n"+gammalib::parformat("Number of parameters"));
599  result.append(gammalib::str(m_npars));
600  result.append("\n"+gammalib::parformat("Number of free parameters"));
601  result.append(gammalib::str(m_nfree));
602  result.append("\n"+gammalib::parformat("Number of iterations"));
603  result.append(gammalib::str(m_iter));
604  result.append("\n"+gammalib::parformat("Lambda"));
605  result.append(gammalib::str(m_lambda));
606 
607  } // endif: chatter was not silent
608 
609  // Return result
610  return result;
611 }
612 
613 
614 /*==========================================================================
615  = =
616  = Private methods =
617  = =
618  ==========================================================================*/
619 
620 /***********************************************************************//**
621  * @brief Initialise class members
622  ***************************************************************************/
624 {
625  // Initialise optimizer parameters
626  m_npars = 0;
627  m_nfree = 0;
628  m_lambda_start = 1.0e-3;
629  m_lambda_inc = 10.0;
630  m_lambda_dec = 0.1;
631  m_eps = 5.0e-3; //!< Changed on 30/10/2014 from 1.0e-6
632  m_accept_dec = 0.0; //!< Do not allow to decrease
633  m_max_iter = 100; //!< Changed on 30/10/2014 from 1000
634  m_max_stall = 10;
635  m_max_hit = 3; //!< Maximum successive boundary hits before freeze
636  m_max_dec = 3; //!< Maximum number of function decreases
637  m_step_adjust = true;
638 
639  // Initialise bookkeeping arrays
640  m_hit_boundary.clear();
641  m_hit_minimum.clear();
642  m_hit_maximum.clear();
643  m_par_freeze.clear();
644  m_par_remove.clear();
645 
646  // Initialise optimizer values
648  m_value = 0.0;
649  m_delta = 0.0;
650  m_status = 0;
651  m_iter = 0;
652  m_num_dec = 0;
653 
654  // Initialise pointer to logger
655  m_logger = NULL;
656 
657  // Return
658  return;
659 }
660 
661 
662 /***********************************************************************//**
663  * @brief Copy class members
664  *
665  * @param[in] opt GOptimizerLM members to be copied.
666  ***************************************************************************/
668 {
669  // Copy attributes
670  m_npars = opt.m_npars;
671  m_nfree = opt.m_nfree;
675  m_eps = opt.m_eps;
677  m_max_iter = opt.m_max_iter;
678  m_max_stall = opt.m_max_stall;
679  m_max_dec = opt.m_max_dec;
686  m_lambda = opt.m_lambda;
687  m_value = opt.m_value;
688  m_delta = opt.m_delta;
689  m_status = opt.m_status;
690  m_iter = opt.m_iter;
691  m_num_dec = opt.m_num_dec;
692  m_logger = opt.m_logger;
693 
694  // Return
695  return;
696 }
697 
698 
699 /***********************************************************************//**
700  * @brief Delete class members
701  ***************************************************************************/
703 {
704  // Return
705  return;
706 }
707 
708 
709 /***********************************************************************//**
710  * @brief Perform one LM iteration
711  *
712  * @param[in] fct Optimizer function.
713  * @param[in] pars Function parameters.
714  * @return Step size for iteration.
715  *
716  * This method performs one LM iteration. Note that the method only acts on
717  * the parameter value, i.e. it does not worry about the true scale of the
718  * parameter. It calls the eval() method of the optimizer function which
719  * is assumed to return the gradient and curvature matrix with respect to the
720  * parameter values (and not the scaled true values).
721  ***************************************************************************/
723 {
724  // Debug option: dump header
725  #if defined(G_DEBUG_ITER)
726  std::cout << "GOptimizerLM::iteration: enter" << std::endl;
727  #endif
728 
729  // Initialise step size
730  double step = 0.0;
731 
732  // Single loop for common exit point
733  do {
734 
735  // Initialise iteration parameters
736  GVector* grad = fct.gradient();
737  GMatrixSparse* curvature = fct.curvature();
738 
739  // Save function value, gradient and curvature matrix
740  double save_value = m_value;
741  GVector save_grad = GVector(*grad);
742  GMatrixSparse save_curvature = GMatrixSparse(*curvature);
743 
744  // Save parameter values in vector
745  GVector save_pars(m_npars);
746  for (int ipar = 0; ipar < m_npars; ++ipar) {
747  save_pars[ipar] = pars[ipar]->factor_value();
748  }
749 
750  // Setup matrix and vector for curvature computation
751  for (int ipar = 0; ipar < m_npars; ++ipar) {
752  (*curvature)(ipar,ipar) *= (1.0 + m_lambda);
753  (*grad)[ipar] = -(*grad)[ipar];
754  }
755 
756  // Debug option: dump gradient and curvature matrix
757  #if defined(G_DEBUG_ITER)
758  std::cout << "Gradient : ";
759  for (int ipar = 0; ipar < m_npars; ++ipar) {
760  std::cout << (*grad)[ipar] << " ";
761  }
762  std::cout << std::endl;
763  std::cout << "Curvature: ";
764  for (int ipar = 0; ipar < m_npars; ++ipar) {
765  for (int jpar = 0; jpar < m_npars; ++jpar) {
766  std::cout << (*curvature)(ipar,jpar) << " ";
767  }
768  }
769  std::cout << std::endl;
770  #endif
771 
772  // Solve: curvature * X = grad. Handle matrix problems
773  try {
774  *grad = curvature->solve(*grad);
775  }
776  catch (GException::runtime_error &e) {
778  if (m_logger != NULL) {
779  *m_logger << "GOptimizerLM::iteration: "
780  << "Curvature matrix not positive definite."
781  << std::endl;
782  }
783  continue;
784  }
785  catch (std::exception &e) {
786  throw;
787  }
788 
789  // Get LM step size
790  step = step_size(*grad, pars);
791 
792  // Debug option: dump step size
793  #if defined(G_DEBUG_ITER)
794  std::cout << "Step size = " << step << std::endl;
795  #endif
796 
797  // Derive new parameter vector
798  for (int ipar = 0; ipar < m_npars; ++ipar) {
799 
800  // Consider only free parameters
801  if (pars[ipar]->is_free()) {
802 
803  // Get actual parameter value and limits
804  double p = pars[ipar]->factor_value();
805  double p_min = pars[ipar]->factor_min();
806  double p_max = pars[ipar]->factor_max();
807 
808  // Compute new parameter value
809  p += (*grad)[ipar] * step;
810 
811  // Debug option: dump new parameter value
812  #if defined(G_DEBUG_ITER)
813  std::cout << "Trial factor value of ";
814  std::cout << pars[ipar]->name() << " = ";
815  std::cout << p << std::endl;
816  #endif
817 
818  // Constrain parameter to within the valid range
819  if (pars[ipar]->has_min() && p < p_min) {
820  if (m_hit_minimum[ipar] >= m_max_hit) {
821  if (m_logger != NULL) {
822  *m_logger << " Parameter \"" << pars[ipar]->name();
823  *m_logger << "\" hits minimum " << p_min << " more than ";
824  *m_logger << m_max_hit << " times.";
825  *m_logger << " Fix parameter at minimum for now." << std::endl;
826  }
827  m_par_freeze[ipar] = true;
828  pars[ipar]->fix();
829  }
830  else {
831  if (m_logger != NULL) {
832  *m_logger << " Parameter \"" << pars[ipar]->name();
833  *m_logger << "\" hits minimum: ";
834  *m_logger << p << " < " << p_min;
835  *m_logger << " (" << m_hit_minimum[ipar]+1 << ")" << std::endl;
836  }
837  }
838  m_hit_minimum[ipar]++;
839  p = p_min;
840  }
841  else if (pars[ipar]->has_max() && p > p_max) {
842  if (m_hit_maximum[ipar] >= m_max_hit) {
843  if (m_logger != NULL) {
844  *m_logger << " Parameter \"" << pars[ipar]->name();
845  *m_logger << "\" hits maximum " << p_max << " more than ";
846  *m_logger << m_max_hit << " times.";
847  *m_logger << " Fix parameter at maximum for now." << std::endl;
848  }
849  m_par_freeze[ipar] = true;
850  pars[ipar]->fix();
851  }
852  else {
853  if (m_logger != NULL) {
854  *m_logger << " Parameter \"" << pars[ipar]->name();
855  *m_logger << "\" hits maximum: ";
856  *m_logger << p << " > " << p_max;
857  *m_logger << " (" << m_hit_maximum[ipar]+1 << ")" << std::endl;
858  }
859  }
860  m_hit_maximum[ipar]++;
861  p = p_max;
862  }
863  else {
864  m_hit_minimum[ipar] = 0;
865  m_hit_maximum[ipar] = 0;
866  }
867 
868  // Set new parameter value
869  pars[ipar]->factor_value(p);
870 
871  // Debug option: dump new parameter value
872  #if defined(G_DEBUG_ITER)
873  std::cout << "New value of ";
874  std::cout << pars[ipar]->name() << " = ";
875  std::cout << pars[ipar]->value() << std::endl;
876  #endif
877 
878  } // endif: Parameter was free
879 
880  } // endfor: computed new parameter vector
881 
882  // Evaluate function at new parameters
883  fct.eval(pars);
884 
885  // If a free parameter has a zero diagonal element in the curvature
886  // matrix then remove this parameter definitely from the fit as it
887  // otherwise will block the fit. The problem appears in the unbinned
888  // fitting where parameter gradients may be zero (due to the truncation
889  // of the PSF), but the Npred gradient is not zero. In principle we
890  // could use the Npred gradient for fitting (I guess), but I still
891  // have to figure out how ... (the diagonal loading was not so
892  // successful as it faked early convergence)
893  for (int ipar = 0; ipar < m_npars; ++ipar) {
894  if (pars[ipar]->is_free()) {
895  if ((*fct.curvature())(ipar,ipar) == 0.0) {
896  if (m_logger != NULL) {
897  *m_logger << " Parameter \"" << pars[ipar]->name();
898  *m_logger << "\" has zero curvature.";
899  *m_logger << " Fix parameter." << std::endl;
900  }
901  m_par_remove[ipar] = true;
902  pars[ipar]->fix();
903  }
904  }
905  }
906 
907  // Fetch new pointers since eval will allocate new memory
908  grad = fct.gradient();
909  curvature = fct.curvature();
910 
911  // Retrieve new function value
912  m_value = fct.value();
913 
914  // Debug option: dump new function value
915  #if defined(G_DEBUG_ITER)
916  std::cout << "Function value " << m_value;
917  std::cout << " (was " << save_value << ")" << std::endl;
918  #endif
919 
920  // Determine how many parameters have changed
921  int par_change = 0;
922  for (int ipar = 0; ipar < m_npars; ++ipar) {
923  if (pars[ipar]->factor_value() != save_pars[ipar]) {
924  par_change++;
925  }
926  }
927 
928  // If the function has decreased then accept the new solution
929  // and decrease lambda ...
930  m_delta = save_value - m_value;
931  if (m_delta > 0.0) {
933  }
934 
935  // ... if function is identical then accept new solution. If the
936  // parameters have changed then increase lambda, otherwise decrease
937  // lambda
938  else if (m_delta == 0.0) {
939  if (par_change) {
941  }
942  else {
944  }
945  }
946 
947  // ... if function worsened by an amount smaller than m_accept_dec
948  // then accept new solution and increase lambda for the next
949  // iteration. This kluge which may help to overcome some local
950  // minimum will only for m_max_dec times at maximum to avoid
951  // keeping oscillating in a loop.
952  else if ((m_num_dec < m_max_dec) && (m_delta >= -std::abs(m_accept_dec))) {
954  m_num_dec++;
955  }
956 
957  // ... otherwise, if the statistics did not improve then use old
958  // parameters and increase lamdba. Restore also the best statistics
959  // value that was reached so far, the gradient vector and the curve
960  // matrix.
961  else {
963  m_value = save_value;
964  *grad = save_grad;
965  *curvature = save_curvature;
966  for (int ipar = 0; ipar < m_npars; ++ipar) {
967  pars[ipar]->factor_value(save_pars[ipar]);
968  }
969  }
970 
971  } while (0); // endwhile: main loop
972 
973  // Debug option: dump trailer
974  #if defined(G_DEBUG_ITER)
975  std::cout << "GOptimizerLM::iteration: exit" << std::endl;
976  #endif
977 
978  // Return step size
979  return step;
980 
981 }
982 
983 
984 /***********************************************************************//**
985  * @brief Return LM step size
986  *
987  * @param[in] grad Function gradient.
988  * @param[in] pars Function parameters.
989  *
990  * Determine the size of the LM step. By default a step size of 1 is taken.
991  * If m_step_adjust=true then the step size will be estimated so that the
992  * next parameter vector should stay within the parameter boundaries
993  * (provided that boundaries exist).
994  ***************************************************************************/
995 double GOptimizerLM::step_size(const GVector& grad, const GOptimizerPars& pars)
996 {
997  // Initialise step size
998  double step = 1.0;
999 
1000  // Check if we should reduce the step size
1001  if (m_step_adjust) {
1002 
1003  // Initialise the parameter index that constrains most the fit
1004  int ipar_bnd = -1;
1005 
1006  // Loop over all parameters
1007  for (int ipar = 0; ipar < pars.size(); ++ipar) {
1008 
1009  // Get parameter attributes
1010  double p = pars[ipar]->factor_value();
1011  double p_min = pars[ipar]->factor_min();
1012  double p_max = pars[ipar]->factor_max();
1013  double delta = grad[ipar];
1014 
1015  // Check if a parameter minimum requires a reduced step size
1016  if (pars[ipar]->has_min()) {
1017  double step_min = (delta < 0.0) ? (p_min - p)/delta : 1.0;
1018  if (step_min > 0.0) {
1019  if (step_min < step) {
1020  if (!m_hit_boundary[ipar]) {
1021  ipar_bnd = ipar;
1022  step = step_min;
1023  }
1024  }
1025  else if (m_hit_boundary[ipar]) {
1026  m_hit_boundary[ipar] = false;
1027  if (m_logger != NULL) {
1028  *m_logger << " Parameter \"";
1029  *m_logger << pars[ipar]->name();
1030  *m_logger << "\" does not drive optimization step anymore.";
1031  *m_logger << std::endl;
1032  }
1033  }
1034  }
1035  }
1036 
1037  // Check if a parameter maximum requires a reduced step size
1038  if (pars[ipar]->has_max()) {
1039  double step_max = (delta > 0.0) ? (p_max - p)/delta : 1.0;
1040  if (step_max > 0.0) {
1041  if (step_max < step) {
1042  if (!m_hit_boundary[ipar]) {
1043  ipar_bnd = ipar;
1044  step = step_max;
1045  }
1046  }
1047  else if (m_hit_boundary[ipar]) {
1048  m_hit_boundary[ipar] = false;
1049  if (m_logger != NULL) {
1050  *m_logger << " Parameter \"";
1051  *m_logger << pars[ipar]->name();
1052  *m_logger << "\" does not drive optimization step anymore.";
1053  *m_logger << std::endl;
1054  }
1055  }
1056  }
1057  }
1058 
1059  } // endfor: looped over all parameters
1060 
1061  // Signal if a parameter is driving the optimization step
1062  if (ipar_bnd != -1) {
1063  m_hit_boundary[ipar_bnd] = true;
1064  if (m_logger != NULL) {
1065  *m_logger << " Parameter \"";
1066  *m_logger << pars[ipar_bnd]->name();
1067  *m_logger << "\" drives optimization step (step=";
1068  *m_logger << step << ")" << std::endl;
1069  }
1070  }
1071 
1072  } // endif: automatic step size adjustment requested
1073 
1074  // Return step size
1075  return step;
1076 }
#define G_LM_SINGULAR
GOptimizerLM & operator=(const GOptimizerLM &opt)
Assignment operator.
const double & lambda_inc(void) const
Return lambda increment value.
double step_size(const GVector &grad, const GOptimizerPars &pars)
Return LM step size.
Optimizer function abstract base class.
void init_members(void)
Initialise class members.
Sparse matrix class interface definition.
#define G_LM_STALLED
GVector abs(const GVector &vector)
Computes absolute of vector elements.
Definition: GVector.cpp:1253
virtual GMatrixSparse * curvature(void)=0
Optimizer parameter container class.
virtual void optimize(GOptimizerFunction &fct, GOptimizerPars &pars)
Optimize function parameters.
double m_lambda_start
Initial start value.
double m_lambda
Actual lambda.
void copy_members(const GOptimizerLM &opt)
Copy class members.
int m_iter
Iteration.
double m_eps
Absolute precision.
std::vector< bool > m_hit_boundary
Bookkeeping array for boundary hits.
virtual void eval(const GOptimizerPars &pars)=0
int m_npars
Number of parameters.
Gammalib tools definition.
#define G_LM_CONVERGED
GOptimizerLM(void)
Void constructor.
virtual GOptimizer & operator=(const GOptimizer &opt)
Assignment operator.
Definition: GOptimizer.cpp:101
std::vector< bool > m_par_remove
Bookkeeping of parameter removal.
int m_max_dec
Maximum number of function decrease.
#define G_LM_BAD_ERRORS
double m_value
Actual function value.
GLog * m_logger
Pointer to optional logger.
Information logger interface definition.
Definition: GLog.hpp:62
virtual GOptimizerLM * clone(void) const
Clone object.
int m_nfree
Number of free parameters.
GVector sqrt(const GVector &vector)
Computes square root of vector elements.
Definition: GVector.cpp:1358
GVector solve(const GVector &vector) const
Solves linear matrix equation.
bool m_step_adjust
Adjust step size to boundaries.
int m_status
Fit status.
virtual void clear(void)
Clear object.
double m_delta
Function improvement.
double m_accept_dec
Acceptable function decrease.
GVector log(const GVector &vector)
Computes natural logarithm of vector elements.
Definition: GVector.cpp:1274
#define G_LM_NOT_POSTIVE_DEFINITE
virtual ~GOptimizerLM(void)
Destructor.
double m_lambda_dec
Lambda decrease.
Abstract optimizer abstract base class.
Definition: GOptimizer.hpp:54
std::vector< bool > m_par_freeze
Bookkeeping of parameter freeze.
GChatter
Definition: GTypemaps.hpp:33
double m_lambda_inc
Lambda increase.
Levenberg Marquardt optimizer class interface definition.
int m_num_dec
Number of function decreases.
void free_members(void)
Delete class members.
const int & npars(void) const
Return number of model parameters.
int m_max_iter
Maximum number of iterations.
int m_max_hit
Maximum number of successive hits.
std::vector< int > m_hit_maximum
Bookkeeping of successive maximum hits.
void init_members(void)
Initialise class members.
Definition: GOptimizer.cpp:137
virtual GVector * gradient(void)=0
GVector cholesky_solver(const GVector &vector, const bool &compress=true) const
Cholesky solver.
int nfree(void) const
Return number of free parameters.
virtual double value(void) const =0
virtual void errors(GOptimizerFunction &fct, GOptimizerPars &pars)
Compute parameter uncertainties.
virtual std::string print(const GChatter &chatter=NORMAL) const
Print optimizer information.
Exception handler interface definition.
Vector class.
Definition: GVector.hpp:46
int size(void) const
Return number of parameters in container.
std::string parformat(const std::string &s, const int &indent=0)
Convert string in parameter format.
Definition: GTools.cpp:1143
std::vector< int > m_hit_minimum
Bookkeeping of successive minimum hits.
GMatrixSparse cholesky_decompose(const bool &compress=true) const
Return Cholesky decomposition.
void free_members(void)
Delete class members.
Definition: GOptimizer.cpp:159
Levenberg Marquardt optimizer class.
int m_max_stall
Maximum number of stalls.
virtual int status(void) const
Return optimizer status.
double iteration(GOptimizerFunction &fct, GOptimizerPars &pars)
Perform one LM iteration.
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
Definition: GTools.cpp:489