Visual Servoing Platform  version 3.0.1
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
vpLevenbergMarquartd.cpp
1 /****************************************************************************
2  *
3  * This file is part of the ViSP software.
4  * Copyright (C) 2005 - 2017 by Inria. All rights reserved.
5  *
6  * This software is free software; you can redistribute it and/or
7  * modify it under the terms of the GNU General Public License
8  * ("GPL") version 2 as published by the Free Software Foundation.
9  * See the file LICENSE.txt at the root directory of this source
10  * distribution for additional information about the GNU GPL.
11  *
12  * For using ViSP with software that can not be combined with the GNU
13  * GPL, please contact Inria about acquiring a ViSP Professional
14  * Edition License.
15  *
16  * See http://visp.inria.fr for more information.
17  *
18  * This software was developed at:
19  * Inria Rennes - Bretagne Atlantique
20  * Campus Universitaire de Beaulieu
21  * 35042 Rennes Cedex
22  * France
23  *
24  * If you have questions regarding the use of this file, please contact
25  * Inria at visp@inria.fr
26  *
27  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
28  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
29  *
30  * Description:
31  * Levenberg Marquartd.
32  *
33  * Authors:
34  * Eric Marchand
35  * Francois Chaumette
36  *
37  *****************************************************************************/
38 
39 #include <cmath> // std::fabs
40 #include <limits> // numeric_limits
41 #include <iostream>
42 #include <algorithm> // std::min
43 
44 #include <visp3/vision/vpLevenbergMarquartd.h>
45 #include <visp3/core/vpMath.h>
46 
47 #define SIGN(x) ((x) < 0 ? -1 : 1)
48 #define SWAP(a,b,c) {(c) = (a); (a) = (b); (b) = (c);}
49 #define MIJ(m,i,j,s) ((m) + ((long) (i) * (long) (s)) + (long) (j))
50 #define TRUE 1
51 #define FALSE 0
52 
53 /*
54  * PROCEDURE : enorm
55  *
56  * ENTREE :
57  *
58  * x Vecteur de taille "n"
59  * n Taille du vecteur "x"
60  *
61  * DESCRIPTION :
62  * La procedure calcule la norme euclidienne d'un vecteur "x" de taille "n"
63  * La norme euclidienne est calculee par accumulation de la somme des carres
64  * dans les trois directions. Les sommes des carres pour les petits et grands
65  * elements sont mis a echelle afin d'eviter les overflows. Des underflows non
66  * destructifs sont autorisee. Les underflows et overflows sont evites dans le
67  * calcul des sommes des carres non encore mis a echelle par les elements
68  * intermediaires. La definition des elements petit, intermediaire et grand
69  * depend de deux constantes : rdwarf et rdiant. Les restrictions principales
70  * sur ces constantes sont rdwarf^2 n'est pas en underflow et rdgiant^2 n'est
71  * pas en overflow. Les constantes donnees ici conviennent pour la plupart des
72  * pc connus.
73  *
74  * RETOUR :
75  * En cas de succes, la valeur retournee est la norme euclidienne du vecteur
76  * Sinon, la valeur -1 est retournee et la variable globale "errno" est
77  * initialisee pour indiquee le type de l'erreur.
78  *
79  */
80 double enorm (const double *x, int n)
81 {
82  const double rdwarf = 3.834e-20;
83  const double rgiant = 1.304e19;
84 
85  int i;
86  double agiant, floatn;
87  double norm_eucl = 0.0;
88  double s1 = 0.0, s2 = 0.0, s3 = 0.0;
89  double x1max = 0.0 , x3max = 0.0;
90 
91  floatn = (double) n;
92  agiant = rgiant / floatn;
93 
94  for (i = 0; i < n; i++)
95  {
96  double xabs = std::fabs(x[i]);
97  if ((xabs > rdwarf) && (xabs < agiant))
98  {
99  /*
100  * somme pour elements intemediaires.
101  */
102  s2 += xabs * xabs;
103  }
104 
105  else if (xabs <= rdwarf)
106  {
107  /*
108  * somme pour elements petits.
109  */
110  if (xabs <= x3max)
111  {
112  //if (xabs != 0.0)
113  if (xabs > std::numeric_limits<double>::epsilon())
114  s3 += (xabs / x3max) * (xabs / x3max);
115  }
116  else
117  {
118  s3 = 1.0 + s3 * (x3max / xabs) * (x3max / xabs);
119  x3max = xabs;
120  }
121  }
122 
123  else
124  {
125  /*
126  * somme pour elements grand.
127  */
128  if (xabs <= x1max)
129  {
130  s1 += (xabs / x1max) * (xabs / x1max);
131  }
132  else
133  {
134  s1 = 1.0 + s1 * (x1max / xabs) * (x1max / xabs);
135  x1max = xabs;
136  }
137  }
138  }
139 
140  /*
141  * calcul de la norme.
142  */
143  //if (s1 == 0.0)
144  if (std::fabs(s1) <= std::numeric_limits<double>::epsilon())
145  {
146  //if (s2 == 0.0)
147  if (std::fabs(s2) <= std::numeric_limits<double>::epsilon())
148  norm_eucl = x3max * sqrt(s3);
149  else if (s2 >= x3max)
150  norm_eucl = sqrt (s2 * (1.0 + ( x3max / s2) * (x3max * s3)));
151  else if (s2 < x3max)
152  norm_eucl = sqrt (x3max * ((s2 / x3max) + (x3max * s3)));
153  }
154  else
155  norm_eucl = x1max * sqrt (s1 + (s2 / x1max) / x1max);
156 
157  return (norm_eucl);
158 }
159 
160 /* PROCEDURE : lmpar
161  *
162  * ENTREE :
163  * n Ordre de la matrice "r".
164  * r Matrice de taille "n" x "n". En entree, la toute la partie
165  * triangulaire superieure doit contenir toute la partie triangulaire
166  * superieure de "r".
167  *
168  * ldr Taille maximum de la matrice "r". "ldr" >= "n".
169  *
170  * ipvt Vecteur de taille "n" qui definit la matrice de permutation "p"
171  * tel que :
172  * a * p = q * r.
173  * La jeme colonne de p la colonne ipvt[j] de la matrice d'identite.
174  *
175  * diag Vecteur de taille "n" contenant les elements diagonaux de la
176  * matrice "d".
177  *
178  * qtb Vecteur de taille "n" contenant les "n" premiers elements du
179  * vecteur (q transpose)*b.
180  *
181  * delta Limite superieure de la norme euclidienne de d * x.
182  *
183  * par Estimee initiale du parametre de Levenberg-Marquardt.
184  * wa1, wa2 Vecteurs de taille "n" de travail.
185  *
186  * DESCRIPTION :
187  * La procedure determine le parametre de Levenberg-Marquardt. Soit une matrice
188  * "a" de taille "m" x "n", une matrice diagonale "d" non singuliere de taille
189  * "n" x "n", un vecteur "b" de taille "m" et un nombre positf delta, le probleme
190  * est le calcul du parametre "par" de telle facon que si "x" resoud le systeme
191  *
192  * a * x = b , sqrt(par) * d * x = 0 ,
193  *
194  * au sens des moindre carre, et dxnorm est la norme euclidienne de d * x
195  * alors "par" vaut 0.0 et (dxnorm - delta) <= 0.1 * delta ,
196  * ou "par" est positif et abs(dxnorm-delta) <= 0.1 * delta.
197  * Cette procedure complete la solution du probleme si on lui fourni les infos
198  * nessaires de la factorisation qr, avec pivotage de colonnes de a.
199  * Donc, si a * p = q * r, ou "p" est une matrice de permutation, les colonnes
200  * de "q" sont orthogonales, et "r" est une matrice triangulaire superieure
201  * avec les elements diagonaux classes par ordre decroissant de leur valeur, lmpar
202  * attend une matrice triangulaire superieure complete, la matrice de permutation
203  * "p" et les "n" premiers elements de (q transpose) * b. En sortie, la procedure
204  * lmpar fournit aussi une matrice triangulaire superieure "s" telle que
205  *
206  * t t t
207  * p * (a * a + par * d * d )* p = s * s .
208  *
209  * "s" est utilise a l'interieure de lmpar et peut etre d'un interet separe.
210  *
211  * Seulement quelques iterations sont necessaire pour la convergence de
212  * l'algorithme. Si neanmoins la limite de 10 iterations est atteinte, la
213  * valeur de sortie "par" aura la derniere meilleure valeur.
214  *
215  * SORTIE :
216  * r En sortie, tout le triangle superieur est inchange, et le
217  * le triangle inferieur contient les elements de la partie
218  * triangulaire superieure (transpose) de la matrice triangulaire
219  * superieure de "s".
220  * par Estimee finale du parametre de Levenberg-Marquardt.
221  * x Vecteur de taille "n" contenant la solution au sens des moindres
222  * carres du systeme a * x = b, sqrt(par) * d * x = 0, pour le
223  * parametre en sortie "par"
224  * sdiag Vecteur de taille "n" contenant les elements diagonaux de la
225  * matrice triangulaire "s".
226  *
227  * RETOUR :
228  * En cas de succes, la valeur 0.0 est retournee.
229  *
230  */
231 int lmpar(int n, double *r, int ldr, int *ipvt, double *diag, double *qtb,
232  double *delta, double *par, double *x, double *sdiag, double *wa1, double *wa2)
233 {
234  const double tol1 = 0.1; /* tolerance a 0.1 */
235 
236  int l;
237  unsigned int iter; /* compteur d'iteration */
238  int nsing; /* nombre de singularite de la matrice */
239  double dxnorm, fp;
240  double temp;
241  double dwarf = DBL_MIN; /* plus petite amplitude positive */
242 
243  /*
244  * calcul et stockage dans "x" de la direction de Gauss-Newton. Si
245  * le jacobien a une rangee de moins, on a une solution au moindre
246  * carres.
247  */
248  nsing = n;
249 
250  for (int i = 0; i < n; i++)
251  {
252  wa1[i] = qtb[i];
253  double *pt = MIJ(r, i, i, ldr);
254  //if (*MIJ(r, i, i, ldr) == 0.0 && nsing == n)
255  if (std::fabs(*pt) <= std::numeric_limits<double>::epsilon() && nsing == n)
256  nsing = i - 1;
257 
258  if (nsing < n)
259  wa1[i] = 0.0;
260  }
261 
262  if (nsing >= 0)
263  {
264  for (int k = 0; k < nsing; k++)
265  {
266  int i = nsing - 1 - k;
267  wa1[i] /= *MIJ(r, i, i, ldr);
268  temp = wa1[i];
269  int jm1 = i - 1;
270 
271  if (jm1 >= 0)
272  {
273  for (unsigned int j = 0; j <= (unsigned int)jm1; j++)
274  wa1[j] -= *MIJ(r, i, j, ldr) * temp;
275  }
276  }
277  }
278 
279  for (int j = 0; j < n; j++)
280  {
281  l = ipvt[j];
282  x[l] = wa1[j];
283  }
284 
285  /*
286  * initialisation du compteur d'iteration.
287  * evaluation de la fonction a l'origine, et test
288  * d'acceptation de la direction de Gauss-Newton.
289  */
290  iter = 0;
291 
292  for (int i = 0; i < n; i++)
293  wa2[i] = diag[i] * x[i];
294 
295  dxnorm = enorm(wa2, n);
296 
297  fp = dxnorm - *delta;
298 
299  if (fp > tol1 * (*delta))
300  {
301  /*
302  * Si le jacobien n'a pas de rangee deficiente,l'etape de
303  * Newton fournit une limite inferieure, parl pour le
304  * zero de la fonction. Sinon cette limite vaut 0.0.
305  */
306  double parl = 0.0;
307 
308  if (nsing >= n)
309  {
310  for (int i = 0; i < n; i++)
311  {
312  l = ipvt[i];
313  wa1[i] = diag[l] * (wa2[l] / dxnorm);
314  }
315 
316  for (int i = 0; i < n; i++)
317  {
318  long im1;
319  double sum = 0.0;
320  im1 = (i - 1L);
321 
322  if (im1 >= 0)
323  {
324  for (unsigned int j = 0; j <= (unsigned int)im1; j++)
325  sum += (*MIJ(r, i, j, ldr) * wa1[j]);
326  }
327  wa1[i] = (wa1[i] - sum) / *MIJ(r, i, i, ldr);
328  }
329 
330  temp = enorm(wa1, n);
331  parl = ((fp / *delta) / temp) / temp;
332  }
333 
334  /*
335  * calcul d'une limite superieure, paru, pour le zero de la
336  * fonction.
337  */
338  for (int i = 0; i < n; i++)
339  {
340  double sum = 0.0;
341 
342  for (int j = 0; j <= i; j++)
343  sum += *MIJ(r, i, j, ldr) * qtb[j];
344 
345  l = ipvt[i];
346  wa1[i] = sum / diag[l];
347  }
348 
349  double gnorm = enorm(wa1, n);
350  double paru = gnorm / *delta;
351 
352  //if (paru == 0.0)
353  if (std::fabs(paru) <= std::numeric_limits<double>::epsilon())
354  paru = dwarf / vpMath::minimum(*delta, tol1);
355 
356  /*
357  * Si "par" en entree tombe hors de l'intervalle (parl,paru),
358  * on le prend proche du point final.
359  */
360 
361  *par = vpMath::maximum(*par, parl);
362  *par = vpMath::maximum(*par, paru);
363 
364  //if (*par == 0.0)
365  if (std::fabs(*par) <= std::numeric_limits<double>::epsilon() )
366  *par = gnorm / dxnorm;
367 
368  /*
369  * debut d'une iteration.
370  */
371  for(;;) //iter >= 0)
372  {
373  iter++;
374 
375  /*
376  * evaluation de la fonction a la valeur courant
377  * de "par".
378  */
379  //if (*par == 0.0)
380  if (std::fabs(*par) <= std::numeric_limits<double>::epsilon() ) {
381  const double tol001 = 0.001; /* tolerance a 0.001 */
382  *par = vpMath::maximum(dwarf, (tol001 * paru));
383  }
384 
385  temp = sqrt(*par);
386 
387  for (int i = 0; i < n; i++)
388  wa1[i] = temp * diag[i];
389 
390  qrsolv(n, r, ldr, ipvt, wa1, qtb, x, sdiag, wa2);
391 
392  for (int i = 0; i < n; i++)
393  wa2[i] = diag[i] * x[i];
394 
395  dxnorm = enorm(wa2, n);
396  temp = fp;
397  fp = dxnorm - *delta;
398 
399  /*
400  * si la fonction est assez petite, acceptation de
401  * la valeur courant de "par". de plus, test des cas
402  * ou parl est nul et ou le nombre d'iteration a
403  * atteint 10.
404  */
405  //if ((std::fabs(fp) <= tol1 * (*delta)) || ((parl == 0.0) && (fp <= temp)
406  // && (temp < 0.0)) || (iter == 10))
407  if ((std::fabs(fp) <= tol1 * (*delta)) || ((std::fabs(parl) <= std::numeric_limits<double>::epsilon()) && (fp <= temp)
408  && (temp < 0.0)) || (iter == 10))
409  {
410  // terminaison.
411 
412  // Remove the two next lines since this is a dead code
413  /* if (iter == 0)
414  *par = 0.0; */
415 
416  return (0);
417  }
418 
419  /*
420  * calcul de la correction de Newton.
421  */
422 
423  for (int i = 0; i < n; i++)
424  {
425  l = ipvt[i];
426  wa1[i] = diag[l] * (wa2[l] / dxnorm);
427  }
428 
429  for (unsigned int i = 0; i < (unsigned int)n; i++)
430  {
431  wa1[i] = wa1[i] / sdiag[i];
432  temp = wa1[i];
433  unsigned int jp1 = i + 1;
434  if ( (unsigned int) n >= jp1)
435  {
436  for (unsigned int j = jp1; j < (unsigned int)n; j++)
437  wa1[j] -= (*MIJ(r, i, j, ldr) * temp);
438  }
439  }
440 
441  temp = enorm(wa1, n);
442  double parc = ((fp / *delta) / temp) / temp;
443 
444  /*
445  * selon le signe de la fonction, mise a jour
446  * de parl ou paru.
447  */
448  if (fp > 0.0)
449  parl = vpMath::maximum(parl, *par);
450 
451  if (fp < 0.0)
452  paru = vpMath::minimum(paru, *par);
453 
454  /*
455  * calcul d'une estimee ameliree de "par".
456  */
457  *par = vpMath::maximum(parl, (*par + parc));
458  }/* fin boucle sur iter */
459  }/* fin fp > tol1 * delta */
460 
461  /*
462  * terminaison.
463  */
464  if (iter == 0)
465  *par = 0.0;
466 
467  return (0);
468 }
469 
470 /*
471  * PROCEDURE : pythag
472  *
473  * ENTREES :
474  * a, b Variables dont on veut la racine carre de leur somme de carre
475  *
476  * DESCRIPTION :
477  * La procedure calcule la racine carre de la somme des carres de deux nombres
478  * en evitant l'overflow ou l'underflow destructif.
479  *
480  * RETOUR :
481  * La procedure retourne la racine carre de a^2 + b^2.
482  *
483  */
484 double pythag (double a, double b)
485 {
486  double pyth, p, r, t;
487 
488  p = vpMath::maximum(std::fabs(a), std::fabs(b));
489 
490  //if (p == 0.0)
491  if (std::fabs(p) <= std::numeric_limits<double>::epsilon() )
492  {
493  pyth = p;
494  return (pyth);
495  }
496 
497  r = (std::min(std::fabs(a), std::fabs(b)) / p) * (std::min(std::fabs(a), std::fabs(b)) / p);
498  t = 4.0 + r;
499 
500  //while (t != 4.0)
501  while (std::fabs(t - 4.0) < std::fabs(vpMath::maximum(t,4.0)) * std::numeric_limits<double>::epsilon() )
502  {
503  double s = r / t;
504  double u = 1.0 + 2.0 * s;
505  p *= u;
506  r *= (s / u) * (s / u);
507  t = 4.0 + r;
508  }
509 
510  pyth = p;
511  return (pyth);
512 }
513 
514 /*
515  * PROCEDURE : qrfac
516  *
517  * ENTREE :
518  * m Nombre de lignes de la matrice "a".
519  * n Nombre de colonne de la matrice "a".
520  * a Matrice de taille "m" x "n". elle contient, en entree la matrice
521  * dont on veut sa factorisation qr.
522  * lda Taille maximale de "a". lda >= m.
523  * pivot Booleen. Si pivot est TRUE, le pivotage de colonnes est realise
524  * Si pivot = FALSE, pas de pivotage.
525  * lipvt Taille du vecteur "ipvt". Si pivot est FALSE, lipvt est de
526  * l'ordre de 1. Sinon lipvt est de l'ordre de "n".
527  * wa Vecteur de travail de taille "n". Si pivot = FALSE "wa"
528  * coincide avec rdiag.
529  *
530  * DESCRIPTION :
531  * La procedure effectue une decomposition de la matrice "a"par la methode qr.
532  * Elle utilise les transformations de householders avec pivotage sur les colonnes
533  * (option) pour calculer la factorisation qr de la matrice "a" de taille "m" x "n".
534  * La procedure determine une matrice orthogonale "q", une matrice de permutation
535  * "p" et une matrice trapesoidale superieure "r" dont les elements diagonaux
536  * sont ordonnes dans l'ordre decroissant de leurs valeurs,tel que a * p = q * r.
537  * La transformation de householder pour la colonne k, k = 1,2,...,min(m,n),
538  * est de la forme
539  * t
540  * i - (1 / u(k)) * u * u
541  *
542  * Ou u a des zeros dans les k-1 premieres positions.
543  *
544  * SORTIE :
545  * a Matrice de taille "m" x "n" dont le trapeze superieur de "a"
546  * contient la partie trapezoidale superieure de "r" et la partie
547  * trapezoidale inferieure de "a" contient une forme factorisee
548  * de "q" (les elements non triviaux du vecteurs "u" sont decrits
549  * ci-dessus).
550  * ipvt Vecteur de taille "n". Il definit la matrice de permutation "p"
551  * tel que a * p = q * r. La jeme colonne de p est la colonne
552  * ipvt[j] de la matrice d'identite. Si pivot = FALSE, ipvt n'est
553  * pas referencee.
554  * rdiag Vecteur de taille "n" contenant les elements diagonaux de la
555  * matrice "r".
556  * acnorm Vecteur de taille "n" contenant les normes des lignes
557  * correspondantes de la matrice "a". Si cette information n'est
558  * pas requise, acnorm coincide avec rdiag.
559  *
560  */
561 int qrfac(int m, int n, double *a, int lda, int *pivot, int *ipvt,
562  int /* lipvt */, double *rdiag, double *acnorm, double *wa)
563 {
564  const double tolerance = 0.05;
565 
566  int i, j, ip1, k, kmax, minmn;
567  double epsmch;
568  double sum, temp, tmp;
569 
570  /*
571  * epsmch est la precision machine.
572  */
573  epsmch = std::numeric_limits<double>::epsilon();
574 
575  /*
576  * calcul des normes initiales des lignes et initialisation
577  * de plusieurs tableaux.
578  */
579  for (i = 0; i < m; i++)
580  {
581  acnorm[i] = enorm(MIJ(a, i, 0, lda), n);
582  rdiag[i] = acnorm[i];
583  wa[i] = rdiag[i];
584 
585  if (pivot)
586  ipvt[i] = i;
587  }
588  /*
589  * reduction de "a" en "r" avec les tranformations de Householder.
590  */
591  minmn = vpMath::minimum(m, n);
592  for (i = 0; i < minmn; i++)
593  {
594  if (pivot)
595  {
596  /*
597  * met la ligne de plus grande norme en position
598  * de pivot.
599  */
600  kmax = i;
601  for (k = i; k < m; k++)
602  {
603  if (rdiag[k] > rdiag[kmax])
604  kmax = k;
605  }
606 
607  if (kmax != i)
608  {
609  for (j = 0; j < n; j++)
610  SWAP(*MIJ(a, i, j, lda),
611  *MIJ(a, kmax, j, lda), tmp);
612 
613  rdiag[kmax] = rdiag[i];
614  wa[kmax] = wa[i];
615 
616  SWAP( ipvt[i], ipvt[kmax], k);
617  }
618  }
619 
620  /*
621  * calcul de al transformationde Householder afin de reduire
622  * la jeme ligne de "a" a un multiple du jeme vecteur unite.
623  */
624  double ajnorm = enorm(MIJ(a, i, i, lda), n - i);
625 
626  //if (ajnorm != 0.0)
627  if (std::fabs(ajnorm) > std::numeric_limits<double>::epsilon() )
628  {
629  if (*MIJ(a, i, i, lda) < 0.0)
630  ajnorm = -ajnorm;
631 
632  for (j = i; j < n; j++)
633  *MIJ(a, i, j, lda) /= ajnorm;
634  *MIJ(a, i, i, lda) += 1.0;
635 
636  /*
637  * application de la tranformation aux lignes
638  * restantes et mise a jour des normes.
639  */
640  ip1 = i + 1;
641 
642  if (m >= ip1)
643  {
644  for (k = ip1; k < m; k++)
645  {
646  sum = 0.0;
647  for (j = i; j < n; j++)
648  sum += *MIJ(a, i, j, lda) * *MIJ(a, k, j, lda);
649 
650  temp = sum / *MIJ(a, i, i, lda);
651 
652  for (j = i; j < n; j++)
653  *MIJ(a, k, j, lda) -= temp * *MIJ(a, i, j, lda);
654 
655  //if (pivot && rdiag[k] != 0.0)
656  if (pivot && (std::fabs(rdiag[k]) > std::numeric_limits<double>::epsilon()) )
657  {
658  temp = *MIJ (a, k, i, lda) / rdiag[k];
659  rdiag[k] *= sqrt(vpMath::maximum(0.0, (1.0 - temp * temp)));
660 
661  if (tolerance * (rdiag[k] / wa[k]) * (rdiag[k] / wa[k]) <= epsmch)
662  {
663  rdiag[k] = enorm(MIJ(a, k, ip1, lda), (n -1 - (int) i));
664  wa[k] = rdiag[k];
665  }
666  }
667  }/* fin boucle for k */
668 
669  }
670 
671  } /* fin if (ajnorm) */
672 
673  rdiag[i] = -ajnorm;
674  } /* fin for (i = 0; i < minmn; i++) */
675  return (0);
676 }
677 
678 /* PROCEDURE : qrsolv
679  *
680  * ENTREE :
681  * n Ordre de la matrice "r".
682  * r Matrice de taille "n" x "n". En entree, la partie triangulaire
683  * complete de "r" doit contenir la partie triangulaire superieure
684  * complete de "r".
685  * ldr Taille maximale de la matrice "r". "ldr" >= n.
686  * ipvt Vecteur de taille "n" definissant la matrice de permutation "p"
687  * La jeme colonne de de "p" est la colonne ipvt[j] de la matrice
688  * identite.
689  * diag Vecteur de taille "n" contenant les elements diagonaux de la
690  * matrice "d".
691  * qtb Vecteur de taille "n" contenant les "n" premiers elements du
692  * vecteur (q transpose) * b.
693  * wa Vecteur de travail de taille "n".
694  *
695  * DESCRIPTION :
696  * La procedure complete la solution du probleme, si on fournit les information
697  * necessaires de la factorisation qr, avec pivotage des colonnes.
698  * Soit une matrice "a" de taille "m" x "n" donnee, une matrice diagonale "d" de
699  * taille "n" x "n" et un vecteur "b" de taille "m", le probleme est la determination
700  * un vecteur "x" qui est solution du systeme
701  *
702  * a*x = b , d*x = 0 ,
703  *
704  * Au sens des moindres carres.
705  *
706  * Soit a * p = q * r, ou p est une matrice de permutation, les colonnes de "q"
707  * sont orthogonales et "r" est une matrice traingulaire superieure dont les
708  * elements diagonaux sont classes de l'ordre decroissant de leur valeur. Cette
709  * procedure attend donc la matrice triangulaire superieure remplie "r", la
710  * matrice de permutaion "p" et les "n" premiers elements de (q transpose) * b.
711  * Le systeme
712  *
713  * a * x = b, d * x = 0, est alors equivalent a
714  *
715  * t t
716  * r * z = q * b , p * d * p * z = 0 ,
717  *
718  * Ou x = p * z. Si ce systeme ne possede pas de rangee pleine, alors une
719  * solution au moindre carre est obtenue. En sortie, la procedure fournit aussi
720  * une matrice triangulaire superieure "s" tel que
721  *
722  * t t t
723  * p * (a * a + d * d) * p = s * s .
724  *
725  * "s" est calculee a l'interieure de qrsolv et peut etre hors interet.
726  *
727  * SORTIE :
728  * r En sortie, le triangle superieur n'est pas altere, et la partie
729  * triangulaire inferieure contient la partie triangulaire superieure
730  * (transpose) de la matrice triangulaire "s".
731  * x Vecteur de taille "n" contenant les solutions au moindres carres du
732  * systeme a * x = b, d * x = 0.
733  * sdiag Vecteur de taille "n" contenant les elements diagonaux de la
734  * matrice triangulaire superieure "s".
735  *
736  */
737 int qrsolv (int n, double *r, int ldr, int *ipvt, double *diag, double *qtb,
738  double *x, double *sdiag, double *wa)
739 {
740  int i, j, k, kp1, l; /* compteur de boucle */
741  int nsing;
742  double cosi, cotg, qtbpj, sinu, tg, temp;
743 
744  /*
745  * copie de r et (q transpose) * b afin de preserver l'entree
746  * et initialisation de "s". En particulier, sauvegarde des elements
747  * diagonaux de "r" dans "x".
748  */
749  for (i = 0; i < n; i++)
750  {
751  for (j = i; j < n; j++)
752  *MIJ(r, i, j, ldr) = *MIJ(r, j, i, ldr);
753 
754  x[i] = *MIJ(r, i, i, ldr);
755  wa[i] = qtb[i];
756  }
757 
758  /*
759  * Elimination de la matrice diagonale "d" en utlisant une rotation
760  * connue.
761  */
762 
763  for (i = 0; i < n; i++)
764  {
765  /*
766  * preparation de la colonne de d a eliminer, reperage de
767  * l'element diagonal par utilisation de p de la
768  * factorisation qr.
769  */
770  l = ipvt[i];
771 
772  //if (diag[l] != 0.0)
773  if (std::fabs(diag[l]) > std::numeric_limits<double>::epsilon())
774  {
775  for (k = i; k < n; k++)
776  sdiag[k] = 0.0;
777 
778  sdiag[i] = diag[l];
779 
780  /*
781  * Les transformations qui eliminent la colonne de d
782  * modifient seulement qu'un seul element de
783  * (q transpose)*b avant les n premiers elements
784  * lesquels sont inialement nuls.
785  */
786 
787  qtbpj = 0.0;
788 
789  for (k = i; k < n; k++)
790  {
791  /*
792  * determination d'une rotation qui elimine
793  * les elements appropriees dans la colonne
794  * courante de d.
795  */
796 
797  //if (sdiag[k] != 0.0)
798  if (std::fabs(sdiag[k]) > std::numeric_limits<double>::epsilon())
799  {
800  if (std::fabs(*MIJ(r, k, k, ldr)) >= std::fabs(sdiag[k]))
801  {
802  tg = sdiag[k] / *MIJ(r, k, k, ldr);
803  cosi = 0.5 / sqrt(0.25 + 0.25 * (tg * tg));
804  sinu = cosi * tg;
805  }
806  else
807  {
808  cotg = *MIJ(r, k, k, ldr) / sdiag[k];
809  sinu = 0.5 / sqrt(0.25 + 0.25 * (cotg * cotg));
810  cosi = sinu * cotg;
811  }
812 
813  /*
814  * calcul des elements de la diagonale modifiee
815  * de r et des elements modifies de
816  * ((q transpose)*b,0).
817  */
818  *MIJ(r, k, k, ldr) = cosi * *MIJ(r, k, k, ldr) + sinu * sdiag[k];
819  temp = cosi * wa[k] + sinu * qtbpj;
820  qtbpj = -sinu * wa[k] + cosi * qtbpj;
821  wa[k] = temp;
822 
823  /*
824  * accumulation des tranformations dans
825  * les lignes de s.
826  */
827 
828  kp1 = k + 1;
829 
830  if ( n >= kp1)
831  {
832  for (j = kp1; j < n; j++)
833  {
834  temp = cosi * *MIJ(r, k, j, ldr) +
835  sinu * sdiag[j];
836  sdiag[j] = - sinu * *MIJ(r, k, j, ldr) +
837  cosi * sdiag[j];
838  *MIJ(r, k, j, ldr) = temp;
839  }
840  }
841  }/* fin if diag[] !=0 */
842  } /* fin boucle for k -> n */
843  }/* fin if diag =0 */
844 
845  /*
846  * stokage de l'element diagonal de s et restauration de
847  * l'element diagonal correspondant de r.
848  */
849  sdiag[i] = *MIJ(r, i, i, ldr);
850  *MIJ(r, i, i, ldr) = x[i];
851  } /* fin boucle for j -> n */
852 
853  /*
854  * resolution du systeme triangulaire pour z. Si le systeme est
855  * singulier, on obtient une solution au moindres carres.
856  */
857  nsing = n;
858 
859  for (i = 0; i < n; i++)
860  {
861  //if (sdiag[i] == 0.0 && nsing == n)
862  if ( (std::fabs(sdiag[i]) <= std::numeric_limits<double>::epsilon()) && nsing == n)
863  nsing = i - 1;
864 
865  if (nsing < n)
866  wa[i] = 0.0;
867  }
868 
869  if (nsing >= 0)
870  {
871  for (k = 0; k < nsing; k++)
872  {
873  i = nsing - 1 - k;
874  double sum = 0.0;
875  int jp1 = i + 1;
876 
877  if (nsing >= jp1)
878  {
879  for (j = jp1; j < nsing; j++)
880  sum += *MIJ(r, i, j, ldr) * wa[j];
881  }
882  wa[i] = (wa[i] - sum) / sdiag[i];
883  }
884  }
885  /*
886  * permutation arriere des composants de z et des componants de x.
887  */
888 
889  for (j = 0; j < n; j++)
890  {
891  l = ipvt[j];
892  x[l] = wa[j];
893 
894  }
895  return (0);
896 }
897 
898 /*
899  * PROCEDURE : lmder
900  *
901  *
902  * ENTREE :
903  * fcn Fonction qui calcule la fonction et le jacobien de la fonction.
904  * m Nombre de fonctions.
905  * n Nombre de variables. n <= m
906  * x Vecteur de taille "n" contenant en entree une estimation
907  * initiale de la solution.
908  * ldfjac Taille dominante de la matrice "fjac". ldfjac >= "m".
909  * ftol Erreur relative desiree dans la somme des carre. La terminaison
910  * survient quand les preductions estimee et vraie de la somme des
911  * carres sont toutes deux au moins egal a ftol.
912  * xtol Erreur relative desiree dans la solution approximee. La
913  * terminaison survient quand l'erreur relative entre deux
914  * iterations consecutives est au moins egal a xtol.
915  * gtol Mesure de l'orthogonalite entre le vecteur des fonctions et les
916  * colonnes du jacobien. La terminaison survient quand le cosinus
917  * de l'angle entre fvec et n'importe quelle colonne du jacobien
918  * est au moins egal a gtol, en valeur absolue.
919  * maxfev Nombre d'appel maximum. La terminaison se produit lorsque le
920  * nombre d'appel a fcn avec iflag = 1 a atteint "maxfev".
921  * diag Vecteur de taille "n". Si mode = 1 (voir ci-apres), diag est
922  * initialisee en interne. Si mode = 2, diag doit contenir les
923  * entree positives qui servent de facteurs d'echelle aux variables.
924  * mode Si mode = 1, les variables seront mis a l'echelle en interne.
925  * Si mode = 2, la mise a l'echelle est specifie par l'entree diag.
926  * Les autres valeurs de mode sont equivalents a mode = 1.
927  * factor Definit la limite de l'etape initial. Cette limite est initialise
928  * au produit de "factor" et de la norme euclidienne de "diag" * "x"
929  * sinon nul. ou a "factor" lui meme. Dans la plupart des cas,
930  * "factor" doit se trouve dans l'intervalle (1, 100); ou 100 est
931  * la valeur recommandee.
932  * nprint Controle de l'impression des iterees (si valeur positive).
933  * Dans ce cas, fcn est appelle avec iflag = 0 au debut de la
934  * premiere iteration et apres chaque nprint iteration, x, fvec,
935  * et fjac sont disponible pour impression, cela avant de quitter
936  * la procedure. Si "nprint" est negatif, aucun appel special de
937  * fcn est faite.
938  * wa1, wa2, wa3 Vecteur de travail de taille "n".
939  * wa4 Vecteur de travail de taille "m".
940  *
941  *
942  * SORTIE :
943  * x Vecteur de taille "n" contenant en sortie l'estimee finale
944  * de la solution.
945  * fvec Vecteur de taille "m" contenant les fonctions evaluee en "x".
946  * fjac Matrice de taille "m" x "n". La sous matrice superieure de
947  * taille "n" x "n" de fjac contient une matrice triangulaire
948  * superieure r dont les elements diagonaux, classe dans le sens
949  * decroissant de leur valeur, sont de la forme :
950  *
951  * T T T
952  * p * (jac * jac) * p = r * r
953  *
954  * Ou p est une matrice de permutation et jac est le jacobien
955  * final calcule.
956  * La colonne j de p est la colonne ipvt (j) (voir ci apres) de
957  * la matrice identite. La partie trapesoidale inferieure de fjac
958  * contient les information genere durant le calcul de r.
959  * info Information de l'execution de la procedure. Lorsque la procedure
960  * a termine son execution, "info" est inialisee a la valeur
961  * (negative) de iflag. sinon elle prend les valeurs suivantes :
962  * info = 0 : parametres en entree non valides.
963  * info = 1 : les reductions relatives reelle et estimee de la
964  * somme des carres sont au moins egales a ftol.
965  * info = 2 : erreur relative entre deux iteres consecutives sont
966  * egaux a xtol.
967  * info = 3 : conditions info = 1 et info = 2 tous deux requis.
968  * info = 4 : le cosinus de l'angle entre fvec et n'importe quelle
969  * colonne du jacobien est au moins egal a gtol, en
970  * valeur absolue.
971  * info = 5 : nombre d'appels a fcn avec iflag = 1 a atteint
972  * maxfev.
973  * info = 6 : ftol est trop petit. Plus moyen de reduire de la
974  * somme des carres.
975  * info = 7 : xtol est trop petit. Plus d'amelioration possible
976  * pour approximer la solution x.
977  * info = 8 : gtol est trop petit. "fvec" est orthogonal aux
978  * colonnes du jacobien a la precision machine pres.
979  * nfev Nombre d'appel a "fcn" avec iflag = 1.
980  * njev Nombre d'appel a "fcn" avec iflag = 2.
981  * ipvt Vecteur de taille "n". Il definit une matrice de permutation p
982  * tel que jac * p = q * p, ou jac est le jacbien final calcule,
983  * q est orthogonal (non socke) et r est triangulaire superieur,
984  * avec les elements diagonaux classes en ordre decroissant de
985  * leur valeur. La colonne j de p est ipvt[j] de la matrice identite.
986  * qtf Vecteur de taille n contenant les n premiers elements du
987  * vecteur qT * fvec.
988  *
989  * DESCRIPTION :
990  * La procedure minimize la somme de carre de m equation non lineaire a n
991  * variables par une modification de l'algorithme de Levenberg - Marquardt.
992  *
993  * REMARQUE :
994  * L'utilisateur doit fournir une procedure "fcn" qui calcule la fonction et
995  * le jacobien.
996  * "fcn" doit etre declare dans une instruction externe a la procedure et doit
997  * etre appele comme suit :
998  * fcn (int m, int n, int ldfjac, double *x, double *fvec, double *fjac, int *iflag)
999  *
1000  * si iflag = 1 calcul de la fonction en x et retour de ce vecteur dans fvec.
1001  * fjac n'est pas modifie.
1002  * si iflag = 2 calcul du jacobien en x et retour de cette matrice dans fjac.
1003  * fvec n'est pas modifie.
1004  *
1005  * RETOUR :
1006  * En cas de succes, la valeur zero est retournee.
1007  * Sinon la valeur -1 est retournee.
1008  */
1009 int lmder (void (*ptr_fcn)(int m, int n, double *xc, double *fvecc,
1010  double *jac, int ldfjac, int iflag), int m, int n, double *x,
1011  double *fvec, double *fjac, int ldfjac, double ftol, double xtol,
1012  double gtol, unsigned int maxfev, double *diag, int mode,
1013  const double factor, int nprint, int *info, unsigned int *nfev,
1014  int *njev, int *ipvt, double *qtf, double *wa1, double *wa2,
1015  double *wa3, double *wa4)
1016 {
1017  const double tol1 = 0.1, tol5 = 0.5, tol25 = 0.25, tol75 = 0.75, tol0001 = 0.0001;
1018  int oncol = TRUE;
1019  int iflag, iter;
1020  int count = 0;
1021  int i, j, l;
1022  double actred, delta, dirder, epsmch, fnorm, fnorm1;
1023  double ratio = std::numeric_limits<double>::epsilon();
1024  double par, pnorm, prered;
1025  double sum, temp, temp1, temp2, xnorm = 0.0;
1026 
1027  /* epsmch est la precision machine. */
1028  epsmch = std::numeric_limits<double>::epsilon();
1029 
1030  *info = 0;
1031  iflag = 0;
1032  *nfev = 0;
1033  *njev = 0;
1034 
1035  /* verification des parametres d'entree. */
1036 
1037  /*if (n <= 0)
1038  return 0;*/
1039  if (m < n)
1040  return 0;
1041  if (ldfjac < m)
1042  return 0;
1043  if (ftol < 0.0)
1044  return 0;
1045  if (xtol < 0.0)
1046  return 0;
1047  if (gtol < 0.0)
1048  return 0;
1049  if (maxfev == 0)
1050  return 0;
1051  if (factor <= 0.0)
1052  return 0;
1053  if ((n <= 0) || (m < n) || (ldfjac < m) || (ftol < 0.0) || (xtol < 0.0)
1054  || (gtol < 0.0) || (maxfev == 0) || (factor <= 0.0))
1055  {
1056  /*
1057  * termination, normal ou imposee par l'utilisateur.
1058  */
1059  if (iflag < 0)
1060  *info = iflag;
1061 
1062  iflag = 0;
1063 
1064  if (nprint > 0)
1065  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1066 
1067  return (count);
1068  }
1069 
1070  if (mode == 2)
1071  {
1072  for (j = 0; j < n; j++)
1073  {
1074  if (diag[j] <= 0.0)
1075  {
1076  /*
1077  * termination, normal ou imposee par l'utilisateur.
1078  */
1079  if (iflag < 0)
1080  *info = iflag;
1081 
1082  iflag = 0;
1083 
1084  if (nprint > 0)
1085  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1086 
1087  return (count);
1088  }
1089  }
1090  }
1091 
1092  /*
1093  * evaluation de la fonction au point de depart
1094  * et calcul de sa norme.
1095  */
1096  iflag = 1;
1097 
1098  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1099 
1100 
1101  *nfev = 1;
1102 
1103  if (iflag < 0)
1104  {
1105  /*
1106  * termination, normal ou imposee par l'utilisateur.
1107  */
1108  if (iflag < 0)
1109  *info = iflag;
1110 
1111  iflag = 0;
1112 
1113  if (nprint > 0)
1114  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1115 
1116  return (count);
1117  }
1118 
1119  fnorm = enorm(fvec, m);
1120 
1121  /*
1122  * initialisation du parametre de Levenberg-Marquardt
1123  * et du conteur d'iteration.
1124  */
1125 
1126  par = 0.0;
1127  iter = 1;
1128 
1129  /*
1130  * debut de la boucle la plus externe.
1131  */
1132  while (count < (int)maxfev)
1133  {
1134  count++;
1135  /*
1136  * calcul de la matrice jacobienne.
1137  */
1138 
1139  iflag = 2;
1140 
1141  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1142 
1143  (*njev) ++;
1144 
1145  if (iflag < 0)
1146  {
1147  /*
1148  * termination, normal ou imposee par l'utilisateur.
1149  */
1150  if (iflag < 0)
1151  *info = iflag;
1152 
1153  iflag = 0;
1154 
1155  if (nprint > 0)
1156  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1157 
1158  return (count);
1159  }
1160 
1161  /*
1162  * si demandee, appel de fcn pour impression des iterees.
1163  */
1164  if (nprint > 0)
1165  {
1166  iflag = 0;
1167  if ((iter-1) % nprint == 0)
1168  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1169 
1170  if (iflag < 0)
1171  {
1172  /*
1173  * termination, normal ou imposee par l'utilisateur.
1174  */
1175  if (iflag < 0)
1176  *info = iflag;
1177 
1178  iflag = 0;
1179 
1180  if (nprint > 0)
1181  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1182 
1183  return (count);
1184  }
1185  }
1186 
1187  /*
1188  * calcul de la factorisation qr du jacobien.
1189  */
1190  qrfac(n, m, fjac, ldfjac, &oncol, ipvt, n, wa1, wa2, wa3);
1191 
1192  /*
1193  * a la premiere iteration et si mode est 1, mise a l'echelle
1194  * en accord avec les normes des colonnes du jacobien initial.
1195  */
1196 
1197  if (iter == 1)
1198  {
1199  if (mode != 2)
1200  {
1201  for (j = 0; j < n; j++)
1202  {
1203  diag[j] = wa2[j];
1204  //if (wa2[j] == 0.0)
1205  if (std::fabs(wa2[j]) <= std::numeric_limits<double>::epsilon())
1206  diag[j] = 1.0;
1207  }
1208  }
1209 
1210  /*
1211  * a la premiere iteration, calcul de la norme de x mis
1212  * a l'echelle et initialisation de la limite delta de
1213  * l'etape.
1214  */
1215 
1216  for (j = 0; j < n; j++)
1217  wa3[j] = diag[j] * x[j];
1218 
1219  xnorm = enorm (wa3, n);
1220  delta = factor * xnorm;
1221 
1222  //if (delta == 0.0)
1223  if (std::fabs(delta) <= std::numeric_limits<double>::epsilon())
1224  delta = factor;
1225  }
1226 
1227  /*
1228  * formation de (q transpose) * fvec et stockage des n premiers
1229  * composants dans qtf.
1230  */
1231  for (i = 0; i < m; i++)
1232  wa4[i] = fvec[i];
1233 
1234  for (i = 0; i < n; i++)
1235  {
1236  double *pt = MIJ(fjac, i, i, ldfjac);
1237  //if (*MIJ(fjac, i, i, ldfjac) != 0.0)
1238  if (std::fabs(*pt) > std::numeric_limits<double>::epsilon() )
1239  {
1240  sum = 0.0;
1241 
1242  for (j = i; j < m; j++)
1243  sum += *MIJ(fjac, i, j, ldfjac) * wa4[j];
1244 
1245  temp = - sum / *MIJ(fjac, i, i, ldfjac);
1246 
1247  for (j = i; j < m; j++)
1248  wa4[j] += *MIJ(fjac, i, j, ldfjac) * temp;
1249 
1250  }
1251 
1252  *MIJ(fjac, i, i, ldfjac) = wa1[i];
1253  qtf[i] = wa4[i];
1254  }
1255 
1256  /*
1257  * calcul de la norme du gradient mis a l'echelle.
1258  */
1259 
1260  double gnorm = 0.0;
1261 
1262  //if (fnorm != 0.0)
1263  if (std::fabs(fnorm) > std::numeric_limits<double>::epsilon())
1264  {
1265  for (i = 0; i < n; i++)
1266  {
1267  l = ipvt[i];
1268  //if (wa2[l] != 0.0)
1269  if (std::fabs(wa2[l]) > std::numeric_limits<double>::epsilon())
1270  {
1271  sum = 0.0;
1272  for (j = 0; j <= i; j++)
1273  sum += *MIJ(fjac, i, j, ldfjac) * (qtf[j] / fnorm);
1274 
1275  gnorm = vpMath::maximum(gnorm, std::fabs(sum / wa2[l]));
1276  }
1277  }
1278  }
1279 
1280  /*
1281  * test pour la convergence de la norme du gradient .
1282  */
1283 
1284  if (gnorm <= gtol)
1285  *info = 4;
1286 
1287  if (*info != 0)
1288  {
1289  /*
1290  * termination, normal ou imposee par l'utilisateur.
1291  */
1292  if (iflag < 0)
1293  *info = iflag;
1294 
1295  iflag = 0;
1296 
1297  if (nprint > 0)
1298  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1299 
1300  return (count);
1301  }
1302 
1303  /*
1304  * remise a l'echelle si necessaire.
1305  */
1306 
1307  if (mode != 2)
1308  {
1309  for (j = 0; j < n; j++)
1310  diag[j] = vpMath::maximum(diag[j], wa2[j]);
1311  }
1312 
1313  /*
1314  * debut de la boucle la plus interne.
1315  */
1316  ratio = 0.0;
1317  while (ratio < tol0001)
1318  {
1319 
1320  /*
1321  * determination du parametre de Levenberg-Marquardt.
1322  */
1323  lmpar(n, fjac, ldfjac, ipvt, diag, qtf, &delta, &par, wa1,
1324  wa2, wa3, wa4);
1325 
1326  /*
1327  * stockage de la direction p et x + p. calcul de la norme de p.
1328  */
1329 
1330  for (j = 0; j < n; j++)
1331  {
1332  wa1[j] = - wa1[j];
1333  wa2[j] = x[j] + wa1[j];
1334  wa3[j] = diag[j] * wa1[j];
1335  }
1336 
1337  pnorm = enorm(wa3, n);
1338 
1339  /*
1340  * a la premiere iteration, ajustement de la premiere limite de
1341  * l'etape.
1342  */
1343 
1344  if (iter == 1)
1345  delta = vpMath::minimum(delta, pnorm);
1346 
1347  /*
1348  * evaluation de la fonction en x + p et calcul de leur norme.
1349  */
1350 
1351  iflag = 1;
1352  (*ptr_fcn)(m, n, wa2, wa4, fjac, ldfjac, iflag);
1353 
1354  (*nfev)++;
1355 
1356  if (iflag < 0)
1357  {
1358  // termination, normal ou imposee par l'utilisateur.
1359  if (iflag < 0)
1360  *info = iflag;
1361 
1362  iflag = 0;
1363 
1364  if (nprint > 0)
1365  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1366 
1367  return (count);
1368  }
1369 
1370  fnorm1 = enorm(wa4, m);
1371 
1372  /*
1373  * calcul de la reduction reelle mise a l'echelle.
1374  */
1375 
1376  actred = - 1.0;
1377 
1378  if ((tol1 * fnorm1) < fnorm)
1379  actred = 1.0 - ((fnorm1 / fnorm) * (fnorm1 / fnorm));
1380 
1381  /*
1382  * calcul de la reduction predite mise a l'echelle et
1383  * de la derivee directionnelle mise a l'echelle.
1384  */
1385 
1386  for (i = 0; i < n; i++)
1387  {
1388  wa3[i] = 0.0;
1389  l = ipvt[i];
1390  temp = wa1[l];
1391  for (j = 0; j <= i; j++)
1392  wa3[j] += *MIJ(fjac, i, j, ldfjac) * temp;
1393  }
1394 
1395  temp1 = enorm(wa3, n) / fnorm;
1396  temp2 = (sqrt(par) * pnorm) / fnorm;
1397  prered = (temp1 * temp1) + (temp2 * temp2) / tol5;
1398  dirder = - ((temp1 * temp1) + (temp2 * temp2));
1399 
1400  /*
1401  * calcul du rapport entre la reduction reel et predit.
1402  */
1403 
1404  ratio = 0.0;
1405 
1406  //if (prered != 0.0)
1407  if (std::fabs(prered) > std::numeric_limits<double>::epsilon())
1408  ratio = actred / prered;
1409 
1410  /*
1411  * mise a jour de la limite de l'etape.
1412  */
1413 
1414  if (ratio > tol25)
1415  {
1416  //if ((par == 0.0) || (ratio <= tol75))
1417  if ((std::fabs(par) <= std::numeric_limits<double>::epsilon()) || (ratio <= tol75))
1418  {
1419  delta = pnorm / tol5;
1420  par *= tol5;
1421  }
1422  }
1423  else
1424  {
1425  if (actred >= 0.0)
1426  temp = tol5;
1427 
1428  else
1429  temp = tol5 * dirder / (dirder + tol5 * actred);
1430 
1431  if ((tol1 * fnorm1 >= fnorm) || (temp < tol1))
1432  temp = tol1;
1433 
1434  delta = temp * vpMath::minimum(delta, (pnorm / tol1));
1435  par /= temp;
1436  }
1437 
1438  /*
1439  * test pour une iteration reussie.
1440  */
1441  if (ratio >= tol0001)
1442  {
1443  /*
1444  * iteration reussie. mise a jour de x, de fvec, et de
1445  * leurs normes.
1446  */
1447 
1448  for (j = 0; j < n; j++)
1449  {
1450  x[j] = wa2[j];
1451  wa2[j] = diag[j] * x[j];
1452  }
1453 
1454  for (i = 0; i < m; i++)
1455  fvec[i] = wa4[i];
1456 
1457  xnorm = enorm(wa2, n);
1458  fnorm = fnorm1;
1459  iter++;
1460  }
1461 
1462  /*
1463  * tests pour convergence.
1464  */
1465 
1466  if ((std::fabs(actred) <= ftol) && (prered <= ftol) && (tol5 * ratio <= 1.0))
1467  *info = 1;
1468 
1469  if (delta <= xtol * xnorm)
1470  *info = 2;
1471 
1472  if ((std::fabs(actred) <= ftol) && (prered <= ftol) && (tol5 * ratio <= 1.0)
1473  && *info == 2)
1474  *info = 3;
1475 
1476  if (*info != 0)
1477  {
1478  /*
1479  * termination, normal ou imposee par l'utilisateur.
1480  */
1481  if (iflag < 0)
1482  *info = iflag;
1483 
1484  iflag = 0;
1485 
1486  if (nprint > 0)
1487  (*ptr_fcn)(m,n,x,fvec,fjac,ldfjac, iflag);
1488 
1489  return (count);
1490  }
1491  /*
1492  * tests pour termination et
1493  * verification des tolerances.
1494  */
1495 
1496  if (*nfev >= maxfev)
1497  *info = 5;
1498 
1499  if ((std::fabs(actred) <= epsmch) && (prered <= epsmch)
1500  && (tol5 * ratio <= 1.0))
1501  *info = 6;
1502 
1503  if (delta <= epsmch * xnorm)
1504  *info = 7;
1505 
1506  if (gnorm <= epsmch)
1507  *info = 8;
1508 
1509  if (*info != 0)
1510  {
1511  /*
1512  * termination, normal ou imposee par l'utilisateur.
1513  */
1514  if (iflag < 0)
1515  *info = iflag;
1516 
1517  iflag = 0;
1518 
1519  if (nprint > 0)
1520  (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1521 
1522  return (count);
1523  }
1524  }/* fin while ratio >=tol0001 */
1525  }/*fin while 1*/
1526 
1527  return 0 ;
1528 }
1529 
1530 
1531 
1532 /*
1533  * PROCEDURE : lmder1
1534  *
1535  * ENTREE :
1536  * fcn Fonction qui calcule la fonction et le jacobien de la fonction.
1537  * m Nombre de fonctions.
1538  * n Nombre de variables (parametres). n <= m
1539  * x Vecteur de taille "n" contenant en entree une estimation
1540  * initiale de la solution.
1541  * ldfjac Taille maximale de la matrice "fjac". ldfjac >= "m".
1542  * tol Tolerance. La terminaison de la procedure survient quand
1543  * l'algorithme estime que l'erreur relative dans la somme des
1544  * carres est au moins egal a tol ou bien que l'erreur relative
1545  * entre x et la solution est au moins egal atol.
1546  * wa Vecteur de travail de taille "lwa".
1547  * lwa Taille du vecteur "wa". wa >= 5 * n + m.
1548  *
1549  *
1550  * SORTIE :
1551  * x Vecteur de taille "n" contenant en sortie l'estimee finale
1552  * de la solution.
1553  * fvec Vecteur de taille "m" contenant les fonctions evaluee en "x".
1554  * fjac Matrice de taille "m" x "n". La sous matrice superieure de
1555  * taille "n" x "n" de fjac contient une matrice triangulaire
1556  * superieure r dont les elements diagonaux, classe dans le sens
1557  * decroissant de leur valeur, sont de la forme :
1558  *
1559  * T T T
1560  * p * (jac * jac) * p = r * r
1561  *
1562  * Ou p est une matrice de permutation et jac est le jacobien
1563  * final calcule.
1564  * La colonne j de p est la colonne ipvt (j) (voir ci apres) de
1565  * la matrice identite. La partie trapesoidale inferieure de fjac
1566  * contient les information genere durant le calcul de r.
1567  * info Information de l'executionde la procedure. Lorsque la procedure
1568  * a termine son execution, "info" est inialisee a la valeur
1569  * (negative) de iflag. sinon elle prend les valeurs suivantes :
1570  * info = 0 : parametres en entre non valides.
1571  * info = 1 : estimation par l'algorithme que l'erreur relative
1572  * de la somme des carre est egal a tol.
1573  * info = 2 : estimation par l'algorithme que l'erreur relative
1574  * entre x et la solution est egal a tol.
1575  * info = 3 : conditions info = 1 et info = 2 tous deux requis.
1576  * info = 4 : fvec est orthogonal aux colonnes du jacobien.
1577  * info = 5 : nombre d'appels a fcn avec iflag = 1 a atteint
1578  * 100 * (n + 1).
1579  * info = 6 : tol est trop petit. Plus moyen de reduire de la
1580  * somme des carres.
1581  * info = 7 : tol est trop petit. Plus d'amelioration possible
1582  * d' approximer la solution x.
1583  * ipvt Vecteur de taille "n". Il definit une matrice de permutation p
1584  * tel que jac * p = q * p, ou jac est le jacbien final calcule,
1585  * q est orthogonal (non socke) et r est triangulaire superieur,
1586  * avec les elements diagonaux classes en ordre decroissant de
1587  * leur valeur. La colonne j de p est ipvt[j] de la matrice identite.
1588  *
1589  * DESCRIPTION :
1590  * La procedure minimize la somme de carre de m equation non lineaire a n
1591  * variables par une modification de l'algorithme de Levenberg - Marquardt.
1592  * Cette procedure appele la procedure generale au moindre carre lmder.
1593  *
1594  * REMARQUE :
1595  * L'utilisateur doit fournir une procedure "fcn" qui calcule la fonction et
1596  * le jacobien.
1597  * "fcn" doit etre declare dans une instruction externe a la procedure et doit
1598  * etre appele comme suit :
1599  * fcn (int m, int n, int ldfjac, double *x, double *fvec, double *fjac, int *iflag)
1600  *
1601  * si iflag = 1 calcul de la fonction en x et retour de ce vecteur dans fvec.
1602  * fjac n'est pas modifie.
1603  * si iflag = 2 calcul du jacobien en x et retour de cette matrice dans fjac.
1604  * fvec n'est pas modifie.
1605  *
1606  * RETOUR :
1607  * En cas de succes, la valeur zero est retournee.
1608  * Sinon, la valeur -1.
1609  *
1610  */
1611 int lmder1 (void (*ptr_fcn)(int m, int n, double *xc, double *fvecc,
1612  double *jac, int ldfjac, int iflag),
1613  int m, int n, double *x, double *fvec, double *fjac,
1614  int ldfjac, double tol, int *info, int *ipvt, int lwa, double *wa)
1615 {
1616  const double factor = 100.0;
1617  unsigned int maxfev, nfev;
1618  int mode, njev, nprint;
1619  double ftol, gtol, xtol;
1620 
1621  *info = 0;
1622 
1623  /* verification des parametres en entree qui causent des erreurs */
1624 
1625 
1626  if (/*(n <= 0) ||*/ (m < n) || (ldfjac < m) || (tol < 0.0) ||
1627  (lwa < (5 * n + m)) )
1628  {
1629  printf("%d %d %d %d \n", (m < n), (ldfjac < m), (tol < 0.0) , (lwa < (5 * n + m))) ;
1630  return (-1);
1631  }
1632 
1633  /* appel a lmder */
1634 
1635  maxfev = (unsigned int)(100 * (n + 1));
1636  ftol = tol;
1637  xtol = tol;
1638  gtol = 0.0;
1639  mode = 1;
1640  nprint = 0;
1641 
1642  lmder (ptr_fcn , m, n, x, fvec, fjac, ldfjac, ftol, xtol, gtol, maxfev, wa,
1643  mode, factor, nprint, info, &nfev, &njev, ipvt, &wa[n], &wa[2 * n],
1644  &wa[3 * n], &wa[4 * n], &wa[5 * n]);
1645 
1646  if (*info == 8)
1647  *info = 4;
1648 
1649  return (0);
1650 }
1651 
1652 #undef TRUE
1653 #undef FALSE
1654 
1655 
1656 /*
1657  * Local variables:
1658  * c-basic-offset: 2
1659  * End:
1660  */
static Type maximum(const Type &a, const Type &b)
Definition: vpMath.h:140
static Type minimum(const Type &a, const Type &b)
Definition: vpMath.h:151