def main():

    s1name = "geometry/BCC_to_FCC/POSCAR_BCC"
    s2name = "geometry/BCC_to_FCC/POSCAR_FCC"
    s1 = Poscar.from_file(s1name).structure
    s2 = Poscar.from_file(s2name).structure
    #    print s1,s2
    smatcher = StructureMatcher(primitive_cell=False)
    res = smatcher.get_transformation(s1, s2)
    print res
Exemplo n.º 2
0
def match_structures(s1, s2, scale_lattice=False, rh_only=True):
    """
    Args
        s1: high sym structure
        s2: low sym structure
        scale_lattice (optional): high_sym_superlcell has same lattice vectors as s2 (no strain)
    Returns
        basis: should be the basis that when applied to s1 makes a supercell of the size and orentation of s2
        origin: any additional translation to best match (applied before applying the basis change to match what isotropy does)
        displacements
        high_sym_supercell
"""
    sm = StructureMatcher(ltol=0.3,
                          stol=0.3,
                          angle_tol=15,
                          scale=True,
                          attempt_supercell=True,
                          primitive_cell=False)
    basis, origin, mapping = sm.get_transformation(s1, s2, rh_only=rh_only)

    struct_hs_supercell = sm.get_s2_like_s1(s1, s2, rh_only=rh_only)

    # change origin from the supercell basis to the high sym basis
    origin = np.round_(frac_vec_convert(origin, struct_hs_supercell.lattice,
                                        s2.lattice),
                       decimals=5)
    if scale_lattice:
        hs_lat = struct_hs_supercell.lattice
        hs_std = pmg.Lattice.from_lengths_and_angles(hs_lat.abc, hs_lat.angles)
        ls_lat = s2.lattice
        ls_std = pmg.Lattice.from_lengths_and_angles(ls_lat.abc, ls_lat.angles)
        for aligned, rot, scale in hs_lat.find_all_mappings(hs_std):
            if (abs(aligned.matrix - hs_lat.matrix) < 1.e-3).all():
                strained_hs_lattice = pmg.Lattice(np.inner(ls_std.matrix, rot))
        struct_hs_supercell = pmg.Structure(
            strained_hs_lattice, struct_hs_supercell.species,
            [site.frac_coords for site in struct_hs_supercell])
    displacements = []
    for s_hs, s_ls in zip(struct_hs_supercell, s1):
        disp = np.round_(smallest_disp(s_hs.frac_coords, s_ls.frac_coords),
                         decimals=5)
        displacements.append(disp)
    return basis, origin, displacements, struct_hs_supercell
