Skip to content

Commit 69ce22f

Browse files
committed
Add example for lmdif1()
1 parent 63350e8 commit 69ce22f

File tree

2 files changed

+63
-0
lines changed

2 files changed

+63
-0
lines changed

examples/CMakeLists.txt

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -2,3 +2,6 @@ include_directories(${PROJECT_BINARY_DIR}/src)
22

33
add_executable(example_lmder1 example_lmder1.f90)
44
target_link_libraries(example_lmder1 minpack)
5+
6+
add_executable(example_lmdif1 example_lmdif1.f90)
7+
target_link_libraries(example_lmdif1 minpack)

examples/example_lmdif1.f90

Lines changed: 60 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,60 @@
1+
module testmod_dif
2+
implicit none
3+
private
4+
public fcn, dp
5+
6+
integer, parameter :: dp=kind(0d0)
7+
8+
contains
9+
10+
subroutine fcn(m, n, x, fvec, iflag)
11+
integer, intent(in) :: m, n, iflag
12+
real(dp), intent(in) :: x(n)
13+
real(dp), intent(out) :: fvec(m)
14+
15+
integer :: i
16+
real(dp) :: tmp1, tmp2, tmp3, y(15)
17+
! Suppress compiler warning:
18+
y(1) = iflag
19+
y = [1.4D-1, 1.8D-1, 2.2D-1, 2.5D-1, 2.9D-1, 3.2D-1, 3.5D-1, 3.9D-1, &
20+
3.7D-1, 5.8D-1, 7.3D-1, 9.6D-1, 1.34D0, 2.1D0, 4.39D0]
21+
22+
do i = 1, 15
23+
tmp1 = i
24+
tmp2 = 16 - i
25+
tmp3 = tmp1
26+
if (i .gt. 8) tmp3 = tmp2
27+
fvec(i) = y(i) - (x(1) + tmp1/(x(2)*tmp2 + x(3)*tmp3))
28+
end do
29+
end subroutine
30+
31+
end module
32+
33+
34+
program example_lmdif1
35+
use minpack, only: enorm, dpmpar, lmdif1
36+
use testmod_dif, only: dp, fcn
37+
implicit none
38+
39+
integer :: info, m, n
40+
real(dp) :: tol, x(3), fvec(15)
41+
integer :: iwa(size(x))
42+
real(dp), allocatable :: wa(:)
43+
44+
! The following starting values provide a rough fit.
45+
x = [1._dp, 1._dp, 1._dp]
46+
47+
! Set tol to the square root of the machine precision. Unless high precision
48+
! solutions are required, this is the recommended setting.
49+
tol = sqrt(dpmpar(1))
50+
51+
m = size(fvec)
52+
n = size(x)
53+
allocate(wa(m*n + 5*n + m))
54+
call lmdif1(fcn, size(fvec), size(x), x, fvec, tol, info, iwa, wa, size(wa))
55+
print 1000, enorm(size(fvec), fvec), info, x
56+
1000 format(5x, 'FINAL L2 NORM OF THE RESIDUALS', d15.7 // &
57+
5x, 'EXIT PARAMETER', 16x, i10 // &
58+
5x, 'FINAL APPROXIMATE SOLUTION' // &
59+
5x, 3d15.7)
60+
end program

0 commit comments

Comments
 (0)