Skip to content

coordinates

Coordinates

Coordinates(coords: ndarray | list[list[float]] = None)

Responsible for handling StructureBase coordinates by storing in a numpy.ndarray with shape (n, 3) where n is the number of atoms in the structure and the 3 dimensions represent x, y, and z coordinates

Parameters:

  • coords (ndarray | list[list[float]], default: None ) –

    The coordinates to store with shape (N, 3). If none are passed an empty container will be generated

Source code in symdesign/structure/coordinates.py
18
19
20
21
22
23
24
25
26
27
28
29
30
def __init__(self, coords: np.ndarray | list[list[float]] = None):
    """Construct the instance

    Args:
        coords: The coordinates to store with shape (N, 3). If none are passed an empty container will be generated
    """
    if coords is None:
        self.coords = np.array([])
    elif not isinstance(coords, (np.ndarray, list)):
        raise TypeError(f"Can't initialize {self.__class__.__name__} with {type(coords).__name__}. Type must be a "
                        f'numpy.ndarray of float with shape (n, 3) or list[list[float]]')
    else:
        self.coords = np.array(coords, np.float_)

delete

delete(indices: Sequence[int])

Delete coordinates from the instance

Parameters:

  • indices (Sequence[int]) –

    The indices to delete from the Coords array

Source code in symdesign/structure/coordinates.py
32
33
34
35
36
37
38
def delete(self, indices: Sequence[int]):
    """Delete coordinates from the instance

    Args:
        indices: The indices to delete from the Coords array
    """
    self.coords = np.delete(self.coords, indices, axis=0)

insert

insert(at: int, new_coords: ndarray | list[list[float]])

Insert additional coordinates into the instance

Parameters:

  • at (int) –

    The index to perform the insert at

  • new_coords (ndarray | list[list[float]]) –

    The coordinate values to insert into Coords

Source code in symdesign/structure/coordinates.py
40
41
42
43
44
45
46
47
def insert(self, at: int, new_coords: np.ndarray | list[list[float]]):
    """Insert additional coordinates into the instance

    Args:
        at: The index to perform the insert at
        new_coords: The coordinate values to insert into Coords
    """
    self.coords = np.concatenate((self.coords[:at], new_coords, self.coords[at:]))

append

append(new_coords: ndarray | list[list[float]])

Append additional coordinates onto the instance

Parameters:

  • new_coords (ndarray | list[list[float]]) –

    The coordinate values to append to Coords

Source code in symdesign/structure/coordinates.py
49
50
51
52
53
54
55
def append(self, new_coords: np.ndarray | list[list[float]]):
    """Append additional coordinates onto the instance

    Args:
        new_coords: The coordinate values to append to Coords
    """
    self.coords = np.concatenate((self.coords, new_coords))

replace

replace(indices: Sequence[int], new_coords: ndarray | list[list[float]])

Replace existing coordinates in the instance with new coordinates

Parameters:

  • indices (Sequence[int]) –

    The indices to replace in the Coords array

  • new_coords (ndarray | list[list[float]]) –

    The coordinate values to replace in Coords

Source code in symdesign/structure/coordinates.py
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
def replace(self, indices: Sequence[int], new_coords: np.ndarray | list[list[float]]):
    """Replace existing coordinates in the instance with new coordinates

    Args:
        indices: The indices to replace in the Coords array
        new_coords: The coordinate values to replace in Coords
    """
    try:
        self.coords[indices] = new_coords
    except ValueError as error:
        # They are probably different lengths or another numpy indexing/setting issue
        if len(self.coords) == 0:  # There are no coords. Use the .set() mechanism
            self.set(new_coords)
        else:
            raise ValueError(
                f"The selected indices aren't the same shape as the 'new_coords': {error}")

set

set(coords: ndarray | list[list[float]])

Set self.coords to the provided coordinates

Parameters:

  • coords (ndarray | list[list[float]]) –

    The coordinate values to set

Source code in symdesign/structure/coordinates.py
74
75
76
77
78
79
80
def set(self, coords: np.ndarray | list[list[float]]):
    """Set self.coords to the provided coordinates

    Args:
        coords: The coordinate values to set
    """
    self.coords = coords

guide_superposition

guide_superposition(fixed_coords: ndarray, moving_coords: ndarray, number_of_points: int = 4) -> tuple[ndarray, ndarray]

TTakes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations should be applied to the "moving_coords" to place them in the setting of the fixed_coords

