public release 2.2.0 - see README.md and CHANGES.md for details

This commit is contained in:
2020-09-04 16:22:42 +02:00
parent fbd2d4fa8c
commit 7c61eb1b41
67 changed files with 2934 additions and 682 deletions

View File

@ -17,7 +17,7 @@ pip install --user periodictable
@author Matthias Muntwiler
@copyright (c) 2015-19 by Paul Scherrer Institut @n
@copyright (c) 2015-20 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
@ -54,14 +54,14 @@ if sys.version_info[0] >= 3:
else:
_SYMBOL_TYPE = 'S2'
## numpy.array datatype of Cluster.data array
DTYPE_CLUSTER_INTERNAL = [('i', 'i4'), ('t', 'i4'), ('s', _SYMBOL_TYPE), ('x', 'f4'), ('y', 'f4'), ('z', 'f4'),
('e', 'u1'), ('q', 'f4'), ('c', 'i4')]
## file format of internal Cluster.data array
## numpy.array datatype of internal Cluster.data array
DTYPE_CLUSTER_INTERNAL = [('i', 'i4'), ('t', 'i4'), ('s', _SYMBOL_TYPE), ('c', 'i4'),
('x', 'f4'), ('y', 'f4'), ('z', 'f4'), ('e', 'u1'), ('q', 'f4')]
## string formatting of native file format
FMT_CLUSTER_INTERNAL = ["%5u", "%2u", "%s", "%5u", "%7.3f", "%7.3f", "%7.3f", "%1u", "%7.3f"]
## field (column) names of internal Cluster.data array
## field (column) names of native file format
FIELDS_CLUSTER_INTERNAL = ['i', 't', 's', 'c', 'x', 'y', 'z', 'e', 'q']
## column names for export
## column names of native file format
NAMES_CLUSTER_INTERNAL = {'i': 'index', 't': 'element', 's': 'symbol', 'c': 'class', 'x': 'x', 'y': 'y', 'z': 'z',
'e': 'emitter', 'q': 'charge'}
@ -178,12 +178,12 @@ class Cluster(object):
# @arg @c 'i' (int) atom index (1-based)
# @arg @c 't' (int) atom type (chemical element number)
# @arg @c 's' (string) chemical element symbol
# @arg @c 'c' (int) scatterer class
# @arg @c 'x' (float32) x coordinate of the atom position
# @arg @c 'y' (float32) t coordinate of the atom position
# @arg @c 'z' (float32) z coordinate of the atom position
# @arg @c 'e' (uint8) 1 = emitter, 0 = regular atom
# @arg @c 'q' (float32) charge/ionicity
# @arg @c 'c' (int) scatterer class
## @var comment (str)
# one-line comment that can be included in some cluster files
@ -227,7 +227,7 @@ class Cluster(object):
"""
self.rmax = r
def build_element(self, index, element_number, x, y, z, emitter, charge=0., scatterer=0):
def build_element(self, index, element_number, x, y, z, emitter, charge=0., scatterer_class=0):
"""
build a tuple in the format of the internal data array.
@ -241,10 +241,10 @@ class Cluster(object):
@param charge: (float) ionicity. default = 0
@param scatterer: (int) scatterer class. default = 0.
@param scatterer_class: (int) scatterer class. default = 0.
"""
symbol = pt.elements[element_number].symbol
element = (index, element_number, symbol, x, y, z, int(emitter), charge, scatterer)
element = (index, element_number, symbol, scatterer_class, x, y, z, int(emitter), charge)
return element
def add_atom(self, atomtype, v_pos, is_emitter=False, charge=0.):
@ -284,17 +284,18 @@ class Cluster(object):
"""
r_great = max(self.rmax, np.linalg.norm(v_pos))
n0 = self.data.shape[0] + 1
n1 = max(int(r_great / np.linalg.norm(v_lat1)) + 1, 3) * 2
n2 = max(int(r_great / np.linalg.norm(v_lat2)) + 1, 3) * 2
nn = 0
buf = np.empty((2 * n1 + 1) * (2 * n2 + 1), dtype=self.dtype)
for i1 in range(-n1, n1 + 1):
for i2 in range(-n2, n2 + 1):
v = v_pos + v_lat1 * i1 + v_lat2 * i2
if np.linalg.norm(v) <= self.rmax:
buf[nn] = self.build_element(nn + n0, atomtype, v[0], v[1], v[2], 0)
nn += 1
buf = np.resize(buf, nn)
n1 = max(int(r_great / np.linalg.norm(v_lat1)) + 1, 4) * 3
n2 = max(int(r_great / np.linalg.norm(v_lat2)) + 1, 4) * 3
idx = np.mgrid[-n1:n1+1, -n2:n2+1]
idx = idx.reshape(idx.shape[0], -1)
lat = np.array([v_lat1, v_lat2])
v = v_pos + np.matmul(idx.T, lat)
rsq = np.sum(np.square(v), axis=-1)
b1 = rsq <= self.rmax**2
sel = b1.nonzero()[0]
buf = np.empty((len(sel)), dtype=self.dtype)
for nn, ii in enumerate(sel):
buf[nn] = self.build_element(nn + n0, atomtype, v[ii, 0], v[ii, 1], v[ii, 2], 0)
self.data = np.append(self.data, buf)
def add_bulk(self, atomtype, v_pos, v_lat1, v_lat2, v_lat3, z_surf=0.0):
@ -322,16 +323,18 @@ class Cluster(object):
n1 = max(int(r_great / np.linalg.norm(v_lat1)) + 1, 4) * 3
n2 = max(int(r_great / np.linalg.norm(v_lat2)) + 1, 4) * 3
n3 = max(int(r_great / np.linalg.norm(v_lat3)) + 1, 4) * 3
nn = 0
buf = np.empty((2 * n1 + 1) * (2 * n2 + 1) * (n3 + 1), dtype=self.dtype)
for i1 in range(-n1, n1 + 1):
for i2 in range(-n2, n2 + 1):
for i3 in range(-n3, n3 + 1):
v = v_pos + v_lat1 * i1 + v_lat2 * i2 + v_lat3 * i3
if np.linalg.norm(v) <= self.rmax and v[2] <= z_surf:
buf[nn] = self.build_element(nn + n0, atomtype, v[0], v[1], v[2], 0)
nn += 1
buf = np.resize(buf, nn)
idx = np.mgrid[-n1:n1+1, -n2:n2+1, -n3:n3+1]
idx = idx.reshape(idx.shape[0], -1)
lat = np.array([v_lat1, v_lat2, v_lat3])
v = v_pos + np.matmul(idx.T, lat)
rsq = np.sum(np.square(v), axis=-1)
b1 = rsq <= self.rmax**2
b2 = v[:, 2] <= z_surf
ba = np.all([b1, b2], axis=0)
sel = ba.nonzero()[0]
buf = np.empty((len(sel)), dtype=self.dtype)
for nn, ii in enumerate(sel):
buf[nn] = self.build_element(nn + n0, atomtype, v[ii, 0], v[ii, 1], v[ii, 2], 0)
self.data = np.append(self.data, buf)
def add_cluster(self, cluster, check_rmax=False, check_unique=False, tol=0.001):
@ -426,15 +429,47 @@ class Cluster(object):
idx = np.where(b_all)
self.data['z'][idx] += z_shift
return idx
return idx[0]
def translate(self, vector, element=0):
def get_center(self, element=None):
"""
get the geometric center of the cluster or a class of atoms.
@param element: chemical element number (int) or symbol (str)
if atoms of a specific element should be considered only.
by default (element == None or 0 or ""),
all atoms are included in the calculation.
@return: (numpy.ndarray) 3-dimensional vector.
"""
if element:
try:
sel = self.data['t'] == int(element)
except ValueError:
sel = self.data['s'] == element
else:
sel = np.ones_like(self.data['t'])
idx = np.where(sel)
center = np.zeros(3)
center[0] = np.mean(self.data['x'][idx])
center[1] = np.mean(self.data['y'][idx])
center[2] = np.mean(self.data['z'][idx])
return center
def translate(self, vector, element=None):
"""
translate the cluster or all atoms of a specified element.
translation shifts each selected atom by the given vector.
@param vector: (numpy.ndarray) 3-dimensional displacement vector.
@param element: (int) chemical element number if atoms of a specific element should be affected.
by default (element = 0), all atoms are moved.
@param element: chemical element number (int) or symbol (str)
if atoms of a specific element should be affected only.
by default (element == None or 0 or ""),
all atoms are translated.
@return: (numpy.ndarray) indices of the atoms that have been shifted.
"""
if element:
@ -449,7 +484,7 @@ class Cluster(object):
self.data['y'][idx] += vector[1]
self.data['z'][idx] += vector[2]
return idx
return idx[0]
def matrix_transform(self, matrix):
"""
@ -461,47 +496,49 @@ class Cluster(object):
@return: None
"""
for atom in self.data:
v = np.matrix([atom['x'], atom['y'], atom['z']])
w = matrix * v.transpose()
atom['x'] = float(w[0])
atom['y'] = float(w[1])
atom['z'] = float(w[2])
pos = np.empty((3, self.data.shape[0]), np.float32)
pos[0, :] = self.data['x']
pos[1, :] = self.data['y']
pos[2, :] = self.data['z']
pos = np.matmul(matrix, pos)
self.data['x'] = pos[0, :]
self.data['y'] = pos[1, :]
self.data['z'] = pos[2, :]
def rotate_x(self, angle):
"""
rotate cluster about the surface normal axis
rotate cluster about the x-axis
@param angle (float) in degrees
"""
angle = math.radians(angle)
s = math.sin(angle)
c = math.cos(angle)
matrix = np.matrix([[1, 0, 0], [0, c, -s], [0, s, c]])
matrix = np.array([[1, 0, 0], [0, c, -s], [0, s, c]])
self.matrix_transform(matrix)
def rotate_y(self, angle):
"""
rotate cluster about the surface normal axis
rotate cluster about the y-axis
@param angle (float) in degrees
"""
angle = math.radians(angle)
s = math.sin(angle)
c = math.cos(angle)
matrix = np.matrix([[c, 0, s], [0, 1, 0], [-s, 0, c]])
matrix = np.array([[c, 0, s], [0, 1, 0], [-s, 0, c]])
self.matrix_transform(matrix)
def rotate_z(self, angle):
"""
rotate cluster about the surface normal axis
rotate cluster about the z-axis (surface normal)
@param angle (float) in degrees
"""
angle = math.radians(angle)
s = math.sin(angle)
c = math.cos(angle)
matrix = np.matrix([[c, -s, 0], [s, c, 0], [0, 0, 1]])
matrix = np.array([[c, -s, 0], [s, c, 0], [0, 0, 1]])
self.matrix_transform(matrix)
def find_positions(self, pos, tol=0.001):
@ -794,6 +831,53 @@ class Cluster(object):
idx = self.data['e'] != 0
return np.sum(idx)
def calc_scattering_angles(self, index_emitter, radius):
"""
calculate forward-scattering angles of the cluster atoms
for each atom within a given radius of the emitter,
the connecting vector between emitter and scatterer is calculated
and returned in cartesian and polar coordinates.
@param index_emitter: atom index of the emitter.
all angles are calculated with respect to this atom.
@param radius: include only atoms within this radius of the emitter.
@note back-scattering angles can be obtained by inverting the angle on the unit sphere:
th' = 180 - th, ph' = -ph.
@return dictionary with results.
each item is a numpy.ndarray of shape (N, M)
where N is the number of scatterers
and M = 3 for dict['xyz'] and M = 1 otherwise.
@arg dict['index']: atom index into the cluster array.
@arg dict['xyz']: connecting vector between the emitter and the atom in cartesian coordinates.
@arg dict['dist']: distance between the emitter and the atom.
@arg dict['polar']: polar angle with respect to the z-axis.
@arg dict['azimuth']: azimuthal angle with respect to the x-axis.
"""
# position of emitter atom
em = self.data[index_emitter]
em = np.asarray((em['x'], em['y'], em['z']))
# relative positions of scattering atoms
xyz = self.get_positions()
xyz -= em
dist = np.linalg.norm(xyz, axis=1)
sel1 = dist <= radius
sel2 = dist > 0.
idx = np.where(np.all([sel1, sel2], axis=0))
xyz = xyz[idx]
dist = dist[idx]
# angles
v1 = np.asarray([0, 0, 1])
v2 = np.transpose(xyz / dist.reshape((dist.shape[0], 1)))
th = np.degrees(np.arccos(np.clip(np.dot(v1, v2), -1., 1.)))
ph = np.degrees(np.arctan2(v2[1], v2[0]))
return {'index': idx[0], 'xyz': xyz, 'dist': dist, 'polar': th, 'azimuth': ph}
def load_from_file(self, f, fmt=FMT_DEFAULT):
"""
load a cluster from a file created by the scattering program.
@ -848,7 +932,7 @@ class Cluster(object):
else:
raise ValueError("unknown file format {}".format(fmt))
data = np.genfromtxt(f, dtype=dtype, skip_header=sh)
data = np.atleast_1d(np.genfromtxt(f, dtype=dtype, skip_header=sh))
if fmt == FMT_PHAGEN_IN and data['t'][-1] < 1:
data = data[:-1]
@ -920,7 +1004,7 @@ class Cluster(object):
or left at the default value 0 in which case PMSCO sets the correct values.
if the scattering factors are loaded from existing files,
the atom class corresponds to the key of the pmsco.project.Params.phase_files dictionary.
the atom class corresponds to the key of the pmsco.project.CalculatorParams.phase_files dictionary.
in this case the meaning of the class value is up to the project,
and the class must be set either by the cluster generator
or the project's after_atomic_scattering hook.
@ -956,7 +1040,7 @@ class Cluster(object):
the other cluster must contain the same atoms (same coordinates) in a possibly random order.
the atoms of this and the other cluster are matched up by sorting them by coordinate.
atomic scattering calculators often change the order of atoms in a cluster based on symmetry,
atomic scattering calculators often change the order of atoms in a cluster based on domain,
and return atom classes versus atomic coordinates.
this method allows to import the atom classes into the original cluster.
@ -1071,9 +1155,9 @@ class ClusterGenerator(object):
def count_emitters(self, model, index):
"""
return the number of emitter configurations for a particular model, scan and symmetry.
return the number of emitter configurations for a particular model, scan and domain.
the number of emitter configurations may depend on the model parameters, scan index and symmetry index.
the number of emitter configurations may depend on the model parameters, scan index and domain index.
by default, the method returns 1, which means that there is only one emitter configuration.
emitter configurations are mainly a way to distribute the calculations to multiple processes
@ -1100,9 +1184,9 @@ class ClusterGenerator(object):
@param index (named tuple CalcID) calculation index.
the method should consider only the following attributes:
@arg @c scan scan index (index into Project.scans)
@arg @c sym symmetry index (index into Project.symmetries)
@arg @c emit emitter index must be -1.
@arg scan scan index (index into Project.scans)
@arg domain domain index (index into Project.domains)
@arg emit emitter index must be -1.
@return number of emitter configurations.
this implementation returns the default value of 1.
@ -1114,23 +1198,23 @@ class ClusterGenerator(object):
create a Cluster object given the model parameters and calculation index.
the generated cluster will typically depend on the model parameters.
depending on the project, it may also depend on the scan index, symmetry index and emitter index.
depending on the project, it may also depend on the scan index, domain index and emitter index.
the scan index can be used to generate a different cluster for different scan geometry,
e.g., if some atoms can be excluded due to a longer mean free path.
if this is not the case for the specific project, the scan index can be ignored.
the symmetry index may select a particular domain that has a different atomic arrangement.
in this case, depending on the value of index.sym, the function must generate a cluster corresponding
to the particular domain/symmetry.
the method can ignore the symmetry index if the project defines only one symmetry,
or if the symmetry does not correspond to a different atomic structure.
the domain index may select a particular domain that has a different atomic arrangement.
in this case, depending on the value of index.domain, the function must generate a cluster corresponding
to the particular domain.
the method can ignore the domain index if the project defines only one domain,
or if the domain does not correspond to a different atomic structure.
the emitter index selects a particular emitter configuration.
depending on the value of the emitter index, the method must react differently:
1. if the value is -1, return the full cluster and mark all inequivalent emitter atoms.
emitters which are reproduced by a symmetry expansion in combine_emitters() should not be marked.
emitters which are reproduced by a domain expansion in combine_emitters() should not be marked.
the full diffraction scan will be calculated in one calculation.
2. if the value is greater or equal to zero, generate the cluster with the emitter configuration
@ -1152,9 +1236,9 @@ class ClusterGenerator(object):
@param index (named tuple CalcID) calculation index.
the method should consider only the following attributes:
@arg @c scan scan index (index into Project.scans)
@arg @c sym symmetry index (index into Project.symmetries)
@arg @c emit emitter index.
@arg scan scan index (index into Project.scans)
@arg domain domain index (index into Project.domains)
@arg emit emitter index.
if -1, generate the full cluster and mark all emitters.
if greater or equal to zero, the value is a zero-based index of the emitter configuration.