1 # calibcant - tools for thermally calibrating AFM cantilevers
3 # Copyright (C) 2008-2011 W. Trevor King <wking@drexel.edu>
5 # This file is part of calibcant.
7 # calibcant is free software: you can redistribute it and/or
8 # modify it under the terms of the GNU Lesser General Public
9 # License as published by the Free Software Foundation, either
10 # version 3 of the License, or (at your option) any later version.
12 # calibcant 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 Lesser General Public License for more details.
17 # You should have received a copy of the GNU Lesser General Public
18 # License along with calibcant. If not, see
19 # <http://www.gnu.org/licenses/>.
21 """Surface bump analysis (measures photodiode sensitivity).
23 Separate the more general `bump_analyze()` from the other `bump_*()`
24 functions in calibcant.
26 The relevant physical quantities are:
28 * `Vzp_out` Output z-piezo voltage (what we generate)
29 * `Vzp` Applied z-piezo voltage (after external amplification)
30 * `Zp` The z-piezo position
31 * `Zcant` The cantilever vertical deflection
32 * `Vphoto` The photodiode vertical deflection voltage (what we measure)
34 Which are related by the parameters:
36 * `zp_gain` Vzp_out / Vzp
37 * `zp_sensitivity` Zp / Vzp
38 * `photo_sensitivity` Vphoto / Zcant
40 `photo_sensitivity` is measured by bumping the cantilever against the
41 surface, where `Zp = Zcant` (see `calibrate.bump_acquire()`). The
42 measured slope `Vphoto/Vout` is converted to `photo_sensitivity` with
46 >>> from pprint import pprint
49 >>> from h5config.storage.hdf5 import pprint_HDF5
50 >>> from pypiezo.config import ChannelConfig, AxisConfig
51 >>> from .config import BumpConfig
53 >>> fd,filename = tempfile.mkstemp(suffix='.h5', prefix='calibcant-')
56 >>> bump_config = BumpConfig()
57 >>> z_channel_config = ChannelConfig()
58 >>> z_channel_config['name'] = 'z'
59 >>> z_channel_config['maxdata'] = 200
60 >>> z_channel_config['conversion-coefficients'] = (0,1)
61 >>> z_channel_config['conversion-origin'] = 0
62 >>> z_axis_config = AxisConfig()
63 >>> z_axis_config['channel'] = z_channel_config
64 >>> deflection_channel_config = ChannelConfig()
65 >>> deflection_channel_config['name'] = 'deflection'
66 >>> deflection_channel_config['maxdata'] = 200
67 >>> deflection_channel_config['conversion-coefficients'] = (0,1)
68 >>> deflection_channel_config['conversion-origin'] = 0
71 ... 'z': numpy.arange(100, dtype=numpy.uint16),
72 ... 'deflection': numpy.arange(100, dtype=numpy.uint16),
74 >>> raw_bump['deflection'][:50] = 50
75 >>> processed_bump = bump_analyze(
76 ... raw_bump, bump_config, z_axis_config, deflection_channel_config)
77 >>> bump_plot(data=raw_bump) # TODO: convert to V and m
78 >>> bump_save(filename=filename, group='/bump/', raw_bump=raw_bump,
79 ... bump_config=bump_config, z_axis_config=z_axis_config,
80 ... deflection_channel_config=deflection_channel_config,
81 ... processed_bump=processed_bump)
83 >>> pprint_HDF5(filename) # doctest: +ELLIPSIS, +REPORT_UDIFF
88 <HDF5 dataset "far-steps": shape (), type "<i4">
90 <HDF5 dataset "initial-position": shape (), type "<f8">
92 <HDF5 dataset "model": shape (), type "|S9">
94 <HDF5 dataset "push-depth": shape (), type "<f8">
96 <HDF5 dataset "push-speed": shape (), type "<f8">
98 <HDF5 dataset "samples": shape (), type "<i4">
100 <HDF5 dataset "setpoint": shape (), type "<f8">
102 /bump/config/deflection
103 /bump/config/deflection/channel
104 <HDF5 dataset "channel": shape (), type "<i4">
106 <HDF5 dataset "conversion-coefficients": shape (2,), type "<i4">
108 <HDF5 dataset "conversion-origin": shape (), type "<i4">
110 <HDF5 dataset "device": shape (), type "|S12">
112 <HDF5 dataset "inverse-conversion-coefficients": shape (), type "|S1">
114 <HDF5 dataset "inverse-conversion-origin": shape (), type "<f8">
116 <HDF5 dataset "maxdata": shape (), type "<i4">
118 <HDF5 dataset "name": shape (), type "|S10">
120 <HDF5 dataset "range": shape (), type "<i4">
122 <HDF5 dataset "subdevice": shape (), type "<i4">
126 /bump/config/z/axis/channel
127 <HDF5 dataset "channel": shape (), type "<i4">
129 <HDF5 dataset "conversion-coefficients": shape (2,), type "<i4">
131 <HDF5 dataset "conversion-origin": shape (), type "<i4">
133 <HDF5 dataset "device": shape (), type "|S12">
135 <HDF5 dataset "inverse-conversion-coefficients": shape (), type "|S1">
137 <HDF5 dataset "inverse-conversion-origin": shape (), type "<f8">
139 <HDF5 dataset "maxdata": shape (), type "<i4">
141 <HDF5 dataset "name": shape (), type "|S1">
143 <HDF5 dataset "range": shape (), type "<i4">
145 <HDF5 dataset "subdevice": shape (), type "<i4">
147 <HDF5 dataset "gain": shape (), type "<f8">
149 <HDF5 dataset "maximum": shape (), type "|S4">
151 <HDF5 dataset "minimum": shape (), type "|S4">
153 <HDF5 dataset "monitor": shape (), type "|S1">
155 <HDF5 dataset "sensitivity": shape (), type "<f8">
157 <HDF5 dataset "processed": shape (), type "<f8">
160 <HDF5 dataset "deflection": shape (100,), type "<u2">
161 [50 50 ... 50 51 52 ... 97 98 99]
162 <HDF5 dataset "z": shape (100,), type "<u2">
163 [ 0 1 2 3 ... 97 98 99]
165 >>> (raw_bump,bump_config,z_axis_config,deflection_channel_config,
166 ... processed_bump) = bump_load(filename=filename, group='/bump/')
168 >>> pprint(raw_bump) # doctest: +ELLIPSIS
169 {'deflection': array([50, 50, ... 51, 52, 53, ..., 97, 98, 99], dtype=uint16),
170 'z': array([ 0, 1, 2, ..., 97, 98, 99], dtype=uint16)}
171 >>> processed_bump # doctest: +ELLIPSIS
174 >>> os.remove(filename)
178 import numpy as _numpy
179 from scipy.optimize import leastsq as _leastsq
182 import matplotlib as _matplotlib
183 import matplotlib.pyplot as _matplotlib_pyplot
184 import time as _time # for timestamping lines on plots
185 except (ImportError, RuntimeError), e:
187 _matplotlib_import_error = e
189 from h5config.storage.hdf5 import HDF5_Storage as _HDF5_Storage
190 from h5config.storage.hdf5 import h5_create_group as _h5_create_group
191 from pypiezo.base import convert_bits_to_volts as _convert_bits_to_volts
192 from pypiezo.base import convert_bits_to_meters as _convert_bits_to_meters
193 from pypiezo.config import ChannelConfig as _ChannelConfig
194 from pypiezo.config import AxisConfig as _AxisConfig
196 from . import LOG as _LOG
197 from . import package_config as _package_config
198 from .config import Linear as _Linear
199 from .config import Quadratic as _Quadratic
200 from .config import BumpConfig as _BumpConfig
203 def bump_analyze(data, bump_config, z_axis_config,
204 deflection_channel_config, plot=False):
205 """Return the slope of the bump.
208 data dictionary of data in DAC/ADC bits
209 bump_config `.config._BumpConfig` instance
210 z_axis_config z `pypiezo.config.AxisConfig` instance
211 deflection_channel_config
212 deflection `pypiezo.config.ChannelConfig` instance
213 plot boolean overriding matplotlib config setting.
215 photo_sensitivity (Vphoto/Zcant) in Volts/m
217 Checks for strong correlation (r-value) and low randomness chance
220 z = _convert_bits_to_meters(z_axis_config, data['z'])
221 deflection = _convert_bits_to_volts(
222 deflection_channel_config, data['deflection'])
223 high_voltage_rail = _convert_bits_to_volts(
224 deflection_channel_config, deflection_channel_config['maxdata'])
225 if bump_config['model'] == _Linear:
227 'param_guesser': limited_linear_param_guess,
228 'model': limited_linear,
229 'sensitivity_from_fit_params': limited_linear_sensitivity,
233 'param_guesser': limited_quadratic_param_guess,
234 'model': limited_quadratic,
235 'sensitivity_from_fit_params': limited_quadratic_sensitivity,
237 photo_sensitivity = bump_fit(
238 z, deflection, high_voltage_rail=high_voltage_rail, plot=plot,
240 return photo_sensitivity
242 def limited_linear(x, params, high_voltage_rail):
245 flat region (off-surface)
246 linear region (in-contact)
247 flat region (high-voltage-rail)
249 x_contact (x value for the surface-contact kink)
250 y_contact (y value for the surface-contact kink)
251 slope (dy/dx at the surface-contact kink)
253 x_contact,y_contact,slope = params
254 off_surface_mask = x <= x_contact
255 on_surface_mask = x > x_contact
256 y = (off_surface_mask * y_contact +
257 on_surface_mask * (y_contact + slope*(x-x_contact)))
258 y = _numpy.clip(y, y_contact, high_voltage_rail)
261 def limited_linear_param_guess(x, y):
263 Guess rough parameters for a limited_linear model. Assumes the
264 bump approaches (raising the deflection as it does so) first.
265 Retracting after the approach is optional. Approximates the contact
266 position and an on-surface (high) position by finding first crossings
267 of thresholds 0.3 and 0.7 of the y value's total range. Not the
268 most efficient algorithm, but it seems fairly robust.
270 y_contact = float(y.min())
271 y_max = float(y.max())
273 y_low = y_contact + 0.3 * (y_max-y_contact)
274 y_high = y_contact + 0.7 * (y_max-y_contact)
281 x_contact = float(x[i_low])
282 x_high = float(x[i_high])
283 if x_high == x_contact:
284 x.tofile('x-bad.dat', sep='\n')
285 y.tofile('y-bad.dat', sep='\n')
286 slope = (y_high - y_contact) / (x_high - x_contact)
287 return (x_contact, y_contact, slope)
289 def limited_linear_sensitivity(params):
291 Return the estimated sensitivity to small deflections according to
292 limited_linear fit parameters.
297 def limited_quadratic(x, params, high_voltage_rail):
300 flat region (off-surface)
301 quadratic region (in-contact)
302 flat region (high-voltage-rail)
304 x_contact (x value for the surface-contact kink)
305 y_contact (y value for the surface-contact kink)
306 slope (dy/dx at the surface-contact kink)
307 quad (d**2 y / dx**2, allow decreasing sensitivity with increased x)
309 x_contact,y_contact,slope,quad = params
310 off_surface_mask = x <= x_contact
311 on_surface_mask = x > x_contact
312 y = (off_surface_mask * y_contact +
314 y_contact + slope*(x-x_contact) + quad*(x-x_contact)**2))
315 y = _numpy.clip(y, y_contact, high_voltage_rail)
318 def limited_quadratic_param_guess(x, y):
320 Guess rough parameters for a limited_quadratic model. Assumes the
321 bump approaches (raising the deflection as it does so) first.
322 Retracting after the approach is optional. Approximates the contact
323 position and an on-surface (high) position by finding first crossings
324 of thresholds 0.3 and 0.7 of the y value's total range. Not the
325 most efficient algorithm, but it seems fairly robust.
327 x_contact,y_contact,linear_slope = limited_linear_param_guess(x,y)
328 contact_depth = x.max() - x_contact
329 slope = linear_slope / 2
330 quad = slope / contact_depth
331 # The above slope and quad were chosen so
333 # x*slope + x**2 * slope == x * linear_slope
334 return (x_contact, y_contact, slope, quad)
336 def limited_quadratic_sensitivity(params):
338 Return the estimated sensitivity to small deflections according to
339 limited_quadratic fit parameters.
344 def bump_fit(z, deflection, high_voltage_rail,
345 param_guesser=limited_quadratic_param_guess,
346 model=limited_quadratic,
347 sensitivity_from_fit_params=limited_quadratic_sensitivity,
349 """Fit a aurface bump and return the photodiode sensitivitiy.
352 z piezo position in meters
353 deflection photodiode deflection in volts
354 param_guesser function that guesses initial model parameters
355 model parametric model of deflection vs. z
356 sensitivity_from_fit_params given fit params, return the sensitivity
357 plot boolean overriding matplotlib config setting.
359 photodiode_sensitivity photodiode volts per piezo meter
361 _LOG.debug('fit bump data with model %s' % model)
362 def residual(p, deflection, z):
363 return model(z, p, high_voltage_rail=high_voltage_rail) - deflection
364 param_guess = param_guesser(z, deflection)
365 p,cov,info,mesg,ier = _leastsq(
366 residual, param_guess, args=(deflection, z), full_output=True,
368 _LOG.debug('fitted params: %s' % p)
369 _LOG.debug('covariance matrix: %s' % cov)
370 #_LOG.debug('info: %s' % info)
371 _LOG.debug('message: %s' % mesg)
373 _LOG.debug('solution converged')
375 _LOG.debug('solution did not converge')
376 if plot or _package_config['matplotlib']:
377 yguess = model(z, param_guess, high_voltage_rail=high_voltage_rail)
378 yfit = model(z, p, high_voltage_rail=high_voltage_rail)
379 bump_plot({'z': z, 'deflection': deflection}, yguess=yguess, yfit=yfit)
380 return sensitivity_from_fit_params(p)
382 def bump_save(filename, group='/', raw_bump=None, bump_config=None,
383 z_axis_config=None, deflection_channel_config=None,
384 processed_bump=None):
385 with _h5py.File(filename, 'a') as f:
386 cwg = _h5_create_group(f, group)
387 if raw_bump is not None:
388 for k in ['z', 'deflection']:
390 del cwg['raw/{}'.format(k)]
393 cwg['raw/z'] = raw_bump['z']
394 cwg['raw/deflection'] = raw_bump['deflection']
395 storage = _HDF5_Storage()
396 for config,key in [(bump_config, 'config/bump'),
397 (z_axis_config, 'config/z/axis'),
398 (deflection_channel_config,
399 'config/deflection/channel')]:
402 config_cwg = _h5_create_group(cwg, key)
403 storage.save(config=config, group=config_cwg)
404 if processed_bump is not None:
409 cwg['processed'] = processed_bump
411 def bump_load(filename, group='/'):
412 assert group.endswith('/')
413 raw_bump = processed_bump = None
415 with _h5py.File(filename, 'a') as f:
418 'z': f[group+'raw/z'][...],
419 'deflection': f[group+'raw/deflection'][...],
423 for Config,key in [(_BumpConfig, 'config/bump'),
424 (_AxisConfig, 'config/z/axis'),
425 (_ChannelConfig, 'config/deflection/channel')]:
426 config = Config(storage=_HDF5_Storage(
427 filename=filename, group=group+key))
428 configs.append(config)
430 processed_bump = float(f[group+'processed'][...])
435 ret.append(processed_bump)
436 for config in configs:
440 def bump_plot(data, yguess=None, yfit=None):
441 "Plot the bump (Vphoto vs Vzp)"
443 raise _matplotlib_import_error
444 figure = _matplotlib_pyplot.figure()
446 axes = figure.add_subplot(1, 1, 1)
448 axes = figure.add_subplot(2, 1, 1)
449 residual_axes = figure.add_subplot(2, 1, 2)
450 timestamp = _time.strftime('%H%M%S')
453 axes.plot(data['z'], data['deflection'], '.', label='bump')
455 axes.plot(data['z'], yguess, 'g-', label='guess')
457 axes.plot(data['z'], yfit, 'r-', label='fit')
459 axes.set_title('bump surface %s' % timestamp)
460 #axes.legend(loc='upper left')
461 axes.set_xlabel('Z piezo (meters)')
462 axes.set_ylabel('Photodiode (Volts)')
464 # second subplot for residual
465 residual_axes.plot(data['z'], data['deflection'] - yfit,
466 'r-', label='residual')
467 #residual_axes.legend(loc='upper right')
468 residual_axes.set_xlabel('Z piezo (meters)')
469 residual_axes.set_ylabel('Photodiode (Volts)')