Commit 20ce530d4f1b4017a2f49a5d325d46fed9ece341

Authored by wdaniau
1 parent e9d7fe24b4

Wrote documentation in source for fvn_lm

git-svn-id: https://lxsd.femto-st.fr/svn/fvn@78 b657c933-2333-4658-acf2-d3c7c2708721

Showing 1 changed file with 122 additions and 141 deletions Side-by-side Diff

fvn_misc/fvn_misc.f90
... ... @@ -7,6 +7,11 @@
7 7 module procedure fvn_z_muller
8 8 end interface fvn_muller
9 9  
  10 +public
  11 +private :: lmdif,lmpar,qrfac,qrsolv,enorm,fdjac2
  12 +! These are made private to not interfere with a
  13 +! possibly linked minpack
  14 +
10 15 contains
11 16 !
12 17 ! Muller
13 18  
14 19  
15 20  
16 21  
... ... @@ -320,31 +325,94 @@
320 325 ! This uses MINPACK Routines (http://www.netlib.org/minpack)
321 326 ! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
322 327 !
  328 +! The purpose of fvn_lm is to minimize the sum of the squares of m nonlinear
  329 +! functions in n variables by a modification of the Levenberg-Marquardt
  330 +! algorithm. This is done by using the more general least-squares
  331 +! solver lmdif. The user must provide a subroutine which calculates the
  332 +! functions. The jacobian is then calculated by a forward-difference
  333 +! approximation.
  334 +!
  335 +! call fvn_lm(fcn,m,n,a,info,tol)
  336 +!
  337 +! fcn : fcn is the user-supplied subroutine which calculates
  338 +! the functions. fcn must follow the following interface that must
  339 +! be declared in the calling subroutine :
  340 +!
  341 +! interface
  342 +! subroutine fcn(m,n,a,fvec,iflag)
  343 +! use fvn_common
  344 +! integer(ip_kind), intent(in) :: m
  345 +! integer(ip_kind), intent(in) :: n
  346 +! real(dp_kind), dimension(:), intent(in) :: a
  347 +! real(dp_kind), dimension(:), intent(inout) :: fvec
  348 +! integer(ip_kind), intent(inout) :: iflag
  349 +! end subroutine
  350 +! end interface
  351 +!
  352 +! This is the function which calculate the differences for which which square sum
  353 +! will be minimized outputing this difference in vector fvec.
  354 +! Parameters of fcn are as follows :
  355 +! m : positive integer input variable set to the number of functions
  356 +! (number of measurement points)
  357 +! n : positive integer input variable set to the number of variables
  358 +! (number of parameters in the function to fit)
  359 +! a : vector of length n containing parameters for which fcn should
  360 +! perform the calculation
  361 +! fvec : vector of length m containing the resulting evaluation
  362 +! iflag : integer normally not used, can be used to exit the
  363 +! the algorithm by setting it to a negative value
  364 +!
  365 +! m : positive integer input variable set to the number of functions
  366 +! (number of measurement points)
  367 +! n : positive integer input variable set to the number of variables
  368 +! (number of parameters in the function to fit)
  369 +! a : vector of length n, on input must contains an initial guess (or unity vector)
  370 +! and on output the solution vector
  371 +! info : is an output positive integer
  372 +! info = 0 improper input parameters.
  373 +! info = 1 algorithm estimates that the relative error
  374 +! in the sum of squares is at most tol.
  375 +! info = 2 algorithm estimates that the relative error
  376 +! between x and the solution is at most tol.
  377 +! info = 3 conditions for info = 1 and info = 2 both hold.
  378 +! info = 4 fvec is orthogonal to the columns of the
  379 +! jacobian to machine precision.
  380 +! info = 5 number of calls to fcn has reached or exceeded 200*(n+1).
  381 +! info = 6 tol is too small. no further reduction in
  382 +! the sum of squares is possible.
  383 +! info = 7 tol is too small. No further improvement in
  384 +! the approximate solution x is possible.
  385 +! tol : is an optional positive value. Termination occurs when the
  386 +! algorithm estimates either that the relative error in the sum of
  387 +! squares is at most tol or that the relative error between x and the
  388 +! solution is at most tol. If not provided default value is :
  389 +! sqrt(epsilon(0.0d0))
  390 +!
  391 +subroutine fvn_lm(fcn,m,n,a,info,tol)
  392 + integer(ip_kind), intent(in) :: m
  393 + integer(ip_kind), intent(in) :: n
  394 + real(dp_kind), dimension(:), intent(inout) :: a
  395 + integer(ip_kind), intent(out) :: info
  396 + real(dp_kind), intent(in), optional :: tol
323 397  
324   -subroutine fvn_lmdif(zef,m,n,a,info,tol)
325   - integer(4), intent(in) :: m
326   - integer(4), intent(in) :: n
327   - real(8), dimension(:), intent(inout) :: a
328   - integer(4), intent(out) :: info
329   - real(8), intent(in), optional :: tol
  398 + real(dp_kind) :: rtol
  399 + real(dp_kind), dimension(:), allocatable :: fvec
  400 + integer(ip_kind), dimension(:), allocatable :: iwa
330 401  
331   - real(8) :: rtol
332   - real(8), dimension(:), allocatable :: fvec
333   - integer(4), dimension(:), allocatable :: iwa
334   -
335 402 interface
336   - subroutine zef(m,n,a,fvec,iflag)
337   - integer(4), intent(in) :: m
338   - integer(4), intent(in) :: n
339   - real(8), dimension(:), intent(in) :: a
340   - real(8), dimension(:), intent(inout) :: fvec
341   - integer(4), intent(inout) :: iflag
  403 + subroutine fcn(m,n,a,fvec,iflag)
  404 + use fvn_common
  405 + integer(ip_kind), intent(in) :: m
  406 + integer(ip_kind), intent(in) :: n
  407 + real(dp_kind), dimension(:), intent(in) :: a
  408 + real(dp_kind), dimension(:), intent(inout) :: fvec
  409 + integer(ip_kind), intent(inout) :: iflag
342 410 end subroutine
343 411 end interface
344 412  
345   - integer(4) :: maxfev, mode, nfev, nprint
346   - real(8) :: epsfcn, ftol, gtol, xtol, fjac(m,n)
347   - real(8), parameter :: factor = 100._8, zero = 0.0_8
  413 + integer(ip_kind) :: maxfev, mode, nfev, nprint
  414 + real(dp_kind) :: epsfcn, ftol, gtol, xtol, fjac(m,n)
  415 + real(dp_kind), parameter :: factor = 100._8, zero = 0.0_8
348 416  
349 417 allocate(fvec(m),iwa(n))
350 418  
... ... @@ -367,7 +435,7 @@
367 435 mode = 1
368 436 nprint = 0
369 437  
370   - call lmdif(zef, m, n, a, fvec, ftol, xtol, gtol, maxfev, epsfcn, &
  438 + call lmdif(fcn, m, n, a, fvec, ftol, xtol, gtol, maxfev, epsfcn, &
371 439 mode, factor, nprint, info, nfev, fjac, iwa)
372 440  
373 441 if (info == 8) info = 4
374 442  
... ... @@ -376,113 +444,14 @@
376 444 end subroutine
377 445  
378 446  
379   -! MINPACK routines which are used by both LMDIF & LMDER
380   -! 25 October 2001:
381   -! Changed INTENT of iflag in several places to IN OUT.
382   -! Changed INTENT of fvec to IN OUT in user routine FCN.
383   -! Removed arguments diag and qtv from LMDIF & LMDER.
384   -! Replaced several DO loops with array operations.
385   -! amiller @ bigpond.net.au
386   -
387   -
388   -! **********
389   -
390   -! subroutine lmdif1
391   -
392   -! The purpose of lmdif1 is to minimize the sum of the squares of m nonlinear
393   -! functions in n variables by a modification of the Levenberg-Marquardt
394   -! algorithm. This is done by using the more general least-squares
395   -! solver lmdif. The user must provide a subroutine which calculates the
396   -! functions. The jacobian is then calculated by a forward-difference
397   -! approximation.
398   -
399   -! the subroutine statement is
400   -
401   -! subroutine lmdif1(fcn, m, n, x, fvec, tol, info, iwa)
402   -
403   -! where
404   -
405   -! fcn is the name of the user-supplied subroutine which calculates
406   -! the functions. fcn must be declared in an external statement in the
407   -! user calling program, and should be written as follows.
408   -
409   -! subroutine fcn(m, n, x, fvec, iflag)
410   -! integer m, n, iflag
411   -! REAL (dp) x(n), fvec(m)
412   -! ----------
413   -! calculate the functions at x and return this vector in fvec.
414   -! ----------
415   -! return
416   -! end
417   -
418   -! the value of iflag should not be changed by fcn unless
419   -! the user wants to terminate execution of lmdif1.
420   -! In this case set iflag to a negative integer.
421   -
422   -! m is a positive integer input variable set to the number of functions.
423   -
424   -! n is a positive integer input variable set to the number of variables.
425   -! n must not exceed m.
426   -
427   -! x is an array of length n. On input x must contain an initial estimate
428   -! of the solution vector. On output x contains the final estimate of
429   -! the solution vector.
430   -
431   -! fvec is an output array of length m which contains
432   -! the functions evaluated at the output x.
433   -
434   -! tol is a nonnegative input variable. Termination occurs when the
435   -! algorithm estimates either that the relative error in the sum of
436   -! squares is at most tol or that the relative error between x and the
437   -! solution is at most tol.
438   -
439   -! info is an integer output variable. If the user has terminated execution,
440   -! info is set to the (negative) value of iflag. See description of fcn.
441   -! Otherwise, info is set as follows.
442   -
443   -! info = 0 improper input parameters.
444   -
445   -! info = 1 algorithm estimates that the relative error
446   -! in the sum of squares is at most tol.
447   -
448   -! info = 2 algorithm estimates that the relative error
449   -! between x and the solution is at most tol.
450   -
451   -! info = 3 conditions for info = 1 and info = 2 both hold.
452   -
453   -! info = 4 fvec is orthogonal to the columns of the
454   -! jacobian to machine precision.
455   -
456   -! info = 5 number of calls to fcn has reached or exceeded 200*(n+1).
457   -
458   -! info = 6 tol is too small. no further reduction in
459   -! the sum of squares is possible.
460   -
461   -! info = 7 tol is too small. No further improvement in
462   -! the approximate solution x is possible.
463   -
464   -! iwa is an integer work array of length n.
465   -
466   -! wa is a work array of length lwa.
467   -
468   -! lwa is a positive integer input variable not less than m*n+5*n+m.
469   -
470   -! subprograms called
471   -
472   -! user-supplied ...... fcn
473   -
474   -! minpack-supplied ... lmdif
475   -
476   -! argonne national laboratory. minpack project. march 1980.
477   -! burton s. garbow, kenneth e. hillstrom, jorge j. more
478   -
479   -! **********
480   -
  447 +! lmdif
  448 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  449 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
481 450 SUBROUTINE lmdif(fcn, m, n, x, fvec, ftol, xtol, gtol, maxfev, epsfcn, &
482 451 mode, factor, nprint, info, nfev, fjac, ipvt)
483 452  
484 453 ! N.B. Arguments LDFJAC, DIAG, QTF, WA1, WA2, WA3 & WA4 have been removed.
485   -INTEGER, PARAMETER :: dp = 8
  454 +INTEGER, PARAMETER :: dp = dp_kind
486 455 INTEGER, INTENT(IN) :: m
487 456 INTEGER, INTENT(IN) :: n
488 457 REAL (dp), INTENT(IN OUT) :: x(:)
... ... @@ -504,10 +473,11 @@
504 473  
505 474 INTERFACE
506 475 SUBROUTINE fcn(m, n, x, fvec, iflag)
507   - INTEGER(4), INTENT(IN) :: m, n
508   - REAL (8), INTENT(IN) :: x(:)
509   - REAL (8), INTENT(IN OUT) :: fvec(:)
510   - INTEGER(4), INTENT(IN OUT) :: iflag
  476 + use fvn_common
  477 + INTEGER(ip_kind), INTENT(IN) :: m, n
  478 + REAL (dp_kind), INTENT(IN) :: x(:)
  479 + REAL (dp_kind), INTENT(IN OUT) :: fvec(:)
  480 + INTEGER(ip_kind), INTENT(IN OUT) :: iflag
511 481 END SUBROUTINE fcn
512 482 END INTERFACE
513 483  
... ... @@ -912,7 +882,9 @@
912 882  
913 883 ! **********
914 884  
915   -
  885 +! lmpar
  886 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  887 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
916 888 SUBROUTINE lmpar(n, r, ipvt, diag, qtb, delta, par, x, sdiag)
917 889  
918 890 ! Code converted using TO_F90 by Alan Miller
... ... @@ -920,7 +892,7 @@
920 892  
921 893 ! N.B. Arguments LDR, WA1 & WA2 have been removed.
922 894  
923   -INTEGER, PARAMETER :: dp = 8
  895 +INTEGER, PARAMETER :: dp = dp_kind
924 896 INTEGER, INTENT(IN) :: n
925 897 REAL (dp), INTENT(IN OUT) :: r(:,:)
926 898 INTEGER, INTENT(IN) :: ipvt(:)
927 899  
... ... @@ -1158,14 +1130,16 @@
1158 1130 END SUBROUTINE lmpar
1159 1131  
1160 1132  
1161   -
  1133 +! qrfac
  1134 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  1135 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
1162 1136 SUBROUTINE qrfac(m, n, a, pivot, ipvt, rdiag, acnorm)
1163 1137  
1164 1138 ! Code converted using TO_F90 by Alan Miller
1165 1139 ! Date: 1999-12-09 Time: 12:46:17
1166 1140  
1167 1141 ! N.B. Arguments LDA, LIPVT & WA have been removed.
1168   -INTEGER, PARAMETER :: dp = 8
  1142 +INTEGER, PARAMETER :: dp = dp_kind
1169 1143 INTEGER, INTENT(IN) :: m
1170 1144 INTEGER, INTENT(IN) :: n
1171 1145 REAL (dp), INTENT(IN OUT) :: a(:,:)
1172 1146  
... ... @@ -1320,11 +1294,13 @@
1320 1294 END SUBROUTINE qrfac
1321 1295  
1322 1296  
1323   -
  1297 +! qrsolv
  1298 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  1299 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
1324 1300 SUBROUTINE qrsolv(n, r, ipvt, diag, qtb, x, sdiag)
1325 1301  
1326 1302 ! N.B. Arguments LDR & WA have been removed.
1327   -INTEGER, PARAMETER :: dp = 8
  1303 +INTEGER, PARAMETER :: dp = dp_kind
1328 1304 INTEGER, INTENT(IN) :: n
1329 1305 REAL (dp), INTENT(IN OUT) :: r(:,:)
1330 1306 INTEGER, INTENT(IN) :: ipvt(:)
1331 1307  
... ... @@ -1507,12 +1483,14 @@
1507 1483 END SUBROUTINE qrsolv
1508 1484  
1509 1485  
1510   -
  1486 +! enorm
  1487 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  1488 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
1511 1489 FUNCTION enorm(n,x) RESULT(fn_val)
1512 1490  
1513 1491 ! Code converted using TO_F90 by Alan Miller
1514 1492 ! Date: 1999-12-09 Time: 12:45:34
1515   -INTEGER, PARAMETER :: dp = 8
  1493 +INTEGER, PARAMETER :: dp = dp_kind
1516 1494 INTEGER, INTENT(IN) :: n
1517 1495 REAL (dp), INTENT(IN) :: x(:)
1518 1496 REAL (dp) :: fn_val
1519 1497  
... ... @@ -1615,14 +1593,16 @@
1615 1593 END FUNCTION enorm
1616 1594  
1617 1595  
1618   -
  1596 +! fdjac2
  1597 +! MINPACK Subroutine (http://www.netlib.org/minpack)
  1598 +! Converted to fortran90 by Alan Miller amiller @ bigpond.net.au
1619 1599 SUBROUTINE fdjac2(fcn, m, n, x, fvec, fjac, iflag, epsfcn)
1620 1600  
1621 1601 ! Code converted using TO_F90 by Alan Miller
1622 1602 ! Date: 1999-12-09 Time: 12:45:44
1623 1603  
1624 1604 ! N.B. Arguments LDFJAC & WA have been removed.
1625   -INTEGER, PARAMETER :: dp = 8
  1605 +INTEGER, PARAMETER :: dp = dp_kind
1626 1606 INTEGER, INTENT(IN) :: m
1627 1607 INTEGER, INTENT(IN) :: n
1628 1608 REAL (dp), INTENT(IN OUT) :: x(n)
... ... @@ -1633,10 +1613,11 @@
1633 1613  
1634 1614 INTERFACE
1635 1615 SUBROUTINE fcn(m, n, x, fvec, iflag)
1636   - INTEGER(4), INTENT(IN) :: m, n
1637   - REAL (8), INTENT(IN) :: x(:)
1638   - REAL (8), INTENT(IN OUT) :: fvec(:)
1639   - INTEGER(4), INTENT(IN OUT) :: iflag
  1616 + use fvn_common
  1617 + INTEGER(ip_kind), INTENT(IN) :: m, n
  1618 + REAL (dp_kind), INTENT(IN) :: x(:)
  1619 + REAL (dp_kind), INTENT(IN OUT) :: fvec(:)
  1620 + INTEGER(ip_kind), INTENT(IN OUT) :: iflag
1640 1621 END SUBROUTINE fcn
1641 1622 END INTERFACE
1642 1623