add files for public distribution

based on internal repository 0a462b6 2017-11-22 14:41:39 +0100
This commit is contained in:
2017-11-22 14:55:20 +01:00
parent 96d206fc7b
commit bbd16d0f94
102 changed files with 230209 additions and 0 deletions

1
tests/__init__.py Normal file
View File

@ -0,0 +1 @@
__author__ = 'muntwiler_m'

321
tests/test_cluster.py Normal file
View File

@ -0,0 +1,321 @@
"""
@package tests.test_cluster
unit tests for pmsco.cluster
the purpose of these tests is to check whether the code runs as expected in a particular environment.
to run the tests, change to the directory which contains the tests directory, and execute =nosetests=.
@pre nose must be installed (python-nose package on Debian).
@author Matthias Muntwiler, matthias.muntwiler@psi.ch
@copyright (c) 2015-17 by Paul Scherrer Institut @n
Licensed under the Apache License, Version 2.0 (the "License"); @n
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
"""
import unittest
import math
import numpy as np
import pmsco.cluster as mc
class TestClusterFunctions(unittest.TestCase):
def setUp(self):
# before each test method
pass
def tearDown(self):
# after each test method
pass
@classmethod
def setup_class(cls):
# before any methods in this class
pass
@classmethod
def teardown_class(cls):
# teardown_class() after any methods in this class
pass
@staticmethod
def create_cube():
"""
create a cluster object with atoms on the corners, faces and body center of the unit cube.
the atom types are unique in an arbitrary sequence.
the emitter is at the origin, atom type 1.
@return: cluster.Cluster object.
"""
clu = mc.Cluster()
clu.add_atom(1, np.asarray([0, 0, 0]), 1)
clu.add_atom(2, np.asarray([1, 0, 0]), 0)
clu.add_atom(3, np.asarray([0, 1, 0]), 0)
clu.add_atom(4, np.asarray([0, 0, 1]), 0)
clu.add_atom(5, np.asarray([-1, 0, 0]), 0)
clu.add_atom(6, np.asarray([0, -1, 0]), 0)
clu.add_atom(7, np.asarray([0, 0, -1]), 0)
clu.add_atom(8, np.asarray([1, 1, 0]), 0)
clu.add_atom(9, np.asarray([0, 1, 1]), 0)
clu.add_atom(10, np.asarray([1, 0, 1]), 0)
clu.add_atom(11, np.asarray([-1, 1, 0]), 0)
clu.add_atom(12, np.asarray([0, -1, 1]), 0)
clu.add_atom(13, np.asarray([1, 0, -1]), 0)
clu.add_atom(14, np.asarray([-1, -1, 0]), 0)
clu.add_atom(15, np.asarray([0, -1, -1]), 0)
clu.add_atom(16, np.asarray([-1, 0, -1]), 0)
clu.add_atom(17, np.asarray([1, -1, 0]), 0)
clu.add_atom(18, np.asarray([0, 1, -1]), 0)
clu.add_atom(19, np.asarray([-1, 0, 1]), 0)
clu.add_atom(20, np.asarray([1, 1, 1]), 0)
clu.add_atom(21, np.asarray([-1, 1, 1]), 0)
clu.add_atom(22, np.asarray([1, -1, 1]), 0)
clu.add_atom(23, np.asarray([1, 1, -1]), 0)
clu.add_atom(24, np.asarray([-1, -1, -1]), 0)
clu.add_atom(25, np.asarray([1, -1, -1]), 0)
clu.add_atom(26, np.asarray([-1, 1, -1]), 0)
clu.add_atom(27, np.asarray([-1, -1, 1]), 0)
return clu
def test_numpy_extract(self):
"""
test array extraction code which should be compatible to numpy versions before and after 1.14.
numpy 1.14 introduces changes to multi-column indexing of structured arrays such as data[['x','y']].
first, it will return a view instead of a copy.
second, it will assign fields by position rather than by name.
the first change affects our cluster code in several places
where we extract XYZ coordinates from the cluster data array.
this test checks whether the new code works with a particular numpy version.
@return: None
"""
clu = self.create_cube()
xy2 = clu.data[['x', 'y']].copy()
xy3 = xy2.view((xy2.dtype[0], len(xy2.dtype.names)))
ctr = np.asarray((1.0, 0.0, 0.0))
dist = np.linalg.norm(xy3 - ctr[0:2], axis=1)
self.assertAlmostEqual(1.0, dist[0])
self.assertAlmostEqual(0.0, dist[1])
clu.clear()
xy2 = clu.data[['x', 'y']].copy()
xy3 = xy2.view((xy2.dtype[0], len(xy2.dtype.names)))
ctr = np.asarray((1.0, 0.0, 0.0))
dist = np.linalg.norm(xy3 - ctr[0:2], axis=1)
self.assertEqual(0, dist.shape[0])
def test_get_positions(self):
"""
check that we get an independent copy of the original data.
@return: None
"""
clu = self.create_cube()
pos = clu.get_positions()
self.assertEqual(clu.data.shape[0], pos.shape[0])
self.assertEqual(3, pos.shape[1])
self.assertEqual(np.float32, pos.dtype)
self.assertEqual(1.0, pos[1, 0])
self.assertEqual(0.0, pos[1, 1])
self.assertEqual(0.0, pos[1, 2])
pos[0, 0] = 15.0
self.assertEqual(0.0, clu.data['x'][0])
# empty cluster
clu.clear()
self.assertEqual(clu.data.shape[0], 0)
self.assertEqual(3, pos.shape[1])
self.assertEqual(np.float32, pos.dtype)
def test_set_positions(self):
clu = mc.Cluster()
clu.data = np.zeros(2, dtype=clu.dtype)
pos = np.array([[1., 2., 3.], [4., 5., 6.]])
clu.set_positions(pos)
self.assertEqual(1., clu.data['x'][0])
self.assertEqual(2., clu.data['y'][0])
self.assertEqual(3., clu.data['z'][0])
self.assertEqual(4., clu.data['x'][1])
self.assertEqual(5., clu.data['y'][1])
self.assertEqual(6., clu.data['z'][1])
def test_get_emitters(self):
clu = self.create_cube()
clu.set_emitter(idx=0)
clu.set_emitter(idx=9)
self.assertEqual(2, clu.get_emitter_count())
result = clu.get_emitters()
expect = [(0., 0., 0., 1), (1., 0., 1., 10)]
self.assertItemsEqual(expect, result)
def test_get_z_layers(self):
clu = mc.Cluster()
clu.add_atom(1, np.asarray([1, 0, 0.1]), 0)
clu.add_atom(2, np.asarray([0, 1, -0.3]), 0)
clu.add_atom(1, np.asarray([0, 1, -0.2]), 0)
clu.add_atom(1, np.asarray([1, 0, 0]), 1)
clu.add_atom(1, np.asarray([0, 1, -0.2001]), 0)
clu.add_atom(2, np.asarray([0, 1, -0.1]), 0)
clu.add_atom(1, np.asarray([0, 1, -0.1999]), 0)
layers = clu.get_z_layers(0.01)
np.testing.assert_allclose(layers, np.asarray([-0.3, -0.2, -0.1, 0.0, +0.1]), atol=0.001)
def test_relax(self):
clu = mc.Cluster()
clu.add_atom(1, np.asarray([1, 0, 1]), 0)
clu.add_atom(1, np.asarray([1, 0, 0]), 1)
clu.add_atom(2, np.asarray([0, 1, -1]), 0)
clu.add_atom(1, np.asarray([0, 1, -2]), 0)
clu.add_atom(2, np.asarray([0, 1, -3]), 0)
idx = clu.relax(-0.3, -0.1, 2)
np.testing.assert_almost_equal(idx, np.asarray([[2, 4]]))
np.testing.assert_allclose(clu.get_position(0), np.asarray([1, 0, 1]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(1), np.asarray([1, 0, 0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(2), np.asarray([0, 1, -1.1]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(3), np.asarray([0, 1, -2.0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(4), np.asarray([0, 1, -3.1]), atol=1e-6)
def test_rotate_x(self):
clu = mc.Cluster()
clu.add_atom(1, np.asarray([1, 0, 0]), 0)
clu.add_atom(1, np.asarray([0, 1, 0]), 0)
clu.add_atom(1, np.asarray([0, 0, 1]), 0)
clu.rotate_x(90)
np.testing.assert_allclose(clu.get_position(0), np.asarray([1, 0, 0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(1), np.asarray([0, 0, 1]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(2), np.asarray([0, -1, 0]), atol=1e-6)
def test_rotate_y(self):
clu = mc.Cluster()
clu.add_atom(1, np.asarray([1, 0, 0]), 0)
clu.add_atom(1, np.asarray([0, 1, 0]), 0)
clu.add_atom(1, np.asarray([0, 0, 1]), 0)
clu.rotate_y(90)
np.testing.assert_allclose(clu.get_position(0), np.asarray([0, 0, -1]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(1), np.asarray([0, 1, 0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(2), np.asarray([1, 0, 0]), atol=1e-6)
def test_rotate_z(self):
clu = mc.Cluster()
clu.add_atom(1, np.asarray([1, 0, 0]), 0)
clu.add_atom(1, np.asarray([0, 1, 0]), 0)
clu.add_atom(1, np.asarray([0, 0, 1]), 0)
clu.rotate_z(90)
np.testing.assert_allclose(clu.get_position(0), np.asarray([0, 1, 0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(1), np.asarray([-1, 0, 0]), atol=1e-6)
np.testing.assert_allclose(clu.get_position(2), np.asarray([0, 0, 1]), atol=1e-6)
def test_add_layer(self):
clu = mc.Cluster()
# from hbncu project
b_surf = 2.50
clu.set_rmax(4.0)
b1 = np.array((b_surf, 0.0, 0.0))
b2 = np.array((b_surf / 2.0, b_surf * math.sqrt(3.0) / 2.0, 0.0))
a1 = -10.0 * b1 - 10.0 * b2
emitter = np.array((0.0, 0.0, 0.0))
clu.add_layer(7, a1, b1, b2)
pos = clu.find_positions(pos=emitter)
self.assertEqual(len(pos), 1)
def test_add_cluster(self):
clu1 = mc.Cluster()
clu1.add_atom(1, np.asarray([0, 0, 0]), 1)
clu1.add_atom(2, np.asarray([1, 0, 0]), 0)
clu1.add_atom(3, np.asarray([0, 1, 0]), 0)
clu1.add_atom(4, np.asarray([0, 0, -1]), 0)
clu1.add_atom(5, np.asarray([0, 0, -2]), 0)
clu2 = mc.Cluster()
clu2.add_atom(3, np.asarray([-0.2, 0, 0]), 0)
clu2.add_atom(4, np.asarray([0, -0.2, 0]), 0)
clu2.add_atom(5, np.asarray([0, 0.05, -1]), 0)
clu2.add_atom(5, np.asarray([0, 0, -1.01]), 0)
clu2.add_atom(6, np.asarray([0, 0, -1.99]), 0)
clu1.set_rmax(1.5)
clu1.add_cluster(clu2, check_rmax=True, check_unique=True, tol=0.1)
self.assertEqual(clu1.get_atom_count(), 5+2)
def test_find_positions(self):
clu = mc.Cluster()
# from hbncu project
b_surf = 2.50
clu.set_rmax(b_surf * 10.0)
b1 = np.array((b_surf, 0.0, 0.0))
b2 = np.array((b_surf / 2.0, b_surf * math.sqrt(3.0) / 2.0, 0.0))
a_N = np.array((0.0, 0.0, 0.0))
a_B = np.array(((b1[0] + b2[0]) / 3.0, (b1[1] + b2[1]) / 3.0, 0.0))
emitter = a_N + b1 * 1 + b2 * 2
clu.add_layer(7, a_N, b1, b2)
clu.add_layer(5, a_B, b1, b2)
pos = clu.find_positions(pos=emitter)
self.assertEqual(len(pos), 1)
self.assertEqual(pos[0], 206)
def test_find_index_cylinder(self):
clu = self.create_cube()
pos = np.array((0.8, 0.8, 0.8))
rxy = 0.5
rz = 1.0
idx = clu.find_index_cylinder(pos, rxy, rz, None)
self.assertEqual(len(idx), 2)
self.assertEqual(clu.get_atomtype(idx[0]), 8)
self.assertEqual(clu.get_atomtype(idx[1]), 20)
idx = clu.find_index_cylinder(pos, rxy, rz, 8)
self.assertEqual(len(idx), 1)
def test_trim_cylinder(self):
clu = mc.Cluster()
clu.set_rmax(10.0)
v_pos = np.asarray([0, 0, 0])
v_lat1 = np.asarray([1, 0, 0])
v_lat2 = np.asarray([0, 1, 0])
v_lat3 = np.asarray([0, 0, 1])
clu.add_bulk(7, v_pos, v_lat1, v_lat2, v_lat3)
clu.set_emitter(pos=v_pos)
clu.trim_cylinder(2.3, 4.2)
self.assertEqual(clu.data.dtype, clu.dtype)
self.assertEqual(clu.data.shape[0], 21 * 5)
self.assertEqual(clu.data[1]['i'], 2)
self.assertEqual(clu.data[1]['s'], 'N')
self.assertEqual(clu.data[1]['t'], 7)
self.assertEqual(clu.get_emitter_count(), 1)
def test_trim_sphere(self):
clu = mc.Cluster()
clu.set_rmax(10.0)
v_pos = np.asarray([0, 0, 0])
v_lat1 = np.asarray([1, 0, 0])
v_lat2 = np.asarray([0, 1, 0])
v_lat3 = np.asarray([0, 0, 1])
clu.add_bulk(7, v_pos, v_lat1, v_lat2, v_lat3)
clu.set_emitter(pos=v_pos)
clu.trim_sphere(2.3)
self.assertEqual(clu.data.dtype, clu.dtype)
self.assertEqual(clu.data.shape[0], 39)
self.assertEqual(clu.data[1]['i'], 2)
self.assertEqual(clu.data[1]['s'], 'N')
self.assertEqual(clu.data[1]['t'], 7)
self.assertEqual(clu.get_emitter_count(), 1)
def test_trim_slab(self):
clu = self.create_cube()
clu.trim_slab('x', 0.5, 1.1)
self.assertEqual(clu.data.dtype, clu.dtype)
self.assertEqual(clu.data.shape[0], 9 * 2)
self.assertEqual(clu.get_emitter_count(), 1)

428
tests/test_data.py Normal file
View File

@ -0,0 +1,428 @@
"""
@package tests.test_data
unit tests for pmsco.data
the purpose of these tests is to mainly to check the syntax, and correct data types,
i.e. anything that could cause a run-time error.
calculation results are sometimes checked for plausibility but not exact values,
depending on the level of debugging required for a specific part of the code.
to run the tests, change to the directory which contains the tests directory, and execute =nosetests=.
@pre nose must be installed (python-nose package on Debian).
@author Matthias Muntwiler, matthias.muntwiler@psi.ch
@copyright (c) 2015 by Paul Scherrer Institut @n
Licensed under the Apache License, Version 2.0 (the "License"); @n
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
"""
import unittest
import math
import numpy as np
import pmsco.data as md
class TestDataFunctions(unittest.TestCase):
def setUp(self):
# before each test method
shape = (10, )
self.e_scan = md.create_data(shape, dtype=md.DTYPE_EI)
self.e_scan['e'] = np.linspace(100.0, 200.0, shape[0])
self.e_scan['i'] = np.linspace(0.0, 10.0, shape[0])
shape = (12, )
self.ea_scan = md.create_data(shape, dtype=md.DTYPE_ETPAI)
self.ea_scan[0] = (1.0, 0.0, 0.0, -1.0, 1.0)
self.ea_scan[1] = (1.0, 0.0, 0.0, 0.0, 1.0)
self.ea_scan[2] = (1.0, 0.0, 0.0, 1.0, 1.0)
self.ea_scan[3] = (1.0, 0.0, 0.0, 2.0, 1.0)
self.ea_scan[4] = (2.0, 0.0, 0.0, -1.0, 2.0)
self.ea_scan[5] = (2.0, 0.0, 0.0, 0.0, 2.0)
self.ea_scan[6] = (2.0, 0.0, 0.0, 1.0, 2.0)
self.ea_scan[7] = (2.0, 0.0, 0.0, 2.0, 2.0)
self.ea_scan[8] = (3.0, 0.0, 0.0, -1.0, 3.0)
self.ea_scan[9] = (3.0, 0.0, 0.0, 0.0, 3.0)
self.ea_scan[10] = (3.0, 0.0, 0.0, 1.0, 3.0)
self.ea_scan[11] = (3.0, 0.0, 0.0, 2.0, 3.0)
shape = (10, )
self.holo_scan = md.create_data(shape, dtype=md.DTYPE_ETPI)
self.holo_scan[0] = (1.0, 0.0, 0.0, 1.0)
self.holo_scan[1] = (1.0, 1.0, 0.0, 2.0)
self.holo_scan[2] = (1.0, 1.0, 120.0, 3.0)
self.holo_scan[3] = (1.0, 1.0, 240.0, 4.0)
self.holo_scan[4] = (1.0, 2.0, 0.0, 5.0)
self.holo_scan[5] = (1.0, 2.0, 60.0, 6.0)
self.holo_scan[6] = (1.0, 2.0, 120.0, 7.0)
self.holo_scan[7] = (1.0, 2.0, 180.0, 8.0)
self.holo_scan[8] = (1.0, 2.0, 240.0, 9.0)
self.holo_scan[9] = (1.0, 2.0, 300.0, 10.0)
def tearDown(self):
# after each test method
pass
@classmethod
def setup_class(cls):
# before any methods in this class
pass
@classmethod
def teardown_class(cls):
# teardown_class() after any methods in this class
pass
def test_create_data(self):
shape = (10, )
data = md.create_data(shape, dtype=md.DTYPE_ETPAIS)
expected_names = ('e', 't', 'p', 'a', 'i', 's')
self.assertItemsEqual(data.dtype.names, expected_names)
self.assertEqual(data.shape, shape)
def test_detect_scan_mode_1d(self):
scan_mode, scan_positions = md.detect_scan_mode(self.e_scan)
expected_mode = ['e']
expected_positions = {}
expected_positions['e'] = np.linspace(0.0, 10.0, 10)
self.assertItemsEqual(scan_mode, expected_mode)
self.assertItemsEqual(scan_positions, expected_positions)
def test_detect_scan_mode_2d(self):
scan_mode, scan_positions = md.detect_scan_mode(self.ea_scan)
expected_mode = ['e', 'a']
expected_positions = {}
expected_positions['e'] = np.asarray((1.0, 2.0, 3.0))
expected_positions['t'] = np.zeros((1))
expected_positions['p'] = np.zeros((1))
expected_positions['a'] = np.asarray((-1.0, 0.0, 1.0, 2.0))
self.assertItemsEqual(scan_mode, expected_mode)
self.assertItemsEqual(scan_positions, expected_positions)
def test_detect_scan_mode_holo(self):
scan_mode, scan_positions = md.detect_scan_mode(self.holo_scan)
expected_mode = ['t', 'p']
expected_positions = {}
expected_positions['e'] = np.ones((1))
expected_positions['t'] = self.holo_scan['t']
expected_positions['p'] = self.holo_scan['p']
self.assertItemsEqual(scan_mode, expected_mode)
self.assertItemsEqual(scan_positions, expected_positions)
def test_detect_scan_mode_theta(self):
scan = self.holo_scan
scan['t'] = np.linspace(1.0, 2.0, scan.shape[0])
scan['p'] = 3.3
scan_mode, scan_positions = md.detect_scan_mode(scan)
expected_mode = ['t']
expected_positions = {}
expected_positions['e'] = np.ones((1))
expected_positions['t'] = np.linspace(1.0, 2.0, scan.shape[0])
expected_positions['p'] = np.ones((1)) * 3.3
self.assertItemsEqual(scan_mode, expected_mode)
self.assertItemsEqual(scan_positions, expected_positions)
def test_calc_modfunc_mean_1d(self):
modf = md.calc_modfunc_mean(self.e_scan)
exp_modf = self.e_scan.copy()
exp_modf['i'] = (self.e_scan['i'] - 5.0) / 5.0
np.testing.assert_allclose(modf['e'], exp_modf['e'])
np.testing.assert_allclose(modf['i'], exp_modf['i'])
def test_calc_modfunc_mean_2d(self):
modf = md.calc_modfunc_mean(self.ea_scan)
exp_modf = self.ea_scan.copy()
exp_modf['i'][0] = -0.5
exp_modf['i'][1] = -0.5
exp_modf['i'][2] = -0.5
exp_modf['i'][3] = -0.5
exp_modf['i'][4] = 0.0
exp_modf['i'][5] = 0.0
exp_modf['i'][6] = 0.0
exp_modf['i'][7] = 0.0
exp_modf['i'][8] = 0.5
exp_modf['i'][9] = 0.5
exp_modf['i'][10] = 0.5
exp_modf['i'][11] = 0.5
np.testing.assert_allclose(modf['e'], exp_modf['e'])
np.testing.assert_allclose(modf['t'], exp_modf['t'])
np.testing.assert_allclose(modf['p'], exp_modf['p'])
np.testing.assert_allclose(modf['a'], exp_modf['a'])
np.testing.assert_allclose(modf['i'], exp_modf['i'])
def test_calc_modfunc_loess_1d(self):
"""
check that the result of msc_data.calc_modfunc_loess() is between -1 and 1.
"""
modf = md.calc_modfunc_loess(self.e_scan)
self.assertEqual(self.e_scan.shape, modf.shape)
exp_modf = self.e_scan.copy()
np.testing.assert_allclose(modf['e'], exp_modf['e'])
exp_modf['i'] = -1.000001
np.testing.assert_array_less(exp_modf['i'], modf['i'])
exp_modf['i'] = +1.000001
np.testing.assert_array_less(modf['i'], exp_modf['i'])
def test_calc_modfunc_loess_1d_nan(self):
"""
check that data.calc_modfunc_loess() ignores NaNs gracefully.
"""
modified_index = 2
self.e_scan['i'][modified_index] = np.nan
modf = md.calc_modfunc_loess(self.e_scan)
exp_modf = self.e_scan.copy()
self.assertEqual(self.e_scan.shape, modf.shape)
np.testing.assert_allclose(modf['e'], exp_modf['e'])
self.assertTrue(math.isnan(modf['i'][modified_index]))
modf['i'][modified_index] = 0.0
exp_modf['i'] = -1.000001
np.testing.assert_array_less(exp_modf['i'], modf['i'])
exp_modf['i'] = +1.000001
np.testing.assert_array_less(modf['i'], exp_modf['i'])
def test_calc_modfunc_loess_2d(self):
"""
check that the msc_data.calc_modfunc_loess() function does approximately what we want for a two-dimensional dataset.
"""
n_e = 10
n_a = 15
shape = (n_e * n_a, )
scan = md.create_data(shape, dtype=md.DTYPE_ETPAI)
e_range = np.linspace(100.0, 200.0, n_e)
a_range = np.linspace(-15.0, 15.0, n_a)
a_grid, e_grid = np.meshgrid(a_range, e_range)
scan['e'] = np.ravel(e_grid)
scan['a'] = np.ravel(a_grid)
scan['t'] = 0.0
scan['p'] = 90.0
scan['i'] = 0.02 * scan['e'] + 0.03 * scan['a'] + np.sin((scan['e'] - 150) / 50 * math.pi) + np.sin(scan['a'] / 180 * math.pi)
modf = md.calc_modfunc_loess(scan)
self.assertEqual(scan.shape, modf.shape)
exp_modf = scan.copy()
np.testing.assert_allclose(modf['e'], exp_modf['e'])
np.testing.assert_allclose(modf['t'], exp_modf['t'])
np.testing.assert_allclose(modf['p'], exp_modf['p'])
np.testing.assert_allclose(modf['a'], exp_modf['a'])
exp_modf['i'] = -1.000001
np.testing.assert_array_less(exp_modf['i'], modf['i'])
exp_modf['i'] = +1.000001
np.testing.assert_array_less(modf['i'], exp_modf['i'])
# this is rough estimate of the result, manually optimized by trial and error in Igor.
# the R factor should be sensitive enough to detect mixed-up axes.
exp_modf['i'] = 0.03 * np.sin((scan['e'] - 150) / 50 * math.pi)
rf = md.rfactor(modf, exp_modf)
print rf
self.assertLessEqual(rf, 0.50)
def test_alpha_mirror_average(self):
n_e = 10
n_a = 15
shape = (n_e * n_a, )
scan = md.create_data(shape, dtype=md.DTYPE_ETPAI)
e_range = np.linspace(100.0, 200.0, n_e)
a_range = np.linspace(-15.0, 15.0, n_a)
a_grid, e_grid = np.meshgrid(a_range, e_range)
scan['e'] = np.ravel(e_grid)
scan['a'] = np.ravel(a_grid)
scan['t'] = 0.0
scan['p'] = 90.0
scan['i'] = 0.02 * scan['e'] + 0.03 * scan['a']
act_result = md.alpha_mirror_average(scan)
exp_result = scan.copy()
exp_result['i'] = 0.02 * scan['e']
np.testing.assert_allclose(act_result['e'], exp_result['e'])
np.testing.assert_allclose(act_result['t'], exp_result['t'])
np.testing.assert_allclose(act_result['p'], exp_result['p'])
np.testing.assert_allclose(act_result['a'], exp_result['a'])
np.testing.assert_allclose(act_result['i'], exp_result['i'])
def test_alpha_average(self):
data = md.create_data((20), dtype=md.DTYPE_ETPAIS)
data['e'][0:10] = 11
data['e'][10:20] = 22
data['a'][0:10] = np.arange(10)
data['a'][10:20] = np.arange(10)
data['i'] = np.arange(20)
data['s'] = np.arange(20) / 10.0
exp_result = md.create_data((2), dtype=md.DTYPE_ETPIS)
exp_result['e'][0] = 11
exp_result['e'][1] = 22
exp_result['i'][0] = 4.5
exp_result['i'][1] = 14.5
exp_result['s'] = exp_result['i'] / 10.0
act_result = md.alpha_average(data)
np.testing.assert_allclose(act_result['e'], exp_result['e'])
np.testing.assert_allclose(act_result['t'], exp_result['t'])
np.testing.assert_allclose(act_result['p'], exp_result['p'])
np.testing.assert_allclose(act_result['i'], exp_result['i'])
np.testing.assert_allclose(act_result['s'], exp_result['s'])
def test_phi_average(self):
data = self.holo_scan.copy()
exp_result = md.create_data((3), dtype=[('e', 'f4'), ('t', 'f4'), ('i', 'f4')])
exp_result['e'] = np.asarray([1, 1, 1])
exp_result['t'] = np.asarray([0, 1, 2])
exp_result['i'] = np.asarray([1, 3, 7.5])
act_result = md.phi_average(data)
np.testing.assert_allclose(act_result['e'], exp_result['e'])
np.testing.assert_allclose(act_result['t'], exp_result['t'])
np.testing.assert_allclose(act_result['i'], exp_result['i'])
def test_filter_tp(self):
data = md.create_data((17), dtype=md.DTYPE_ETPAIS)
data['e'] = 100.0
data['t'] = np.array([0.0, 0.0, 0.0, 0.0, 1.0, 1.0, 1.0, 1.0, 2.0, 2.0, 2.0, 2.0, 3.0, 3.0, 3.0, 3.0, 3.0])
data['p'] = np.array([0.0, 1.0, 2.0, 3.0, 0.0, 1.01, 2.0, 3.0, 0.0, 1.0, 2.0, 3.0, 0.0, 1.0, 2.0, 3.0, 3.1])
data['a'] = 0.0
data['i'] = data['t'] * 10.0 + data['p']
data['s'] = np.sqrt(data['i'])
# duplicate
data['p'][10] = 0.0
filter = md.create_data((10), dtype=md.DTYPE_ETPI)
filter['e'] = 100.0
filter['t'] = np.array([0.0, 1.0, 1.0, 2.0, 2.0, 2.0, 3.0, 3.0, 3.0, 3.0])
filter['p'] = np.array([0.0, 0.0, 1.0, 0.0, 1.0, 2.0, 0.0, 1.0, 2.0, 3.0])
filter['i'] = -99.999
exp_result = md.create_data((10), dtype=md.DTYPE_ETPAIS)
exp_result['e'] = 100.0
exp_result['t'] = np.array([0.0, 1.0, 1.00, 2.0, 2.0, 2.0, 3.0, 3.0, 3.0, 3.0])
exp_result['p'] = np.array([0.0, 0.0, 1.01, 0.0, 1.0, 2.0, 0.0, 1.0, 2.0, 3.0])
exp_result['a'] = 0.0
exp_result['i'] = exp_result['t'] * 10.0 + exp_result['p']
exp_result['s'] = np.sqrt(exp_result['i'])
exp_result['p'][5] = 0.0
act_result = md.filter_tp(data, filter)
np.testing.assert_allclose(act_result['e'], exp_result['e'])
np.testing.assert_allclose(act_result['t'], exp_result['t'])
np.testing.assert_allclose(act_result['p'], exp_result['p'])
np.testing.assert_allclose(act_result['a'], exp_result['a'])
np.testing.assert_allclose(act_result['i'], exp_result['i'])
np.testing.assert_allclose(act_result['s'], exp_result['s'])
def test_interpolate_hemi_scan(self):
n_t = 91
n_p = 11
shape = (n_t * n_p, )
calc_data = md.create_data(shape, dtype=md.DTYPE_TPI)
t_range = np.linspace(0.0, 90.0, n_t)
p_range = np.linspace(0.0, 360.0, n_p)
t_grid, p_grid = np.meshgrid(t_range, p_range)
calc_data['t'] = np.ravel(t_grid)
calc_data['p'] = np.ravel(p_grid)
calc_data['i'] = 100.0 * calc_data['t'] + calc_data['p']
scan_data = md.create_data((10), dtype=md.DTYPE_TPI)
scan_data['t'] = np.array([0.0, 1.0, 1.0, 2.0, 2.0, 2.0, 3.0, 3.0, 3.0, 3.0])
scan_data['p'] = np.array([0.0, 0.0, 1.0, 0.5, 1.7, 2.8, 0.0, 90.0, 180.0, 360.0])
scan_data['i'] = -99.999
exp_result = scan_data.copy()
exp_result['i'] = 100.0 * scan_data['t'] + scan_data['p']
act_result = md.interpolate_hemi_scan(calc_data, scan_data)
np.testing.assert_allclose(act_result['t'], exp_result['t'])
np.testing.assert_allclose(act_result['p'], exp_result['p'])
np.testing.assert_allclose(act_result['i'], exp_result['i'])
def test_scaled_rfactor(self):
n = 20
calc_modf = md.create_etpi((n,), False)
calc_modf['e'] = 0.0
calc_modf['t'] = 0.0
calc_modf['p'] = np.linspace(0, np.pi, n)
calc_modf['i'] = np.sin(calc_modf['p'])
exp_modf = calc_modf.copy()
exp_modf['i'] = 0.6 * np.sin(exp_modf['p'] + np.pi / 100.0)
weights = np.ones_like(exp_modf['i'])
r = md.scaled_rfactor(1.4, exp_modf['i'], weights, calc_modf['i'])
self.assertGreater(r, 0.0)
self.assertLess(r, 0.05)
def test_rfactor(self):
n = 20
calc_modf = md.create_data((n,), dtype=md.DTYPE_ETPI)
calc_modf['e'] = 0.0
calc_modf['t'] = 0.0
calc_modf['p'] = np.linspace(0, np.pi, n)
calc_modf['i'] = np.sin(calc_modf['p'])
exp_modf = md.create_data((n,), dtype=md.DTYPE_ETPIS)
exp_modf['i'] = 0.6 * np.sin(exp_modf['p'] + np.pi / 100.0)
exp_modf['s'] = np.sqrt(np.abs(exp_modf['i']))
r1 = md.rfactor(exp_modf, calc_modf)
self.assertAlmostEqual(r1, 0.95, delta=0.02)
# one nan should not make a big difference
calc_modf['i'][3] = np.nan
r2 = md.rfactor(exp_modf, calc_modf)
self.assertAlmostEqual(r1, r2, delta=0.02)
# all values nan should raise an exception
with self.assertRaises(ValueError):
calc_modf['i'] = np.nan
md.rfactor(exp_modf, calc_modf)
def test_optimize_rfactor(self):
n = 20
calc_modf = md.create_data((n,), dtype=md.DTYPE_ETPI)
calc_modf['e'] = 0.0
calc_modf['t'] = 0.0
calc_modf['p'] = np.linspace(0, np.pi, n)
calc_modf['i'] = np.sin(calc_modf['p'])
exp_modf = md.create_data((n,), dtype=md.DTYPE_ETPIS)
exp_modf['i'] = 0.6 * np.sin(exp_modf['p'] + np.pi / 100.0)
exp_modf['s'] = np.sqrt(np.abs(exp_modf['i']))
r1 = md.optimize_rfactor(exp_modf, calc_modf)
self.assertAlmostEqual(r1, 0.55, delta=0.02)
# one nan should not make a big difference
calc_modf['i'][3] = np.nan
r2 = md.optimize_rfactor(exp_modf, calc_modf)
self.assertAlmostEqual(r1, r2, delta=0.02)
# all values nan should raise an exception
with self.assertRaises(ValueError):
calc_modf['i'] = np.nan
md.optimize_rfactor(exp_modf, calc_modf)
if __name__ == '__main__':
unittest.main()

104
tests/test_files.py Normal file
View File

@ -0,0 +1,104 @@
"""
@package tests.test_files
unit tests for pmsco.files
the purpose of these tests is to help debugging the code.
to run the tests, change to the directory which contains the tests directory, and execute =nosetests=.
@pre nose must be installed (python-nose package on Debian).
@author Matthias Muntwiler, matthias.muntwiler@psi.ch
@copyright (c) 2015 by Paul Scherrer Institut @n
Licensed under the Apache License, Version 2.0 (the "License"); @n
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
"""
import unittest
import mock
import pmsco.files as files
class TestFileTracker(unittest.TestCase):
def setUp(self):
# before each test method
self.files = files.FileTracker()
self.files.keep_rfac = 1
self.files._os_delete_file = mock.Mock(return_value=None)
self.files.add_file("model 1 file 1 cluster K", 1, 'cluster')
self.files.add_file("model 1 file 2 output D", 1, 'output')
self.files.add_file("model 2 file 1 cluster K", 2, 'cluster')
self.files.add_file("model 2 file 2 output D", 2, 'output')
self.files.add_file("model 3 file 1 cluster K", 3, 'cluster')
self.files.add_file("model 3 file 2 output D", 3, 'output')
self.files.add_file("model 4 file 1 cluster K", 4, 'cluster')
self.files.add_file("model 4 file 2 output D", 4, 'output')
self.files.add_file("model 5 file 1 cluster K", 5, 'cluster')
self.files.add_file("model 5 file 2 output D", 5, 'output')
self.files.update_model_rfac(2, 0.0)
self.files.update_model_rfac(3, 0.1)
self.files.update_model_rfac(4, 0.3)
self.files.update_model_rfac(1, 0.5)
self.files.update_model_rfac(5, 0.6)
self.files.set_model_complete(1, True)
self.files.set_model_complete(2, True)
self.files.set_model_complete(3, False)
self.files.set_model_complete(5, True)
def tearDown(self):
# after each test method
pass
@classmethod
def setup_class(cls):
# before any methods in this class
pass
@classmethod
def teardown_class(cls):
# teardown_class() after any methods in this class
pass
def test_add_file(self):
pass
def test_rename_file(self):
pass
def test_remove_file(self):
pass
def test_update_model_rfac(self):
pass
def test_delete_files(self):
self.files.keep_rfac = 10
self.files.delete_files()
self.files._os_delete_file.assert_any_call("model 1 file 2 output D")
self.files._os_delete_file.assert_any_call("model 2 file 2 output D")
self.files._os_delete_file.assert_any_call("model 5 file 2 output D")
self.assertEqual(len(self.files._id_by_path), 5+2)
self.assertEqual(len(self.files._path_by_id), 5+2)
self.assertEqual(len(self.files._model_by_id), 5+2)
self.assertEqual(len(self.files._category_by_id), 5+2)
def test_delete_file(self):
pass
def test_delete_bad_rfac(self):
self.files.delete_bad_rfac(keep=2, force_delete=True)
self.files._os_delete_file.assert_any_call("model 1 file 1 cluster K")
self.files._os_delete_file.assert_any_call("model 5 file 1 cluster K")
self.assertEqual(len(self.files._id_by_path), 6)
self.assertEqual(len(self.files._path_by_id), 6)
self.assertEqual(len(self.files._model_by_id), 6)
self.assertEqual(len(self.files._category_by_id), 6)
def test_delete_category(self):
pass

48
tests/test_hbncu.py Normal file
View File

@ -0,0 +1,48 @@
"""
@package tests.test_hbncu
unit tests for projects.hbncu
the purpose of these tests is to help debugging the code.
to run the tests, change to the directory which contains the tests directory, and execute =nosetests=.
@pre nose must be installed (python-nose package on Debian).
@author Matthias Muntwiler, matthias.muntwiler@psi.ch
@copyright (c) 2015 by Paul Scherrer Institut @n
Licensed under the Apache License, Version 2.0 (the "License"); @n
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
"""
import unittest
import os.path
import tempfile
import shutil
import numpy as np
import projects.hbncu.hbncu as hbncu
import pmsco.data as data
import pmsco.dispatch as dispatch
class TestHbncuProject(unittest.TestCase):
def setUp(self):
self.test_dir = tempfile.mkdtemp()
self.project = hbncu.HbncuProject()
def tearDown(self):
# after each test method
self.project = None
shutil.rmtree(self.test_dir)
@classmethod
def setup_class(cls):
# before any methods in this class
pass
@classmethod
def teardown_class(cls):
# teardown_class() after any methods in this class
pass

364
tests/test_swarm.py Normal file
View File

@ -0,0 +1,364 @@
"""
@package tests.test_swarm
unit tests for pmsco.swarm
the purpose of these tests is to help debugging the code.
to run the tests, change to the directory which contains the tests directory, and execute =nosetests=.
@pre nose must be installed (python-nose package on Debian).
@author Matthias Muntwiler, matthias.muntwiler@psi.ch
@copyright (c) 2015 by Paul Scherrer Institut @n
Licensed under the Apache License, Version 2.0 (the "License"); @n
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
"""
import unittest
import os
import os.path
import tempfile
import shutil
import numpy as np
import pmsco.swarm as mo
import pmsco.project as mp
POP_SIZE = 5
class TestPopulation(unittest.TestCase):
def setUp(self):
self.test_dir = tempfile.mkdtemp()
self.domain = mp.Domain()
self.domain.add_param('A', 1.5, 1.0, 2.0, 0.1)
self.domain.add_param('B', 2.5, 2.0, 3.0, 0.1)
self.domain.add_param('C', 3.5, 3.0, 4.0, 0.1)
self.expected_names = ('A', 'B', 'C', '_particle', '_gen', '_model', '_rfac')
self.size = POP_SIZE
self.pop = mo.Population()
self.optimum1 = {'A': 1.045351, 'B': 2.346212, 'C': 3.873627}
def tearDown(self):
# after each test method
self.pop = None
shutil.rmtree(self.test_dir)
@classmethod
def setup_class(cls):
# before any methods in this class
pass
@classmethod
def teardown_class(cls):
# teardown_class() after any methods in this class
pass
def rfactor1(self, pos):
r = (pos['A'] - self.optimum1['A'])**2 \
+ (pos['B'] - self.optimum1['B'])**2 \
+ (pos['C'] - self.optimum1['C'])**2
r /= 3.0
return r
def test_setup(self):
self.pop.setup(self.size, self.domain)
self.assertItemsEqual(self.pop.pos.dtype.names, self.expected_names)
self.assertEqual(self.pop.pos.shape, (POP_SIZE,))
self.assertItemsEqual(np.arange(POP_SIZE), self.pop.pos['_particle'])
self.assertItemsEqual(np.zeros((POP_SIZE)), self.pop.pos['_gen'])
self.assertItemsEqual(np.arange(POP_SIZE), self.pop.pos['_model'])
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][0], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][1], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][2], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][3], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][4], 3)
self.assertEqual(0, self.pop.generation)
self.assertEqual(POP_SIZE, self.pop.model_count)
def test_setup_with_results(self):
data_dir = os.path.dirname(os.path.abspath(__file__))
data_file = os.path.join(data_dir, "test_swarm.setup_with_results.1.dat")
self.pop.setup(self.size, self.domain, data_file, False)
self.assertItemsEqual(self.pop.pos.dtype.names, self.expected_names)
self.assertEqual(self.pop.pos.shape, (POP_SIZE,))
self.assertEqual(0, self.pop.generation)
self.assertEqual(3, self.pop.model_count)
self.assertItemsEqual(np.arange(POP_SIZE), self.pop.pos['_particle'])
self.assertItemsEqual([-1, -1, 0, 0, 0], self.pop.pos['_gen'])
self.assertItemsEqual([-1, -2, 0, 1, 2], self.pop.pos['_model'])
self.assertAlmostEqual(0.3, self.pop.pos['_rfac'][0], 3)
self.assertAlmostEqual(0.6, self.pop.pos['_rfac'][1], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][2], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][3], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][4], 3)
self.assertAlmostEqual(1.3, self.pop.pos['A'][0], 3)
self.assertAlmostEqual(1.1, self.pop.pos['A'][1], 3)
self.assertAlmostEqual(1.5, self.pop.pos['A'][4], 3)
self.assertAlmostEqual(2.3, self.pop.pos['B'][0], 3)
self.assertAlmostEqual(2.1, self.pop.pos['B'][1], 3)
self.assertAlmostEqual(2.5, self.pop.pos['B'][4], 3)
self.assertAlmostEqual(3.3, self.pop.pos['C'][0], 3)
self.assertAlmostEqual(3.1, self.pop.pos['C'][1], 3)
self.assertAlmostEqual(3.5, self.pop.pos['C'][4], 3)
def test_setup_with_results_recalc(self):
data_dir = os.path.dirname(os.path.abspath(__file__))
data_file = os.path.join(data_dir, "test_swarm.setup_with_results.1.dat")
self.pop.setup(self.size, self.domain, data_file, True)
self.assertItemsEqual(self.pop.pos.dtype.names, self.expected_names)
self.assertEqual(self.pop.pos.shape, (POP_SIZE,))
self.assertEqual(self.pop.generation, 0)
self.assertEqual(self.pop.model_count, POP_SIZE)
self.assertItemsEqual(self.pop.pos['_particle'], np.arange(POP_SIZE))
self.assertItemsEqual(self.pop.pos['_gen'], [0, 0, 0, 0, 0])
self.assertItemsEqual(self.pop.pos['_model'], np.arange(POP_SIZE))
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][0], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][1], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][2], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][3], 3)
self.assertAlmostEqual(2.1, self.pop.pos['_rfac'][4], 3)
self.assertAlmostEqual(1.3, self.pop.pos['A'][0], 3)
self.assertAlmostEqual(1.1, self.pop.pos['A'][1], 3)
self.assertAlmostEqual(1.5, self.pop.pos['A'][4], 3)
self.assertAlmostEqual(2.3, self.pop.pos['B'][0], 3)
self.assertAlmostEqual(2.1, self.pop.pos['B'][1], 3)
self.assertAlmostEqual(2.5, self.pop.pos['B'][4], 3)
self.assertAlmostEqual(3.3, self.pop.pos['C'][0], 3)
self.assertAlmostEqual(3.1, self.pop.pos['C'][1], 3)
self.assertAlmostEqual(3.5, self.pop.pos['C'][4], 3)
def test_pos_gen(self):
self.pop.setup(self.size, self.domain)
for index, item in enumerate(self.pop.pos_gen()):
self.assertIsInstance(item, dict)
self.assertItemsEqual(item.keys(), self.expected_names)
self.assertEqual(item['_particle'], index)
def test_randomize(self):
self.pop.setup(self.size, self.domain)
self.pop.randomize()
m = np.mean(self.pop.pos['A'])
self.assertGreaterEqual(m, self.domain.min['A'])
self.assertLessEqual(m, self.domain.max['A'])
def test_seed(self):
self.pop.setup(self.size, self.domain)
self.pop.seed(self.domain.start)
self.assertAlmostEqual(self.pop.pos['A'][0], self.domain.start['A'], delta=0.001)
def test_best_friend(self):
self.pop.setup(self.size, self.domain)
self.pop.best['_rfac'] = np.arange(self.size)
friend = self.pop.best_friend(0)
self.assertNotIsInstance(friend, np.ndarray)
self.assertItemsEqual(friend.dtype.names, self.expected_names)
def test_advance_particle(self):
self.pop.setup(self.size, self.domain)
self.pop.pos['A'] = np.linspace(1.0, 2.0, POP_SIZE)
self.pop.pos['B'] = np.linspace(2.0, 3.0, POP_SIZE)
self.pop.pos['C'] = np.linspace(3.0, 4.0, POP_SIZE)
self.pop.pos['_rfac'] = np.linspace(2.0, 1.0, POP_SIZE)
self.pop.vel['A'] = np.linspace(-0.1, 0.1, POP_SIZE)
self.pop.vel['B'] = np.linspace(-0.1, 0.1, POP_SIZE)
self.pop.vel['C'] = np.linspace(-0.1, 0.1, POP_SIZE)
pos0 = self.pop.pos['A'][0]
self.pop.advance_particle(0)
pos1 = self.pop.pos['A'][0]
self.assertNotAlmostEqual(pos0, pos1, delta=0.001)
for key in ['A','B','C']:
for pos in self.pop.pos[key]:
self.assertGreaterEqual(pos, self.domain.min[key])
self.assertLessEqual(pos, self.domain.max[key])
def test_add_result(self):
self.pop.setup(self.size, self.domain)
i_sample = 1
i_result = 0
result = self.pop.pos[i_sample]
self.pop.add_result(result, 0.0)
self.assertEqual(self.pop.results.shape[0], 1)
self.assertItemsEqual(self.pop.results[i_result], result)
self.assertItemsEqual(self.pop.best[i_sample], result)
def test_is_converged(self):
self.pop.setup(self.size, self.domain)
self.assertFalse(self.pop.is_converged())
i_sample = 0
result = self.pop.pos[i_sample]
for i in range(POP_SIZE):
rfac = 1.0 - float(i)/POP_SIZE
self.pop.add_result(result, rfac)
self.assertFalse(self.pop.is_converged())
for i in range(POP_SIZE):
rfac = (1.0 - float(i)/POP_SIZE) / 1000.0
self.pop.add_result(result, rfac)
self.assertTrue(self.pop.is_converged())
def test_save_population(self):
self.pop.setup(self.size, self.domain)
filename = os.path.join(self.test_dir, "test_save_population.pop")
self.pop.save_population(filename)
def test_save_results(self):
self.pop.setup(self.size, self.domain)
i_sample = 1
result = self.pop.pos[i_sample]
self.pop.add_result(result, 1.0)
filename = os.path.join(self.test_dir, "test_save_results.dat")
self.pop.save_results(filename)
def test_save_array(self):
self.pop.setup(self.size, self.domain)
filename = os.path.join(self.test_dir, "test_save_array.pos")
self.pop.save_array(filename, self.pop.pos)
def test_load_array(self):
n = 3
filename = os.path.join(self.test_dir, "test_load_array")
self.pop.setup(self.size, self.domain)
# expected array
dt_exp = self.pop.get_model_dtype(self.domain.start)
a_exp = np.zeros((n,), dtype=dt_exp)
a_exp['A'] = np.linspace(0, 1, n)
a_exp['B'] = np.linspace(1, 2, n)
a_exp['C'] = np.linspace(3, 4, n)
a_exp['_rfac'] = np.linspace(5, 6, n)
a_exp['_gen'] = np.array([3, 4, 7])
a_exp['_particle'] = np.array([1, 0, 2])
a_exp['_model'] = np.array([3, 6, 1])
# test array is a expected array with different column order
dt_test = [('A', 'f4'), ('_particle', 'i4'), ('_rfac', 'f4'), ('C', 'f4'), ('_gen', 'i4'), ('B', 'f4'),
('_model', 'i4')]
names_test = [a[0] for a in dt_test]
a_test = np.zeros((n,), dtype=dt_test)
for name in names_test:
a_test[name] = a_exp[name]
header = " ".join(names_test)
np.savetxt(filename, a_test, fmt='%g', header=header)
result = np.zeros((n,), dtype=dt_exp)
result = self.pop.load_array(filename, result)
self.assertItemsEqual(result.dtype.names, a_exp.dtype.names)
for name in a_exp.dtype.names:
np.testing.assert_almost_equal(result[name], a_exp[name], err_msg=name)
def test_constrain_position(self):
# upper
pos1 = 11.0
vel1 = 5.0
min1 = 0.0
max1 = 10.0
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 're-enter')
self.assertAlmostEqual(pos2, 1.0)
self.assertAlmostEqual(vel2, vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'bounce')
self.assertAlmostEqual(pos2, 9.0)
self.assertAlmostEqual(vel2, -vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'scatter')
self.assertGreaterEqual(pos2, 6.0)
self.assertLessEqual(pos2, 10.0)
self.assertGreaterEqual(vel2, 0.0)
self.assertLessEqual(vel2, vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'stick')
self.assertAlmostEqual(pos2, max1)
self.assertAlmostEqual(vel2, max1 - (pos1 - vel1))
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'random')
self.assertGreaterEqual(pos2, 0.0)
self.assertLessEqual(pos2, 10.0)
self.assertGreaterEqual(vel2, -(max1 - min1))
self.assertLessEqual(vel2, (max1 - min1))
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'expand')
self.assertAlmostEqual(pos2, pos1)
self.assertAlmostEqual(vel2, vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max2, pos1)
self.assertRaises(ValueError, mo.Population.constrain_position, pos1, vel1, min1, max1, 'undefined')
# lower
pos1 = -1.0
vel1 = -5.0
min1 = 0.0
max1 = 10.0
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 're-enter')
self.assertAlmostEqual(pos2, 9.0)
self.assertAlmostEqual(vel2, vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'bounce')
self.assertAlmostEqual(pos2, 1.0)
self.assertAlmostEqual(vel2, -vel1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'scatter')
self.assertGreaterEqual(pos2, 0.0)
self.assertLessEqual(pos2, 4.0)
self.assertGreaterEqual(vel2, vel1)
self.assertLessEqual(vel2, 0.0)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'stick')
self.assertAlmostEqual(pos2, min1)
self.assertAlmostEqual(vel2, min1 - (pos1 - vel1))
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'random')
self.assertGreaterEqual(pos2, 0.0)
self.assertLessEqual(pos2, 10.0)
self.assertGreaterEqual(vel2, -(max1 - min1))
self.assertLessEqual(vel2, max1 - min1)
self.assertAlmostEqual(min1, min2)
self.assertAlmostEqual(max1, max2)
pos2, vel2, min2, max2 = mo.Population.constrain_position(pos1, vel1, min1, max1, 'expand')
self.assertAlmostEqual(pos2, pos1)
self.assertAlmostEqual(vel2, vel1)
self.assertAlmostEqual(min2, pos1)
self.assertAlmostEqual(max2, max1)
self.assertRaises(ValueError, mo.Population.constrain_position, pos1, vel1, min1, max1, 'undefined')
def test_convergence_1(self):
self.pop.setup(self.size, self.domain)
self.pop.pos['A'] = np.linspace(1.0, 2.0, POP_SIZE)
self.pop.pos['B'] = np.linspace(2.0, 3.0, POP_SIZE)
self.pop.pos['C'] = np.linspace(3.0, 4.0, POP_SIZE)
self.pop.pos['_rfac'] = np.linspace(2.0, 1.0, POP_SIZE)
self.pop.vel['A'] = np.linspace(-0.1, 0.1, POP_SIZE)
self.pop.vel['B'] = np.linspace(-0.1, 0.1, POP_SIZE)
self.pop.vel['C'] = np.linspace(-0.1, 0.1, POP_SIZE)
for i in range(10):
self.pop.advance_population()
for pos in self.pop.pos:
self.pop.add_result(pos, self.rfactor1(pos))
for pos in self.pop.pos:
self.assertLess(pos['_rfac'], 0.2)
if __name__ == '__main__':
unittest.main()

View File

@ -0,0 +1,5 @@
# A B C _gen _model _rfac
1.1 2.1 3.1 1 1111 0.6
1.2 2.2 3.2 2 2222 1.5
1.3 2.3 3.3 3 3333 0.3
1.4 2.4 3.4 4 4444 1.0