File: example_fit_jacobian_benchmark.py

package info (click to toggle)
lmfit-py 1.3.3-4
  • links: PTS, VCS
  • area: main
  • in suites: sid, trixie
  • size: 2,332 kB
  • sloc: python: 13,071; makefile: 130; sh: 30
file content (178 lines) | stat: -rw-r--r-- 4,345 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
"""
Benchmarks of methods with and without computing the Jacobian analytically
==========================================================================

Providing a function that calculates the Jacobian matrix analytically can
reduce the time spent finding a solution. The results from benchmarks comparing
two methods (``leastsq`` and ``least_squares``) with and without a function to
calculate the Jacobian matrix analytically are presented below.

First we define the model function, the residual function, and the appropriate
Jacobian functions:
"""
from timeit import timeit
from types import SimpleNamespace

import matplotlib.pyplot as plt
import numpy as np

from lmfit import Parameters, minimize

NUM_JACOBIAN_CALLS = 0


def func(var, x):
    return var[0] * np.exp(-var[1]*x) + var[2]


def residual(pars, x, data):
    a, b, c = pars['a'], pars['b'], pars['c']
    model = func((a, b, c), x)
    return model - data


def dfunc(pars, x, data):
    global NUM_JACOBIAN_CALLS
    NUM_JACOBIAN_CALLS += 1

    a, b = pars['a'], pars['b']
    v = np.exp(-b*x)
    return np.array([v, -a*x*v, np.ones(len(x))])


def jacfunc(pars, x, data):
    global NUM_JACOBIAN_CALLS
    NUM_JACOBIAN_CALLS += 1

    a, b = pars['a'], pars['b']
    v = np.exp(-b*x)
    jac = np.ones((len(x), 3), dtype=np.float64)
    jac[:, 0] = v
    jac[:, 1] = -a * x * v
    return jac


a, b, c = 2.5, 1.3, 0.8

x = np.linspace(0, 4, 50)
y = func([a, b, c], x)

data = y + 0.15*np.random.RandomState(seed=2021).normal(size=x.size)


###############################################################################
# Then we define the different cases to benchmark (i.e., different methods with
# and without a function to calculate the Jacobian analytically) and the number
# of repetitions per case:
cases = (
    dict(
        method='leastsq',
    ),
    dict(
        method='leastsq',
        Dfun=dfunc,
        col_deriv=1,
    ),
    dict(
        method='least_squares',
    ),
    dict(
        method='least_squares',
        jac=jacfunc,
    ),
)

num_repeats = 100
results = []

for kwargs in cases:
    params = Parameters()
    params.add('a', value=10)
    params.add('b', value=10)
    params.add('c', value=10)

    wrapper = lambda: minimize(
        residual,
        params,
        args=(x,),
        kws={'data': data},
        **kwargs,
    )
    time = timeit(wrapper, number=num_repeats) / num_repeats

    NUM_JACOBIAN_CALLS = 0
    fit = wrapper()

    results.append(SimpleNamespace(
        time=time,
        num_jacobian_calls=NUM_JACOBIAN_CALLS,
        fit=fit,
        kwargs=kwargs,
    ))


###############################################################################
# Finally, we present the results:
labels = []

for result in results:
    label = result.kwargs['method']
    if result.num_jacobian_calls > 0:
        label += ' with Jac.'

    labels.append(label)

label_width = max(map(len, labels))
lines = [
    '| '
    + ' | '.join([
        'Method'.ljust(label_width),
        'Avg. time (ms)',
        '# func. (+ Jac.) calls',
        'Chi-squared',
        'a'.ljust(5),
        'b'.ljust(5),
        'c'.ljust(6),
    ])
    + '|'
]

print(f'The "true" parameters are: a = {a:.3f}, b = {b:.3f}, c = {c:.3f}\n')
fig, ax = plt.subplots()
ax.plot(x, data, marker='.', linestyle='none', label='data')

for (result, label) in zip(results, labels):
    linestyle = '-'
    if result.num_jacobian_calls > 0:
        linestyle = '--'

    a = result.fit.params['a'].value
    b = result.fit.params['b'].value
    c = result.fit.params['c'].value
    y = func([a, b, c], x)
    ax.plot(x, y, label=label, alpha=0.5, linestyle=linestyle)

    columns = [
        label.ljust(label_width),
        f'{result.time * 1000:.2f}'.ljust(14),
        (
            f'{result.fit.nfev}'
            + (
                f' (+{result.num_jacobian_calls})'
                if result.num_jacobian_calls > 0 else
                ''
            )
        ).ljust(22),
        f'{result.fit.chisqr:.3f}'.ljust(11),
        f'{a:.3f}'.ljust(5, '0'),
        f'{b:.3f}'.ljust(5, '0'),
        f'{c:.3f}'.ljust(5, '0'),
    ]
    lines.append('| ' + ' | '.join(columns) + ' |')

lines.insert(1, '|-' + '-|-'.join('-' * len(col) for col in columns) + '-|')
print('\n'.join(lines))

ax.set_xlabel('x')
ax.set_ylabel('y')
ax.legend()