This function implements a more general variant of the method from: R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216 This version has been augmented slightly. The version in the original paper only considers rotation and translation and does not allow the coordinates of either object to be rescaled (multiplication by a scalar). (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries form a vector (of length sin(θ/2)), pointing along the axis of rotation. Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

MIT License. Copyright (c) 2016, Andrew Jewett Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

Parameters:

  • fixed_coords (ndarray) –

    The coordinates for the 'frozen' object

  • moving_coords (ndarray) –

    The coordinates for the 'mobile' object

  • number_of_points (int, default: 4 ) –

    The number of points included in the coordinate sets

Returns: rotation, translation_vector

Source code in symdesign/structure/coordinates.py
 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
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
@jit(nopython=True)  # , cache=True)
def guide_superposition(fixed_coords: np.ndarray, moving_coords: np.ndarray, number_of_points: int = 4) \
        -> tuple[np.ndarray, np.ndarray]:
    """TTakes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation
    operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations
    should be applied to the "moving_coords" to place them in the setting of the fixed_coords

    This function implements a more general variant of the method from:
    R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216
    This version has been augmented slightly. The version in the original paper only considers rotation and translation
    and does not allow the coordinates of either object to be rescaled (multiplication by a scalar).
    (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

    The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries
    form a vector (of length sin(θ/2)), pointing along the axis of rotation.
    Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

    MIT License. Copyright (c) 2016, Andrew Jewett
    Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated
    documentation files (the "Software"), to deal in the Software without restriction, including without limitation the
    rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to
    permit persons to whom the Software is furnished to do so, subject to the following conditions:

    The above copyright notice and this permission notice shall be included in all copies or substantial portions of the
    Software.

    THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE
    WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS
    OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
    OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

    Args:
        fixed_coords: The coordinates for the 'frozen' object
        moving_coords: The coordinates for the 'mobile' object
        number_of_points: The number of points included in the coordinate sets
    Returns:
        rotation, translation_vector
    """
    # number_of_points = fixed_coords.shape[0]
    # if number_of_points != moving_coords.shape[0]:
    #     raise ValueError(f'{guide_superposition.__name__}: Inputs should have the same size. '
    #                      f'Input 1={number_of_points}, 2={moving_coords.shape[0]}')

    # Find the center of mass of each object:
    # center_of_mass_fixed = fixed_coords.mean(axis=0)
    # center_of_mass_moving = moving_coords.mean(axis=0)
    center_of_mass_fixed = fixed_coords.sum(axis=0)
    center_of_mass_moving = moving_coords.sum(axis=0)
    center_of_mass_fixed /= number_of_points
    center_of_mass_moving /= number_of_points

    # Subtract the centers-of-mass from the original coordinates for each object
    # Translate the center of mass to the origin
    fixed_coords_at_origin = fixed_coords - center_of_mass_fixed
    moving_coords_at_origin = moving_coords - center_of_mass_moving

    # Calculate the "m" array from the Diamond paper (equation 16)
    m = moving_coords_at_origin.T @ fixed_coords_at_origin

    # Calculate "v" (equation 18)
    # v = np.empty(3)
    # v[0] = m[1, 2] - m[2, 1]
    # v[1] = m[2, 0] - m[0, 2]
    # v[2] = m[0, 1] - m[1, 0]
    v = [m[1, 2] - m[2, 1], m[2, 0] - m[0, 2], m[0, 1] - m[1, 0]]

    # Calculate "P" (equation 22)
    matrix_p = np.zeros((4, 4))
    # Calculate "q" (equation 17)
    # q = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[:3, :3] = m + m.T - 2*identity_matrix*np.trace(m)
    # matrix_p[:3, :3] = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[3, :3] = v
    matrix_p[:3, 3] = v
    # [[ q[0, 0] q[0, 1] q[0, 2] v[0] ]
    #  [ q[1, 0] q[1, 1] q[1, 2] v[1] ]
    #  [ q[2, 0] q[2, 1] q[2, 2] v[2] ]
    #  [ v[0]    v[1]    v[2]    0    ]]

    # Calculate "p" - optimal_quat
    # "p" contains the optimal rotation (in backwards-quaternion format)
    # (Note: A discussion of various quaternion conventions is included below)
    # try:
    # The eigenvalues/eigenvector are returned as 1D array in ascending order; largest is last
    a_eigenvals, aa_eigenvects = np.linalg.eigh(matrix_p)
    # except np.linalg.LinAlgError:
    #     singular = True  # I have never seen this happen
    # Pull out the largest magnitude
    optimal_quat = aa_eigenvects[:, -1]
    # normalize the vector
    # (It should be normalized already, but just in case it is not, do it again)
    # optimal_quat /= np.linalg.norm(optimal_quat)

    # Calculate the rotation matrix corresponding to "optimal_quat" which is in scipy quaternion format
    # """
    rotation_matrix = np.empty((3, 3))
    quat0, quat1, quat2, quat3 = optimal_quat
    quat2_0, quat2_1, quat2_2, quat2_3 = optimal_quat**2
    rotation_matrix[0, 1] = quat0*quat1 - quat2*quat3  # 2*( )
    rotation_matrix[1, 0] = quat0*quat1 + quat2*quat3  # 2*( )
    rotation_matrix[1, 2] = quat1*quat2 - quat0*quat3  # 2*( )
    rotation_matrix[2, 1] = quat1*quat2 + quat0*quat3  # 2*( )
    rotation_matrix[0, 2] = quat0*quat2 + quat1*quat3  # 2*( )
    rotation_matrix[2, 0] = quat0*quat2 - quat1*quat3  # 2*( )
    rotation_matrix *= 2
    rotation_matrix[0, 0] = quat2_0 - quat2_1 - quat2_2 + quat2_3
    rotation_matrix[1, 1] = -quat2_0 + quat2_1 - quat2_2 + quat2_3
    rotation_matrix[2, 2] = -quat2_0 - quat2_1 + quat2_2 + quat2_3
    # rotation_matrix[0, 0] = quat0*quat0 - quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[1, 1] = -quat0*quat0 + quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[2, 2] = -quat0*quat0 - quat1*quat1 + quat2*quat2 + quat3*quat3
    """
    # Alternatively, in modern python versions, this code also works:
    rotation_matrix = Rotation.from_quat(optimal_quat).as_matrix()
    """
    # input(f'{rotation_matrix_} {rotation_matrix}')
    # input(rotation_matrix_ == rotation_matrix)
    # """
    # Finally compute the RMSD between the two coordinate sets:
    # First compute E0 from equation 24 of the paper
    # e0 = np.sum((fixed_coords_at_origin - moving_coords_at_origin) ** 2)
    # sum_sqr_dist = max(0, ((fixed_coords_at_origin-moving_coords_at_origin) ** 2).sum() - 2.*pPp)

    # Lastly, calculate the translational offset:
    # Recall that:
    # RMSD=sqrt((Σ_i  w_i * |X_i - (Σ_j c*R_ij*x_j + T_i))|^2) / (Σ_j w_j))
    #    =sqrt((Σ_i  w_i * |X_i - x_i'|^2) / (Σ_j w_j))
    #  where
    # x_i' = Σ_j c*R_ij*x_j + T_i
    #      = Xcm_i + c*R_ij*(x_j - xcm_j)
    #  and Xcm and xcm = center_of_mass for the frozen and mobile point clouds
    #                  = center_of_mass_fixed[]       and       center_of_mass_moving[],  respectively
    # Hence:
    #  T_i = Xcm_i - Σ_j c*R_ij*xcm_j  =  a_translate[i]

    # a_translate = center_of_mass_fixed - np.matmul(c * rotation_matrix, center_of_mass_moving).T.reshape(3,)

    # Calculate the translation
    translation = center_of_mass_fixed - rotation_matrix@center_of_mass_moving

    return rotation_matrix, translation

superposition3d

superposition3d(fixed_coords: ndarray, moving_coords: ndarray) -> tuple[float, ndarray, ndarray]

Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations should be applied to the "moving_coords" to place them in the setting of the fixed_coords

This function implements a more general variant of the method from: R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216 This version has been augmented slightly. The version in the original paper only considers rotation and translation and does not allow the coordinates of either object to be rescaled (multiplication by a scalar). (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries form a vector (of length sin(θ/2)), pointing along the axis of rotation. Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

MIT License. Copyright (c) 2016, Andrew Jewett Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

Parameters:

  • fixed_coords (ndarray) –

    The coordinates for the 'frozen' object

  • moving_coords (ndarray) –

    The coordinates for the 'mobile' object

Raises: ValueError: If coordinates are not the same length Returns: rmsd, rotation, translation_vector

Source code in symdesign/structure/coordinates.py
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
@jit(nopython=True)  # , cache=True)
def superposition3d(fixed_coords: np.ndarray, moving_coords: np.ndarray) -> tuple[float, np.ndarray, np.ndarray]:
    """Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation
    operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations
    should be applied to the "moving_coords" to place them in the setting of the fixed_coords

    This function implements a more general variant of the method from:
    R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216
    This version has been augmented slightly. The version in the original paper only considers rotation and translation
    and does not allow the coordinates of either object to be rescaled (multiplication by a scalar).
    (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

    The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries
    form a vector (of length sin(θ/2)), pointing along the axis of rotation.
    Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

    MIT License. Copyright (c) 2016, Andrew Jewett
    Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated
    documentation files (the "Software"), to deal in the Software without restriction, including without limitation the
    rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to
    permit persons to whom the Software is furnished to do so, subject to the following conditions:

    The above copyright notice and this permission notice shall be included in all copies or substantial portions of the
    Software.

    THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE
    WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS
    OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
    OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

    Args:
        fixed_coords: The coordinates for the 'frozen' object
        moving_coords: The coordinates for the 'mobile' object
    Raises:
        ValueError: If coordinates are not the same length
    Returns:
        rmsd, rotation, translation_vector
    """
    number_of_points = fixed_coords.shape[0]
    if number_of_points != moving_coords.shape[0]:
        raise ValueError('superposition3d: Inputs should have the same size')  # . '
        #                  f'Input 1={number_of_points}, 2={moving_coords.shape[0]}')

    # convert weights into array
    # if a_weights is None or len(a_weights) == 0:
    # a_weights = np.full((number_of_points, 1), 1.)
    # sum_weights = float(number_of_points)
    # else:  # reshape a_eights so multiplications are done column-wise
    #     a_weights = np.array(a_weights).reshape(number_of_points, 1)
    #     sum_weights = np.sum(a_weights, axis=0)

    # Find the center of mass of each object:
    center_of_mass_fixed = fixed_coords.sum(axis=0)
    center_of_mass_moving = moving_coords.sum(axis=0)

    # Subtract the centers-of-mass from the original coordinates for each object
    # if sum_weights != 0:
    # try:
    center_of_mass_fixed /= number_of_points
    center_of_mass_moving /= number_of_points
    # except ZeroDivisionError:
    #     pass  # The weights are a total of zero which is allowed algorithmically, but not possible

    # Translate the center of mass to the origin
    fixed_coords_at_origin = fixed_coords - center_of_mass_fixed
    moving_coords_at_origin = moving_coords - center_of_mass_moving

    # Calculate the "m" array from the Diamond paper (equation 16)
    m = moving_coords_at_origin.T @ fixed_coords_at_origin

    # Calculate "v" (equation 18)
    # v = np.empty(3)
    # v[0] = m[1, 2] - m[2, 1]
    # v[1] = m[2, 0] - m[0, 2]
    # v[2] = m[0, 1] - m[1, 0]
    v = [m[1, 2] - m[2, 1], m[2, 0] - m[0, 2], m[0, 1] - m[1, 0]]

    # Calculate "P" (equation 22)
    matrix_p = np.zeros((4, 4))
    # Calculate "q" (equation 17)
    # q = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[:3, :3] = m + m.T - 2*identity_matrix*np.trace(m)
    # matrix_p[:3, :3] = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[3, :3] = v
    matrix_p[:3, 3] = v
    # [[ q[0, 0] q[0, 1] q[0, 2] v[0] ]
    #  [ q[1, 0] q[1, 1] q[1, 2] v[1] ]
    #  [ q[2, 0] q[2, 1] q[2, 2] v[2] ]
    #  [ v[0]    v[1]    v[2]    0    ]]

    # Calculate "p" - optimal_quat
    # "p" contains the optimal rotation (in backwards-quaternion format)
    # (Note: A discussion of various quaternion conventions is included below)
    if number_of_points < 2:
        # Specify the default values for p, pPp
        optimal_quat = np.array([0., 0., 0., 1.])  # p = [0,0,0,1]    default value
        pPp = 0.  # = p^T * P * p    (zero by default)
    else:
        # try:
        # The eigenvalues/eigenvector are returned as 1D array in ascending order; largest is last
        a_eigenvals, aa_eigenvects = np.linalg.eigh(matrix_p)
        # except np.linalg.LinAlgError:
        #     singular = True  # I have never seen this happen
        # Pull out the largest magnitude
        pPp = a_eigenvals[-1]
        optimal_quat = aa_eigenvects[:, -1]
        # normalize the vector
        # (It should be normalized already, but just in case it is not, do it again)
        # optimal_quat /= np.linalg.norm(optimal_quat)

    # Calculate the rotation matrix corresponding to "optimal_quat" which is in scipy quaternion format
    # """
    rotation_matrix = np.empty((3, 3))
    quat0, quat1, quat2, quat3 = optimal_quat
    quat2_0, quat2_1, quat2_2, quat2_3 = optimal_quat**2
    rotation_matrix[0, 1] = quat0*quat1 - quat2*quat3  # 2*( )
    rotation_matrix[1, 0] = quat0*quat1 + quat2*quat3  # 2*( )
    rotation_matrix[1, 2] = quat1*quat2 - quat0*quat3  # 2*( )
    rotation_matrix[2, 1] = quat1*quat2 + quat0*quat3  # 2*( )
    rotation_matrix[0, 2] = quat0*quat2 + quat1*quat3  # 2*( )
    rotation_matrix[2, 0] = quat0*quat2 - quat1*quat3  # 2*( )
    rotation_matrix *= 2
    rotation_matrix[0, 0] = quat2_0 - quat2_1 - quat2_2 + quat2_3
    rotation_matrix[1, 1] = -quat2_0 + quat2_1 - quat2_2 + quat2_3
    rotation_matrix[2, 2] = -quat2_0 - quat2_1 + quat2_2 + quat2_3
    # rotation_matrix[0, 0] = quat0*quat0 - quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[1, 1] = -quat0*quat0 + quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[2, 2] = -quat0*quat0 - quat1*quat1 + quat2*quat2 + quat3*quat3
    """
    # Alternatively, in modern python versions, this code also works:
    rotation_matrix = Rotation.from_quat(optimal_quat).as_matrix()
    ""
    input(f'{rotation_matrix_.as_matrix()} {rotation_matrix}')
    input(np.allclose(rotation_matrix_.as_matrix(), rotation_matrix))
    # """
    # Finally compute the RMSD between the two coordinate sets:
    # First compute E0 from equation 24 of the paper
    # e0 = np.sum((fixed_coords_at_origin - moving_coords_at_origin) ** 2)
    # sum_sqr_dist = max(0, ((fixed_coords_at_origin-moving_coords_at_origin) ** 2).sum() - 2.*pPp)

    # if sum_weights != 0.:
    # try:
    rmsd = np.sqrt(max(0, ((fixed_coords_at_origin-moving_coords_at_origin)**2).sum() - 2.*pPp) / number_of_points)
    # except ZeroDivisionError:
    #     rmsd = 0.  # The weights are a total of zero which is allowed algorithmically, but not possible

    # Lastly, calculate the translational offset:
    # Recall that:
    # RMSD=sqrt((Σ_i  w_i * |X_i - (Σ_j c*R_ij*x_j + T_i))|^2) / (Σ_j w_j))
    #    =sqrt((Σ_i  w_i * |X_i - x_i'|^2) / (Σ_j w_j))
    #  where
    # x_i' = Σ_j c*R_ij*x_j + T_i
    #      = Xcm_i + c*R_ij*(x_j - xcm_j)
    #  and Xcm and xcm = center_of_mass for the frozen and mobile point clouds
    #                  = center_of_mass_fixed[]       and       center_of_mass_moving[],  respectively
    # Hence:
    #  T_i = Xcm_i - Σ_j c*R_ij*xcm_j  =  a_translate[i]

    # a_translate = center_of_mass_fixed - np.matmul(c * rotation_matrix, center_of_mass_moving).T.reshape(3,)

    # Calculate the translation
    translation = center_of_mass_fixed - rotation_matrix@center_of_mass_moving

    return rmsd, rotation_matrix, translation

superposition3d_quat

superposition3d_quat(fixed_coords: ndarray, moving_coords: ndarray) -> tuple[float, ndarray, ndarray]

Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations should be applied to the "moving_coords" to place them in the setting of the fixed_coords

This function implements a more general variant of the method from: R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216 This version has been augmented slightly. The version in the original paper only considers rotation and translation and does not allow the coordinates of either object to be rescaled (multiplication by a scalar). (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries form a vector (of length sin(θ/2)), pointing along the axis of rotation. Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

MIT License. Copyright (c) 2016, Andrew Jewett Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

Parameters:

  • fixed_coords (ndarray) –

    The coordinates for the 'frozen' object

  • moving_coords (ndarray) –

    The coordinates for the 'mobile' object

Raises: ValueError: If coordinates are not the same length Returns: rmsd, quaternion, translation_vector

Source code in symdesign/structure/coordinates.py
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
@jit(nopython=True)  # , cache=True)
def superposition3d_quat(fixed_coords: np.ndarray, moving_coords: np.ndarray) -> tuple[float, np.ndarray, np.ndarray]:
    """Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotation and translation
    operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations
    should be applied to the "moving_coords" to place them in the setting of the fixed_coords

    This function implements a more general variant of the method from:
    R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216
    This version has been augmented slightly. The version in the original paper only considers rotation and translation
    and does not allow the coordinates of either object to be rescaled (multiplication by a scalar).
    (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

    The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries
    form a vector (of length sin(θ/2)), pointing along the axis of rotation.
    Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

    MIT License. Copyright (c) 2016, Andrew Jewett
    Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated
    documentation files (the "Software"), to deal in the Software without restriction, including without limitation the
    rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to
    permit persons to whom the Software is furnished to do so, subject to the following conditions:

    The above copyright notice and this permission notice shall be included in all copies or substantial portions of the
    Software.

    THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE
    WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS
    OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
    OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

    Args:
        fixed_coords: The coordinates for the 'frozen' object
        moving_coords: The coordinates for the 'mobile' object
    Raises:
        ValueError: If coordinates are not the same length
    Returns:
        rmsd, quaternion, translation_vector
    """
    number_of_points = fixed_coords.shape[0]
    if number_of_points != moving_coords.shape[0]:
        raise ValueError('superposition3d: Inputs should have the same size')  # . '
        #                  f'Input 1={number_of_points}, 2={moving_coords.shape[0]}')

    # convert weights into array
    # if a_weights is None or len(a_weights) == 0:
    # a_weights = np.full((number_of_points, 1), 1.)
    # sum_weights = float(number_of_points)
    # else:  # reshape a_eights so multiplications are done column-wise
    #     a_weights = np.array(a_weights).reshape(number_of_points, 1)
    #     sum_weights = np.sum(a_weights, axis=0)

    # Find the center of mass of each object:
    center_of_mass_fixed = fixed_coords.sum(axis=0)
    center_of_mass_moving = moving_coords.sum(axis=0)

    # Subtract the centers-of-mass from the original coordinates for each object
    # if sum_weights != 0:
    # try:
    center_of_mass_fixed /= number_of_points
    center_of_mass_moving /= number_of_points
    # except ZeroDivisionError:
    #     pass  # The weights are a total of zero which is allowed algorithmically, but not possible

    # Translate the center of mass to the origin
    fixed_coords_at_origin = fixed_coords - center_of_mass_fixed
    moving_coords_at_origin = moving_coords - center_of_mass_moving

    # Calculate the "m" array from the Diamond paper (equation 16)
    m = moving_coords_at_origin.T @ fixed_coords_at_origin

    # Calculate "v" (equation 18)
    # v = np.empty(3)
    # v[0] = m[1, 2] - m[2, 1]
    # v[1] = m[2, 0] - m[0, 2]
    # v[2] = m[0, 1] - m[1, 0]
    v = [m[1, 2] - m[2, 1], m[2, 0] - m[0, 2], m[0, 1] - m[1, 0]]

    # Calculate "P" (equation 22)
    matrix_p = np.zeros((4, 4))
    # Calculate "q" (equation 17)
    # q = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[:3, :3] = m + m.T - 2*identity_matrix*np.trace(m)
    # matrix_p[:3, :3] = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[3, :3] = v
    matrix_p[:3, 3] = v
    # [[ q[0, 0] q[0, 1] q[0, 2] v[0] ]
    #  [ q[1, 0] q[1, 1] q[1, 2] v[1] ]
    #  [ q[2, 0] q[2, 1] q[2, 2] v[2] ]
    #  [ v[0]    v[1]    v[2]    0    ]]

    # Calculate "p" - optimal_quat
    # "p" contains the optimal rotation (in backwards-quaternion format)
    # (Note: A discussion of various quaternion conventions is included below)
    if number_of_points < 2:
        # Specify the default values for p, pPp
        optimal_quat = np.array([0., 0., 0., 1.])  # p = [0,0,0,1]    default value
        pPp = 0.  # = p^T * P * p    (zero by default)
    else:
        # try:
        # The eigenvalues/eigenvector are returned as 1D array in ascending order; largest is last
        a_eigenvals, aa_eigenvects = np.linalg.eigh(matrix_p)
        # except np.linalg.LinAlgError:
        #     singular = True  # I have never seen this happen
        # Pull out the largest magnitude
        pPp = a_eigenvals[-1]
        optimal_quat = aa_eigenvects[:, -1]
        # normalize the vector
        # (It should be normalized already, but just in case it is not, do it again)
        # optimal_quat /= np.linalg.norm(optimal_quat)

    # Calculate the rotation matrix corresponding to "optimal_quat" which is in scipy quaternion format
    # """
    rotation_matrix = np.empty((3, 3))
    quat0, quat1, quat2, quat3 = optimal_quat
    quat2_0, quat2_1, quat2_2, quat2_3 = optimal_quat**2
    rotation_matrix[0, 1] = quat0*quat1 - quat2*quat3  # 2*( )
    rotation_matrix[1, 0] = quat0*quat1 + quat2*quat3  # 2*( )
    rotation_matrix[1, 2] = quat1*quat2 - quat0*quat3  # 2*( )
    rotation_matrix[2, 1] = quat1*quat2 + quat0*quat3  # 2*( )
    rotation_matrix[0, 2] = quat0*quat2 + quat1*quat3  # 2*( )
    rotation_matrix[2, 0] = quat0*quat2 - quat1*quat3  # 2*( )
    rotation_matrix *= 2
    rotation_matrix[0, 0] = quat2_0 - quat2_1 - quat2_2 + quat2_3
    rotation_matrix[1, 1] = -quat2_0 + quat2_1 - quat2_2 + quat2_3
    rotation_matrix[2, 2] = -quat2_0 - quat2_1 + quat2_2 + quat2_3
    # rotation_matrix[0, 0] = quat0*quat0 - quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[1, 1] = -quat0*quat0 + quat1*quat1 - quat2*quat2 + quat3*quat3
    # rotation_matrix[2, 2] = -quat0*quat0 - quat1*quat1 + quat2*quat2 + quat3*quat3
    """
    # Alternatively, in modern python versions, this code also works:
    rotation_matrix = Rotation.from_quat(optimal_quat).as_matrix()
    ""
    input(f'{rotation_matrix_.as_matrix()} {rotation_matrix}')
    input(np.allclose(rotation_matrix_.as_matrix(), rotation_matrix))
    # """
    # Finally compute the RMSD between the two coordinate sets:
    # First compute E0 from equation 24 of the paper
    # e0 = np.sum((fixed_coords_at_origin - moving_coords_at_origin) ** 2)
    # sum_sqr_dist = max(0, ((fixed_coords_at_origin-moving_coords_at_origin) ** 2).sum() - 2.*pPp)

    # if sum_weights != 0.:
    # try:
    rmsd = np.sqrt(max(0, ((fixed_coords_at_origin-moving_coords_at_origin)**2).sum() - 2.*pPp) / number_of_points)
    # except ZeroDivisionError:
    #     rmsd = 0.  # The weights are a total of zero which is allowed algorithmically, but not possible

    # Lastly, calculate the translational offset:
    # Recall that:
    # RMSD=sqrt((Σ_i  w_i * |X_i - (Σ_j c*R_ij*x_j + T_i))|^2) / (Σ_j w_j))
    #    =sqrt((Σ_i  w_i * |X_i - x_i'|^2) / (Σ_j w_j))
    #  where
    # x_i' = Σ_j c*R_ij*x_j + T_i
    #      = Xcm_i + c*R_ij*(x_j - xcm_j)
    #  and Xcm and xcm = center_of_mass for the frozen and mobile point clouds
    #                  = center_of_mass_fixed[]       and       center_of_mass_moving[],  respectively
    # Hence:
    #  T_i = Xcm_i - Σ_j c*R_ij*xcm_j  =  a_translate[i]

    # a_translate = center_of_mass_fixed - np.matmul(c * rotation_matrix, center_of_mass_moving).T.reshape(3,)

    # Calculate the translation
    translation = center_of_mass_fixed - rotation_matrix@center_of_mass_moving
    # The p array is a quaternion that uses this convention:
    # https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.transform.Rotation.from_quat.html
    # However it seems that the following convention is much more popular:
    # https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation
    # https://mathworld.wolfram.com/Quaternion.html
    # So I return "q" (a version of "p" using the more popular convention).
    # rotation_matrix = np.array([p[3], p[0], p[1], p[2]])
    # KM: Disregard above, I am using the scipy version for python continuity which returns X, Y, Z, W
    return rmsd, optimal_quat, translation

superposition3d_weighted

superposition3d_weighted(fixed_coords: ndarray, moving_coords: ndarray, a_weights: ndarray = None, quaternion: bool = False) -> tuple[float, ndarray, ndarray]

Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotations, translations, and (optionally) rescale operations to minimize the root mean squared distance (RMSD) between them. The found transformation operations should be applied to the "moving_coords" to place them in the setting of the fixed_coords

This function implements a more general variant of the method from: R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216 This version has been augmented slightly. The version in the original paper only considers rotation and translation and does not allow the coordinates of either object to be rescaled (multiplication by a scalar). (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries form a vector (of length sin(θ/2)), pointing along the axis of rotation. Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

MIT License. Copyright (c) 2016, Andrew Jewett Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

Parameters:

  • fixed_coords (ndarray) –

    The coordinates for the 'frozen' object

  • moving_coords (ndarray) –

    The coordinates for the 'mobile' object

  • a_weights (ndarray, default: None ) –

    Weights for the calculation of RMSD

  • quaternion (bool, default: False ) –

    Whether to report the rotation angle and axis in Scipy.Rotation quaternion format

Raises: ValueError: If coordinates are not the same length Returns: rmsd, rotation/quaternion_matrix, translation_vector

Source code in symdesign/structure/coordinates.py
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
def superposition3d_weighted(fixed_coords: np.ndarray, moving_coords: np.ndarray, a_weights: np.ndarray = None,
                             quaternion: bool = False) -> tuple[float, np.ndarray, np.ndarray]:
    """Takes two xyz coordinate sets (same length), and attempts to superimpose them using rotations, translations,
    and (optionally) rescale operations to minimize the root mean squared distance (RMSD) between them. The found
    transformation operations should be applied to the "moving_coords" to place them in the setting of the fixed_coords

    This function implements a more general variant of the method from:
    R. Diamond, (1988) "A Note on the Rotational Superposition Problem", Acta Cryst. A44, pp. 211-216
    This version has been augmented slightly. The version in the original paper only considers rotation and translation
    and does not allow the coordinates of either object to be rescaled (multiplication by a scalar).
    (Additional documentation can be found at https://pypi.org/project/superpose3d/ )

    The quaternion_matrix has the last entry storing cos(θ/2) (where θ is the rotation angle). The first 3 entries
    form a vector (of length sin(θ/2)), pointing along the axis of rotation.
    Details: https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation

    MIT License. Copyright (c) 2016, Andrew Jewett
    Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated
    documentation files (the "Software"), to deal in the Software without restriction, including without limitation the
    rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to
    permit persons to whom the Software is furnished to do so, subject to the following conditions:

    The above copyright notice and this permission notice shall be included in all copies or substantial portions of the
    Software.

    THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE
    WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS
    OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
    OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.

    Args:
        fixed_coords: The coordinates for the 'frozen' object
        moving_coords: The coordinates for the 'mobile' object
        a_weights: Weights for the calculation of RMSD
        quaternion: Whether to report the rotation angle and axis in Scipy.Rotation quaternion format
    Raises:
        ValueError: If coordinates are not the same length
    Returns:
        rmsd, rotation/quaternion_matrix, translation_vector
    """
    number_of_points = len(fixed_coords)
    if number_of_points != len(moving_coords):
        raise ValueError(f'{superposition3d.__name__}: Inputs should have the same size. '
                         f'Input 1={number_of_points}, 2={moving_coords.shape[0]}')

    # convert weights into array
    if a_weights is None or len(a_weights) == 0:
        a_weights = np.full((number_of_points, 1), 1.)
        sum_weights = float(number_of_points)
    else:  # reshape a_eights so multiplications are done column-wise
        a_weights = np.array(a_weights).reshape(number_of_points, 1)
        sum_weights = np.sum(a_weights, axis=0)

    # Find the center of mass of each object:
    center_of_mass_fixed = np.sum(fixed_coords * a_weights, axis=0)
    center_of_mass_moving = np.sum(moving_coords * a_weights, axis=0)

    # Subtract the centers-of-mass from the original coordinates for each object
    # if sum_weights != 0:
    try:
        center_of_mass_fixed /= sum_weights
        center_of_mass_moving /= sum_weights
    except ZeroDivisionError:
        pass  # the weights are a total of zero which is allowed algorithmically, but not possible

    aa_xf = fixed_coords - center_of_mass_fixed
    aa_xm = moving_coords - center_of_mass_moving

    # Calculate the "m" array from the Diamond paper (equation 16)
    m = np.matmul(aa_xm.T, (aa_xf * a_weights))

    # Calculate "v" (equation 18)
    # v = np.empty(3)
    # v[0] = m[1, 2] - m[2, 1]
    # v[1] = m[2, 0] - m[0, 2]
    # v[2] = m[0, 1] - m[1, 0]
    v = [m[1, 2] - m[2, 1], m[2, 0] - m[0, 2], m[0, 1] - m[1, 0]]

    # Calculate "P" (equation 22)
    matrix_p = np.zeros((4, 4))
    # Calculate "q" (equation 17)
    # q = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[:3, :3] = m + m.T - 2*identity_matrix*np.trace(m)
    # matrix_p[:3, :3] = m + m.T - 2*utils.symmetry.identity_matrix*np.trace(m)
    matrix_p[3, :3] = v
    matrix_p[:3, 3] = v
    # [[ q[0, 0] q[0, 1] q[0, 2] v[0] ]
    #  [ q[1, 0] q[1, 1] q[1, 2] v[1] ]
    #  [ q[2, 0] q[2, 1] q[2, 2] v[2] ]
    #  [ v[0]    v[1]    v[2]    0    ]]

    # Calculate "p" - optimal_quat
    # "p" contains the optimal rotation (in backwards-quaternion format)
    # (Note: A discussion of various quaternion conventions is included below)
    if number_of_points < 2:
        # Specify the default values for p, pPp
        optimal_quat = np.array([0., 0., 0., 1.])  # p = [0,0,0,1]    default value
        pPp = 0.  # = p^T * P * p    (zero by default)
    else:
        # try:
        a_eigenvals, aa_eigenvects = np.linalg.eigh(matrix_p)
        # except np.linalg.LinAlgError:
        #     singular = True  # I have never seen this happen
        pPp = np.max(a_eigenvals)
        optimal_quat = aa_eigenvects[:, np.argmax(a_eigenvals)]  # pull out the largest magnitude eigenvector
        # normalize the vector
        # (It should be normalized already, but just in case it is not, do it again)
        optimal_quat /= np.linalg.norm(optimal_quat)

    # Calculate the rotation matrix corresponding to "optimal_quat" which is in scipy quaternion format
    """
    rotation_matrix = np.empty((3, 3))
    rotation_matrix[0][0] = (quat0*quat0)-(quat1*quat1)
                     -(quat2*quat2)+(quat3*quat3)
    rotation_matrix[1][1] = -(quat0*quat0)+(quat1*quat1)
                      -(quat2*quat2)+(quat3*quat3)
    rotation_matrix[2][2] = -(quat0*quat0)-(quat1*quat1)
                      +(quat2*quat2)+(quat3*quat3)
    rotation_matrix[0][1] = 2*(quat0*quat1 - quat2*quat3)
    rotation_matrix[1][0] = 2*(quat0*quat1 + quat2*quat3)
    rotation_matrix[1][2] = 2*(quat1*quat2 - quat0*quat3)
    rotation_matrix[2][1] = 2*(quat1*quat2 + quat0*quat3)
    rotation_matrix[0][2] = 2*(quat0*quat2 + quat1*quat3)
    rotation_matrix[2][0] = 2*(quat0*quat2 - quat1*quat3)
    """
    # Alternatively, in modern python versions, this code also works:
    rotation_matrix = Rotation.from_quat(optimal_quat).as_matrix()

    # Finally compute the RMSD between the two coordinate sets:
    # First compute E0 from equation 24 of the paper
    # e0 = np.sum((aa_xf - aa_xm) ** 2)
    # sum_sqr_dist = max(0, ((aa_xf-aa_xm) ** 2).sum() - 2.*pPp)

    # if sum_weights != 0.:
    try:
        rmsd = np.sqrt(max(0, ((aa_xf-aa_xm) ** 2).sum() - 2.*pPp) / sum_weights)
    except ZeroDivisionError:
        rmsd = 0.  # the weights are a total of zero which is allowed algorithmically, but not possible

    # Lastly, calculate the translational offset:
    # Recall that:
    # RMSD=sqrt((Σ_i  w_i * |X_i - (Σ_j c*R_ij*x_j + T_i))|^2) / (Σ_j w_j))
    #    =sqrt((Σ_i  w_i * |X_i - x_i'|^2) / (Σ_j w_j))
    #  where
    # x_i' = Σ_j c*R_ij*x_j + T_i
    #      = Xcm_i + c*R_ij*(x_j - xcm_j)
    #  and Xcm and xcm = center_of_mass for the frozen and mobile point clouds
    #                  = center_of_mass_fixed[]       and       center_of_mass_moving[],  respectively
    # Hence:
    #  T_i = Xcm_i - Σ_j c*R_ij*xcm_j  =  a_translate[i]

    # a_translate = center_of_mass_fixed - np.matmul(c * aa_rotate, center_of_mass_moving).T.reshape(3,)

    # return rmsd, aa_rotate, center_of_mass_fixed - np.matmul(aa_rotate, center_of_mass_moving).T.reshape(3,)
    # Calculate the translation
    translation = center_of_mass_fixed - np.matmul(rotation_matrix, center_of_mass_moving)
    if quaternion:  # does the caller want the quaternion?
        # The p array is a quaternion that uses this convention:
        # https://docs.scipy.org/doc/scipy/reference/generated/scipy.spatial.transform.Rotation.from_quat.html
        # However it seems that the following convention is much more popular:
        # https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation
        # https://mathworld.wolfram.com/Quaternion.html
        # So I return "q" (a version of "p" using the more popular convention).
        # rotation_matrix = np.array([p[3], p[0], p[1], p[2]])
        # KM: Disregard above, I am using the scipy version for python continuity
        return rmsd, optimal_quat, translation
    else:
        return rmsd, rotation_matrix, translation

transform_coordinates

transform_coordinates(coords: ndarray | Iterable, rotation: ndarray | Iterable = None, translation: ndarray | Iterable | int | float = None, rotation2: ndarray | Iterable = None, translation2: ndarray | Iterable | int | float = None) -> ndarray

Take a set of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the order of operations as: rotation, translation, rotation2, translation2

Parameters:

  • coords (ndarray | Iterable) –

    The coordinates to transform, can be shape (number of coordinates, 3)

  • rotation (ndarray | Iterable, default: None ) –

    The first rotation to apply, expected general rotation matrix shape (3, 3)

  • translation (ndarray | Iterable | int | float, default: None ) –

    The first translation to apply, expected shape (3)

  • rotation2 (ndarray | Iterable, default: None ) –

    The second rotation to apply, expected general rotation matrix shape (3, 3)

  • translation2 (ndarray | Iterable | int | float, default: None ) –

    The second translation to apply, expected shape (3)

Returns: The transformed coordinate set with the same shape as the original

Source code in symdesign/structure/coordinates.py
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
def transform_coordinates(coords: np.ndarray | Iterable, rotation: np.ndarray | Iterable = None,
                          translation: np.ndarray | Iterable | int | float = None,
                          rotation2: np.ndarray | Iterable = None,
                          translation2: np.ndarray | Iterable | int | float = None) -> np.ndarray:
    """Take a set of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the order of
    operations as: rotation, translation, rotation2, translation2

    Args:
        coords: The coordinates to transform, can be shape (number of coordinates, 3)
        rotation: The first rotation to apply, expected general rotation matrix shape (3, 3)
        translation: The first translation to apply, expected shape (3)
        rotation2: The second rotation to apply, expected general rotation matrix shape (3, 3)
        translation2: The second translation to apply, expected shape (3)
    Returns:
        The transformed coordinate set with the same shape as the original
    """
    new_coords = coords.copy()

    if rotation is not None:
        np.matmul(new_coords, np.transpose(rotation), out=new_coords)

    if translation is not None:
        new_coords += translation  # No array allocation, sets in place

    if rotation2 is not None:
        np.matmul(new_coords, np.transpose(rotation2), out=new_coords)

    if translation2 is not None:
        new_coords += translation2

    return coords

transform_coordinate_sets_with_broadcast

transform_coordinate_sets_with_broadcast(coord_sets: ndarray, rotation: ndarray = None, translation: ndarray | Iterable | int | float = None, rotation2: ndarray = None, translation2: ndarray | Iterable | int | float = None) -> ndarray

Take stacked sets of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the order of operations as: rotation, translation, rotation2, translation2. Non-efficient memory use

Parameters:

  • coord_sets (ndarray) –

    The coordinates to transform, can be shape (number of sets, number of coordinates, 3)

  • rotation (ndarray, default: None ) –

    The first rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)

  • translation (ndarray | Iterable | int | float, default: None ) –

    The first translation to apply, expected shape (number of sets, 3)

  • rotation2 (ndarray, default: None ) –

    The second rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)

  • translation2 (ndarray | Iterable | int | float, default: None ) –

    The second translation to apply, expected shape (number of sets, 3)

Returns: The transformed coordinate set with the same shape as the original

Source code in symdesign/structure/coordinates.py
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
def transform_coordinate_sets_with_broadcast(coord_sets: np.ndarray,
                                             rotation: np.ndarray = None,
                                             translation: np.ndarray | Iterable | int | float = None,
                                             rotation2: np.ndarray = None,
                                             translation2: np.ndarray | Iterable | int | float = None) \
        -> np.ndarray:
    """Take stacked sets of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the
    order of operations as: rotation, translation, rotation2, translation2. Non-efficient memory use

    Args:
        coord_sets: The coordinates to transform, can be shape (number of sets, number of coordinates, 3)
        rotation: The first rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)
        translation: The first translation to apply, expected shape (number of sets, 3)
        rotation2: The second rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)
        translation2: The second translation to apply, expected shape (number of sets, 3)
    Returns:
        The transformed coordinate set with the same shape as the original
    """
    # in general, the np.tensordot module accomplishes this coordinate set multiplication without stacking
    # np.tensordot(a, b, axes=1)  <-- axes=1 performs the correct multiplication with a 3d (3,3,N) by 2d (3,3) matrix
    # np.matmul solves as well due to broadcasting
    set_shape = getattr(coord_sets, 'shape', None)
    if set_shape is None or set_shape[0] < 1:
        return coord_sets
    # else:  # Create a new array for the result
    #     new_coord_sets = coord_sets.copy()

    if rotation is not None:
        coord_sets = np.matmul(coord_sets, rotation.swapaxes(-2, -1))

    if translation is not None:
        coord_sets += translation  # No array allocation, sets in place

    if rotation2 is not None:
        coord_sets = np.matmul(coord_sets, rotation2.swapaxes(-2, -1))

    if translation2 is not None:
        coord_sets += translation2

    return coord_sets

transform_coordinate_sets

transform_coordinate_sets(coord_sets: ndarray, rotation: ndarray = None, translation: ndarray | Iterable | int | float = None, rotation2: ndarray = None, translation2: ndarray | Iterable | int | float = None) -> ndarray

Take stacked sets of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the order of operations as: rotation, translation, rotation2, translation2. If transformation uses broadcasting, for efficient memory use, the returned array will be the size of the coord_sets multiplied by rotation. Additional broadcasting is not allowed. If that behavior is desired, use "transform_coordinate_sets_with_broadcast()" instead

Parameters:

  • coord_sets (ndarray) –

    The coordinates to transform, can be shape (number of sets, number of coordinates, 3)

  • rotation (ndarray, default: None ) –

    The first rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)

  • translation (ndarray | Iterable | int | float, default: None ) –

    The first translation to apply, expected shape (number of sets, 3)

  • rotation2 (ndarray, default: None ) –

    The second rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)

  • translation2 (ndarray | Iterable | int | float, default: None ) –

    The second translation to apply, expected shape (number of sets, 3)

