GammaLib 2.1.0.dev
Loading...
Searching...
No Matches
GOptimizerLM.cpp
Go to the documentation of this file.
1/***************************************************************************
2 * GOptimizerLM.cpp - Levenberg Marquardt optimizer *
3 * ----------------------------------------------------------------------- *
4 * copyright (C) 2009-2025 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
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
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
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();
166
167 // Initialise 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
231 for (int ipar = 0; ipar < m_npars; ++ipar) {
232 if (pars[ipar]->is_free()) {
233 if ((*fct.curvature())(ipar,ipar) == 0.0) {
234 if (m_logger != NULL) {
235 *m_logger << " Parameter \"" << pars[ipar]->name();
236 *m_logger << "\" = " << pars[ipar]->value();
237 *m_logger << " with gradient " << pars[ipar]->gradient();
238 *m_logger << " has zero curvature.";
239 *m_logger << " Fix parameter." << std::endl;
240 }
241 m_par_remove[ipar] = true;
242 pars[ipar]->fix();
243 }
244 }
245 }
246
247 // Save function value
248 m_value = fct.value();
249
250 // Save initial statistics and lambda values
251 double lambda_old = m_lambda;
252 int lambda_inc = 0;
253
254 // Optionally write initial iteration into logger
255 if (m_logger != NULL) {
256 (*m_logger)(">Iteration %3d: -logL=%.3f, Lambda=%.1e",
257 0, m_value, m_lambda);
258 }
259 #if defined(G_DEBUG_OPT)
260 std::cout << "Initial iteration: func=" << m_value << ", Lambda="
261 << m_lambda << std::endl;
262 #endif
263 #if defined(G_DEBUG_SHOW_GRAD_COVAR)
264 if (m_logger != NULL) {
265 *m_logger << *fct.gradient() << std::endl;
266 *m_logger << *fct.curvature() << std::endl;
267 }
268 #endif
269
270 // Initialise iteration flags
271 bool check_for_freeze = true;
272
273 // Iterative fitting
274 for (m_iter = 1; m_iter <= m_max_iter; ++m_iter) {
275
276 // Perform one iteration
277 double step = iteration(fct, pars);
278
279 // Determine maximum (scaled) gradient
280 double grad_max = 0.0;
281 int grad_imax = -1;
282 for (int ipar = 0; ipar < m_npars; ++ipar) {
283 if (pars[ipar]->is_free()) {
284 double grad = pars[ipar]->factor_gradient();
285 if (std::abs(grad) > std::abs(grad_max)) {
286 grad_max = grad;
287 grad_imax = ipar;
288 }
289 /*
290 if (grad == 0.0) {
291 if (m_logger != NULL) {
292 *m_logger << "Parameter " << ipar;
293 *m_logger << " (" << pars[ipar]->name() << ")";
294 *m_logger << " has a zero gradient." << std::endl;
295 }
296 }
297 */
298 }
299 }
300
301 // Optionally write iteration results into logger
302 if (m_logger != NULL) {
303 std::string stalled = "";
304 std::string status = "";
305 if (m_lambda > lambda_old) {
306 status = " ";
307 stalled = " (stalled)";
308 }
309 else {
310 status = ">";
311 }
312 std::string parname = "";
313 if (grad_imax != -1) {
314 parname = " [" + pars[grad_imax]->name() + ":" +
315 gammalib::str(grad_imax) + "]";
316 }
317 (*m_logger)("%sIteration %3d: -logL=%.3f, Lambda=%.1e,"
318 " delta=%.3f, step=%.1e, max(|grad|)=%f%s%s",
319 status.c_str(), m_iter, m_value, lambda_old,
320 m_delta, step, grad_max,
321 parname.c_str(), stalled.c_str());
322 }
323 #if defined(G_DEBUG_OPT)
324 std::cout << "Iteration " << m_iter << ": func="
325 << m_value << ", Lambda=" << lambda_old
326 << ", delta=" << m_delta << std::endl;
327 #endif
328 #if defined(G_DEBUG_SHOW_GRAD_COVAR)
329 if (m_logger != NULL) {
330 *m_logger << *fct.gradient() << std::endl;
331 *m_logger << *fct.curvature() << std::endl;
332 }
333 #endif
334
335 // Reset lambda increment if we had success
336 if (m_lambda < lambda_old) {
337 lambda_inc = 0;
338 }
339
340 // Stop if convergence was reached. Before stopping, check
341 // if some parameters were frozen, and if this was the case,
342 // free them now and continue. We do this only once, i.e.
343 // the next time a parameter is frozen and convergence is
344 // reached we really stop.
345 if ((m_lambda <= lambda_old) &&
346 (m_delta >= 0.0) &&
347 (m_delta < step*m_eps)) {
348
349 // Check for frozen parameters, and if some exist, free
350 // them and start over
351 if (check_for_freeze) {
352
353 // Signal that we won't check frozen parameters
354 // anymore
355 check_for_freeze = false;
356
357 // Free frozen parameters and determine how many
358 // have been frozen
359 int nfrozen = 0;
360 for (int ipar = 0; ipar < m_npars; ++ipar) {
361 if (m_par_freeze[ipar]) {
362 nfrozen++;
363 pars[ipar]->free();
364 if (m_logger != NULL) {
365 *m_logger << " Free parameter \""
366 << pars[ipar]->name()
367 << "\" after convergence was"
368 " reached with frozen"
369 " parameter." << std::endl;
370 }
371 }
372 }
373
374 // If there were frozen parameters then start over
375 // again (initialise optimizer)
376 if (nfrozen > 0) {
378 lambda_old = m_lambda;
379 lambda_inc = 0;
380 continue;
381 }
382 }
383
384 // ... otherwise, convergence was reached and we can stop
385 // now
386 break;
387
388 } // endif: convergence check
389
390 // Monitor the number of subsequent increases of lambda and
391 // stop if the number of increases exceeds threshold
392 lambda_inc = (m_lambda > lambda_old) ? lambda_inc + 1 : 0;
393 if (lambda_inc >= m_max_stall) {
395 break;
396 }
397
398 // Bookkeeping of actual result (we always store the last
399 // lambda to detect turn arounds in the lambda tendency; however
400 // we always keep the best function value)
401 lambda_old = m_lambda;
402
403 } // endfor: iterations
404
405 // Free now all temporarily frozen parameters so that the resulting
406 // model has the same attributes as the initial model
407 for (int ipar = 0; ipar < m_npars; ++ipar) {
408 if (m_par_freeze[ipar] || m_par_remove[ipar]) {
409 pars[ipar]->free();
410 if (m_logger != NULL) {
411 *m_logger << " Free parameter \""
412 << pars[ipar]->name()
413 << "\" after convergence was"
414 " reached with frozen"
415 " parameter." << std::endl;
416 }
417 }
418 }
419
420 } // endif: there were free parameters to fit
421
422 // ... otherwise just evaluate function and store the parameters
423 else {
424
425 // Evaluate function
426 fct.eval(pars);
427
428 // Save function value
429 m_value = fct.value();
430
431 // Optionally write initial iteration into logger
432 if (m_logger != NULL) {
433 (*m_logger)(">Iteration %3d: -logL=%.3f, Lambda=%.1e (no free parameters)",
434 0, m_value, m_lambda);
435 }
436
437 }
438
439 // Return
440 return;
441}
442
443
444/***********************************************************************//**
445 * @brief Compute parameter uncertainties
446 *
447 * @param[in] fct Optimizer function.
448 * @param[in] pars Function parameters.
449 *
450 * Compute parameter uncertainties from the diagonal elements of the
451 * curvature matrix.
452 ***************************************************************************/
454{
455 // Get number of parameters
456 int npars = pars.size();
457
458 // Perform final parameter evaluation
459 fct.eval(pars);
460
461 // Fetch sparse matrix pointer. We have to do this after the eval()
462 // method since eval() will allocate new memory for the curvature
463 // matrix!
464 GMatrixSparse* curvature = fct.curvature();
465
466 // Save best fitting value
467 m_value = fct.value();
468
469 // Save curvature matrix
470 GMatrixSparse save_curvature = GMatrixSparse(*curvature);
471
472 // Signal no diagonal element loading
473 bool diag_loaded = false;
474
475 // Loop over error computation (maximum 2 turns)
476 for (int i = 0; i < 2; ++i) {
477
478 // Solve: curvature * X = unit
479 try {
480 GMatrixSparse decomposition = curvature->cholesky_decompose(true);
481 GVector unit(npars);
482 for (int ipar = 0; ipar < npars; ++ipar) {
483 unit[ipar] = 1.0;
484 GVector x = decomposition.cholesky_solver(unit, true);
485 if (x[ipar] >= 0.0) {
486 pars[ipar]->factor_error(sqrt(x[ipar]));
487 }
488 else {
489 pars[ipar]->factor_error(0.0);
491 }
492 unit[ipar] = 0.0;
493 }
494 }
495 catch (GException::runtime_error &e) {
496
497 // Load diagonal if this has not yet been tried
498 if (!diag_loaded) {
499
500 // Flag errors as inaccurate
502 if (m_logger != NULL) {
503 *m_logger << "Non-Positive definite curvature matrix encountered."
504 << std::endl;
505 *m_logger << "Load diagonal elements with 1e-10."
506 << " Fit errors may be inaccurate."
507 << std::endl;
508 }
509
510 // Try now with diagonal loaded matrix
511 *curvature = save_curvature;
512 for (int ipar = 0; ipar < npars; ++ipar) {
513 if ((*curvature)(ipar,ipar) != 0.0) {
514 (*curvature)(ipar,ipar) += 1.0e-10;
515 }
516 }
517
518 // Signal loading
519 diag_loaded = true;
520
521 // Try again
522 continue;
523
524 } // endif: diagonal has not yet been loaded
525
526 // ... otherwise signal an error
527 else {
529 if (m_logger != NULL) {
530 *m_logger << "Non-Positive definite curvature matrix encountered,"
531 << " even after diagonal loading." << std::endl;
532 }
533 break;
534 }
535 }
536 catch (std::exception &e) {
537 throw;
538 }
539
540 // If no error occured then break now
541 break;
542
543 } // endfor: looped over error computation
544
545 // Return
546 return;
547}
548
549
550/***********************************************************************//**
551 * @brief Set fit status string
552 *
553 * @return Fit status string.
554 ***************************************************************************/
555std::string GOptimizerLM::status_string(void) const
556{
557 // Initialise status string
558 std::string status;
559
560 // Set status string
561 switch (m_status) {
562 case G_LM_CONVERGED:
563 status = "converged";
564 break;
565 case G_LM_STALLED:
566 status = "stalled";
567 break;
568 case G_LM_SINGULAR:
569 status = "singular curvature matrix encountered";
570 break;
572 status = "curvature matrix not positive definite";
573 break;
574 case G_LM_BAD_ERRORS:
575 status = "errors are inaccurate";
576 break;
577 default:
578 status = "unknown";
579 break;
580 }
581
582 // Return status
583 return status;
584}
585
586
587/***********************************************************************//**
588 * @brief Print optimizer information
589 *
590 * @param[in] chatter Chattiness.
591 * @return String containing optimizer information.
592 ***************************************************************************/
593std::string GOptimizerLM::print(const GChatter& chatter) const
594{
595 // Initialise result string
596 std::string result;
597
598 // Continue only if chatter is not silent
599 if (chatter != SILENT) {
600
601 // Append header
602 result.append("=== GOptimizerLM ===");
603
604 // Append information
605 result.append("\n"+gammalib::parformat("Optimized function value"));
606 result.append(gammalib::str(m_value, 3));
607 result.append("\n"+gammalib::parformat("Absolute precision"));
608 result.append(gammalib::str(m_eps));
609 result.append("\n"+gammalib::parformat("Acceptable value decrease"));
610 result.append(gammalib::str(std::abs(m_accept_dec)));
611
612 // Append status
613 result.append("\n"+gammalib::parformat("Optimization status"));
614 result.append(status_string());
615
616 // Append further information
617 result.append("\n"+gammalib::parformat("Number of parameters"));
618 result.append(gammalib::str(m_npars));
619 result.append("\n"+gammalib::parformat("Number of free parameters"));
620 result.append(gammalib::str(m_nfree));
621 result.append("\n"+gammalib::parformat("Number of iterations"));
622 result.append(gammalib::str(m_iter));
623 result.append("\n"+gammalib::parformat("Lambda"));
624 result.append(gammalib::str(m_lambda));
625
626 } // endif: chatter was not silent
627
628 // Return result
629 return result;
630}
631
632
633/*==========================================================================
634 = =
635 = Private methods =
636 = =
637 ==========================================================================*/
638
639/***********************************************************************//**
640 * @brief Initialise class members
641 ***************************************************************************/
643{
644 // Initialise optimizer parameters
645 m_npars = 0;
646 m_nfree = 0;
647 m_lambda_start = 1.0e-3;
648 m_lambda_inc = 10.0;
649 m_lambda_dec = 0.1;
650 m_eps = 5.0e-3; //!< Changed on 30/10/2014 from 1.0e-6
651 m_accept_dec = 0.0; //!< Do not allow to decrease
652 m_max_iter = 100; //!< Changed on 30/10/2014 from 1000
653 m_max_stall = 10;
654 m_max_hit = 3; //!< Maximum successive boundary hits before freeze
655 m_max_dec = 3; //!< Maximum number of function decreases
656 m_step_adjust = true;
657
658 // Initialise bookkeeping arrays
659 m_hit_boundary.clear();
660 m_hit_minimum.clear();
661 m_hit_maximum.clear();
662 m_par_freeze.clear();
663 m_par_remove.clear();
664
665 // Initialise optimizer values
667 m_value = 0.0;
668 m_delta = 0.0;
669 m_status = 0;
670 m_iter = 0;
671 m_num_dec = 0;
672
673 // Initialise pointer to logger
674 m_logger = NULL;
675
676 // Return
677 return;
678}
679
680
681/***********************************************************************//**
682 * @brief Copy class members
683 *
684 * @param[in] opt GOptimizerLM members to be copied.
685 ***************************************************************************/
687{
688 // Copy attributes
689 m_npars = opt.m_npars;
690 m_nfree = opt.m_nfree;
694 m_eps = opt.m_eps;
698 m_max_dec = opt.m_max_dec;
705 m_lambda = opt.m_lambda;
706 m_value = opt.m_value;
707 m_delta = opt.m_delta;
708 m_status = opt.m_status;
709 m_iter = opt.m_iter;
710 m_num_dec = opt.m_num_dec;
711 m_logger = opt.m_logger;
712
713 // Return
714 return;
715}
716
717
718/***********************************************************************//**
719 * @brief Delete class members
720 ***************************************************************************/
722{
723 // Return
724 return;
725}
726
727
728/***********************************************************************//**
729 * @brief Perform one LM iteration
730 *
731 * @param[in] fct Optimizer function.
732 * @param[in] pars Function parameters.
733 * @return Step size for iteration.
734 *
735 * This method performs one LM iteration. Note that the method only acts on
736 * the parameter value, i.e. it does not worry about the true scale of the
737 * parameter. It calls the eval() method of the optimizer function which
738 * is assumed to return the gradient and curvature matrix with respect to the
739 * parameter values (and not the scaled true values).
740 ***************************************************************************/
742{
743 // Debug option: dump header
744 #if defined(G_DEBUG_ITER)
745 std::cout << "GOptimizerLM::iteration: enter" << std::endl;
746 #endif
747
748 // Initialise step size
749 double step = 0.0;
750
751 // Single loop for common exit point
752 do {
753
754 // Initialise iteration parameters
755 GVector* grad = fct.gradient();
756 GMatrixSparse* curvature = fct.curvature();
757
758 // Save function value, gradient and curvature matrix
759 double save_value = m_value;
760 GVector save_grad = GVector(*grad);
761 GMatrixSparse save_curvature = GMatrixSparse(*curvature);
762
763 // Save parameter values in vector
764 GVector save_pars(m_npars);
765 for (int ipar = 0; ipar < m_npars; ++ipar) {
766 save_pars[ipar] = pars[ipar]->factor_value();
767 }
768
769 // Setup matrix and vector for curvature computation
770 for (int ipar = 0; ipar < m_npars; ++ipar) {
771 (*curvature)(ipar,ipar) *= (1.0 + m_lambda);
772 (*grad)[ipar] = -(*grad)[ipar];
773 }
774
775 // Debug option: dump gradient and curvature matrix
776 #if defined(G_DEBUG_ITER)
777 std::cout << "Gradient : ";
778 for (int ipar = 0; ipar < m_npars; ++ipar) {
779 std::cout << (*grad)[ipar] << " ";
780 }
781 std::cout << std::endl;
782 std::cout << "Curvature: ";
783 for (int ipar = 0; ipar < m_npars; ++ipar) {
784 for (int jpar = 0; jpar < m_npars; ++jpar) {
785 std::cout << (*curvature)(ipar,jpar) << " ";
786 }
787 }
788 std::cout << std::endl;
789 #endif
790
791 // Solve: curvature * X = grad. Handle matrix problems
792 try {
793 *grad = curvature->solve(*grad);
794 }
795 catch (GException::runtime_error &e) {
797 if (m_logger != NULL) {
798 *m_logger << "GOptimizerLM::iteration: "
799 << "Curvature matrix not positive definite."
800 << std::endl;
801 }
802 continue;
803 }
804 catch (std::exception &e) {
805 throw;
806 }
807
808 // Get LM step size
809 step = step_size(*grad, pars);
810
811 // Debug option: dump step size
812 #if defined(G_DEBUG_ITER)
813 std::cout << "Step size = " << step << std::endl;
814 #endif
815
816 // Derive new parameter vector
817 for (int ipar = 0; ipar < m_npars; ++ipar) {
818
819 // Consider only free parameters
820 if (pars[ipar]->is_free()) {
821
822 // Get actual parameter value and limits
823 double p = pars[ipar]->factor_value();
824 double p_min = pars[ipar]->factor_min();
825 double p_max = pars[ipar]->factor_max();
826
827 // Compute new parameter value
828 p += (*grad)[ipar] * step;
829
830 // Debug option: dump new parameter value
831 #if defined(G_DEBUG_ITER)
832 std::cout << "Trial factor value of ";
833 std::cout << pars[ipar]->name() << " = ";
834 std::cout << p << std::endl;
835 #endif
836
837 // Constrain parameter to within the valid range
838 if (pars[ipar]->has_min() && p < p_min) {
839 if (m_hit_minimum[ipar] >= m_max_hit) {
840 if (m_logger != NULL) {
841 *m_logger << " Parameter \"" << pars[ipar]->name();
842 *m_logger << "\" hits minimum " << p_min << " more than ";
843 *m_logger << m_max_hit << " times.";
844 *m_logger << " Fix parameter at minimum for now." << std::endl;
845 }
846 m_par_freeze[ipar] = true;
847 pars[ipar]->fix();
848 }
849 else {
850 if (m_logger != NULL) {
851 *m_logger << " Parameter \"" << pars[ipar]->name();
852 *m_logger << "\" hits minimum: ";
853 *m_logger << p << " < " << p_min;
854 *m_logger << " (" << m_hit_minimum[ipar]+1 << ")" << std::endl;
855 }
856 }
857 m_hit_minimum[ipar]++;
858 p = p_min;
859 }
860 else if (pars[ipar]->has_max() && p > p_max) {
861 if (m_hit_maximum[ipar] >= m_max_hit) {
862 if (m_logger != NULL) {
863 *m_logger << " Parameter \"" << pars[ipar]->name();
864 *m_logger << "\" hits maximum " << p_max << " more than ";
865 *m_logger << m_max_hit << " times.";
866 *m_logger << " Fix parameter at maximum for now." << std::endl;
867 }
868 m_par_freeze[ipar] = true;
869 pars[ipar]->fix();
870 }
871 else {
872 if (m_logger != NULL) {
873 *m_logger << " Parameter \"" << pars[ipar]->name();
874 *m_logger << "\" hits maximum: ";
875 *m_logger << p << " > " << p_max;
876 *m_logger << " (" << m_hit_maximum[ipar]+1 << ")" << std::endl;
877 }
878 }
879 m_hit_maximum[ipar]++;
880 p = p_max;
881 }
882 else {
883 m_hit_minimum[ipar] = 0;
884 m_hit_maximum[ipar] = 0;
885 }
886
887 // Set new parameter value
888 pars[ipar]->factor_value(p);
889
890 // Debug option: dump new parameter value
891 #if defined(G_DEBUG_ITER)
892 std::cout << "New value of ";
893 std::cout << pars[ipar]->name() << " = ";
894 std::cout << pars[ipar]->value() << std::endl;
895 #endif
896
897 } // endif: Parameter was free
898
899 } // endfor: computed new parameter vector
900
901 // Evaluate function at new parameters
902 fct.eval(pars);
903
904 // If a free parameter has a zero diagonal element in the curvature
905 // matrix then remove this parameter definitely from the fit as it
906 // otherwise will block the fit
907 for (int ipar = 0; ipar < m_npars; ++ipar) {
908 if (pars[ipar]->is_free()) {
909 if ((*fct.curvature())(ipar,ipar) == 0.0) {
910 if (m_logger != NULL) {
911 *m_logger << " Parameter \"" << pars[ipar]->name();
912 *m_logger << "\" = " << pars[ipar]->value();
913 *m_logger << " with gradient " << pars[ipar]->gradient();
914 *m_logger << " has zero curvature.";
915 *m_logger << " Fix parameter." << std::endl;
916 }
917 m_par_remove[ipar] = true;
918 pars[ipar]->fix();
919 }
920 }
921 }
922
923 // Fetch new pointers since eval will allocate new memory
924 grad = fct.gradient();
925 curvature = fct.curvature();
926
927 // Retrieve new function value
928 m_value = fct.value();
929
930 // Debug option: dump new function value
931 #if defined(G_DEBUG_ITER)
932 std::cout << "Function value " << m_value;
933 std::cout << " (was " << save_value << ")" << std::endl;
934 #endif
935
936 // Determine how many parameters have changed
937 int par_change = 0;
938 for (int ipar = 0; ipar < m_npars; ++ipar) {
939 if (pars[ipar]->factor_value() != save_pars[ipar]) {
940 par_change++;
941 }
942 }
943
944 // If the function has decreased then accept the new solution
945 // and decrease lambda ...
946 m_delta = save_value - m_value;
947 if (m_delta > 0.0) {
949 }
950
951 // ... if function is identical then accept new solution. If the
952 // parameters have changed then increase lambda, otherwise decrease
953 // lambda
954 else if (m_delta == 0.0) {
955 if (par_change) {
957 }
958 else {
960 }
961 }
962
963 // ... if function worsened by an amount smaller than m_accept_dec
964 // then accept new solution and increase lambda for the next
965 // iteration. This kluge which may help to overcome some local
966 // minimum will only for m_max_dec times at maximum to avoid
967 // keeping oscillating in a loop.
968 else if ((m_num_dec < m_max_dec) && (m_delta >= -std::abs(m_accept_dec))) {
970 m_num_dec++;
971 }
972
973 // ... otherwise, if the statistics did not improve then use old
974 // parameters and increase lamdba. Restore also the best statistics
975 // value that was reached so far, the gradient vector and the curve
976 // matrix.
977 else {
979 m_value = save_value;
980 *grad = save_grad;
981 *curvature = save_curvature;
982 for (int ipar = 0; ipar < m_npars; ++ipar) {
983 pars[ipar]->factor_value(save_pars[ipar]);
984 }
985 }
986
987 } while (0); // endwhile: main loop
988
989 // Debug option: dump trailer
990 #if defined(G_DEBUG_ITER)
991 std::cout << "GOptimizerLM::iteration: exit" << std::endl;
992 #endif
993
994 // Return step size
995 return step;
996
997}
998
999
1000/***********************************************************************//**
1001 * @brief Return LM step size
1002 *
1003 * @param[in] grad Function gradient.
1004 * @param[in] pars Function parameters.
1005 *
1006 * Determine the size of the LM step. By default a step size of 1 is taken.
1007 * If m_step_adjust=true then the step size will be estimated so that the
1008 * next parameter vector should stay within the parameter boundaries
1009 * (provided that boundaries exist).
1010 ***************************************************************************/
1011double GOptimizerLM::step_size(const GVector& grad, const GOptimizerPars& pars)
1012{
1013 // Initialise step size
1014 double step = 1.0;
1015
1016 // Check if we should reduce the step size
1017 if (m_step_adjust) {
1018
1019 // Initialise the parameter index that constrains most the fit
1020 int ipar_bnd = -1;
1021
1022 // Loop over all parameters
1023 for (int ipar = 0; ipar < pars.size(); ++ipar) {
1024
1025 // Get parameter attributes
1026 double p = pars[ipar]->factor_value();
1027 double p_min = pars[ipar]->factor_min();
1028 double p_max = pars[ipar]->factor_max();
1029 double delta = grad[ipar];
1030
1031 // Check if a parameter minimum requires a reduced step size
1032 if (pars[ipar]->has_min()) {
1033 double step_min = (delta < 0.0) ? (p_min - p)/delta : 1.0;
1034 if (step_min > 0.0) {
1035 if (step_min < step) {
1036 if (!m_hit_boundary[ipar]) {
1037 ipar_bnd = ipar;
1038 step = step_min;
1039 }
1040 }
1041 else if (m_hit_boundary[ipar]) {
1042 m_hit_boundary[ipar] = false;
1043 if (m_logger != NULL) {
1044 *m_logger << " Parameter \"";
1045 *m_logger << pars[ipar]->name();
1046 *m_logger << "\" does not drive optimization step anymore.";
1047 *m_logger << std::endl;
1048 }
1049 }
1050 }
1051 }
1052
1053 // Check if a parameter maximum requires a reduced step size
1054 if (pars[ipar]->has_max()) {
1055 double step_max = (delta > 0.0) ? (p_max - p)/delta : 1.0;
1056 if (step_max > 0.0) {
1057 if (step_max < step) {
1058 if (!m_hit_boundary[ipar]) {
1059 ipar_bnd = ipar;
1060 step = step_max;
1061 }
1062 }
1063 else if (m_hit_boundary[ipar]) {
1064 m_hit_boundary[ipar] = false;
1065 if (m_logger != NULL) {
1066 *m_logger << " Parameter \"";
1067 *m_logger << pars[ipar]->name();
1068 *m_logger << "\" does not drive optimization step anymore.";
1069 *m_logger << std::endl;
1070 }
1071 }
1072 }
1073 }
1074
1075 } // endfor: looped over all parameters
1076
1077 // Signal if a parameter is driving the optimization step
1078 if (ipar_bnd != -1) {
1079 m_hit_boundary[ipar_bnd] = true;
1080 if (m_logger != NULL) {
1081 *m_logger << " Parameter \"";
1082 *m_logger << pars[ipar_bnd]->name();
1083 *m_logger << "\" drives optimization step (step=";
1084 *m_logger << step << ")" << std::endl;
1085 }
1086 }
1087
1088 } // endif: automatic step size adjustment requested
1089
1090 // Return step size
1091 return step;
1092}
Exception handler interface definition.
Levenberg Marquardt optimizer class interface definition.
#define G_LM_STALLED
#define G_LM_NOT_POSTIVE_DEFINITE
#define G_LM_BAD_ERRORS
#define G_LM_SINGULAR
#define G_LM_CONVERGED
Gammalib tools definition.
GChatter
Definition GTypemaps.hpp:33
@ SILENT
Definition GTypemaps.hpp:34
GVector log(const GVector &vector)
Computes natural logarithm of vector elements.
Definition GVector.cpp:1342
GVector sqrt(const GVector &vector)
Computes square root of vector elements.
Definition GVector.cpp:1426
Information logger interface definition.
Definition GLog.hpp:62
Sparse matrix class interface definition.
GMatrixSparse cholesky_decompose(const bool &compress=true) const
Return Cholesky decomposition.
GVector solve(const GVector &vector) const
Solves linear matrix equation.
GVector cholesky_solver(const GVector &vector, const bool &compress=true) const
Cholesky solver.
Optimizer function abstract base class.
virtual double value(void) const =0
virtual GMatrixSparse * curvature(void)=0
virtual GVector * gradient(void)=0
virtual void eval(const GOptimizerPars &pars)=0
Levenberg Marquardt optimizer class.
bool m_step_adjust
Adjust step size to boundaries.
int m_npars
Number of parameters.
int m_max_iter
Maximum number of iterations.
GLog * m_logger
Pointer to optional logger.
double m_lambda_inc
Lambda increase.
std::string status_string(void) const
Set fit status string.
int m_max_dec
Maximum number of function decrease.
std::vector< int > m_hit_minimum
Bookkeeping of successive minimum hits.
int m_num_dec
Number of function decreases.
std::vector< bool > m_hit_boundary
Bookkeeping array for boundary hits.
virtual GOptimizerLM * clone(void) const
Clone object.
void init_members(void)
Initialise class members.
GOptimizerLM & operator=(const GOptimizerLM &opt)
Assignment operator.
void free_members(void)
Delete class members.
const double & lambda_inc(void) const
Return lambda increment value.
double m_value
Actual function value.
int m_max_stall
Maximum number of stalls.
double m_lambda_dec
Lambda decrease.
std::vector< bool > m_par_freeze
Bookkeeping of parameter freeze.
GOptimizerLM(void)
Void constructor.
double step_size(const GVector &grad, const GOptimizerPars &pars)
Return LM step size.
virtual void errors(GOptimizerFunction &fct, GOptimizerPars &pars)
Compute parameter uncertainties.
double iteration(GOptimizerFunction &fct, GOptimizerPars &pars)
Perform one LM iteration.
void copy_members(const GOptimizerLM &opt)
Copy class members.
double m_eps
Absolute precision.
virtual std::string print(const GChatter &chatter=NORMAL) const
Print optimizer information.
virtual void optimize(GOptimizerFunction &fct, GOptimizerPars &pars)
Optimize function parameters.
std::vector< int > m_hit_maximum
Bookkeeping of successive maximum hits.
int m_status
Fit status.
std::vector< bool > m_par_remove
Bookkeeping of parameter removal.
double m_lambda_start
Initial start value.
int m_nfree
Number of free parameters.
double m_delta
Function improvement.
virtual void clear(void)
Clear object.
int m_iter
Iteration.
int m_max_hit
Maximum number of successive hits.
double m_accept_dec
Acceptable function decrease.
virtual ~GOptimizerLM(void)
Destructor.
double m_lambda
Actual lambda.
virtual int status(void) const
Return optimizer status.
const int & npars(void) const
Return number of model parameters.
Optimizer parameter container class.
int size(void) const
Return number of parameters in container.
int nfree(void) const
Return number of free parameters.
Abstract optimizer abstract base class.
void free_members(void)
Delete class members.
virtual GOptimizer & operator=(const GOptimizer &opt)
Assignment operator.
void init_members(void)
Initialise class members.
Vector class.
Definition GVector.hpp:46
std::string parformat(const std::string &s, const int &indent=0)
Convert string in parameter format.
Definition GTools.cpp:1162
std::string str(const unsigned short int &value)
Convert unsigned short integer value into string.
Definition GTools.cpp:508