1# coding: utf-8 2# Copyright (c) Pymatgen Development Team. 3# Distributed under the terms of the MIT License. 4 5""" 6Defines the classes relating to 3D lattices. 7""" 8 9import collections 10import itertools 11import math 12import warnings 13from fractions import Fraction 14from functools import reduce 15from typing import Dict, Iterator, List, Optional, Sequence, Tuple, Union 16 17import numpy as np 18from monty.dev import deprecated 19from monty.json import MSONable 20from numpy import dot, pi, transpose 21from numpy.linalg import inv 22 23from pymatgen.util.coord import pbc_shortest_vectors 24from pymatgen.util.num import abs_cap 25from pymatgen.util.typing import ArrayLike 26 27__author__ = "Shyue Ping Ong, Michael Kocher" 28__copyright__ = "Copyright 2011, The Materials Project" 29__maintainer__ = "Shyue Ping Ong" 30__email__ = "shyuep@gmail.com" 31 32 33class Lattice(MSONable): 34 """ 35 A lattice object. Essentially a matrix with conversion matrices. In 36 general, it is assumed that length units are in Angstroms and angles are in 37 degrees unless otherwise stated. 38 """ 39 40 # Properties lazily generated for efficiency. 41 42 def __init__(self, matrix: ArrayLike): 43 """ 44 Create a lattice from any sequence of 9 numbers. Note that the sequence 45 is assumed to be read one row at a time. Each row represents one 46 lattice vector. 47 48 Args: 49 matrix: Sequence of numbers in any form. Examples of acceptable 50 input. 51 i) An actual numpy array. 52 ii) [[1, 0, 0], [0, 1, 0], [0, 0, 1]] 53 iii) [1, 0, 0 , 0, 1, 0, 0, 0, 1] 54 iv) (1, 0, 0, 0, 1, 0, 0, 0, 1) 55 Each row should correspond to a lattice vector. 56 E.g., [[10, 0, 0], [20, 10, 0], [0, 0, 30]] specifies a lattice 57 with lattice vectors [10, 0, 0], [20, 10, 0] and [0, 0, 30]. 58 """ 59 m = np.array(matrix, dtype=np.float64).reshape((3, 3)) 60 m.setflags(write=False) 61 self._matrix = m # type: np.ndarray 62 self._inv_matrix = None # type: Optional[np.ndarray] 63 self._diags = None 64 self._lll_matrix_mappings = {} # type: Dict[float, Tuple[np.ndarray, np.ndarray]] 65 self._lll_inverse = None 66 67 @property 68 def lengths(self) -> Tuple[float, float, float]: 69 """ 70 :return: The lengths (a, b, c) of the lattice. 71 """ 72 return tuple(np.sqrt(np.sum(self._matrix ** 2, axis=1)).tolist()) # type: ignore 73 74 @property 75 def angles(self) -> Tuple[float, float, float]: 76 """ 77 Returns the angles (alpha, beta, gamma) of the lattice. 78 """ 79 m = self._matrix 80 lengths = self.lengths 81 angles = np.zeros(3) 82 for i in range(3): 83 j = (i + 1) % 3 84 k = (i + 2) % 3 85 angles[i] = abs_cap(dot(m[j], m[k]) / (lengths[j] * lengths[k])) 86 angles = np.arccos(angles) * 180.0 / pi 87 return tuple(angles.tolist()) # type: ignore 88 89 @property 90 def is_orthogonal(self) -> bool: 91 """ 92 :return: Whether all angles are 90 degrees. 93 """ 94 return all(abs(a - 90) < 1e-5 for a in self.angles) 95 96 def __format__(self, fmt_spec=""): 97 """ 98 Support format printing. Supported formats are: 99 100 1. "l" for a list format that can be easily copied and pasted, e.g., 101 ".3fl" prints something like 102 "[[10.000, 0.000, 0.000], [0.000, 10.000, 0.000], [0.000, 0.000, 10.000]]" 103 2. "p" for lattice parameters ".1fp" prints something like 104 "{10.0, 10.0, 10.0, 90.0, 90.0, 90.0}" 105 3. Default will simply print a 3x3 matrix form. E.g., 106 10.000 0.000 0.000 107 0.000 10.000 0.000 108 0.000 0.000 10.000 109 """ 110 m = self._matrix.tolist() 111 if fmt_spec.endswith("l"): 112 fmt = "[[{}, {}, {}], [{}, {}, {}], [{}, {}, {}]]" 113 fmt_spec = fmt_spec[:-1] 114 elif fmt_spec.endswith("p"): 115 fmt = "{{{}, {}, {}, {}, {}, {}}}" 116 fmt_spec = fmt_spec[:-1] 117 m = (self.lengths, self.angles) 118 else: 119 fmt = "{} {} {}\n{} {} {}\n{} {} {}" 120 return fmt.format(*[format(c, fmt_spec) for row in m for c in row]) 121 122 def copy(self): 123 """Deep copy of self.""" 124 return self.__class__(self.matrix.copy()) 125 126 @property 127 def matrix(self) -> np.ndarray: 128 """Copy of matrix representing the Lattice""" 129 return self._matrix 130 131 @property 132 def inv_matrix(self) -> np.ndarray: 133 """ 134 Inverse of lattice matrix. 135 """ 136 if self._inv_matrix is None: 137 self._inv_matrix = inv(self._matrix) 138 self._inv_matrix.setflags(write=False) 139 return self._inv_matrix 140 141 @property 142 def metric_tensor(self) -> np.ndarray: 143 """ 144 The metric tensor of the lattice. 145 """ 146 return dot(self._matrix, self._matrix.T) 147 148 def get_cartesian_coords(self, fractional_coords: ArrayLike) -> np.ndarray: 149 """ 150 Returns the cartesian coordinates given fractional coordinates. 151 152 Args: 153 fractional_coords (3x1 array): Fractional coords. 154 155 Returns: 156 Cartesian coordinates 157 """ 158 return dot(fractional_coords, self._matrix) 159 160 def get_fractional_coords(self, cart_coords: ArrayLike) -> np.ndarray: 161 """ 162 Returns the fractional coordinates given cartesian coordinates. 163 164 Args: 165 cart_coords (3x1 array): Cartesian coords. 166 167 Returns: 168 Fractional coordinates. 169 """ 170 return dot(cart_coords, self.inv_matrix) 171 172 def get_vector_along_lattice_directions(self, cart_coords: ArrayLike) -> np.ndarray: 173 """ 174 Returns the coordinates along lattice directions given cartesian coordinates. 175 176 Note, this is different than a projection of the cartesian vector along the 177 lattice parameters. It is simply the fractional coordinates multiplied by the 178 lattice vector magnitudes. 179 180 For example, this method is helpful when analyzing the dipole moment (in 181 units of electron Angstroms) of a ferroelectric crystal. See the `Polarization` 182 class in `pymatgen.analysis.ferroelectricity.polarization`. 183 184 Args: 185 cart_coords (3x1 array): Cartesian coords. 186 187 Returns: 188 Lattice coordinates. 189 """ 190 return self.lengths * self.get_fractional_coords(cart_coords) 191 192 def d_hkl(self, miller_index: ArrayLike) -> float: 193 """ 194 Returns the distance between the hkl plane and the origin 195 196 Args: 197 miller_index ([h,k,l]): Miller index of plane 198 199 Returns: 200 d_hkl (float) 201 """ 202 203 gstar = self.reciprocal_lattice_crystallographic.metric_tensor 204 hkl = np.array(miller_index) 205 return 1 / ((dot(dot(hkl, gstar), hkl.T)) ** (1 / 2)) 206 207 @staticmethod 208 def cubic(a: float) -> "Lattice": 209 """ 210 Convenience constructor for a cubic lattice. 211 212 Args: 213 a (float): The *a* lattice parameter of the cubic cell. 214 215 Returns: 216 Cubic lattice of dimensions a x a x a. 217 """ 218 return Lattice([[a, 0.0, 0.0], [0.0, a, 0.0], [0.0, 0.0, a]]) 219 220 @staticmethod 221 def tetragonal(a: float, c: float) -> "Lattice": 222 """ 223 Convenience constructor for a tetragonal lattice. 224 225 Args: 226 a (float): *a* lattice parameter of the tetragonal cell. 227 c (float): *c* lattice parameter of the tetragonal cell. 228 229 Returns: 230 Tetragonal lattice of dimensions a x a x c. 231 """ 232 return Lattice.from_parameters(a, a, c, 90, 90, 90) 233 234 @staticmethod 235 def orthorhombic(a: float, b: float, c: float) -> "Lattice": 236 """ 237 Convenience constructor for an orthorhombic lattice. 238 239 Args: 240 a (float): *a* lattice parameter of the orthorhombic cell. 241 b (float): *b* lattice parameter of the orthorhombic cell. 242 c (float): *c* lattice parameter of the orthorhombic cell. 243 244 Returns: 245 Orthorhombic lattice of dimensions a x b x c. 246 """ 247 return Lattice.from_parameters(a, b, c, 90, 90, 90) 248 249 @staticmethod 250 def monoclinic(a: float, b: float, c: float, beta: float) -> "Lattice": 251 """ 252 Convenience constructor for a monoclinic lattice. 253 254 Args: 255 a (float): *a* lattice parameter of the monoclinc cell. 256 b (float): *b* lattice parameter of the monoclinc cell. 257 c (float): *c* lattice parameter of the monoclinc cell. 258 beta (float): *beta* angle between lattice vectors b and c in 259 degrees. 260 261 Returns: 262 Monoclinic lattice of dimensions a x b x c with non right-angle 263 beta between lattice vectors a and c. 264 """ 265 return Lattice.from_parameters(a, b, c, 90, beta, 90) 266 267 @staticmethod 268 def hexagonal(a: float, c: float) -> "Lattice": 269 """ 270 Convenience constructor for a hexagonal lattice. 271 272 Args: 273 a (float): *a* lattice parameter of the hexagonal cell. 274 c (float): *c* lattice parameter of the hexagonal cell. 275 276 Returns: 277 Hexagonal lattice of dimensions a x a x c. 278 """ 279 return Lattice.from_parameters(a, a, c, 90, 90, 120) 280 281 @staticmethod 282 def rhombohedral(a: float, alpha: float) -> "Lattice": 283 """ 284 Convenience constructor for a rhombohedral lattice. 285 286 Args: 287 a (float): *a* lattice parameter of the rhombohedral cell. 288 alpha (float): Angle for the rhombohedral lattice in degrees. 289 290 Returns: 291 Rhombohedral lattice of dimensions a x a x a. 292 """ 293 return Lattice.from_parameters(a, a, a, alpha, alpha, alpha) 294 295 @staticmethod 296 @deprecated(message="Use Lattice.from_parameters instead. This will be removed in v2020.*") 297 def from_lengths_and_angles(abc: Sequence[float], ang: Sequence[float]): 298 """ 299 Create a Lattice using unit cell lengths and angles (in degrees). 300 301 Args: 302 abc (3x1 array): Lattice parameters, e.g. (4, 4, 5). 303 ang (3x1 array): Lattice angles in degrees, e.g., (90,90,120). 304 305 Returns: 306 A Lattice with the specified lattice parameters. 307 """ 308 return Lattice.from_parameters(abc[0], abc[1], abc[2], ang[0], ang[1], ang[2]) 309 310 @classmethod 311 def from_parameters( 312 cls, 313 a: float, 314 b: float, 315 c: float, 316 alpha: float, 317 beta: float, 318 gamma: float, 319 vesta: bool = False, 320 ): 321 """ 322 Create a Lattice using unit cell lengths and angles (in degrees). 323 324 Args: 325 a (float): *a* lattice parameter. 326 b (float): *b* lattice parameter. 327 c (float): *c* lattice parameter. 328 alpha (float): *alpha* angle in degrees. 329 beta (float): *beta* angle in degrees. 330 gamma (float): *gamma* angle in degrees. 331 vesta: True if you import Cartesian coordinates from VESTA. 332 333 Returns: 334 Lattice with the specified lattice parameters. 335 """ 336 337 angles_r = np.radians([alpha, beta, gamma]) 338 cos_alpha, cos_beta, cos_gamma = np.cos(angles_r) 339 sin_alpha, sin_beta, sin_gamma = np.sin(angles_r) 340 341 if vesta: 342 c1 = c * cos_beta 343 c2 = (c * (cos_alpha - (cos_beta * cos_gamma))) / sin_gamma 344 345 vector_a = [float(a), 0.0, 0.0] 346 vector_b = [b * cos_gamma, b * sin_gamma, 0] 347 vector_c = [c1, c2, math.sqrt(c ** 2 - c1 ** 2 - c2 ** 2)] 348 349 else: 350 val = (cos_alpha * cos_beta - cos_gamma) / (sin_alpha * sin_beta) 351 # Sometimes rounding errors result in values slightly > 1. 352 val = abs_cap(val) 353 gamma_star = np.arccos(val) 354 355 vector_a = [a * sin_beta, 0.0, a * cos_beta] 356 vector_b = [ 357 -b * sin_alpha * np.cos(gamma_star), 358 b * sin_alpha * np.sin(gamma_star), 359 b * cos_alpha, 360 ] 361 vector_c = [0.0, 0.0, float(c)] 362 363 return Lattice([vector_a, vector_b, vector_c]) 364 365 @classmethod 366 def from_dict(cls, d: Dict, fmt: str = None, **kwargs): 367 """ 368 Create a Lattice from a dictionary containing the a, b, c, alpha, beta, 369 and gamma parameters if fmt is None. 370 371 If fmt == "abivars", the function build a `Lattice` object from a 372 dictionary with the Abinit variables `acell` and `rprim` in Bohr. 373 If acell is not given, the Abinit default is used i.e. [1,1,1] Bohr 374 375 Example: 376 377 Lattice.from_dict(fmt="abivars", acell=3*[10], rprim=np.eye(3)) 378 """ 379 if fmt == "abivars": 380 # pylint: disable=C0415 381 from pymatgen.io.abinit.abiobjects import lattice_from_abivars 382 383 kwargs.update(d) 384 return lattice_from_abivars(cls=cls, **kwargs) 385 386 if "matrix" in d: 387 return cls(d["matrix"]) 388 return cls.from_parameters(d["a"], d["b"], d["c"], d["alpha"], d["beta"], d["gamma"]) 389 390 @property 391 def a(self) -> float: 392 """ 393 *a* lattice parameter. 394 """ 395 return self.lengths[0] 396 397 @property 398 def b(self) -> float: 399 """ 400 *b* lattice parameter. 401 """ 402 return self.lengths[1] 403 404 @property 405 def c(self) -> float: 406 """ 407 *c* lattice parameter. 408 """ 409 return self.lengths[2] 410 411 @property 412 def abc(self) -> Tuple[float, float, float]: 413 """ 414 Lengths of the lattice vectors, i.e. (a, b, c) 415 """ 416 return self.lengths 417 418 @property 419 def alpha(self) -> float: 420 """ 421 Angle alpha of lattice in degrees. 422 """ 423 return self.angles[0] 424 425 @property 426 def beta(self) -> float: 427 """ 428 Angle beta of lattice in degrees. 429 """ 430 return self.angles[1] 431 432 @property 433 def gamma(self) -> float: 434 """ 435 Angle gamma of lattice in degrees. 436 """ 437 return self.angles[2] 438 439 @property 440 def volume(self) -> float: 441 """ 442 Volume of the unit cell. 443 """ 444 m = self._matrix 445 return float(abs(dot(np.cross(m[0], m[1]), m[2]))) 446 447 @property 448 def parameters(self) -> Tuple[float, float, float, float, float, float]: 449 """ 450 Returns: (a, b, c, alpha, beta, gamma). 451 """ 452 return (*self.lengths, *self.angles) 453 454 @property # type: ignore 455 @deprecated(message="Use Lattice.parameters instead. This will be removed in v2020.*") 456 def lengths_and_angles( 457 self, 458 ) -> Tuple[Tuple[float, float, float], Tuple[float, float, float]]: 459 """ 460 Returns (lattice lengths, lattice angles). 461 """ 462 return self.lengths, self.angles 463 464 @property 465 def reciprocal_lattice(self) -> "Lattice": 466 """ 467 Return the reciprocal lattice. Note that this is the standard 468 reciprocal lattice used for solid state physics with a factor of 2 * 469 pi. If you are looking for the crystallographic reciprocal lattice, 470 use the reciprocal_lattice_crystallographic property. 471 The property is lazily generated for efficiency. 472 """ 473 v = np.linalg.inv(self._matrix).T 474 return Lattice(v * 2 * np.pi) 475 476 @property 477 def reciprocal_lattice_crystallographic(self) -> "Lattice": 478 """ 479 Returns the *crystallographic* reciprocal lattice, i.e., no factor of 480 2 * pi. 481 """ 482 return Lattice(self.reciprocal_lattice.matrix / (2 * np.pi)) 483 484 @property 485 def lll_matrix(self) -> np.ndarray: 486 """ 487 :return: The matrix for LLL reduction 488 """ 489 if 0.75 not in self._lll_matrix_mappings: 490 self._lll_matrix_mappings[0.75] = self._calculate_lll() 491 return self._lll_matrix_mappings[0.75][0] 492 493 @property 494 def lll_mapping(self) -> np.ndarray: 495 """ 496 :return: The mapping between the LLL reduced lattice and the original 497 lattice. 498 """ 499 if 0.75 not in self._lll_matrix_mappings: 500 self._lll_matrix_mappings[0.75] = self._calculate_lll() 501 return self._lll_matrix_mappings[0.75][1] 502 503 @property 504 def lll_inverse(self) -> np.ndarray: 505 """ 506 :return: Inverse of self.lll_mapping. 507 """ 508 return np.linalg.inv(self.lll_mapping) 509 510 @property 511 def selling_vector(self) -> np.ndarray: 512 """ 513 Returns the (1,6) array of Selling Scalars. 514 """ 515 a, b, c = self.matrix 516 d = -(a + b + c) 517 tol = 1e-10 518 519 selling_vector = np.array( 520 [ 521 np.dot(b, c), 522 np.dot(a, c), 523 np.dot(a, b), 524 np.dot(a, d), 525 np.dot(b, d), 526 np.dot(c, d), 527 ] 528 ) 529 selling_vector = np.array([s if abs(s) > tol else 0 for s in selling_vector]) 530 531 reduction_matrices = np.array( 532 [ 533 np.array( 534 [ 535 [-1, 0, 0, 0, 0, 0], 536 [1, 1, 0, 0, 0, 0], 537 [1, 0, 0, 0, 1, 0], 538 [-1, 0, 0, 1, 0, 0], 539 [1, 0, 1, 0, 0, 0], 540 [1, 0, 0, 0, 0, 1], 541 ] 542 ), 543 np.array( 544 [ 545 [1, 1, 0, 0, 0, 0], 546 [0, -1, 0, 0, 0, 0], 547 [0, 1, 0, 1, 0, 0], 548 [0, 1, 1, 0, 0, 0], 549 [0, -1, 0, 0, 1, 0], 550 [0, 1, 0, 0, 0, 1], 551 ] 552 ), 553 np.array( 554 [ 555 [1, 0, 1, 0, 0, 0], 556 [0, 0, 1, 1, 0, 0], 557 [0, 0, -1, 0, 0, 0], 558 [0, 1, 1, 0, 0, 0], 559 [0, 0, 1, 0, 1, 0], 560 [0, 0, -1, 0, 0, 1], 561 ] 562 ), 563 np.array( 564 [ 565 [1, 0, 0, -1, 0, 0], 566 [0, 0, 1, 1, 0, 0], 567 [0, 1, 0, 1, 0, 0], 568 [0, 0, 0, -1, 0, 0], 569 [0, 0, 0, 1, 1, 0], 570 [0, 0, 0, 1, 0, 1], 571 ] 572 ), 573 np.array( 574 [ 575 [0, 0, 1, 0, 1, 0], 576 [0, 1, 0, 0, -1, 0], 577 [1, 0, 0, 0, 1, 0], 578 [0, 0, 0, 1, 1, 0], 579 [0, 0, 0, 0, -1, 0], 580 [0, 0, 0, 0, 1, 1], 581 ] 582 ), 583 np.array( 584 [ 585 [0, 1, 0, 0, 0, 1], 586 [1, 0, 0, 0, 0, 1], 587 [0, 0, 1, 0, 0, -1], 588 [0, 0, 0, 1, 0, 1], 589 [0, 0, 0, 0, 1, 1], 590 [0, 0, 0, 0, 0, -1], 591 ] 592 ), 593 ] 594 ) 595 while np.greater(np.max(selling_vector), 0): 596 max_index = selling_vector.argmax() 597 selling_vector = np.dot(reduction_matrices[max_index], selling_vector) 598 599 return selling_vector 600 601 def selling_dist(self, other): 602 """ 603 Returns the minimum Selling distance between two lattices. 604 """ 605 vcp_matrices = [ 606 np.array( 607 [ 608 [-1, 0, 0, 0, 0, 0], 609 [0, 1, 0, 0, 0, 0], 610 [0, 0, 0, 0, 1, 0], 611 [0, 0, 0, 1, 0, 0], 612 [0, 0, 1, 0, 0, 0], 613 [0, 0, 0, 0, 0, 1], 614 ] 615 ), 616 np.array( 617 [ 618 [1, 0, 0, 0, 0, 0], 619 [0, -1, 0, 0, 0, 0], 620 [0, 0, 0, 1, 0, 0], 621 [0, 0, 1, 0, 0, 0], 622 [0, 0, 0, 0, 1, 0], 623 [0, 0, 0, 0, 0, 1], 624 ] 625 ), 626 np.array( 627 [ 628 [1, 0, 0, 0, 0, 0], 629 [0, 0, 0, 1, 0, 0], 630 [0, 0, -1, 0, 0, 0], 631 [0, 1, 0, 0, 0, 0], 632 [0, 0, 0, 0, 1, 0], 633 [0, 0, 0, 0, 0, 1], 634 ] 635 ), 636 np.array( 637 [ 638 [1, 0, 0, 0, 0, 0], 639 [0, 0, 1, 0, 0, 0], 640 [0, 1, 0, 0, 0, 0], 641 [0, 0, 0, -1, 0, 0], 642 [0, 0, 0, 0, 1, 0], 643 [0, 0, 0, 0, 0, 1], 644 ] 645 ), 646 np.array( 647 [ 648 [0, 0, 1, 0, 0, 0], 649 [0, 1, 0, 0, 0, 0], 650 [1, 0, 0, 0, 0, 0], 651 [0, 0, 0, 1, 0, 0], 652 [0, 0, 0, 0, -1, 0], 653 [0, 0, 0, 0, 0, 1], 654 ] 655 ), 656 np.array( 657 [ 658 [0, 1, 0, 0, 0, 0], 659 [1, 0, 0, 0, 0, 0], 660 [0, 0, 1, 0, 0, 0], 661 [0, 0, 0, 1, 0, 0], 662 [0, 0, 0, 0, 1, 0], 663 [0, 0, 0, 0, 0, -1], 664 ] 665 ), 666 ] 667 668 reflection_matrices = [ 669 np.array( 670 [ 671 [1, 0, 0, 0, 0, 0], 672 [0, 1, 0, 0, 0, 0], 673 [0, 0, 1, 0, 0, 0], 674 [0, 0, 0, 1, 0, 0], 675 [0, 0, 0, 0, 1, 0], 676 [0, 0, 0, 0, 0, 1], 677 ] 678 ), 679 np.array( 680 [ 681 [1, 0, 0, 0, 0, 0], 682 [0, 0, 1, 0, 0, 0], 683 [0, 1, 0, 0, 0, 0], 684 [0, 0, 0, 1, 0, 0], 685 [0, 0, 0, 0, 0, 1], 686 [0, 0, 0, 0, 1, 0], 687 ] 688 ), 689 np.array( 690 [ 691 [1, 0, 0, 0, 0, 0], 692 [0, 0, 0, 0, 1, 0], 693 [0, 0, 0, 0, 0, 1], 694 [0, 0, 0, 1, 0, 0], 695 [0, 1, 0, 0, 0, 0], 696 [0, 0, 1, 0, 0, 0], 697 ] 698 ), 699 np.array( 700 [ 701 [1, 0, 0, 0, 0, 0], 702 [0, 0, 0, 0, 0, 1], 703 [0, 0, 0, 0, 1, 0], 704 [0, 0, 0, 1, 0, 0], 705 [0, 0, 1, 0, 0, 0], 706 [0, 1, 0, 0, 0, 0], 707 ] 708 ), 709 np.array( 710 [ 711 [0, 1, 0, 0, 0, 0], 712 [1, 0, 0, 0, 0, 0], 713 [0, 0, 1, 0, 0, 0], 714 [0, 0, 0, 0, 1, 0], 715 [0, 0, 0, 1, 0, 0], 716 [0, 0, 0, 0, 0, 1], 717 ] 718 ), 719 np.array( 720 [ 721 [0, 1, 0, 0, 0, 0], 722 [0, 0, 1, 0, 0, 0], 723 [1, 0, 0, 0, 0, 0], 724 [0, 0, 0, 0, 1, 0], 725 [0, 0, 0, 0, 0, 1], 726 [0, 0, 0, 1, 0, 0], 727 ] 728 ), 729 np.array( 730 [ 731 [0, 1, 0, 0, 0, 0], 732 [0, 0, 0, 1, 0, 0], 733 [0, 0, 0, 0, 0, 1], 734 [0, 0, 0, 0, 1, 0], 735 [1, 0, 0, 0, 0, 0], 736 [0, 0, 1, 0, 0, 0], 737 ] 738 ), 739 np.array( 740 [ 741 [0, 1, 0, 0, 0, 0], 742 [0, 0, 0, 0, 0, 1], 743 [0, 0, 0, 1, 0, 0], 744 [0, 0, 0, 0, 1, 0], 745 [0, 0, 1, 0, 0, 0], 746 [1, 0, 0, 0, 0, 0], 747 ] 748 ), 749 np.array( 750 [ 751 [0, 0, 1, 0, 0, 0], 752 [1, 0, 0, 0, 0, 0], 753 [0, 1, 0, 0, 0, 0], 754 [0, 0, 0, 0, 0, 1], 755 [0, 0, 0, 1, 0, 0], 756 [0, 0, 0, 0, 1, 0], 757 ] 758 ), 759 np.array( 760 [ 761 [0, 0, 1, 0, 0, 0], 762 [0, 1, 0, 0, 0, 0], 763 [1, 0, 0, 0, 0, 0], 764 [0, 0, 0, 0, 0, 1], 765 [0, 0, 0, 0, 1, 0], 766 [0, 0, 0, 1, 0, 0], 767 ] 768 ), 769 np.array( 770 [ 771 [0, 0, 1, 0, 0, 0], 772 [0, 0, 0, 1, 0, 0], 773 [0, 0, 0, 0, 1, 0], 774 [0, 0, 0, 0, 0, 1], 775 [1, 0, 0, 0, 0, 0], 776 [0, 1, 0, 0, 0, 0], 777 ] 778 ), 779 np.array( 780 [ 781 [0, 0, 1, 0, 0, 0], 782 [0, 0, 0, 0, 1, 0], 783 [0, 0, 0, 1, 0, 0], 784 [0, 0, 0, 0, 0, 1], 785 [0, 1, 0, 0, 0, 0], 786 [1, 0, 0, 0, 0, 0], 787 ] 788 ), 789 np.array( 790 [ 791 [0, 0, 0, 1, 0, 0], 792 [0, 1, 0, 0, 0, 0], 793 [0, 0, 0, 0, 0, 1], 794 [1, 0, 0, 0, 0, 0], 795 [0, 0, 0, 0, 1, 0], 796 [0, 0, 1, 0, 0, 0], 797 ] 798 ), 799 np.array( 800 [ 801 [0, 0, 0, 1, 0, 0], 802 [0, 0, 1, 0, 0, 0], 803 [0, 0, 0, 0, 1, 0], 804 [1, 0, 0, 0, 0, 0], 805 [0, 0, 0, 0, 0, 1], 806 [0, 1, 0, 0, 0, 0], 807 ] 808 ), 809 np.array( 810 [ 811 [0, 0, 0, 1, 0, 0], 812 [0, 0, 0, 0, 1, 0], 813 [0, 0, 1, 0, 0, 0], 814 [1, 0, 0, 0, 0, 0], 815 [0, 1, 0, 0, 0, 0], 816 [0, 0, 0, 0, 0, 1], 817 ] 818 ), 819 np.array( 820 [ 821 [0, 0, 0, 1, 0, 0], 822 [0, 0, 0, 0, 0, 1], 823 [0, 1, 0, 0, 0, 0], 824 [1, 0, 0, 0, 0, 0], 825 [0, 0, 1, 0, 0, 0], 826 [0, 0, 0, 0, 1, 0], 827 ] 828 ), 829 np.array( 830 [ 831 [0, 0, 0, 0, 1, 0], 832 [1, 0, 0, 0, 0, 0], 833 [0, 0, 0, 0, 0, 1], 834 [0, 1, 0, 0, 0, 0], 835 [0, 0, 0, 1, 0, 0], 836 [0, 0, 1, 0, 0, 0], 837 ] 838 ), 839 np.array( 840 [ 841 [0, 0, 0, 0, 1, 0], 842 [0, 0, 1, 0, 0, 0], 843 [0, 0, 0, 1, 0, 0], 844 [0, 1, 0, 0, 0, 0], 845 [0, 0, 0, 0, 0, 1], 846 [1, 0, 0, 0, 0, 0], 847 ] 848 ), 849 np.array( 850 [ 851 [0, 0, 0, 0, 1, 0], 852 [0, 0, 0, 1, 0, 0], 853 [0, 0, 1, 0, 0, 0], 854 [0, 1, 0, 0, 0, 0], 855 [1, 0, 0, 0, 0, 0], 856 [0, 0, 0, 0, 0, 1], 857 ] 858 ), 859 np.array( 860 [ 861 [0, 0, 0, 0, 1, 0], 862 [0, 0, 0, 0, 0, 1], 863 [1, 0, 0, 0, 0, 0], 864 [0, 1, 0, 0, 0, 0], 865 [0, 0, 1, 0, 0, 0], 866 [0, 0, 0, 1, 0, 0], 867 ] 868 ), 869 np.array( 870 [ 871 [0, 0, 0, 0, 0, 1], 872 [1, 0, 0, 0, 0, 0], 873 [0, 0, 0, 0, 1, 0], 874 [0, 0, 1, 0, 0, 0], 875 [0, 0, 0, 1, 0, 0], 876 [0, 1, 0, 0, 0, 0], 877 ] 878 ), 879 np.array( 880 [ 881 [0, 0, 0, 0, 0, 1], 882 [0, 1, 0, 0, 0, 0], 883 [0, 0, 0, 1, 0, 0], 884 [0, 0, 1, 0, 0, 0], 885 [0, 0, 0, 0, 1, 0], 886 [1, 0, 0, 0, 0, 0], 887 ] 888 ), 889 np.array( 890 [ 891 [0, 0, 0, 0, 0, 1], 892 [0, 0, 0, 1, 0, 0], 893 [0, 1, 0, 0, 0, 0], 894 [0, 0, 1, 0, 0, 0], 895 [1, 0, 0, 0, 0, 0], 896 [0, 0, 0, 0, 1, 0], 897 ] 898 ), 899 np.array( 900 [ 901 [0, 0, 0, 0, 0, 1], 902 [0, 0, 0, 0, 1, 0], 903 [1, 0, 0, 0, 0, 0], 904 [0, 0, 1, 0, 0, 0], 905 [0, 1, 0, 0, 0, 0], 906 [0, 0, 0, 1, 0, 0], 907 ] 908 ), 909 ] 910 911 selling1 = self.selling_vector 912 selling2 = other.selling_vector 913 914 vcps = np.dot(selling1, vcp_matrices)[0] 915 916 all_reflections = [] 917 for vcp in vcps: 918 for reflection_matrix in reflection_matrices: 919 all_reflections.append(np.dot(vcp, reflection_matrix)) 920 921 for reflection_matrix in reflection_matrices: 922 all_reflections.append(np.dot(selling1, reflection_matrix)) 923 924 return min([np.linalg.norm(reflection - selling2) for reflection in all_reflections]) 925 926 def __repr__(self): 927 outs = [ 928 "Lattice", 929 " abc : " + " ".join(map(repr, self.lengths)), 930 " angles : " + " ".join(map(repr, self.angles)), 931 " volume : " + repr(self.volume), 932 " A : " + " ".join(map(repr, self._matrix[0])), 933 " B : " + " ".join(map(repr, self._matrix[1])), 934 " C : " + " ".join(map(repr, self._matrix[2])), 935 ] 936 return "\n".join(outs) 937 938 def __eq__(self, other): 939 """ 940 A lattice is considered to be equal to another if the internal matrix 941 representation satisfies np.allclose(matrix1, matrix2) to be True. 942 """ 943 if other is None: 944 return False 945 # shortcut the np.allclose if the memory addresses are the same 946 # (very common in Structure.from_sites) 947 return self is other or np.allclose(self.matrix, other.matrix) 948 949 def __ne__(self, other): 950 return not self.__eq__(other) 951 952 def __hash__(self): 953 return 7 954 955 def __str__(self): 956 return "\n".join([" ".join(["%.6f" % i for i in row]) for row in self._matrix]) 957 958 def as_dict(self, verbosity: int = 0) -> Dict: 959 """ 960 Json-serialization dict representation of the Lattice. 961 962 Args: 963 verbosity (int): Verbosity level. Default of 0 only includes the 964 matrix representation. Set to 1 for more details. 965 """ 966 967 d = { 968 "@module": self.__class__.__module__, 969 "@class": self.__class__.__name__, 970 "matrix": self._matrix.tolist(), 971 } 972 a, b, c, alpha, beta, gamma = self.parameters 973 if verbosity > 0: 974 d.update( 975 { 976 "a": a, 977 "b": b, 978 "c": c, 979 "alpha": alpha, 980 "beta": beta, 981 "gamma": gamma, 982 "volume": self.volume, 983 } 984 ) 985 986 return d 987 988 def find_all_mappings( 989 self, 990 other_lattice: "Lattice", 991 ltol: float = 1e-5, 992 atol: float = 1, 993 skip_rotation_matrix: bool = False, 994 ) -> Iterator[Tuple["Lattice", Optional[np.ndarray], np.ndarray]]: 995 """ 996 Finds all mappings between current lattice and another lattice. 997 998 Args: 999 other_lattice (Lattice): Another lattice that is equivalent to 1000 this one. 1001 ltol (float): Tolerance for matching lengths. Defaults to 1e-5. 1002 atol (float): Tolerance for matching angles. Defaults to 1. 1003 skip_rotation_matrix (bool): Whether to skip calculation of the 1004 rotation matrix 1005 1006 Yields: 1007 (aligned_lattice, rotation_matrix, scale_matrix) if a mapping is 1008 found. aligned_lattice is a rotated version of other_lattice that 1009 has the same lattice parameters, but which is aligned in the 1010 coordinate system of this lattice so that translational points 1011 match up in 3D. rotation_matrix is the rotation that has to be 1012 applied to other_lattice to obtain aligned_lattice, i.e., 1013 aligned_matrix = np.inner(other_lattice, rotation_matrix) and 1014 op = SymmOp.from_rotation_and_translation(rotation_matrix) 1015 aligned_matrix = op.operate_multi(latt.matrix) 1016 Finally, scale_matrix is the integer matrix that expresses 1017 aligned_matrix as a linear combination of this 1018 lattice, i.e., aligned_matrix = np.dot(scale_matrix, self.matrix) 1019 1020 None is returned if no matches are found. 1021 """ 1022 lengths = other_lattice.lengths 1023 (alpha, beta, gamma) = other_lattice.angles 1024 1025 frac, dist, _, _ = self.get_points_in_sphere( 1026 [[0, 0, 0]], [0, 0, 0], max(lengths) * (1 + ltol), zip_results=False 1027 ) 1028 cart = self.get_cartesian_coords(frac) # type: ignore 1029 # this can't be broadcast because they're different lengths 1030 inds = [np.logical_and(dist / l < 1 + ltol, dist / l > 1 / (1 + ltol)) for l in lengths] # type: ignore 1031 c_a, c_b, c_c = (cart[i] for i in inds) 1032 f_a, f_b, f_c = (frac[i] for i in inds) 1033 l_a, l_b, l_c = (np.sum(c ** 2, axis=-1) ** 0.5 for c in (c_a, c_b, c_c)) 1034 1035 def get_angles(v1, v2, l1, l2): 1036 x = np.inner(v1, v2) / l1[:, None] / l2 1037 x[x > 1] = 1 1038 x[x < -1] = -1 1039 angles = np.arccos(x) * 180.0 / pi 1040 return angles 1041 1042 alphab = np.abs(get_angles(c_b, c_c, l_b, l_c) - alpha) < atol 1043 betab = np.abs(get_angles(c_a, c_c, l_a, l_c) - beta) < atol 1044 gammab = np.abs(get_angles(c_a, c_b, l_a, l_b) - gamma) < atol 1045 1046 for i, all_j in enumerate(gammab): 1047 inds = np.logical_and(all_j[:, None], np.logical_and(alphab, betab[i][None, :])) 1048 for j, k in np.argwhere(inds): 1049 scale_m = np.array((f_a[i], f_b[j], f_c[k]), dtype=np.int_) # type: ignore 1050 if abs(np.linalg.det(scale_m)) < 1e-8: 1051 continue 1052 1053 aligned_m = np.array((c_a[i], c_b[j], c_c[k])) 1054 1055 if skip_rotation_matrix: 1056 rotation_m = None 1057 else: 1058 rotation_m = np.linalg.solve(aligned_m, other_lattice.matrix) 1059 1060 yield Lattice(aligned_m), rotation_m, scale_m 1061 1062 def find_mapping( 1063 self, 1064 other_lattice: "Lattice", 1065 ltol: float = 1e-5, 1066 atol: float = 1, 1067 skip_rotation_matrix: bool = False, 1068 ) -> Optional[Tuple["Lattice", Optional[np.ndarray], np.ndarray]]: 1069 """ 1070 Finds a mapping between current lattice and another lattice. There 1071 are an infinite number of choices of basis vectors for two entirely 1072 equivalent lattices. This method returns a mapping that maps 1073 other_lattice to this lattice. 1074 1075 Args: 1076 other_lattice (Lattice): Another lattice that is equivalent to 1077 this one. 1078 ltol (float): Tolerance for matching lengths. Defaults to 1e-5. 1079 atol (float): Tolerance for matching angles. Defaults to 1. 1080 1081 Returns: 1082 (aligned_lattice, rotation_matrix, scale_matrix) if a mapping is 1083 found. aligned_lattice is a rotated version of other_lattice that 1084 has the same lattice parameters, but which is aligned in the 1085 coordinate system of this lattice so that translational points 1086 match up in 3D. rotation_matrix is the rotation that has to be 1087 applied to other_lattice to obtain aligned_lattice, i.e., 1088 aligned_matrix = np.inner(other_lattice, rotation_matrix) and 1089 op = SymmOp.from_rotation_and_translation(rotation_matrix) 1090 aligned_matrix = op.operate_multi(latt.matrix) 1091 Finally, scale_matrix is the integer matrix that expresses 1092 aligned_matrix as a linear combination of this 1093 lattice, i.e., aligned_matrix = np.dot(scale_matrix, self.matrix) 1094 1095 None is returned if no matches are found. 1096 """ 1097 for x in self.find_all_mappings(other_lattice, ltol, atol, skip_rotation_matrix=skip_rotation_matrix): 1098 return x 1099 return None 1100 1101 def get_lll_reduced_lattice(self, delta: float = 0.75) -> "Lattice": 1102 """ 1103 :param delta: Delta parameter. 1104 :return: LLL reduced Lattice. 1105 """ 1106 if delta not in self._lll_matrix_mappings: 1107 self._lll_matrix_mappings[delta] = self._calculate_lll() 1108 return Lattice(self._lll_matrix_mappings[delta][0]) 1109 1110 def _calculate_lll(self, delta: float = 0.75) -> Tuple[np.ndarray, np.ndarray]: 1111 """ 1112 Performs a Lenstra-Lenstra-Lovasz lattice basis reduction to obtain a 1113 c-reduced basis. This method returns a basis which is as "good" as 1114 possible, with "good" defined by orthongonality of the lattice vectors. 1115 1116 This basis is used for all the periodic boundary condition calculations. 1117 1118 Args: 1119 delta (float): Reduction parameter. Default of 0.75 is usually 1120 fine. 1121 1122 Returns: 1123 Reduced lattice matrix, mapping to get to that lattice. 1124 """ 1125 # Transpose the lattice matrix first so that basis vectors are columns. 1126 # Makes life easier. 1127 # pylint: disable=E1136,E1137 1128 a = self._matrix.copy().T 1129 1130 b = np.zeros((3, 3)) # Vectors after the Gram-Schmidt process 1131 u = np.zeros((3, 3)) # Gram-Schmidt coefficients 1132 m = np.zeros(3) # These are the norm squared of each vec. 1133 1134 b[:, 0] = a[:, 0] 1135 m[0] = dot(b[:, 0], b[:, 0]) 1136 for i in range(1, 3): 1137 u[i, 0:i] = dot(a[:, i].T, b[:, 0:i]) / m[0:i] 1138 b[:, i] = a[:, i] - dot(b[:, 0:i], u[i, 0:i].T) 1139 m[i] = dot(b[:, i], b[:, i]) 1140 1141 k = 2 1142 1143 mapping = np.identity(3, dtype=np.double) 1144 while k <= 3: 1145 # Size reduction. 1146 for i in range(k - 1, 0, -1): 1147 q = round(u[k - 1, i - 1]) 1148 if q != 0: 1149 # Reduce the k-th basis vector. 1150 a[:, k - 1] = a[:, k - 1] - q * a[:, i - 1] 1151 mapping[:, k - 1] = mapping[:, k - 1] - q * mapping[:, i - 1] 1152 uu = list(u[i - 1, 0 : (i - 1)]) 1153 uu.append(1) 1154 # Update the GS coefficients. 1155 u[k - 1, 0:i] = u[k - 1, 0:i] - q * np.array(uu) 1156 1157 # Check the Lovasz condition. 1158 if dot(b[:, k - 1], b[:, k - 1]) >= (delta - abs(u[k - 1, k - 2]) ** 2) * dot(b[:, (k - 2)], b[:, (k - 2)]): 1159 # Increment k if the Lovasz condition holds. 1160 k += 1 1161 else: 1162 # If the Lovasz condition fails, 1163 # swap the k-th and (k-1)-th basis vector 1164 v = a[:, k - 1].copy() 1165 a[:, k - 1] = a[:, k - 2].copy() 1166 a[:, k - 2] = v 1167 1168 v_m = mapping[:, k - 1].copy() 1169 mapping[:, k - 1] = mapping[:, k - 2].copy() 1170 mapping[:, k - 2] = v_m 1171 1172 # Update the Gram-Schmidt coefficients 1173 for s in range(k - 1, k + 1): 1174 u[s - 1, 0 : (s - 1)] = dot(a[:, s - 1].T, b[:, 0 : (s - 1)]) / m[0 : (s - 1)] 1175 b[:, s - 1] = a[:, s - 1] - dot(b[:, 0 : (s - 1)], u[s - 1, 0 : (s - 1)].T) 1176 m[s - 1] = dot(b[:, s - 1], b[:, s - 1]) 1177 1178 if k > 2: 1179 k -= 1 1180 else: 1181 # We have to do p/q, so do lstsq(q.T, p.T).T instead. 1182 p = dot(a[:, k:3].T, b[:, (k - 2) : k]) 1183 q = np.diag(m[(k - 2) : k]) 1184 result = np.linalg.lstsq(q.T, p.T, rcond=None)[0].T # type: ignore 1185 u[k:3, (k - 2) : k] = result 1186 1187 return a.T, mapping.T 1188 1189 def get_lll_frac_coords(self, frac_coords: ArrayLike) -> np.ndarray: 1190 """ 1191 Given fractional coordinates in the lattice basis, returns corresponding 1192 fractional coordinates in the lll basis. 1193 """ 1194 return dot(frac_coords, self.lll_inverse) 1195 1196 def get_frac_coords_from_lll(self, lll_frac_coords: ArrayLike) -> np.ndarray: 1197 """ 1198 Given fractional coordinates in the lll basis, returns corresponding 1199 fractional coordinates in the lattice basis. 1200 """ 1201 return dot(lll_frac_coords, self.lll_mapping) 1202 1203 def get_niggli_reduced_lattice(self, tol: float = 1e-5) -> "Lattice": 1204 """ 1205 Get the Niggli reduced lattice using the numerically stable algo 1206 proposed by R. W. Grosse-Kunstleve, N. K. Sauter, & P. D. Adams, 1207 Acta Crystallographica Section A Foundations of Crystallography, 2003, 1208 60(1), 1-6. doi:10.1107/S010876730302186X 1209 1210 Args: 1211 tol (float): The numerical tolerance. The default of 1e-5 should 1212 result in stable behavior for most cases. 1213 1214 Returns: 1215 Niggli-reduced lattice. 1216 """ 1217 # lll reduction is more stable for skewed cells 1218 matrix = self.lll_matrix 1219 e = tol * self.volume ** (1 / 3) 1220 1221 # Define metric tensor 1222 G = np.dot(matrix, matrix.T) 1223 1224 # This sets an upper limit on the number of iterations. 1225 for count in range(100): 1226 # The steps are labelled as Ax as per the labelling scheme in the 1227 # paper. 1228 (A, B, C, E, N, Y) = ( 1229 G[0, 0], 1230 G[1, 1], 1231 G[2, 2], 1232 2 * G[1, 2], 1233 2 * G[0, 2], 1234 2 * G[0, 1], 1235 ) 1236 1237 if A > B + e or (abs(A - B) < e and abs(E) > abs(N) + e): 1238 # A1 1239 M = [[0, -1, 0], [-1, 0, 0], [0, 0, -1]] 1240 G = dot(transpose(M), dot(G, M)) 1241 if (B > C + e) or (abs(B - C) < e and abs(N) > abs(Y) + e): 1242 # A2 1243 M = [[-1, 0, 0], [0, 0, -1], [0, -1, 0]] 1244 G = dot(transpose(M), dot(G, M)) 1245 continue 1246 1247 l = 0 if abs(E) < e else E / abs(E) 1248 m = 0 if abs(N) < e else N / abs(N) 1249 n = 0 if abs(Y) < e else Y / abs(Y) 1250 if l * m * n == 1: 1251 # A3 1252 i = -1 if l == -1 else 1 1253 j = -1 if m == -1 else 1 1254 k = -1 if n == -1 else 1 1255 M = [[i, 0, 0], [0, j, 0], [0, 0, k]] 1256 G = dot(transpose(M), dot(G, M)) 1257 elif l * m * n == 0 or l * m * n == -1: 1258 # A4 1259 i = -1 if l == 1 else 1 1260 j = -1 if m == 1 else 1 1261 k = -1 if n == 1 else 1 1262 1263 if i * j * k == -1: 1264 if n == 0: 1265 k = -1 1266 elif m == 0: 1267 j = -1 1268 elif l == 0: 1269 i = -1 1270 M = [[i, 0, 0], [0, j, 0], [0, 0, k]] 1271 G = dot(transpose(M), dot(G, M)) 1272 1273 (A, B, C, E, N, Y) = ( 1274 G[0, 0], 1275 G[1, 1], 1276 G[2, 2], 1277 2 * G[1, 2], 1278 2 * G[0, 2], 1279 2 * G[0, 1], 1280 ) 1281 1282 # A5 1283 if abs(E) > B + e or (abs(E - B) < e and 2 * N < Y - e) or (abs(E + B) < e and Y < -e): 1284 M = [[1, 0, 0], [0, 1, -E / abs(E)], [0, 0, 1]] 1285 G = dot(transpose(M), dot(G, M)) 1286 continue 1287 1288 # A6 1289 if abs(N) > A + e or (abs(A - N) < e and 2 * E < Y - e) or (abs(A + N) < e and Y < -e): 1290 M = [[1, 0, -N / abs(N)], [0, 1, 0], [0, 0, 1]] 1291 G = dot(transpose(M), dot(G, M)) 1292 continue 1293 1294 # A7 1295 if abs(Y) > A + e or (abs(A - Y) < e and 2 * E < N - e) or (abs(A + Y) < e and N < -e): 1296 M = [[1, -Y / abs(Y), 0], [0, 1, 0], [0, 0, 1]] 1297 G = dot(transpose(M), dot(G, M)) 1298 continue 1299 1300 # A8 1301 if E + N + Y + A + B < -e or (abs(E + N + Y + A + B) < e < Y + (A + N) * 2): 1302 M = [[1, 0, 1], [0, 1, 1], [0, 0, 1]] 1303 G = dot(transpose(M), dot(G, M)) 1304 continue 1305 1306 break 1307 1308 A = G[0, 0] 1309 B = G[1, 1] 1310 C = G[2, 2] 1311 E = 2 * G[1, 2] 1312 N = 2 * G[0, 2] 1313 Y = 2 * G[0, 1] 1314 a = math.sqrt(A) 1315 b = math.sqrt(B) 1316 c = math.sqrt(C) 1317 alpha = math.acos(E / 2 / b / c) / math.pi * 180 1318 beta = math.acos(N / 2 / a / c) / math.pi * 180 1319 gamma = math.acos(Y / 2 / a / b) / math.pi * 180 1320 1321 latt = Lattice.from_parameters(a, b, c, alpha, beta, gamma) 1322 1323 mapped = self.find_mapping(latt, e, skip_rotation_matrix=True) 1324 if mapped is not None: 1325 if np.linalg.det(mapped[0].matrix) > 0: 1326 return mapped[0] 1327 return Lattice(-mapped[0].matrix) 1328 1329 raise ValueError("can't find niggli") 1330 1331 def scale(self, new_volume: float) -> "Lattice": 1332 """ 1333 Return a new Lattice with volume new_volume by performing a 1334 scaling of the lattice vectors so that length proportions and angles 1335 are preserved. 1336 1337 Args: 1338 new_volume: 1339 New volume to scale to. 1340 1341 Returns: 1342 New lattice with desired volume. 1343 """ 1344 versors = self.matrix / self.abc 1345 1346 geo_factor = abs(dot(np.cross(versors[0], versors[1]), versors[2])) 1347 1348 ratios = np.array(self.abc) / self.c 1349 1350 new_c = (new_volume / (geo_factor * np.prod(ratios))) ** (1 / 3.0) 1351 1352 return Lattice(versors * (new_c * ratios)) 1353 1354 def get_wigner_seitz_cell(self) -> List[List[np.ndarray]]: 1355 """ 1356 Returns the Wigner-Seitz cell for the given lattice. 1357 1358 Returns: 1359 A list of list of coordinates. 1360 Each element in the list is a "facet" of the boundary of the 1361 Wigner Seitz cell. For instance, a list of four coordinates will 1362 represent a square facet. 1363 """ 1364 vec1 = self._matrix[0] 1365 vec2 = self._matrix[1] 1366 vec3 = self._matrix[2] 1367 1368 list_k_points = [] 1369 for i, j, k in itertools.product([-1, 0, 1], [-1, 0, 1], [-1, 0, 1]): 1370 list_k_points.append(i * vec1 + j * vec2 + k * vec3) 1371 # pylint: disable=C0415 1372 from scipy.spatial import Voronoi 1373 1374 tess = Voronoi(list_k_points) 1375 to_return = [] 1376 for r in tess.ridge_dict: 1377 if r[0] == 13 or r[1] == 13: 1378 to_return.append([tess.vertices[i] for i in tess.ridge_dict[r]]) 1379 1380 return to_return 1381 1382 def get_brillouin_zone(self) -> List[List[np.ndarray]]: 1383 """ 1384 Returns the Wigner-Seitz cell for the reciprocal lattice, aka the 1385 Brillouin Zone. 1386 1387 Returns: 1388 A list of list of coordinates. 1389 Each element in the list is a "facet" of the boundary of the 1390 Brillouin Zone. For instance, a list of four coordinates will 1391 represent a square facet. 1392 """ 1393 return self.reciprocal_lattice.get_wigner_seitz_cell() 1394 1395 def dot(self, coords_a: ArrayLike, coords_b: ArrayLike, frac_coords: bool = False) -> np.ndarray: 1396 """ 1397 Compute the scalar product of vector(s). 1398 1399 Args: 1400 coords_a, coords_b: Array-like objects with the coordinates. 1401 frac_coords (bool): Boolean stating whether the vector 1402 corresponds to fractional or cartesian coordinates. 1403 1404 Returns: 1405 one-dimensional `numpy` array. 1406 """ 1407 coords_a, coords_b = ( 1408 np.reshape(coords_a, (-1, 3)), 1409 np.reshape(coords_b, (-1, 3)), 1410 ) 1411 1412 if len(coords_a) != len(coords_b): 1413 raise ValueError("") 1414 1415 if np.iscomplexobj(coords_a) or np.iscomplexobj(coords_b): 1416 raise TypeError("Complex array!") 1417 1418 if not frac_coords: 1419 cart_a, cart_b = coords_a, coords_b 1420 else: 1421 cart_a = np.reshape([self.get_cartesian_coords(vec) for vec in coords_a], (-1, 3)) # type: ignore 1422 cart_b = np.reshape([self.get_cartesian_coords(vec) for vec in coords_b], (-1, 3)) # type: ignore 1423 1424 return np.array([dot(a, b) for a, b in zip(cart_a, cart_b)]) 1425 1426 def norm(self, coords: ArrayLike, frac_coords: bool = True) -> float: 1427 """ 1428 Compute the norm of vector(s). 1429 1430 Args: 1431 coords: 1432 Array-like object with the coordinates. 1433 frac_coords: 1434 Boolean stating whether the vector corresponds to fractional or 1435 cartesian coordinates. 1436 1437 Returns: 1438 one-dimensional `numpy` array. 1439 """ 1440 return np.sqrt(self.dot(coords, coords, frac_coords=frac_coords)) 1441 1442 def get_points_in_sphere( 1443 self, 1444 frac_points: ArrayLike, 1445 center: ArrayLike, 1446 r: float, 1447 zip_results=True, 1448 ) -> Union[List[Tuple[np.ndarray, float, int, np.ndarray]], List[np.ndarray], List]: 1449 """ 1450 Find all points within a sphere from the point taking into account 1451 periodic boundary conditions. This includes sites in other periodic 1452 images. 1453 1454 Algorithm: 1455 1456 1. place sphere of radius r in crystal and determine minimum supercell 1457 (parallelpiped) which would contain a sphere of radius r. for this 1458 we need the projection of a_1 on a unit vector perpendicular 1459 to a_2 & a_3 (i.e. the unit vector in the direction b_1) to 1460 determine how many a_1"s it will take to contain the sphere. 1461 1462 Nxmax = r * length_of_b_1 / (2 Pi) 1463 1464 2. keep points falling within r. 1465 1466 Args: 1467 frac_points: All points in the lattice in fractional coordinates. 1468 center: Cartesian coordinates of center of sphere. 1469 r: radius of sphere. 1470 zip_results (bool): Whether to zip the results together to group by 1471 point, or return the raw fcoord, dist, index arrays 1472 1473 Returns: 1474 if zip_results: 1475 [(fcoord, dist, index, supercell_image) ...] since most of the time, subsequent 1476 processing requires the distance, index number of the atom, or index of the image 1477 else: 1478 fcoords, dists, inds, image 1479 """ 1480 try: 1481 # pylint: disable=C0415 1482 from pymatgen.optimization.neighbors import ( 1483 find_points_in_spheres, # type: ignore 1484 ) 1485 except ImportError: 1486 return self.get_points_in_sphere_py(frac_points=frac_points, center=center, r=r, zip_results=zip_results) 1487 else: 1488 frac_points = np.ascontiguousarray(frac_points, dtype=np.float_) # type: ignore 1489 r = float(r) 1490 lattice_matrix = np.array(self.matrix) 1491 lattice_matrix = np.ascontiguousarray(lattice_matrix) 1492 cart_coords = self.get_cartesian_coords(frac_points) 1493 _, indices, images, distances = find_points_in_spheres( 1494 all_coords=cart_coords, 1495 center_coords=np.ascontiguousarray([center], dtype=float), 1496 r=r, 1497 pbc=np.array([1, 1, 1]), 1498 lattice=lattice_matrix, 1499 tol=1e-8, 1500 ) 1501 if len(indices) < 1: 1502 return [] if zip_results else [()] * 4 1503 fcoords = frac_points[indices] + images 1504 if zip_results: 1505 return list( 1506 zip( 1507 fcoords, 1508 distances, 1509 indices, 1510 images, 1511 ) 1512 ) 1513 return [ 1514 fcoords, 1515 distances, 1516 indices, 1517 images, 1518 ] 1519 1520 def get_points_in_sphere_py( 1521 self, 1522 frac_points: ArrayLike, 1523 center: ArrayLike, 1524 r: float, 1525 zip_results=True, 1526 ) -> Union[List[Tuple[np.ndarray, float, int, np.ndarray]], List[np.ndarray],]: 1527 """ 1528 Find all points within a sphere from the point taking into account 1529 periodic boundary conditions. This includes sites in other periodic 1530 images. 1531 1532 Algorithm: 1533 1534 1. place sphere of radius r in crystal and determine minimum supercell 1535 (parallelpiped) which would contain a sphere of radius r. for this 1536 we need the projection of a_1 on a unit vector perpendicular 1537 to a_2 & a_3 (i.e. the unit vector in the direction b_1) to 1538 determine how many a_1"s it will take to contain the sphere. 1539 1540 Nxmax = r * length_of_b_1 / (2 Pi) 1541 1542 2. keep points falling within r. 1543 1544 Args: 1545 frac_points: All points in the lattice in fractional coordinates. 1546 center: Cartesian coordinates of center of sphere. 1547 r: radius of sphere. 1548 zip_results (bool): Whether to zip the results together to group by 1549 point, or return the raw fcoord, dist, index arrays 1550 1551 Returns: 1552 if zip_results: 1553 [(fcoord, dist, index, supercell_image) ...] since most of the time, subsequent 1554 processing requires the distance, index number of the atom, or index of the image 1555 else: 1556 fcoords, dists, inds, image 1557 """ 1558 cart_coords = self.get_cartesian_coords(frac_points) 1559 neighbors = get_points_in_spheres( 1560 all_coords=cart_coords, 1561 center_coords=np.array([center]), 1562 r=r, 1563 pbc=True, 1564 numerical_tol=1e-8, 1565 lattice=self, 1566 return_fcoords=True, 1567 )[0] 1568 if len(neighbors) < 1: 1569 return [] if zip_results else [()] * 4 # type: ignore 1570 if zip_results: 1571 return neighbors 1572 return [np.array(i) for i in list(zip(*neighbors))] 1573 1574 @deprecated(get_points_in_sphere, "This is retained purely for checking purposes.") 1575 def get_points_in_sphere_old( 1576 self, 1577 frac_points: ArrayLike, 1578 center: ArrayLike, 1579 r: float, 1580 zip_results=True, 1581 ) -> Union[ 1582 List[Tuple[np.ndarray, float, int, np.ndarray]], 1583 Tuple[List[np.ndarray], List[float], List[int], List[np.ndarray]], 1584 ]: 1585 """ 1586 Find all points within a sphere from the point taking into account 1587 periodic boundary conditions. This includes sites in other periodic 1588 images. 1589 1590 Algorithm: 1591 1592 1. place sphere of radius r in crystal and determine minimum supercell 1593 (parallelpiped) which would contain a sphere of radius r. for this 1594 we need the projection of a_1 on a unit vector perpendicular 1595 to a_2 & a_3 (i.e. the unit vector in the direction b_1) to 1596 determine how many a_1"s it will take to contain the sphere. 1597 1598 Nxmax = r * length_of_b_1 / (2 Pi) 1599 1600 2. keep points falling within r. 1601 1602 Args: 1603 frac_points: All points in the lattice in fractional coordinates. 1604 center: Cartesian coordinates of center of sphere. 1605 r: radius of sphere. 1606 zip_results (bool): Whether to zip the results together to group by 1607 point, or return the raw fcoord, dist, index arrays 1608 1609 Returns: 1610 if zip_results: 1611 [(fcoord, dist, index, supercell_image) ...] since most of the time, subsequent 1612 processing requires the distance, index number of the atom, or index of the image 1613 else: 1614 fcoords, dists, inds, image 1615 """ 1616 # TODO: refactor to use lll matrix (nmax will be smaller) 1617 # Determine the maximum number of supercells in each direction 1618 # required to contain a sphere of radius n 1619 recp_len = np.array(self.reciprocal_lattice.abc) / (2 * pi) 1620 nmax = float(r) * recp_len + 0.01 1621 1622 # Get the fractional coordinates of the center of the sphere 1623 pcoords = self.get_fractional_coords(center) 1624 center = np.array(center) 1625 1626 # Prepare the list of output atoms 1627 n = len(frac_points) # type: ignore 1628 fcoords = np.array(frac_points) % 1 1629 indices = np.arange(n) 1630 1631 # Generate all possible images that could be within `r` of `center` 1632 mins = np.floor(pcoords - nmax) 1633 maxes = np.ceil(pcoords + nmax) 1634 arange = np.arange(start=mins[0], stop=maxes[0], dtype=np.int_) 1635 brange = np.arange(start=mins[1], stop=maxes[1], dtype=np.int_) 1636 crange = np.arange(start=mins[2], stop=maxes[2], dtype=np.int_) 1637 arange = arange[:, None] * np.array([1, 0, 0], dtype=np.int_)[None, :] 1638 brange = brange[:, None] * np.array([0, 1, 0], dtype=np.int_)[None, :] 1639 crange = crange[:, None] * np.array([0, 0, 1], dtype=np.int_)[None, :] 1640 images = arange[:, None, None] + brange[None, :, None] + crange[None, None, :] 1641 1642 # Generate the coordinates of all atoms within these images 1643 shifted_coords = fcoords[:, None, None, None, :] + images[None, :, :, :, :] 1644 1645 # Determine distance from `center` 1646 cart_coords = self.get_cartesian_coords(fcoords) 1647 cart_images = self.get_cartesian_coords(images) 1648 coords = cart_coords[:, None, None, None, :] + cart_images[None, :, :, :, :] 1649 coords -= center[None, None, None, None, :] 1650 coords **= 2 1651 d_2 = np.sum(coords, axis=4) 1652 1653 # Determine which points are within `r` of `center` 1654 within_r = np.where(d_2 <= r ** 2) 1655 # `within_r` now contains the coordinates of each image that is 1656 # inside of the cutoff distance. It has 4 coordinates: 1657 # 0 - index of the image within `frac_points` 1658 # 1,2,3 - index of the supercell which holds the images in the x, y, z directions 1659 1660 if zip_results: 1661 return list( 1662 zip( 1663 shifted_coords[within_r], 1664 np.sqrt(d_2[within_r]), 1665 indices[within_r[0]], 1666 images[within_r[1:]], 1667 ) 1668 ) 1669 return ( 1670 shifted_coords[within_r], 1671 np.sqrt(d_2[within_r]), 1672 indices[within_r[0]], 1673 images[within_r[1:]], 1674 ) 1675 1676 def get_all_distances( 1677 self, 1678 fcoords1: ArrayLike, 1679 fcoords2: ArrayLike, 1680 ) -> np.ndarray: 1681 """ 1682 Returns the distances between two lists of coordinates taking into 1683 account periodic boundary conditions and the lattice. Note that this 1684 computes an MxN array of distances (i.e. the distance between each 1685 point in fcoords1 and every coordinate in fcoords2). This is 1686 different functionality from pbc_diff. 1687 1688 Args: 1689 fcoords1: First set of fractional coordinates. e.g., [0.5, 0.6, 1690 0.7] or [[1.1, 1.2, 4.3], [0.5, 0.6, 0.7]]. It can be a single 1691 coord or any array of coords. 1692 fcoords2: Second set of fractional coordinates. 1693 1694 Returns: 1695 2d array of cartesian distances. E.g the distance between 1696 fcoords1[i] and fcoords2[j] is distances[i,j] 1697 """ 1698 v, d2 = pbc_shortest_vectors(self, fcoords1, fcoords2, return_d2=True) 1699 return np.sqrt(d2) 1700 1701 def is_hexagonal(self, hex_angle_tol: float = 5, hex_length_tol: float = 0.01) -> bool: 1702 """ 1703 :param hex_angle_tol: Angle tolerance 1704 :param hex_length_tol: Length tolerance 1705 :return: Whether lattice corresponds to hexagonal lattice. 1706 """ 1707 lengths = self.lengths 1708 angles = self.angles 1709 right_angles = [i for i in range(3) if abs(angles[i] - 90) < hex_angle_tol] 1710 hex_angles = [ 1711 i for i in range(3) if abs(angles[i] - 60) < hex_angle_tol or abs(angles[i] - 120) < hex_angle_tol 1712 ] 1713 1714 return ( 1715 len(right_angles) == 2 1716 and len(hex_angles) == 1 1717 and abs(lengths[right_angles[0]] - lengths[right_angles[1]]) < hex_length_tol 1718 ) 1719 1720 def get_distance_and_image( 1721 self, 1722 frac_coords1: ArrayLike, 1723 frac_coords2: ArrayLike, 1724 jimage: Optional[ArrayLike] = None, 1725 ) -> Tuple[float, np.ndarray]: 1726 """ 1727 Gets distance between two frac_coords assuming periodic boundary 1728 conditions. If the index jimage is not specified it selects the j 1729 image nearest to the i atom and returns the distance and jimage 1730 indices in terms of lattice vector translations. If the index jimage 1731 is specified it returns the distance between the frac_coords1 and 1732 the specified jimage of frac_coords2, and the given jimage is also 1733 returned. 1734 1735 Args: 1736 frac_coords1 (3x1 array): Reference fcoords to get distance from. 1737 frac_coords2 (3x1 array): fcoords to get distance from. 1738 jimage (3x1 array): Specific periodic image in terms of 1739 lattice translations, e.g., [1,0,0] implies to take periodic 1740 image that is one a-lattice vector away. If jimage is None, 1741 the image that is nearest to the site is found. 1742 1743 Returns: 1744 (distance, jimage): distance and periodic lattice translations 1745 of the other site for which the distance applies. This means that 1746 the distance between frac_coords1 and (jimage + frac_coords2) is 1747 equal to distance. 1748 """ 1749 if jimage is None: 1750 v, d2 = pbc_shortest_vectors(self, frac_coords1, frac_coords2, return_d2=True) 1751 fc = self.get_fractional_coords(v[0][0]) + frac_coords1 - frac_coords2 1752 fc = np.array(np.round(fc), dtype=np.int_) 1753 return np.sqrt(d2[0, 0]), fc 1754 1755 jimage = np.array(jimage) 1756 mapped_vec = self.get_cartesian_coords(jimage + frac_coords2 - frac_coords1) 1757 return np.linalg.norm(mapped_vec), jimage 1758 1759 def get_miller_index_from_coords( 1760 self, 1761 coords: ArrayLike, 1762 coords_are_cartesian: bool = True, 1763 round_dp: int = 4, 1764 verbose: bool = True, 1765 ) -> Tuple[int, int, int]: 1766 """ 1767 Get the Miller index of a plane from a list of site coordinates. 1768 1769 A minimum of 3 sets of coordinates are required. If more than 3 sets of 1770 coordinates are given, the best plane that minimises the distance to all 1771 points will be calculated. 1772 1773 Args: 1774 coords (iterable): A list or numpy array of coordinates. Can be 1775 cartesian or fractional coordinates. If more than three sets of 1776 coordinates are provided, the best plane that minimises the 1777 distance to all sites will be calculated. 1778 coords_are_cartesian (bool, optional): Whether the coordinates are 1779 in cartesian space. If using fractional coordinates set to 1780 False. 1781 round_dp (int, optional): The number of decimal places to round the 1782 miller index to. 1783 verbose (bool, optional): Whether to print warnings. 1784 1785 Returns: 1786 (tuple): The Miller index. 1787 """ 1788 if coords_are_cartesian: 1789 coords = [self.get_fractional_coords(c) for c in coords] # type: ignore 1790 1791 coords = np.asarray(coords) 1792 g = coords.sum(axis=0) / coords.shape[0] 1793 1794 # run singular value decomposition 1795 _, _, vh = np.linalg.svd(coords - g) 1796 1797 # get unitary normal vector 1798 u_norm = vh[2, :] 1799 return get_integer_index(u_norm, round_dp=round_dp, verbose=verbose) 1800 1801 def get_recp_symmetry_operation(self, symprec: float = 0.01) -> List: 1802 """ 1803 Find the symmetric operations of the reciprocal lattice, 1804 to be used for hkl transformations 1805 Args: 1806 symprec: default is 0.001 1807 """ 1808 recp_lattice = self.reciprocal_lattice_crystallographic 1809 # get symmetry operations from input conventional unit cell 1810 # Need to make sure recp lattice is big enough, otherwise symmetry 1811 # determination will fail. We set the overall volume to 1. 1812 recp_lattice = recp_lattice.scale(1) 1813 # need a localized import of structure to build a 1814 # pseudo empty lattice for SpacegroupAnalyzer 1815 # pylint: disable=C0415 1816 from pymatgen.core.structure import Structure 1817 from pymatgen.symmetry.analyzer import SpacegroupAnalyzer 1818 1819 recp = Structure(recp_lattice, ["H"], [[0, 0, 0]]) 1820 # Creates a function that uses the symmetry operations in the 1821 # structure to find Miller indices that might give repetitive slabs 1822 analyzer = SpacegroupAnalyzer(recp, symprec=symprec) 1823 recp_symmops = analyzer.get_symmetry_operations() 1824 1825 return recp_symmops 1826 1827 1828def get_integer_index(miller_index: Sequence[float], round_dp: int = 4, verbose: bool = True) -> Tuple[int, int, int]: 1829 """ 1830 Attempt to convert a vector of floats to whole numbers. 1831 1832 Args: 1833 miller_index (list of float): A list miller indexes. 1834 round_dp (int, optional): The number of decimal places to round the 1835 miller index to. 1836 verbose (bool, optional): Whether to print warnings. 1837 1838 Returns: 1839 (tuple): The Miller index. 1840 """ 1841 mi = np.asarray(miller_index) 1842 # deal with the case we have small irregular floats 1843 # that are all equal or factors of each other 1844 mi /= min([m for m in mi if m != 0]) 1845 mi /= np.max(np.abs(mi)) 1846 1847 # deal with the case we have nice fractions 1848 md = [Fraction(n).limit_denominator(12).denominator for n in mi] 1849 mi *= reduce(lambda x, y: x * y, md) 1850 int_miller_index = np.int_(np.round(mi, 1)) 1851 mi /= np.abs(reduce(math.gcd, int_miller_index)) 1852 1853 # round to a reasonable precision 1854 mi = np.array([round(h, round_dp) for h in mi]) 1855 1856 # need to recalculate this after rounding as values may have changed 1857 int_miller_index = np.int_(np.round(mi, 1)) 1858 if np.any(np.abs(mi - int_miller_index) > 1e-6) and verbose: 1859 warnings.warn("Non-integer encountered in Miller index") 1860 else: 1861 mi = int_miller_index 1862 1863 # minimise the number of negative indexes 1864 mi += 0 # converts -0 to 0 1865 1866 def n_minus(index): 1867 return len([h for h in index if h < 0]) 1868 1869 if n_minus(mi) > n_minus(mi * -1): 1870 mi *= -1 1871 1872 # if only one index is negative, make sure it is the smallest 1873 # e.g. (-2 1 0) -> (2 -1 0) 1874 if sum(mi != 0) == 2 and n_minus(mi) == 1 and abs(min(mi)) > max(mi): 1875 mi *= -1 1876 1877 return tuple(mi) # type: ignore 1878 1879 1880def get_points_in_spheres( 1881 all_coords: np.ndarray, 1882 center_coords: np.ndarray, 1883 r: float, 1884 pbc: Union[bool, List[bool]] = True, 1885 numerical_tol: float = 1e-8, 1886 lattice: Lattice = None, 1887 return_fcoords: bool = False, 1888) -> List[List[Tuple[np.ndarray, float, int, np.ndarray]]]: 1889 """ 1890 For each point in `center_coords`, get all the neighboring points in `all_coords` that are within the 1891 cutoff radius `r`. 1892 1893 Args: 1894 all_coords: (list of cartesian coordinates) all available points 1895 center_coords: (list of cartesian coordinates) all centering points 1896 r: (float) cutoff radius 1897 pbc: (bool or a list of bool) whether to set periodic boundaries 1898 numerical_tol: (float) numerical tolerance 1899 lattice: (Lattice) lattice to consider when PBC is enabled 1900 return_fcoords: (bool) whether to return fractional coords when pbc is set. 1901 Returns: 1902 List[List[Tuple[coords, distance, index, image]]] 1903 """ 1904 if isinstance(pbc, bool): 1905 pbc = [pbc] * 3 1906 pbc = np.array(pbc, dtype=np.bool_) # type: ignore 1907 if return_fcoords and lattice is None: 1908 raise ValueError("Lattice needs to be supplied to compute fractional coordinates") 1909 center_coords_min = np.min(center_coords, axis=0) 1910 center_coords_max = np.max(center_coords, axis=0) 1911 # The lower bound of all considered atom coords 1912 global_min = center_coords_min - r - numerical_tol 1913 global_max = center_coords_max + r + numerical_tol 1914 if np.any(pbc): 1915 if lattice is None: 1916 raise ValueError("Lattice needs to be supplied when considering periodic boundary") 1917 recp_len = np.array(lattice.reciprocal_lattice.abc) 1918 maxr = np.ceil((r + 0.15) * recp_len / (2 * math.pi)) 1919 frac_coords = lattice.get_fractional_coords(center_coords) 1920 nmin_temp = np.floor(np.min(frac_coords, axis=0)) - maxr 1921 nmax_temp = np.ceil(np.max(frac_coords, axis=0)) + maxr 1922 nmin = np.zeros_like(nmin_temp) 1923 nmin[pbc] = nmin_temp[pbc] 1924 nmax = np.ones_like(nmax_temp) 1925 nmax[pbc] = nmax_temp[pbc] 1926 all_ranges = [np.arange(x, y, dtype="int64") for x, y in zip(nmin, nmax)] 1927 matrix = lattice.matrix 1928 # temporarily hold the fractional coordinates 1929 image_offsets = lattice.get_fractional_coords(all_coords) 1930 all_fcoords = [] 1931 # only wrap periodic boundary 1932 for k in range(3): 1933 if pbc[k]: # type: ignore 1934 all_fcoords.append(np.mod(image_offsets[:, k : k + 1], 1)) 1935 else: 1936 all_fcoords.append(image_offsets[:, k : k + 1]) 1937 all_fcoords = np.concatenate(all_fcoords, axis=1) 1938 image_offsets = image_offsets - all_fcoords 1939 coords_in_cell = np.dot(all_fcoords, matrix) 1940 # Filter out those beyond max range 1941 valid_coords = [] 1942 valid_images = [] 1943 valid_indices = [] 1944 for image in itertools.product(*all_ranges): 1945 coords = np.dot(image, matrix) + coords_in_cell 1946 valid_index_bool = np.all( 1947 np.bitwise_and(coords > global_min[None, :], coords < global_max[None, :]), 1948 axis=1, 1949 ) 1950 ind = np.arange(len(all_coords)) 1951 if np.any(valid_index_bool): 1952 valid_coords.append(coords[valid_index_bool]) 1953 valid_images.append(np.tile(image, [np.sum(valid_index_bool), 1]) - image_offsets[valid_index_bool]) 1954 valid_indices.extend([k for k in ind if valid_index_bool[k]]) 1955 if len(valid_coords) < 1: 1956 return [[]] * len(center_coords) 1957 valid_coords = np.concatenate(valid_coords, axis=0) 1958 valid_images = np.concatenate(valid_images, axis=0) 1959 1960 else: 1961 valid_coords = all_coords # type: ignore 1962 valid_images = [[0, 0, 0]] * len(valid_coords) 1963 valid_indices = np.arange(len(valid_coords)) 1964 1965 # Divide the valid 3D space into cubes and compute the cube ids 1966 all_cube_index = _compute_cube_index(valid_coords, global_min, r) # type: ignore 1967 nx, ny, nz = _compute_cube_index(global_max, global_min, r) + 1 1968 all_cube_index = _three_to_one(all_cube_index, ny, nz) 1969 site_cube_index = _three_to_one(_compute_cube_index(center_coords, global_min, r), ny, nz) 1970 # create cube index to coordinates, images, and indices map 1971 cube_to_coords = collections.defaultdict(list) # type: Dict[int, List] 1972 cube_to_images = collections.defaultdict(list) # type: Dict[int, List] 1973 cube_to_indices = collections.defaultdict(list) # type: Dict[int, List] 1974 for i, j, k, l in zip(all_cube_index.ravel(), valid_coords, valid_images, valid_indices): 1975 cube_to_coords[i].append(j) 1976 cube_to_images[i].append(k) 1977 cube_to_indices[i].append(l) 1978 1979 # find all neighboring cubes for each atom in the lattice cell 1980 site_neighbors = find_neighbors(site_cube_index, nx, ny, nz) 1981 neighbors = [] # type: List[List[Tuple[np.ndarray, float, int, np.ndarray]]] 1982 1983 for i, j in zip(center_coords, site_neighbors): 1984 l1 = np.array(_three_to_one(j, ny, nz), dtype=int).ravel() 1985 # use the cube index map to find the all the neighboring 1986 # coords, images, and indices 1987 ks = [k for k in l1 if k in cube_to_coords] 1988 if not ks: 1989 neighbors.append([]) 1990 continue 1991 nn_coords = np.concatenate([cube_to_coords[k] for k in ks], axis=0) 1992 nn_images = itertools.chain(*[cube_to_images[k] for k in ks]) 1993 nn_indices = itertools.chain(*[cube_to_indices[k] for k in ks]) 1994 dist = np.linalg.norm(nn_coords - i[None, :], axis=1) 1995 nns: List[Tuple[np.ndarray, float, int, np.ndarray]] = [] 1996 for coord, index, image, d in zip(nn_coords, nn_indices, nn_images, dist): 1997 # filtering out all sites that are beyond the cutoff 1998 # Here there is no filtering of overlapping sites 1999 if d < r + numerical_tol: 2000 if return_fcoords and (lattice is not None): 2001 coord = np.round(lattice.get_fractional_coords(coord), 10) 2002 nn = (coord, float(d), int(index), image) 2003 nns.append(nn) 2004 neighbors.append(nns) 2005 return neighbors 2006 2007 2008# The following internal methods are used in the get_points_in_sphere method. 2009def _compute_cube_index(coords: np.ndarray, global_min: float, radius: float) -> np.ndarray: 2010 """ 2011 Compute the cube index from coordinates 2012 Args: 2013 coords: (nx3 array) atom coordinates 2014 global_min: (float) lower boundary of coordinates 2015 radius: (float) cutoff radius 2016 2017 Returns: (nx3 array) int indices 2018 2019 """ 2020 return np.array(np.floor((coords - global_min) / radius), dtype=int) 2021 2022 2023def _one_to_three(label1d: np.ndarray, ny: int, nz: int) -> np.ndarray: 2024 """ 2025 Convert a 1D index array to 3D index array 2026 2027 Args: 2028 label1d: (array) 1D index array 2029 ny: (int) number of cells in y direction 2030 nz: (int) number of cells in z direction 2031 2032 Returns: (nx3) int array of index 2033 2034 """ 2035 last = np.mod(label1d, nz) 2036 second = np.mod((label1d - last) / nz, ny) 2037 first = (label1d - last - second * nz) / (ny * nz) 2038 return np.concatenate([first, second, last], axis=1) 2039 2040 2041def _three_to_one(label3d: np.ndarray, ny: int, nz: int) -> np.ndarray: 2042 """ 2043 The reverse of _one_to_three 2044 """ 2045 return np.array(label3d[:, 0] * ny * nz + label3d[:, 1] * nz + label3d[:, 2]).reshape((-1, 1)) 2046 2047 2048def find_neighbors(label: np.ndarray, nx: int, ny: int, nz: int) -> List[np.ndarray]: 2049 """ 2050 Given a cube index, find the neighbor cube indices 2051 2052 Args: 2053 label: (array) (n,) or (n x 3) indice array 2054 nx: (int) number of cells in y direction 2055 ny: (int) number of cells in y direction 2056 nz: (int) number of cells in z direction 2057 2058 Returns: neighbor cell indices 2059 2060 """ 2061 2062 array = [[-1, 0, 1]] * 3 2063 neighbor_vectors = np.array(list(itertools.product(*array)), dtype=int) 2064 if np.shape(label)[1] == 1: 2065 label3d = _one_to_three(label, ny, nz) 2066 else: 2067 label3d = label 2068 all_labels = label3d[:, None, :] - neighbor_vectors[None, :, :] 2069 filtered_labels = [] 2070 # filter out out-of-bound labels i.e., label < 0 2071 for labels in all_labels: 2072 ind = (labels[:, 0] < nx) * (labels[:, 1] < ny) * (labels[:, 2] < nz) * np.all(labels > -1e-5, axis=1) 2073 filtered_labels.append(labels[ind]) 2074 return filtered_labels 2075