Exemplo n.º 3
0
def find_common_transformation(asecell,
                               layer_indices,
                               ltol=0.05,
                               stol=0.05,
                               angle_tol=2.0):  # pylint: disable=too-many-arguments,too-many-locals,too-many-branches,too-many-statements
    """
    Given an input structure, in ASE format, and the list with 
    the indices of atoms belonging to each layer, determine
    if there exists a common transformation that brings one
    layer to the next.

    :note: it MUST receive the layers already rearranged so that they are ordered
        along the stacking axis, and with all atoms nearby (i.e., when removing the PBC
        along the third axis, I should still see the layer and not it broken in two or more parts).
        This is performed by the function `find_layers`.

    :param asecell: ASE structure of the bulk, where the first two 
              lattice vectors have been re-oriented to lie 
              in the plane of the layers
    :param layer_indices: list of lists containing the indices of the
                    atoms belonging to each layer
    :param ltol: tolerance on cell length 
    :param stol: tolerance on atomic site positions
    :param angle_tol: tolerance on cell angles
    :return: a tuple of length three: either (rot, transl, None) if there is a common transformation,
        or (None, None, message) if a common transformation could not be found.
        IMPORTANT: the code will try to find a transformation matrix so that rot[2,2] > 0, if it exists.
        If it does not find it, and it returns rot[2,2] < 0, it means it's a category III material (see
        definition in the text). Otherwise, it's either category I or II depending on the monolayer symmetry).
    """
    # instance of the adaptor to convert ASE structures to pymatgen format
    adaptor = AseAtomsAdaptor()

    # Create an instance of Structure Matcher with the given tolerances
    sm = StructureMatcher(ltol, stol, angle_tol, primitive_cell=False)

    # Number of layers
    num_layers = len(layer_indices)
    # If there is only one layer, the transformation is
    # simply a translation along the third axis
    if num_layers == 1:
        return np.eye(3), asecell.cell[2], np.zeros(3), None

    # If we are here, there are at least two layers
    # First check that all layers are identical
    layers = [asecell[layer] for layer in layer_indices]
    if not layers_match(layers, ltol=ltol, stol=stol, angle_tol=angle_tol):
        return None, None, None, "Layers are not identical"

    # Layers are already ordered according to their
    # projection along the stacking direction
    # we start by comparing the first and second layer
    # As we did in function `layers_match`, we put the third vector orthogonal to the plane and
    # extend the third direction slightly
    layer0 = layers[0].copy()
    layer0.cell[2] = [0, 0, 2 * asecell.cell[2, 2]]
    layer1 = layers[1].copy()
    layer1.cell[2] = [0, 0, 2 * asecell.cell[2, 2]]
    str0 = adaptor.get_structure(layer0)
    str1 = adaptor.get_structure(layer1)
    # This is the transformation that brings layer0 onto layer1
    # it is a tuple with a supercell matrix, a translation vector, and the indices
    # of how atoms are rearranged
    transformation01 = sm.get_transformation(str1, str0)

    # There are still cases that, even if layers_match was OK, here we get `None` for transformation.
    # Probably the difference with the code above is in the thresholds and/or in the way the pymatgen
    # code internally implements it. Note: probably we could skip the call to `layers_match` above
    # and just keep the logic below!
    # To avoid crashes, we cope with this here
    if transformation01 is None:
        return None, None, None, "No transformation to match first and second layer"

    # define the common lattice of both layers (which is the one of the bulk)
    # we use twice the same "common lattice" because they are identical in our case
    # in general we should have str0.lattice and str1.lattice (not necessarily in this order)
    common_lattice = str0.lattice
    # find the expression of the rotation in the common lattice
    rot01 = np.linalg.solve(np.dot(transformation01[0], common_lattice.matrix),
                            common_lattice.matrix)
    # translation vector in Cartesian coordinates
    tr01 = np.matmul(common_lattice.matrix.T, transformation01[1])
    # this is the "symmetry" operation that brings layer0 into layer1
    op01 = SymmOp.from_rotation_and_translation(rot01, tr01)

    # While on the x-y plane we are OK if we find a periodic image, on z we want to
    # drop periodic boundary conditions, so finding a periodic image (using the periodicity
    # of the bulk) is not OK.
    # We therefore check if we can identify the same translation along z for all atoms
    # Note that I only take the z coordinate (last column); this is in angstrom
    required_z_shift_per_atom = (
        layers[1].positions -
        op01.operate_multi(layers[0].positions)[transformation01[2]])[:, 2]
    # Note however that if the layers are not properly defined with all atoms closeby,
    # The values of this array might not be all identical.
    # E.g. if atoms of the same layer are not closeby unless one considers PBC.

    z_shift = required_z_shift_per_atom[0]
    for atom_idx in range(1, len(required_z_shift_per_atom)):
        # NOTE: This test is going to work because the function called before this
        # already reorganised layers to have all atoms close to each other
        assert np.abs(z_shift - required_z_shift_per_atom[atom_idx]) < SYMPREC

    # Get the spacegroup of the monolayer (useful below)
    spg_mono = SpacegroupAnalyzer(adaptor.get_structure(layer0),
                                  symprec=SYMPREC)
    # If the transformation involves a flip in the z-direction
    # we check if, by combining it with a symmetry of the monolayer,
    # we get a transformation that does NOT flip z
    # As soon as I find one, I replace tr01, rot01 and op01 (the combination of tr01 and rot01)
    # with those we found, so that the axis is not flipped and transformation01[0][2, 2] > 0
    if transformation01[0][2, 2] < 0:
        for op in spg_mono.get_symmetry_operations(cartesian=True):
            affine_prod = np.dot(op01.affine_matrix, op.affine_matrix)
            if affine_prod[2, 2] > 0:
                tr01 = affine_prod[0:3][:, 3]
                rot01 = affine_prod[0:3][:, 0:3]
                op01 = SymmOp.from_rotation_and_translation(rot01, tr01)
                break

    # Here, there are now two options:
    # - I couldn't find any operation: then I still have transformation01[0][2, 2] < 0
    #   and I will categorize it as category III
    # - I could find it, then transformation01[0][2, 2] > 0 and it's either category I or II
    #   (depending on the monolayer symmetry)
    #   It could also be category III in weird cases like identical layers that are invariant
    #   under z -> -z but with two alternating interlayer distances, i.e. "dimerized"

    # I now use op01 to check if, with op01, I can bring each layer onto the next,
    # and layer num_layers onto num_layers+1, i.e. the first one + the third lattice vector
    # If op01 does not work we might need to combine it with symmetry operations of the monolayer
    # before concluding that the system is not a MDO polytype
    # print(op01)

    for symop in sorted(
            spg_mono.get_symmetry_operations(cartesian=True),
            key=lambda op: -op.affine_matrix[2, 2],
    ):
        # We sort the symmetry operations of the monolayer so that the ones that
        # do NOT flip z are considered first.
        # Indeed in principles operations that flip z should be possible
        # only in category I, but would result in a coincidence operation
        # that flip z, which is not necessary in this case (category I), as in this case
        # we'll find another one that does the same job without flipping, so we could even skip these operations.
        # (also because we want to give the guarantee that, if symop.affine_matrix[2, 2] < 0, then we are
        # in category III)
        # Still, operations that flip z could be relevant in some weird cases of category III, called "dimerized" above
        # So that in the end we will find a coincidence operation that flips z even if op01 did not
        #
        # Combine the coincidence operation with the symmetry operation of the monolayer
        affine_prod = np.dot(op01.affine_matrix, symop.affine_matrix)
        this_rot = affine_prod[0:3][:, 0:3]
        this_tr = affine_prod[0:3][:, 3]
        # The coincidence operation, elevated to the number of layers,
        # should be a symmetry operation of the monolayer, and we need the
        # corresponding fractional translation
        frac_tr = get_fractional_translation(this_rot, num_layers, spg_mono)
        if frac_tr is None:
            found_common = False
            continue

        # Check if the coincidence operation satisfies a complex condition, implemented
        # in 'check_transformation'. If not, there is no need to continue
        # with this transformation
        if not check_transformation(this_rot, this_tr, frac_tr, asecell,
                                    num_layers):
            found_common = False
            continue

        # If the condition is instead satisfied, we can continue and decompose the translation
        # vector this_tr, into a component which is just associated with the fact that
        # the coincidence operation has to be performed around a point this_tr0 different from the origin
        # and true traslation, which is either purely vertical or has a component invariant
        # under the coincidence transormation (e.g. a component parallel to a mirror plane)
        this_tr0 = get_transformation_center(this_rot, this_tr, asecell,
                                             num_layers)

        # Once we have identified possible component associated with a non-trivial center of the
        # coincidence operation, this_tr0, we can subtract it from this_tr to obtain the true translation
        this_tr -= np.dot(np.identity(3) - this_rot, this_tr0)

        # The coincidence operation is then defined in terms of this_tr only
        this_op = SymmOp.from_rotation_and_translation(this_rot, this_tr)

        # Bring back the inplane components of the vector this_tr0 identifying the center of
        # the coincidence operation inside the unit cell
        this_tr0[:2] = np.dot(
            np.dot(this_tr0, np.linalg.inv(asecell.cell))[:2] % 1,
            asecell.cell[:2])[:2]

        # It is also useful to define a vector along the stacking direction
        vec = np.array([0.0, 0.0, asecell.cell[2, 2] / num_layers])

        found_common = True
        # NOTE: here we start from layer ZERO! The reason is that we want to
        # double check that now that we made a combination with a symmetry operation,
        # we still send layer 1 into layer 2.
        # So we want to check that case as well.
        for il in range(0, num_layers):
            # We need to copy the layers as we'll change them in place
            layer0 = layers[il].copy()
            layer1 = layers[(il + 1) % num_layers].copy()
            # translate back the two layers to put at the origin the
            # center around which we want to perform the coincidence operation
            # if layer1 is the layer num_layer + 1 we need to translate it
            # by a full lattice vector
            layer0.translate(-il * vec - this_tr0)
            layer1.translate((-il * vec - this_tr0 + np.floor(
                (il + 1.0) / num_layers) * asecell.cell[2]))
            # the transformed positions of the atoms in the first layer
            pos0 = this_op.operate_multi(layer0.positions)
            # that should be identical to the ones of the second layer
            pos1 = layer1.positions
            # we already know from above that the species in each layer are identical
            # so we take the set of the ones in the first layer
            for an in set(layer0.get_chemical_symbols()):
                # check, species by species, that the atomic positions are identical
                # up to a threshold
                # The second variable would be the mapping, but I'm not using it
                distance, _ = are_same_except_order(
                    np.array([
                        _[0] for _ in zip(pos0, layer0.get_chemical_symbols())
                        if _[1] == an
                    ]),
                    np.array([
                        _[0] for _ in zip(pos1, layer1.get_chemical_symbols())
                        if _[1] == an
                    ]),
                    common_lattice.matrix,
                )
                # if the distance for any of the atoms of the given species
                # is larger than the threshold the transformation is not the same
                # between all consecutive layers
                found_common = found_common and distance.max() < stol
                # if this transformation does not work it is useless con continue with
                # other species of this layer
                if not found_common:
                    break
            # if this transformation does not work it is useless con continue with
            # subsequent layers
            if not found_common:
                break
        # if the transformation works, no need to test additional transformations
        if found_common:
            break
    if not found_common:
        return (
            None,
            None,
            None,
            "The transformation between consecutive layers is not always the same",
        )
    return this_rot, this_tr, this_tr0, None