Returns: The transformed coordinate set with the same shape as the original

Source code in symdesign/structure/coordinates.py
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
def transform_coordinate_sets(coord_sets: np.ndarray,
                              rotation: np.ndarray = None, translation: np.ndarray | Iterable | int | float = None,
                              rotation2: np.ndarray = None, translation2: np.ndarray | Iterable | int | float = None) \
        -> np.ndarray:
    """Take stacked sets of x,y,z coordinates and transform. Transformation proceeds by matrix multiplication with the
    order of operations as: rotation, translation, rotation2, translation2. If transformation uses broadcasting, for
    efficient memory use, the returned array will be the size of the coord_sets multiplied by rotation. Additional
    broadcasting is not allowed. If that behavior is desired, use "transform_coordinate_sets_with_broadcast()" instead

    Args:
        coord_sets: The coordinates to transform, can be shape (number of sets, number of coordinates, 3)
        rotation: The first rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)
        translation: The first translation to apply, expected shape (number of sets, 3)
        rotation2: The second rotation to apply, expected general rotation matrix shape (number of sets, 3, 3)
        translation2: The second translation to apply, expected shape (number of sets, 3)
    Returns:
        The transformed coordinate set with the same shape as the original
    """
    # in general, the np.tensordot module accomplishes this coordinate set multiplication without stacking
    # np.tensordot(a, b, axes=1)  <-- axes=1 performs the correct multiplication with a 3d (3,3,N) by 2d (3,3) matrix
    # np.matmul solves as well due to broadcasting
    set_shape = getattr(coord_sets, 'shape', None)
    if set_shape is None or set_shape[0] < 1:
        return coord_sets

    if rotation is not None:
        new_coord_sets = np.matmul(coord_sets, rotation.swapaxes(-2, -1))
    else:  # Create a new array for the result
        new_coord_sets = coord_sets.copy()

    if translation is not None:
        new_coord_sets += translation  # No array allocation, sets in place

    if rotation2 is not None:
        np.matmul(new_coord_sets, rotation2.swapaxes(-2, -1), out=new_coord_sets)
        # new_coord_sets[:] = np.matmul(new_coord_sets, rotation2.swapaxes(-2, -1))
        # new_coord_sets = np.matmul(new_coord_sets, rotation2.swapaxes(-2, -1))

    if translation2 is not None:
        new_coord_sets += translation2

    return new_coord_sets