1# coding: utf-8
2"""
3This module defines basic objects representing the crystalline structure.
4"""
5import sys
6import os
7import collections
8import tempfile
9import numpy as np
10import pickle
11import pymatgen.core.units as pmg_units
12
13from pprint import pformat
14from collections import OrderedDict
15from monty.collections import AttrDict, dict2namedtuple
16from monty.functools import lazy_property
17from monty.string import is_string, marquee, list_strings
18from monty.termcolor import cprint
19from pymatgen.core.structure import Structure as pmg_Structure
20from pymatgen.core.sites import PeriodicSite
21from pymatgen.core.lattice import Lattice
22from pymatgen.symmetry.analyzer import SpacegroupAnalyzer
23from abipy.tools.plotting import add_fig_kwargs, get_ax_fig_plt, get_axarray_fig_plt
24from abipy.flowtk import PseudoTable
25from abipy.core.mixins import NotebookWriter
26from abipy.core.symmetries import AbinitSpaceGroup
27from abipy.iotools import as_etsfreader, Visualizer
28from abipy.flowtk.abiobjects import structure_from_abivars, structure_to_abivars, species_by_znucl
29
30
31__all__ = [
32    "mp_match_structure",
33    "mp_search",
34    "cod_search",
35    "Structure",
36    "dataframes_from_structures",
37]
38
39
40def mp_match_structure(obj, api_key=None, endpoint=None, final=True):
41    """
42    Finds matching structures on the Materials Project database.
43
44    Args:
45        obj: filename or |Structure| object.
46        api_key (str): A String API key for accessing the MaterialsProject REST interface.
47        endpoint (str): Url of endpoint to access the MaterialsProject REST interface.
48        final (bool): Whether to get the final structure, or the initial
49            (pre-relaxation) structure. Defaults to True.
50
51    Returns:
52        :class:`MpStructures` object with
53            structures: List of matching structures and list of Materials Project identifier.
54    """
55    structure = Structure.as_structure(obj)
56    # Must use pymatgen structure else server does not know how to handle the JSON doc.
57    structure.__class__ = pmg_Structure
58
59    from abipy.core import restapi
60    structures = []
61    with restapi.get_mprester(api_key=api_key, endpoint=endpoint) as rest:
62        try:
63            mpids = rest.find_structure(structure)
64            if mpids:
65                structures = [Structure.from_mpid(mid, final=final, api_key=api_key, endpoint=endpoint)
66                              for mid in mpids]
67
68        except rest.Error as exc:
69            cprint(str(exc), "red")
70
71        finally:
72            # Back to abipy structure
73            structure = Structure.as_structure(structure)
74            structures.insert(0, structure)
75            mpids.insert(0, "this")
76            return restapi.MpStructures(structures=structures, ids=mpids)
77
78
79def mp_search(chemsys_formula_id, api_key=None, endpoint=None):
80    """
81    Connect to the materials project database.
82    Get a list of structures corresponding to a chemical system, formula, or materials_id.
83
84    Args:
85        chemsys_formula_id (str): A chemical system (e.g., Li-Fe-O),
86            or formula (e.g., Fe2O3) or materials_id (e.g., mp-1234).
87        api_key (str): A String API key for accessing the MaterialsProject REST interface.
88            If this is None, the code will check if there is a `PMG_MAPI_KEY` in your .pmgrc.yaml.
89        endpoint (str): Url of endpoint to access the MaterialsProject REST interface.
90
91    Returns:
92        :class:`MpStructures` object with
93            List of Structure objects, Materials project ids associated to structures.
94            and List of dictionaries with MP data (same order as structures).
95
96        Note that the attributes evalute to False if no match is found
97    """
98    chemsys_formula_id = chemsys_formula_id.replace(" ", "")
99
100    structures, mpids, data = [], [], None
101    from abipy.core import restapi
102    with restapi.get_mprester(api_key=api_key, endpoint=endpoint) as rest:
103        try:
104            data = rest.get_data(chemsys_formula_id, prop="")
105            if data:
106                structures = [Structure.from_str(d["cif"], fmt="cif", primitive=False, sort=False)
107                              for d in data]
108                mpids = [d["material_id"] for d in data]
109                # Want AbiPy structure.
110                structures = list(map(Structure.as_structure, structures))
111
112        except rest.Error as exc:
113            cprint(str(exc), "magenta")
114
115        return restapi.MpStructures(structures, mpids, data=data)
116
117
118def cod_search(formula, primitive=False):
119    """
120    Connect to the COD_ database. Get list of structures corresponding to a chemical formula
121
122    Args:
123        formula (str): Chemical formula (e.g., Fe2O3)
124        primitive (bool): True if primitive structures are wanted. Note that many COD structures are not primitive.
125
126    Returns:
127        :class:`CodStructures` object with
128            List of Structure objects, COD ids associated to structures.
129            and List of dictionaries with COD data (same order as structures).
130
131        Note that the attributes evalute to False if no match is found
132    """
133    from pymatgen.ext.cod import COD
134    data = COD().get_structure_by_formula(formula)
135
136    cod_ids = [e.pop("cod_id") for e in data]
137    # Want AbiPy structure.
138    structures = list(map(Structure.as_structure, [e.pop("structure") for e in data]))
139    if primitive:
140        structures = [s.get_primitive_structure() for s in structures]
141
142    from abipy.core import restapi
143    return restapi.CodStructures(structures, cod_ids, data=data)
144
145
146class Structure(pmg_Structure, NotebookWriter):
147    """
148    Extends :class:`pymatgen.core.structure.Structure` with Abinit-specific methods.
149
150    A jupyter_ notebook documenting the usage of this object is available at
151    <https://nbviewer.jupyter.org/github/abinit/abitutorials/blob/master/abitutorials/structure.ipynb>
152
153    For the pymatgen project see :cite:`Ong2013`.
154
155    .. rubric:: Inheritance Diagram
156    .. inheritance-diagram:: Structure
157    """
158    @classmethod
159    def as_structure(cls, obj):
160        """
161        Convert obj into a |Structure|. Accepts:
162
163            - Structure object.
164            - Filename
165            - Dictionaries (JSON_ format or dictionaries with abinit variables).
166            - Objects with a ``structure`` attribute.
167        """
168        if isinstance(obj, cls): return obj
169        if isinstance(obj, pmg_Structure):
170            obj.__class__ = cls
171            return obj
172
173        if is_string(obj):
174            return cls.from_file(obj)
175
176        if isinstance(obj, collections.abc.Mapping):
177            if "@module" in obj:
178                return Structure.from_dict(obj)
179            else:
180                return Structure.from_abivars(obj)
181
182        if hasattr(obj, "structure"):
183            return cls.as_structure(obj.structure)
184        elif hasattr(obj, "final_structure"):
185            # This for HIST.nc file
186            return cls.as_structure(obj.final_structure)
187
188        raise TypeError("Don't know how to convert %s into a structure" % type(obj))
189
190    @classmethod
191    def from_file(cls, filepath, primitive=False, sort=False):
192        """
193        Reads a structure from a file. For example, anything ending in
194        a "cif" is assumed to be a Crystallographic Information Format file.
195        Supported formats include CIF_, POSCAR/CONTCAR, CHGCAR, LOCPOT,
196        vasprun.xml, CSSR, Netcdf and pymatgen's JSON serialized structures.
197
198        Netcdf files supported:
199            All files produced by ABINIT with info of the crystalline geometry
200            HIST.nc, in this case the last structure of the history is returned.
201
202        Args:
203            filename (str): The filename to read from.
204            primitive (bool): Whether to convert to a primitive cell
205                Only available for cifs, POSCAR, CSSR, JSON, YAML
206                Defaults to True.
207            sort (bool): Whether to sort sites. Default to False.
208
209        Returns: |Structure| object
210        """
211        #zipped_exts = (".bz2", ".gz", ".z"):
212        root, ext = os.path.splitext(filepath)
213
214        if filepath.endswith("_HIST.nc"):
215            # Abinit history file. In this case we return the last structure!
216            # Note that HIST does not follow the etsf-io conventions.
217            from abipy.dynamics.hist import HistFile
218            with HistFile(filepath) as hist:
219                return hist.structures[-1]
220
221        elif filepath.endswith(".nc"):
222            # Generic netcdf file.
223            ncfile, closeit = as_etsfreader(filepath)
224
225            new = ncfile.read_structure(cls=cls)
226            new.set_abi_spacegroup(AbinitSpaceGroup.from_ncreader(ncfile))
227
228            # Try to read indsym table from file (added in 8.9.x)
229            indsym = ncfile.read_value("indsym", default=None)
230            if indsym is not None:
231                # Fortran --> C convention
232                indsym[:, :, 3] -= 1
233                new.indsym = indsym
234
235            if closeit: ncfile.close()
236
237        elif filepath.endswith(".abi") or filepath.endswith(".in"):
238            # Abinit input file.
239            # Here I assume that the input file contains a single structure.
240            from abipy.abio.abivars import AbinitInputFile
241            return AbinitInputFile.from_file(filepath).structure
242
243        elif filepath.endswith(".abo") or filepath.endswith(".out"):
244            # Abinit output file. We can have multi-datasets and multiple initial/final structures!
245            # By desing, we return the last structure if out is completed else the initial one.
246            # None is returned if the structures are different.
247            from abipy.abio.outputs import AbinitOutputFile
248            with AbinitOutputFile(filepath) as out:
249                #print("initial_structures:\n", out.initial_structures, "\nfinal_structures:\n", out.final_structures)
250                if out.final_structures: return out.final_structure
251                if out.initial_structures: return out.initial_structure
252            raise ValueError("Cannot find structure in Abinit output file `%s`" % filepath)
253
254        elif filepath.endswith("_DDB") or root.endswith("_DDB"):
255            # DDB file.
256            from abipy.abilab import abiopen
257            with abiopen(filepath) as abifile:
258                return abifile.structure
259
260        elif filepath.endswith(".pickle"):
261            # From pickle.
262            with open(filepath, "rb") as fh:
263                new = pickle.load(fh)
264                if not isinstance(new, pmg_Structure):
265                    # Is it a object with a structure property?
266                    if hasattr(new, "structure"): new = new.structure
267
268                if not isinstance(new, pmg_Structure):
269                    raise TypeError("Don't know how to extract a Structure from file %s, received type %s" %
270                        (filepath, type(new)))
271
272                if new.__class__ != cls: new.__class__ = cls
273
274        else:
275            # Invoke pymatgen and change class
276            # Note that AbinitSpacegroup is missing here.
277            new = super().from_file(filepath, primitive=primitive, sort=sort)
278            if new.__class__ != cls: new.__class__ = cls
279
280        return new
281
282    @classmethod
283    def from_mpid(cls, material_id, final=True, api_key=None, endpoint=None):
284        """
285        Get a Structure corresponding to a material_id.
286
287        Args:
288            material_id (str): Materials Project material_id (a string, e.g., mp-1234).
289            final (bool): Whether to get the final structure, or the initial
290                (pre-relaxation) structure. Defaults to True.
291            api_key (str): A String API key for accessing the MaterialsProject
292                REST interface. Please apply on the Materials Project website for one.
293                If this is None, the code will check if there is a ``PMG_MAPI_KEY`` in your .pmgrc.yaml.
294                If so, it will use that environment
295                This makes easier for heavy users to simply add this environment variable
296                to their setups and MPRester can then be called without any arguments.
297            endpoint (str): Url of endpoint to access the MaterialsProject REST interface.
298                Defaults to the standard Materials Project REST address, but
299                can be changed to other urls implementing a similar interface.
300
301        Returns: |Structure| object.
302        """
303        # Get pytmatgen structure and convert it to abipy structure
304        from abipy.core import restapi
305        with restapi.get_mprester(api_key=api_key, endpoint=endpoint) as rest:
306            new = rest.get_structure_by_material_id(material_id, final=final)
307            return cls.as_structure(new)
308
309    @classmethod
310    def from_cod_id(cls, cod_id, primitive=False, **kwargs):
311        """
312        Queries the COD_ for a structure by id. Returns |Structure| object.
313
314        Args:
315            cod_id (int): COD id.
316            primitive (bool): True if primitive structures are wanted. Note that many COD structures are not primitive.
317            kwargs: Arguments passed to ``get_structure_by_id``
318
319        Returns: |Structure| object.
320        """
321        from pymatgen.ext.cod import COD
322        new = COD().get_structure_by_id(cod_id, **kwargs)
323        if primitive: new = new.get_primitive_structure()
324        return cls.as_structure(new)
325
326    @classmethod
327    def from_ase_atoms(cls, atoms):
328        """
329        Returns structure from ASE Atoms.
330
331        Args:
332            atoms: ASE Atoms object
333
334        Returns:
335            Equivalent Structure
336        """
337        import pymatgen.io.ase as aio
338        return aio.AseAtomsAdaptor.get_structure(atoms, cls=cls)
339
340    def to_ase_atoms(self):
341        """
342        Returns ASE_ Atoms object from structure.
343        """
344        import pymatgen.io.ase as aio
345        return aio.AseAtomsAdaptor.get_atoms(self)
346
347    @classmethod
348    def boxed_molecule(cls, pseudos, cart_coords, acell=3*(10,)):
349        """
350        Creates a molecule in a periodic box of lengths acell [Bohr]
351
352        Args:
353            pseudos: List of pseudopotentials
354            cart_coords: Cartesian coordinates
355            acell: Lengths of the box in *Bohr*
356        """
357        from pymatgen.core.structure import Molecule
358        cart_coords = np.atleast_2d(cart_coords)
359        molecule = Molecule([p.symbol for p in pseudos], cart_coords)
360        l = pmg_units.ArrayWithUnit(acell, "bohr").to("ang")
361
362        new = molecule.get_boxed_structure(l[0], l[1], l[2])
363        return cls.as_structure(new)
364
365    @classmethod
366    def boxed_atom(cls, pseudo, cart_coords=3*(0,), acell=3*(10,)):
367        """
368        Creates an atom in a periodic box of lengths acell [Bohr]
369
370        Args:
371            pseudo: Pseudopotential object.
372            cart_coords: Cartesian coordinates in Angstrom
373            acell: Lengths of the box in *Bohr* (Abinit input variable)
374        """
375        return cls.boxed_molecule([pseudo], cart_coords, acell=acell)
376
377    @classmethod
378    def bcc(cls, a, species, primitive=True, units="ang", **kwargs):
379        """
380        Build a primitive or a conventional bcc crystal structure.
381
382        Args:
383            a: Lattice parameter (Angstrom if units is not given)
384            species: Chemical species. See __init__ method of |pymatgen-Structure|
385            primitive: if True a primitive cell will be produced, otherwise a conventional one
386            units: Units of input lattice parameters e.g. "bohr", "pm"
387            kwargs: All keyword arguments accepted by |pymatgen-Structure|.
388        """
389        a = pmg_units.Length(a, units).to("ang")
390        if primitive:
391            lattice = 0.5 * float(a) * np.array([
392                -1,  1,  1,
393                 1, -1,  1,
394                 1,  1, -1])
395
396            coords = [[0, 0, 0]]
397
398        else:
399            lattice = float(a) * np.eye(3)
400            coords = [[0, 0, 0],
401                      [0.5, 0.5, 0.5]]
402            species = np.repeat(species, 2)
403
404        return cls(lattice, species, coords=coords,  **kwargs)
405
406    @classmethod
407    def fcc(cls, a, species, primitive=True, units="ang", **kwargs):
408        """
409        Build a primitive or a conventional fcc crystal structure.
410
411        Args:
412            a: Lattice parameter (Angstrom if units is not given)
413            species: Chemical species. See __init__ method of :class:`pymatgen.Structure`
414            primitive: if True a primitive cell will be produced, otherwise a conventional one
415            units: Units of input lattice parameters e.g. "bohr", "pm"
416            kwargs: All keyword arguments accepted by :class:`pymatgen.Structure`
417        """
418        a = pmg_units.Length(a, units).to("ang")
419        if primitive:
420            lattice = 0.5 * float(a) * np.array([
421                0,  1,  1,
422                1,  0,  1,
423                1,  1,  0])
424            coords = [[0, 0, 0]]
425        else:
426            lattice = float(a) * np.eye(3)
427            species = np.repeat(species, 4)
428            coords = [[0, 0, 0],
429                      [0.5, 0.5, 0],
430                      [0.5, 0, 0.5],
431                      [0, 0.5, 0.5]]
432
433        return cls(lattice, species, coords=coords, **kwargs)
434
435    @classmethod
436    def zincblende(cls, a, species, units="ang", **kwargs):
437        """
438        Build a primitive zincblende crystal structure.
439
440        Args:
441            a: Lattice parameter (Angstrom if units is not given)
442            species: Chemical species. See __init__ method of :class:`pymatgen.Structure`
443            units: Units of input lattice parameters e.g. "bohr", "pm"
444            kwargs: All keyword arguments accepted by :class:`pymatgen.Structure`
445
446        Example::
447
448            Structure.zincblende(a, ["Zn", "S"])
449
450        """
451        a = pmg_units.Length(a, units).to("ang")
452        lattice = 0.5 * float(a) * np.array([
453            0,  1,  1,
454            1,  0,  1,
455            1,  1,  0])
456
457        frac_coords = np.reshape([0, 0, 0, 0.25, 0.25, 0.25], (2, 3))
458        return cls(lattice, species, frac_coords, coords_are_cartesian=False, **kwargs)
459
460    @classmethod
461    def rocksalt(cls, a, species, units="ang", **kwargs):
462        """
463        Build a primitive fcc crystal structure.
464
465        Args:
466            a: Lattice parameter (Angstrom if units is not given)
467            units: Units of input lattice parameters e.g. "bohr", "pm"
468            species: Chemical species. See __init__ method of :class:`pymatgen.Structure`
469            kwargs: All keyword arguments accepted by :class:`pymatgen.Structure`
470
471        Example::
472
473            Structure.rocksalt(a, ["Na", "Cl"])
474
475        """
476        a = pmg_units.Length(a, units).to("ang")
477        lattice = 0.5 * float(a) * np.array([
478            0,  1,  1,
479            1,  0,  1,
480            1,  1,  0])
481
482        frac_coords = np.reshape([0, 0, 0, 0.5, 0.5, 0.5], (2, 3))
483        return cls(lattice, species, frac_coords, coords_are_cartesian=False, **kwargs)
484
485    @classmethod
486    def ABO3(cls, a, species, units="ang", **kwargs):
487        """
488        Peroviskite structures.
489
490        Args:
491            a: Lattice parameter (Angstrom if units is not given)
492            species: Chemical species. See __init__ method of :class:`pymatgen.Structure`
493            units: Units of input lattice parameters e.g. "bohr", "pm"
494            kwargs: All keyword arguments accepted by :class:`pymatgen.Structure`
495        """
496        a = pmg_units.Length(a, units).to("ang")
497        lattice = float(a) * np.eye(3)
498        frac_coords = np.reshape([
499            0,     0,   0,  # A (2a)
500            0.5, 0.5, 0.5,  # B (2a)
501            0.5, 0.5, 0.0,  # O (6b)
502            0.5, 0.0, 0.5,  # O (6b)
503            0.0, 0.5, 0.5,  # O (6b)
504        ], (5, 3))
505
506        return cls(lattice, species, frac_coords, coords_are_cartesian=False, **kwargs)
507
508    @classmethod
509    def from_abistring(cls, string):
510        """Initialize Structure from string with Abinit input variables."""
511        from abipy.abio.abivars import AbinitInputFile
512        return AbinitInputFile.from_string(string).structure
513
514    @classmethod
515    def from_abivars(cls, *args, **kwargs):
516        """
517        Build a |Structure| object from a dictionary with ABINIT variables.
518
519        Example::
520
521            al_structure = Structure.from_abivars(
522                acell=3*[7.5],
523                rprim=[0.0, 0.5, 0.5,
524                       0.5, 0.0, 0.5,
525                       0.5, 0.5, 0.0],
526                typat=1,
527                xred=[0.0, 0.0, 0.0],
528                ntypat=1,
529                znucl=13,
530            )
531
532        ``xred`` can be replaced with ``xcart`` or ``xangst``.
533        """
534        return structure_from_abivars(cls, *args, **kwargs)
535
536    @property
537    def species_by_znucl(self):
538        """
539        Return list of unique specie found in the structure **ordered according to sites**.
540
541        Example:
542
543            Site0: 0.5 0 0 O
544            Site1: 0   0 0 Si
545
546        produces [Specie_O, Specie_Si] and not set([Specie_O, Specie_Si]) as in `types_of_specie`.
547
548        Important:: We call this method `species_by_znucl` but this does not mean that the list can automagically
549        reproduce the value of `znucl(ntypat)` specified in an **arbitrary** ABINIT input file created by the user.
550        This array is ordered as the znucl list produced by AbiPy when writing the structure to the input file.
551        """
552        return species_by_znucl(self)
553
554    def __str__(self):
555        return self.to_string()
556
557    def to_string(self, title=None, verbose=0):
558        """String representation."""
559        lines = []; app = lines.append
560        if title is not None: app(marquee(title, mark="="))
561        if verbose:
562            app(self.spget_summary(verbose=verbose))
563        else:
564            app(super().__str__())
565
566        if self.abi_spacegroup is not None:
567            app("\nAbinit Spacegroup: %s" % self.abi_spacegroup.to_string(verbose=verbose))
568
569        return "\n".join(lines)
570
571    def to(self, fmt=None, filename=None, **kwargs):
572        __doc__ = pmg_Structure.to.__doc__ + \
573            "\n Accepts also fmt=`abinit` or `abivars` or `.abi` as Abinit input file extension"
574
575        filename = filename or ""
576        fmt = "" if fmt is None else fmt.lower()
577        fname = os.path.basename(filename)
578
579        if fmt in ("abi", "abivars", "abinit") or fname.endswith(".abi"):
580            if filename:
581                with open(filename, "wt") as f:
582                    f.write(self.abi_string)
583            else:
584                return self.abi_string
585        else:
586            return super().to(fmt=fmt, filename=filename, **kwargs)
587
588    def write_cif_with_spglib_symms(self, filename, symprec=1e-3, angle_tolerance=5.0, significant_figures=8,
589                                    ret_string=False):
590        """
591        Args:
592            symprec (float): If not none, finds the symmetry of the structure
593                and writes the cif with symmetry information. Passes symprec
594                to the SpacegroupAnalyzer.
595            significant_figures (int): Specifies precision for formatting of floats.
596                Defaults to 8.
597            angle_tolerance (float): Angle tolerance for symmetry finding. Passes
598                angle_tolerance to the SpacegroupAnalyzer. Used only if symprec
599                is not None.
600        """
601        from pymatgen.io.cif import CifWriter
602        cif_str = str(CifWriter(self,
603                      symprec=symprec, significant_figures=significant_figures, angle_tolerance=angle_tolerance,
604                      refine_struct=False))
605
606        if not ret_string:
607            with open(filename, "wt") as fh:
608                fh.write(cif_str)
609        else:
610            return cif_str
611
612    def __mul__(self, scaling_matrix):
613        """
614        Makes a supercell. Allowing to have sites outside the unit cell
615        See pymatgen for docs.
616
617        Wraps __mul__ operator of pymatgen structure to return abipy structure
618        """
619        new = super().__mul__(scaling_matrix)
620        return self.__class__.as_structure(new)
621
622    __rmul__ = __mul__
623
624    def to_abivars(self, enforce_znucl=None, enforce_typat=None, **kwargs):
625        """
626        Returns a dictionary with the ABINIT variables.
627
628        enforce_znucl[ntypat] =
629        enforce_typat[natom] = Fortran conventions. Start to count from 1.
630        """
631        return structure_to_abivars(self, enforce_znucl=enforce_znucl, enforce_typat=enforce_typat, **kwargs)
632
633    @property
634    def latex_formula(self):
635        """LaTeX formatted formula. E.g., Fe2O3 is transformed to Fe$_{2}$O$_{3}$."""
636        from pymatgen.util.string import latexify
637        return latexify(self.formula)
638
639    @property
640    def abi_string(self):
641        """Return a string with the ABINIT input associated to this structure."""
642        from abipy.abio.variable import InputVariable
643        lines = []
644        app = lines.append
645        abivars = self.to_abivars()
646        for varname, value in abivars.items():
647            app(str(InputVariable(varname, value)))
648
649        return("\n".join(lines))
650
651    def get_panel(self):
652        """Build panel with widgets to interact with the structure either in a notebook or in a bokeh app"""
653        from abipy.panels.structure import StructurePanel
654        return StructurePanel(self).get_panel()
655
656    def get_conventional_standard_structure(self, international_monoclinic=True,
657                                           symprec=1e-3, angle_tolerance=5):
658        """
659        Gives a structure with a conventional cell according to certain
660        standards. The standards are defined in :cite:`Setyawan2010`
661        They basically enforce as much as possible norm(a1) < norm(a2) < norm(a3)
662
663        Returns: The structure in a conventional standardized cell
664        """
665        spga = SpacegroupAnalyzer(self, symprec=symprec, angle_tolerance=angle_tolerance)
666        new = spga.get_conventional_standard_structure(international_monoclinic=international_monoclinic)
667        return self.__class__.as_structure(new)
668
669    def abi_primitive(self, symprec=1e-3, angle_tolerance=5, no_idealize=0):
670        #TODO: this should be moved to pymatgen in the get_refined_structure or so ...
671        # to be considered in February 2016
672        import spglib
673        from pymatgen.io.ase import AseAtomsAdaptor
674        try:
675            from ase.atoms import Atoms
676        except ImportError:
677            raise ImportError('Could not import Atoms from ase. Install it with `conda install ase` or pip')
678
679        s = self.get_sorted_structure()
680        ase_adaptor = AseAtomsAdaptor()
681        ase_atoms = ase_adaptor.get_atoms(structure=s)
682        standardized = spglib.standardize_cell(ase_atoms, to_primitive=1, no_idealize=no_idealize,
683                                               symprec=symprec, angle_tolerance=angle_tolerance)
684        standardized_ase_atoms = Atoms(scaled_positions=standardized[1], numbers=standardized[2], cell=standardized[0])
685        standardized_structure = ase_adaptor.get_structure(standardized_ase_atoms)
686
687        return self.__class__.as_structure(standardized_structure)
688
689    def refine(self, symprec=1e-3, angle_tolerance=5):
690        """
691        Get the refined structure based on detected symmetry. The refined
692        structure is a *conventional* cell setting with atoms moved to the
693        expected symmetry positions.
694
695        Returns: Refined structure.
696        """
697        sym_finder = SpacegroupAnalyzer(structure=self, symprec=symprec, angle_tolerance=angle_tolerance)
698        new = sym_finder.get_refined_structure()
699        return self.__class__.as_structure(new)
700
701    def abi_sanitize(self, symprec=1e-3, angle_tolerance=5, primitive=True, primitive_standard=False):
702        """
703        Returns a new structure in which:
704
705            * Structure is refined.
706            * Reduced to primitive settings.
707            * Lattice vectors are exchanged if the triple product is negative
708
709        Args:
710            symprec (float): Symmetry precision used to refine the structure.
711            angle_tolerance (float): Tolerance on angles.
712                if ``symprec`` is None and `angle_tolerance` is None, no structure refinement is peformed.
713            primitive (bool): Whether to convert to a primitive cell following :cite:`Setyawan2010`
714            primitive_standard (bool): Returns most primitive structure found.
715        """
716        from pymatgen.transformations.standard_transformations import PrimitiveCellTransformation, SupercellTransformation
717        structure = self.__class__.from_sites(self)
718
719        # Refine structure
720        if symprec is not None and angle_tolerance is not None:
721            structure = structure.refine(symprec=symprec, angle_tolerance=angle_tolerance)
722
723        # Convert to primitive structure.
724        #primitive_standard = False
725        if primitive:
726            if primitive_standard:
727                # Setyawan, W., & Curtarolo, S.
728                sym_finder_prim = SpacegroupAnalyzer(structure=structure, symprec=symprec, angle_tolerance=angle_tolerance)
729                structure = sym_finder_prim.get_primitive_standard_structure(international_monoclinic=False)
730            else:
731                # Find most primitive structure.
732                get_prim = PrimitiveCellTransformation()
733                structure = get_prim.apply_transformation(structure)
734
735        # Exchange last two lattice vectors if triple product is negative.
736        m = structure.lattice.matrix
737        x_prod = np.dot(np.cross(m[0], m[1]), m[2])
738        if x_prod < 0:
739            #print("Negative triple product --> exchanging last two lattice vectors.")
740            trans = SupercellTransformation(((1, 0, 0), (0, 0, 1), (0, 1, 0)))
741            structure = trans.apply_transformation(structure)
742            m = structure.lattice.matrix
743            x_prod = np.dot(np.cross(m[0], m[1]), m[2])
744            if x_prod < 0: raise RuntimeError("x_prod is still negative!")
745
746        return self.__class__.as_structure(structure)
747
748    def get_oxi_state_decorated(self, **kwargs):
749        """
750        Use :class:`pymatgen.analysis.bond_valence.BVAnalyzer` to estimate oxidation states
751        Return oxidation state decorated structure.
752        This currently works only for ordered structures only.
753
754        Args:
755            kwargs: Arguments passed to BVAnalyzer
756
757        Returns:
758            A modified structure that is oxidation state decorated.
759        """
760        from pymatgen.analysis.bond_valence import BVAnalyzer
761        new = BVAnalyzer(**kwargs).get_oxi_state_decorated_structure(self)
762        return self.__class__.as_structure(new)
763
764    @property
765    def reciprocal_lattice(self):
766        """
767        Reciprocal lattice of the structure. Note that this is the standard
768        reciprocal lattice used for solid state physics with a factor of 2 * pi
769        i.e.  a_j . b_j = 2pi delta_ij
770
771        If you are looking for the crystallographic reciprocal lattice,
772        use the reciprocal_lattice_crystallographic property.
773        """
774        return self._lattice.reciprocal_lattice
775
776    def lattice_vectors(self, space="r"):
777        """
778        Returns the vectors of the unit cell in Angstrom.
779
780        Args:
781            space: "r" for real space vectors, "g" for reciprocal space basis vectors.
782        """
783        if space.lower() == "r":
784            return self.lattice.matrix
785        if space.lower() == "g":
786            return self.lattice.reciprocal_lattice.matrix
787        raise ValueError("Wrong value for space: %s " % str(space))
788
789    def spget_lattice_type(self, symprec=1e-3, angle_tolerance=5):
790        """
791        Call spglib to get the lattice for the structure, e.g., (triclinic,
792        orthorhombic, cubic, etc.).This is the same than the
793        crystal system with the exception of the hexagonal/rhombohedral lattice
794
795        Args:
796            symprec (float): Symmetry precision for distance
797            angle_tolerance (float): Tolerance on angles.
798
799        Returns:
800            (str): Lattice type for structure or None if type cannot be detected.
801        """
802        spgan = SpacegroupAnalyzer(self, symprec=symprec, angle_tolerance=angle_tolerance)
803        return spgan.get_lattice_type()
804
805    def spget_equivalent_atoms(self, symprec=1e-3, angle_tolerance=5, printout=False):
806        """
807        Call spglib_ to find the inequivalent atoms and build symmetry tables.
808
809        Args:
810            symprec (float): Symmetry precision for distance.
811            angle_tolerance (float): Tolerance on angles.
812            printout (bool): True to print symmetry tables.
813
814        Returns:
815            ``namedtuple`` (irred_pos, eqmap, spgdata) with the following attributes::
816
817                * irred_pos: array giving the position of the i-th irred atom in the structure.
818                    The number of irred atoms is len(irred_pos).
819                *   eqmap: Mapping irred atom position --> list with positions of symmetrical atoms.
820                *   wyckoffs: Wyckoff letters.
821                *   wyck_mult: Array with Wyckoff multiplicity.
822                *   wyck_labels: List of labels with Wyckoff multiplicity and letter e.g. 3a
823                *   site_labels: Labels for each site in computed from element symbol and wyckoff positions e.g Si2a
824                *   spgdata: spglib dataset with additional data reported by spglib_.
825
826         :Example:
827
828            for irr_pos in irred_pos:
829                eqmap[irr_pos]   # List of symmetrical positions associated to the irr_pos atom.
830        """
831        natom = len(self)
832        spgan = SpacegroupAnalyzer(self, symprec=symprec, angle_tolerance=angle_tolerance)
833        spgdata = spgan.get_symmetry_dataset()
834        equivalent_atoms = spgdata["equivalent_atoms"]
835        wyckoffs = np.array(spgdata["wyckoffs"])
836
837        wyck_mult = [np.count_nonzero(equivalent_atoms == equivalent_atoms[i]) for i in range(natom)]
838        wyck_mult = np.array(wyck_mult, dtype=int)
839
840        irred_pos = []
841        eqmap = collections.defaultdict(list)
842        for pos, eqpos in enumerate(equivalent_atoms):
843            eqmap[eqpos].append(pos)
844            # Add it to irred_pos if it's irreducible.
845            if pos == eqpos: irred_pos.append(pos)
846
847        # Convert to numpy arrays
848        irred_pos = np.array(irred_pos)
849        for eqpos in eqmap:
850            eqmap[eqpos] = np.array(eqmap[eqpos], dtype=int)
851
852        if printout:
853            print("Found %d inequivalent position(s):" % len(irred_pos))
854            for i, irr_pos in enumerate(sorted(eqmap.keys())):
855                print("Wyckoff position: (%s%s)" % (wyck_mult[irr_pos], wyckoffs[irr_pos]))
856                print("\t[%d]: %s" % (irr_pos, repr(self[irr_pos])))
857                for eqind in eqmap[irr_pos]:
858                    if eqind == irr_pos: continue
859                    print("\t[%d]: %s" % (eqind, repr(self[eqind])))
860            print("")
861
862        # Build list of labels from multiplicity and name: e.g. 3a
863        wyck_labels = np.array(["%s%s" % (wmul, wsymb) for wsymb, wmul in zip(wyckoffs, wyck_mult)])
864
865        # Build labels for sites with chemical element.
866        site_labels = []
867        for i, (site, wsymb, wmul) in enumerate(zip(self, wyckoffs, wyck_mult)):
868            site_labels.append("%s%d (%s%s)" % (site.specie.symbol, i, wmul, wsymb))
869
870        return dict2namedtuple(irred_pos=irred_pos, eqmap=eqmap, wyckoffs=wyckoffs,
871                               wyck_mult=wyck_mult, wyck_labels=wyck_labels,
872                               site_labels=np.array(site_labels), spgdata=spgdata)
873
874    def spget_summary(self, symprec=1e-3, angle_tolerance=5, site_symmetry=False, verbose=0):
875        """
876        Return string with full information about crystalline structure i.e.
877        space group, point group, wyckoff positions, equivalent sites.
878
879        Args:
880            symprec (float): Symmetry precision for distance.
881            angle_tolerance (float): Tolerance on angles.
882            site_symmetry: True to show site symmetries i.e. the point group operations that leave the site invariant.
883            verbose (int): Verbosity level.
884        """
885        spgan = SpacegroupAnalyzer(self, symprec=symprec, angle_tolerance=angle_tolerance)
886        spgdata = spgan.get_symmetry_dataset()
887        # Get spacegroup number computed by Abinit if available.
888        abispg_number = None if self.abi_spacegroup is None else self.abi_spacegroup.spgid
889
890        # Print lattice info
891        outs = ["Full Formula ({s})".format(s=self.composition.formula),
892                "Reduced Formula: {}".format(self.composition.reduced_formula)]
893        app = outs.append
894        to_s = lambda x: "%0.6f" % x
895        outs.append("abc   : " + " ".join([to_s(i).rjust(10)
896                                           for i in self.lattice.abc]))
897        outs.append("angles: " + " ".join([to_s(i).rjust(10)
898                                           for i in self.lattice.angles]))
899        app("")
900        app("Spglib space group info (magnetic symmetries not taken into account).")
901        app("Spacegroup: %s (%s), Hall: %s, Abinit spg_number: %s" % (
902             spgan.get_space_group_symbol(), spgan.get_space_group_number(), spgan.get_hall(), str(abispg_number)))
903        app("Crystal_system: %s, Lattice_type: %s, Point_group: %s" % (
904            spgan.get_crystal_system(), spgan.get_lattice_type(), spgan.get_point_group_symbol()))
905        app("")
906
907        wickoffs, equivalent_atoms = spgdata["wyckoffs"], spgdata["equivalent_atoms"]
908        header = ["Idx", "Symbol", "Reduced_Coords", "Wyckoff", "EqIdx"]
909
910        if site_symmetry:
911            header.append("site_symmetry")
912            sitesym_labels = self.spget_site_symmetries()
913
914        table = [header]
915        for i, site in enumerate(self):
916            mult = np.count_nonzero(equivalent_atoms == equivalent_atoms[i])
917            row = [
918                i,
919                site.species_string,
920                "%+.5f %+.5f %+.5f" % tuple(site.frac_coords),
921                "(%s%s)" % (mult, wickoffs[i]),
922                "%d" % equivalent_atoms[i],
923            ]
924            if site_symmetry: row.append(sitesym_labels[i])
925
926            table.append(row)
927
928        from tabulate import tabulate
929        app(tabulate(table, headers="firstrow"))
930
931        # Print entire dataset.
932        if verbose > 1:
933            app("\nSpglib dataset:")
934            app(pformat(spgdata, indent=4))
935
936        return "\n".join(outs)
937
938    @property
939    def abi_spacegroup(self):
940        """
941        :class:`AbinitSpaceGroup` instance with Abinit symmetries read from the netcd file.
942        None if abinit symmetries are not available e.g. if the structure has been created
943        from a CIF file.
944        """
945        try:
946            return self._abi_spacegroup
947        except AttributeError:
948            return None
949
950    def set_abi_spacegroup(self, spacegroup):
951        """``AbinitSpaceGroup`` setter."""
952        self._abi_spacegroup = spacegroup
953
954    @property
955    def has_abi_spacegroup(self):
956        """True is the structure contains info on the spacegroup."""
957        return self.abi_spacegroup is not None
958
959    def spgset_abi_spacegroup(self, has_timerev, overwrite=False):
960        """
961        Call spglib to find the spacegroup of the crystal, create new
962        :class:`AbinitSpaceGroup` object and store it in ``self.abi_spacegroup``.
963
964        Args:
965            has_timerev (bool): True if time-reversal can be used.
966            overwrite (bool): By default, the method raises `ValueError` if the object
967                already has the list of symmetries found by Abinit.
968
969        Returns: :class:`AbinitSpaceGroup`
970
971        .. warning:
972
973            This method should be called only if the Abipy structure does not have
974            spacegroup symmetries e.g. if we are reading a CIF file or if the structure
975            is initialized from an output file produced by another code.
976        """
977        if self.has_abi_spacegroup and not overwrite:
978            raise ValueError(("Structure object already has an Abinit spacegroup object.\n"
979                              "Use `overwrite=True` to allow modification."))
980
981        msg = ("Structure object does not have symmetry operations computed from Abinit.\n"
982               "Calling spglib to get symmetry operations.")
983        cprint(msg, "magenta")
984
985        spglib_data = SpacegroupAnalyzer(self).get_symmetry_dataset()
986        spgid = spglib_data["number"]
987        symrel, tnons = spglib_data["rotations"], spglib_data["translations"]
988        # TODO: Anti-ferromagnetic symmetries are not supported by spglib
989        symafm = [1] * len(symrel)
990
991        abispg = AbinitSpaceGroup(spgid, symrel, tnons, symafm, has_timerev, inord="C")
992        self.set_abi_spacegroup(abispg)
993
994        return abispg
995
996    @property
997    def indsym(self):
998        """
999        Compute indsym (natom, nsym, 4) array.
1000
1001        For each isym,iatom, the fourth element is label of atom into
1002        which iatom is sent by INVERSE of symmetry operation isym;
1003        first three elements are the primitive translations which must be
1004        subtracted after the transformation to get back to the original unit cell (see symatm.F90).
1005        """
1006        if getattr(self, "_indsym", None) is not None: return self._indsym
1007        if not self.has_abi_spacegroup:
1008            self.spgset_abi_spacegroup(has_timerev=True, overwrite=False)
1009
1010        from abipy.core.symmetries import indsym_from_symrel
1011        self._indsym = indsym_from_symrel(self.abi_spacegroup.symrel, self.abi_spacegroup.tnons, self, tolsym=1e-8)
1012        return self._indsym
1013
1014    @indsym.setter
1015    def indsym(self, indsym):
1016        """Set indsym array."""
1017        if getattr(self, "_indsym", None) is not None:
1018            cprint("structure.indsym is already set!", "yellow")
1019        self._indsym = indsym
1020
1021    @lazy_property
1022    def site_symmetries(self):
1023        """Object with SiteSymmetries."""
1024        from abipy.core.site_symmetries import SiteSymmetries
1025        return SiteSymmetries(self)
1026
1027    # TODO: site_symmetry or spget_site_symmetries?
1028    def spget_site_symmetries(self):
1029        import spglib
1030        indsym = self.indsym
1031        symrel, symafm = self.abi_spacegroup.symrel, self.abi_spacegroup.symafm
1032        nsym = len(symrel)
1033        sitesym_labels = []
1034        for iatom, site in enumerate(self):
1035            rotations = [symrel[isym] for isym in range(nsym) if
1036                         indsym[iatom, isym, 3] == iatom and symafm[isym] == +1]
1037            # Passing a 0-length rotations list to spglib can segfault.
1038            herm_symbol, ptg_num = "1", 1
1039            if len(rotations) != 0:
1040                herm_symbol, ptg_num, trans_mat = spglib.get_pointgroup(rotations)
1041
1042            sitesym_labels.append("%s (#%d,nsym:%d)" % (herm_symbol.strip(), ptg_num, len(rotations)))
1043
1044        return sitesym_labels
1045
1046    def abiget_spginfo(self, tolsym=None, pre=None):
1047        """
1048        Call Abinit to get spacegroup information.
1049        Return dictionary with e.g.
1050        {'bravais': 'Bravais cF (face-center cubic)', 'spg_number': 227, 'spg_symbol': 'Fd-3m'}.
1051
1052        Args:
1053            tolsym: Abinit tolsym input variable. None correspondes to the default value.
1054            pre: Keywords in dictionary are prepended with this string
1055        """
1056        from abipy.data.hgh_pseudos import HGH_TABLE
1057        from abipy.abio import factories
1058        gsinp = factories.gs_input(self, HGH_TABLE, spin_mode="unpolarized")
1059        gsinp["chkprim"] = 0
1060        d = gsinp.abiget_spacegroup(tolsym=tolsym, retdict=True)
1061        if pre: d = {pre + k: v for k, v in d.items()}
1062        return d
1063
1064    def print_neighbors(self, radius=2.0):
1065        """
1066        Get neighbors for each atom in the unit cell, out to a distance ``radius`` in Angstrom
1067        Print results.
1068        """
1069        print(" ")
1070        print("Finding neighbors for each atom in the unit cell, out to a distance %s (Angstrom)" % radius)
1071        print(" ")
1072
1073        ns = self.get_all_neighbors_old(radius, include_index=False)
1074        for i, (site, sited_list) in enumerate(zip(self, ns)):
1075            print("[%s] site %s has %s neighbors:" % (i, repr(site), len(sited_list)))
1076            for s, dist in sorted(sited_list, key=lambda t: t[1]):
1077                print("\t", repr(s), " at distance", dist)
1078            print("")
1079
1080    @lazy_property
1081    def hsym_kpath(self):
1082        """
1083        Returns an instance of :class:`pymatgen.symmetry.bandstructure.HighSymmKpath`.
1084        (Database of high symmetry k-points and high symmetry lines).
1085        """
1086        from pymatgen.symmetry.bandstructure import HighSymmKpath
1087        return HighSymmKpath(self)
1088
1089    @lazy_property
1090    def hsym_kpoints(self):
1091        """|KpointList| object with the high-symmetry K-points."""
1092        # Get mapping name --> frac_coords for the special k-points in the database.
1093        name2frac_coords = self.hsym_kpath.kpath["kpoints"]
1094        kpath = self.hsym_kpath.kpath["path"]
1095
1096        frac_coords, names = [], []
1097        for segment in kpath:
1098            for name in segment:
1099                fc = name2frac_coords[name]
1100                frac_coords.append(fc)
1101                names.append(name)
1102
1103        # Build KpointList instance.
1104        from .kpoints import KpointList
1105        return KpointList(self.reciprocal_lattice, frac_coords, weights=None, names=names)
1106
1107    def get_kcoords_from_names(self, knames, cart_coords=False):
1108        """
1109        Return numpy array with the fractional coordinates of the high-symmetry k-points listed in `knames`.
1110
1111        Args:
1112            knames: List of strings with the k-point labels.
1113            cart_coords: True if the ``coords`` dataframe should contain Cartesian cordinates
1114                instead of Reduced coordinates.
1115        """
1116        kname2frac = {k.name: k.frac_coords for k in self.hsym_kpoints}
1117        # Add aliases for Gamma.
1118        if r"$\Gamma$" in kname2frac:
1119            kname2frac["G"] = kname2frac[r"$\Gamma$"]
1120            kname2frac["Gamma"] = kname2frac[r"$\Gamma$"]
1121
1122        try:
1123            kcoords = np.reshape([kname2frac[name] for name in list_strings(knames)], (-1, 3))
1124        except KeyError:
1125            cprint("Internal list of high-symmetry k-points:\n%s" % str(self.hsym_kpoints))
1126            raise
1127
1128        if cart_coords:
1129            kcoords = self.reciprocal_lattice.get_cartesian_coords(kcoords)
1130
1131        return kcoords
1132
1133    @lazy_property
1134    def hsym_stars(self):
1135        """
1136        List of |KpointStar| objects. Each star is associated to one of the special k-points
1137        present in the pymatgen database.
1138        """
1139        # Construct the stars.
1140        return [kpoint.compute_star(self.abi_spacegroup.fm_symmops) for kpoint in self.hsym_kpoints]
1141
1142    # TODO
1143    #def get_star_kpoint(self, kpoint):
1144
1145    #    # Call spglib to get spacegroup if Abinit spacegroup is not available.
1146    #    if self.abi_spacegroup is None:
1147    #        self.spgset_abi_spacegroup(has_timerev=not options.no_time_reversal)
1148
1149    #    kpoint = Kpoint(options.kpoint, self.reciprocal_lattice)
1150    #    kstar = kpoint.compute_star(self.abi_spacegroup, wrap_tows=True)
1151    #    return kstar
1152    #    #print("Found %s points in the star of %s\n" % (len(kstar), repr(kpoint)))
1153    #    #for k in kstar:
1154    #    #    print(4 * " ", repr(k))
1155
1156    def get_sorted_structure_z(self):
1157        """Order the structure according to increasing Z of the elements"""
1158        return self.__class__.from_sites(sorted(self.sites, key=lambda site: site.specie.Z))
1159
1160    def findname_in_hsym_stars(self, kpoint):
1161        """
1162        Returns the name of the special k-point, None if kpoint is unknown.
1163        """
1164        if self.abi_spacegroup is None: return None
1165
1166        from .kpoints import Kpoint
1167        kpoint = Kpoint.as_kpoint(kpoint, self.reciprocal_lattice)
1168
1169        # Try to find kpoint in hsym_stars without taking into accout symmetry operation (compare with base_point)
1170        # Important if there are symmetry equivalent k-points in hsym_kpoints e.g. K and U in FCC lattice
1171        # as U should not be mapped onto K as done in the second loop below.
1172        from .kpoints import issamek
1173        for star in self.hsym_stars:
1174            if issamek(kpoint.frac_coords, star.base_point.frac_coords):
1175                return star.name
1176
1177        # Now check if kpoint is in one of the stars.
1178        for star in self.hsym_stars:
1179            i = star.find(kpoint)
1180            if i != -1:
1181                #print("input kpt:", kpoint, "star image", star[i], star[i].name)
1182                return star.name
1183        else:
1184            return None
1185
1186    def get_symbol2indices(self):
1187        """
1188        Return a dictionary mapping chemical symbols to numpy array with the position of the atoms.
1189
1190        Example:
1191
1192            MgB2 --> {Mg: [0], B: [1, 2]}
1193        """
1194        return {symbol: np.array(self.indices_from_symbol(symbol)) for symbol in self.symbol_set}
1195
1196    def get_symbol2coords(self):
1197        """
1198        Return a dictionary mapping chemical symbols to a [ntype_symbol, 3] numpy array
1199        with the fractional coordinates.
1200        """
1201        # TODO:
1202        #use structure.frac_coords but add reshape in pymatgen.
1203        #fcoords = np.reshape([s.frac_coords for s in self], (-1, 3))
1204        coords = {}
1205        for symbol in self.symbol_set:
1206            coords[symbol] = np.reshape(
1207                [site.frac_coords for site in self if site.specie.symbol == symbol], (-1, 3))
1208
1209        return coords
1210
1211    def dot(self, coords_a, coords_b, space="r", frac_coords=False):
1212        """
1213        Compute the scalar product of vector(s) either in real space or
1214        reciprocal space.
1215
1216        Args:
1217            coords (3x1 array): Array-like object with the coordinates.
1218            space (str): "r" for real space, "g" for reciprocal space.
1219            frac_coords (bool): Whether the vector corresponds to fractional or
1220                cartesian coordinates.
1221
1222        Returns:
1223            one-dimensional `numpy` array.
1224        """
1225        lattice = {"r": self.lattice,
1226                   "g": self.reciprocal_lattice}[space.lower()]
1227        return lattice.dot(coords_a, coords_b, frac_coords=frac_coords)
1228
1229    def norm(self, coords, space="r", frac_coords=True):
1230        """
1231        Compute the norm of vector(s) either in real space or reciprocal space.
1232
1233        Args:
1234            coords (3x1 array): Array-like object with the coordinates.
1235            space (str): "r" for real space, "g" for reciprocal space.
1236            frac_coords (bool): Whether the vector corresponds to fractional or
1237                cartesian coordinates.
1238
1239        Returns:
1240            one-dimensional `numpy` array.
1241        """
1242        return np.sqrt(self.dot(coords, coords, space=space, frac_coords=frac_coords))
1243
1244    def scale_lattice(self, new_volume):
1245        """
1246        Return a new |Structure| with volume new_volume by performing a
1247        scaling of the lattice vectors so that length proportions and angles are preserved.
1248        """
1249        new_lattice = self.lattice.scale(new_volume)
1250        return self.__class__(new_lattice, self.species, self.frac_coords)
1251
1252    def get_dict4pandas(self, symprec=1e-2, angle_tolerance=5.0, with_spglib=True):
1253        """
1254        Return a :class:`OrderedDict` with the most important structural parameters:
1255
1256            - Chemical formula and number of atoms.
1257            - Lattice lengths, angles and volume.
1258            - The spacegroup number computed by Abinit (set to None if not available).
1259            - The spacegroup number and symbol computed by spglib (if `with_spglib`).
1260
1261        Useful to construct pandas DataFrames
1262
1263        Args:
1264            with_spglib (bool): If True, spglib is invoked to get the spacegroup symbol and number
1265            symprec (float): Symmetry precision used to refine the structure.
1266            angle_tolerance (float): Tolerance on angles.
1267        """
1268        abc, angles = self.lattice.abc, self.lattice.angles
1269
1270        # Get spacegroup info from spglib.
1271        spglib_symbol, spglib_number, spglib_lattice_type = None, None, None
1272        if with_spglib:
1273            try:
1274                spglib_symbol, spglib_number = self.get_space_group_info(symprec=symprec,
1275                                                                         angle_tolerance=angle_tolerance)
1276                spglib_lattice_type = self.spget_lattice_type(symprec=symprec, angle_tolerance=angle_tolerance)
1277            except Exception as exc:
1278                cprint("Spglib couldn't find space group symbol and number for composition: `%s`" %
1279                        str(self.composition), "red")
1280                print("Exception:\n", exc)
1281
1282        # Get spacegroup number computed by Abinit if available.
1283        abispg_number = None if self.abi_spacegroup is None else self.abi_spacegroup.spgid
1284
1285        od = OrderedDict([
1286            ("formula", self.formula), ("natom", self.num_sites),
1287            ("alpha", angles[0]), ("beta", angles[1]), ("gamma", angles[2]),
1288            ("a", abc[0]), ("b", abc[1]), ("c", abc[2]), ("volume", self.volume),
1289            ("abispg_num", abispg_number),
1290        ])
1291        if with_spglib:
1292            od["spglib_symb"] = spglib_symbol
1293            od["spglib_num"] = spglib_number
1294            od["spglib_lattice_type"] = spglib_lattice_type
1295
1296        return od
1297
1298    @add_fig_kwargs
1299    def plot(self, **kwargs):
1300        """
1301        Plot structure in 3D with matplotlib. Return matplotlib Figure
1302        See plot_structure for kwargs
1303        """
1304        from abipy.tools.plotting import plot_structure
1305        return plot_structure(self, **kwargs)
1306
1307    @add_fig_kwargs
1308    def plot_bz(self, ax=None, pmg_path=True, with_labels=True, **kwargs):
1309        """
1310        Use matplotlib to plot the symmetry line path in the Brillouin Zone.
1311
1312        Args:
1313            ax: matplotlib :class:`Axes` or None if a new figure should be created.
1314            pmg_path (bool): True if the default path used in pymatgen should be show.
1315            with_labels (bool): True to plot k-point labels.
1316
1317        Returns: |matplotlib-Figure|.
1318        """
1319        from pymatgen.electronic_structure.plotter import plot_brillouin_zone, plot_brillouin_zone_from_kpath
1320        labels = None if not with_labels else self.hsym_kpath.kpath["kpoints"]
1321        if pmg_path:
1322            return plot_brillouin_zone_from_kpath(self.hsym_kpath, ax=ax, show=False, **kwargs)
1323        else:
1324            return plot_brillouin_zone(self.reciprocal_lattice, ax=ax, labels=labels, show=False, **kwargs)
1325
1326    @add_fig_kwargs
1327    def plot_xrd(self, wavelength="CuKa", symprec=0, debye_waller_factors=None,
1328                 two_theta_range=(0, 90), annotate_peaks=True, ax=None, **kwargs):
1329        """
1330        Use pymatgen :class:`XRDCalculator` to show the XRD plot.
1331
1332        Args:
1333            wavelength (str/float): The wavelength can be specified as either a
1334                float or a string. If it is a string, it must be one of the
1335                supported definitions in the AVAILABLE_RADIATION class
1336                variable, which provides useful commonly used wavelengths.
1337                If it is a float, it is interpreted as a wavelength in
1338                angstroms. Defaults to "CuKa", i.e, Cu K_alpha radiation.
1339            symprec (float): Symmetry precision for structure refinement. If
1340                set to 0, no refinement is done. Otherwise, refinement is
1341                performed using spglib_ with provided precision.
1342            debye_waller_factors ({element symbol: float}): Allows the
1343                specification of Debye-Waller factors. Note that these
1344                factors are temperature dependent.
1345            two_theta_range ([float of length 2]): Tuple for range of
1346                two_thetas to calculate in degrees. Defaults to (0, 90). Set to
1347                None if you want all diffracted beams within the limiting
1348                sphere of radius 2 / wavelength.
1349            annotate_peaks (bool): Whether to annotate the peaks with plane information.
1350            ax: matplotlib :class:`Axes` or None if a new figure should be created.
1351
1352        Returns: |matplotlib-Figure|
1353        """
1354        ax, fig, plt = get_ax_fig_plt(ax=ax)
1355        from pymatgen.analysis.diffraction.xrd import XRDCalculator
1356        xrd = XRDCalculator(wavelength=wavelength, symprec=symprec, debye_waller_factors=debye_waller_factors)
1357        xrd.get_plot(self, two_theta_range=two_theta_range, annotate_peaks=annotate_peaks, ax=ax)
1358
1359        return fig
1360
1361    def yield_figs(self, **kwargs):  # pragma: no cover
1362        """
1363        This function *generates* a predefined list of matplotlib figures with minimal input from the user.
1364        """
1365        yield self.plot(show=False)
1366        yield self.plot_bz(show=False)
1367
1368    def export(self, filename, visu=None, verbose=1):
1369        """
1370        Export the crystalline structure to file ``filename``.
1371
1372        Args:
1373            filename (str): String specifying the file path and the file format.
1374                The format is defined by the file extension. filename="prefix.xsf", for example,
1375                will produce a file in XSF format. An *empty* prefix, e.g. ".xsf" makes the code use a temporary file.
1376            visu: |Visualizer| subclass. By default, this method returns the first available
1377                visualizer that supports the given file format. If visu is not None, an
1378                instance of visu is returned. See |Visualizer| for the list of applications and formats supported.
1379            verbose: Verbosity level
1380
1381        Returns: ``Visulizer`` instance.
1382        """
1383        if "." not in filename:
1384            raise ValueError("Cannot detect extension in filename %s:" % filename)
1385
1386        tokens = filename.strip().split(".")
1387        ext = tokens[-1]
1388        #print("tokens", tokens, "ext", ext)
1389        #if ext == "POSCAR":
1390
1391        if not tokens[0]:
1392            # filename == ".ext" ==> Create temporary file.
1393            # nbworkdir in cwd is needed when we invoke the method from a notebook.
1394            from abipy.core.globals import abinb_mkstemp
1395            _, rpath = abinb_mkstemp(force_abinb_workdir=False, use_relpath=False,
1396                                     suffix="." + ext, text=True)
1397            #if abilab.in_notebook():
1398            #    _, filename = tempfile.mkstemp(suffix="." + ext, dir=abilab.get_abipy_nbworkdir(), text=True)
1399            #else:
1400            #    _, filename = tempfile.mkstemp(suffix="." + ext, text=True)
1401
1402        if ext.lower() in ("xsf", "poscar", "cif"):
1403            if verbose:
1404                print("Writing data to:", filename, "with fmt:", ext.lower())
1405            s = self.to(fmt=ext)
1406            with open(filename, "wt") as fh:
1407                fh.write(s)
1408
1409        if visu is None:
1410            return Visualizer.from_file(filename)
1411        else:
1412            return visu(filename)
1413
1414    def get_chemview(self, **kwargs): # pragma: no cover
1415        """
1416        Visualize structure inside the jupyter notebook using chemview package.
1417        """
1418        from pymatgen.vis.structure_chemview import quick_view
1419        return quick_view(self, **kwargs)
1420
1421    def plot_vtk(self, show=True, **kwargs):
1422        """
1423        Visualize structure with VTK. Requires vVTK python bindings.
1424
1425        Args:
1426            show: True to show structure immediately.
1427            kwargs: keyword arguments passed to :class:`StructureVis`.
1428
1429        Return: StructureVis object.
1430        """
1431        from pymatgen.vis.structure_vtk import StructureVis
1432        vis = StructureVis(**kwargs)
1433        vis.set_structure(self, to_unit_cell=True)
1434        if show: vis.show()
1435        return vis
1436
1437    def plot_mayaview(self, figure=None, show=True, **kwargs):
1438        """Visualize structure with mayavi."""
1439        from abipy.display import mvtk
1440        return mvtk.plot_structure(self, figure=figure, show=show, **kwargs)
1441
1442    @add_fig_kwargs
1443    def plot_atoms(self, rotations="default", **kwargs):
1444        """
1445        Plot 2d representation with matplotlib using ASE `plot_atoms` function.
1446
1447        Args:
1448            rotations: String or List of strings.
1449                Each string defines a rotation (in degrees) in the form '10x,20y,30z'
1450                Note that the order of rotation matters, i.e. '50x,40z' is different from '40z,50x'.
1451            kwargs: extra kwargs passed to plot_atoms ASE function.
1452
1453        Returns: |matplotlib-Figure|
1454        """
1455        if rotations == "default":
1456            rotations = [
1457                "", "90x", "90y",
1458                "45x,45y", "45y,45z", "45x,45z",
1459            ]
1460        else:
1461            rotations = list_strings(rotations)
1462
1463        nrows, ncols, num_plots = 1, 1, len(rotations)
1464        if num_plots > 1:
1465            ncols = 3
1466            nrows = num_plots // ncols + num_plots % ncols
1467
1468        ax_mat, fig, plt = get_axarray_fig_plt(None, nrows=nrows, ncols=ncols,
1469                                                sharex=False, sharey=True, squeeze=False)
1470
1471        # don't show the last ax if num_plots is odd.
1472        if num_plots % ncols != 0: ax_mat[-1, -1].axis("off")
1473
1474        from ase.visualize.plot import plot_atoms
1475        atoms = self.to_ase_atoms()
1476        for rotation, ax in zip(rotations, ax_mat.flat):
1477            plot_atoms(atoms, ax=ax, rotation=rotation, **kwargs)
1478            ax.set_axis_off()
1479            if rotation:
1480                ax.set_title("rotation: %s" % str(rotation), fontsize=6)
1481
1482        return fig
1483
1484    def get_ngl_view(self): # pragma: no cover
1485        """
1486        Visualize the structure with nglview inside the jupyter notebook.
1487        """
1488        try:
1489            import nglview as nv
1490        except ImportError:
1491            raise ImportError("nglview is not installed. See https://github.com/arose/nglview")
1492
1493        view = nv.show_pymatgen(self)
1494        view.add_unitcell()
1495        return view
1496
1497    def get_crystaltk_view(self): # pragma: no cover
1498        """
1499        Visualize the structure with crystal_toolkit inside the jupyter notebook.
1500        """
1501        try:
1502            from crystal_toolkit import view
1503        except ImportError:
1504            raise ImportError("crystal_toolkit is not installed. See https://docs.crystaltoolkit.org/jupyter")
1505
1506        return view(self)
1507
1508    def get_jsmol_view(self, symprec=None, verbose=0, **kwargs): # pragma: no cover
1509        """
1510        Visualize the structure with jsmol inside the jupyter notebook.
1511
1512        Args:
1513            symprec (float): If not none, finds the symmetry of the structure
1514                and writes the CIF with symmetry information.
1515                Passes symprec to the spglib SpacegroupAnalyzer.
1516            verbose: Verbosity level.
1517        """
1518        try:
1519            from jupyter_jsmol import JsmolView
1520        except ImportError:
1521            raise ImportError("jupyter_jsmol is not installed. See https://github.com/fekad/jupyter-jsmol")
1522
1523        cif_str = self.write_cif_with_spglib_symms(None, symprec=symprec, ret_string=True)
1524        print("cif_str:\n", cif_str)
1525        #return JsmolView.from_str(cif_str)
1526
1527        from IPython.display import display, HTML
1528        # FIXME TEMPORARY HACK TO LOAD JSMOL.js
1529        # See discussion at
1530        #   https://stackoverflow.com/questions/16852885/ipython-adding-javascript-scripts-to-ipython-notebook
1531        display(HTML('<script type="text/javascript" src="/nbextensions/jupyter-jsmol/jsmol/JSmol.min.js"></script>'))
1532
1533        jsmol = JsmolView(color='white')
1534        display(jsmol)
1535        cmd = 'load inline "%s" {1 1 1}' % cif_str
1536        if verbose: print("executing cmd:", cmd)
1537        jsmol.script(cmd)
1538
1539        return jsmol
1540
1541    def visualize(self, appname="vesta"):
1542        """
1543        Visualize the crystalline structure with visualizer.
1544        See |Visualizer| for the list of applications and formats supported.
1545        """
1546        if appname in ("mpl", "matplotlib"): return self.plot()
1547        if appname == "vtk": return self.plot_vtk()
1548        if appname == "mayavi": return self.plot_mayaview()
1549
1550        # Get the Visualizer subclass from the string.
1551        visu = Visualizer.from_name(appname)
1552
1553        # Try to export data to one of the formats supported by the visualizer
1554        # Use a temporary file (note "." + ext)
1555        for ext in visu.supported_extensions():
1556            ext = "." + ext
1557            try:
1558                return self.export(ext, visu=visu)()
1559            except visu.Error as exc:
1560                cprint(str(exc), color="red")
1561                pass
1562        else:
1563            raise visu.Error("Don't know how to export data for %s" % appname)
1564
1565    def convert(self, fmt="cif", **kwargs):
1566        """
1567        Return string with the structure in the given format `fmt`
1568        Options include "abivars", "cif", "xsf", "poscar", "siesta", "wannier90", "cssr", "json".
1569        """
1570        if fmt in ("abivars", "abinit"):
1571            return self.abi_string
1572        elif fmt == "abipython":
1573            return pformat(self.to_abivars(), indent=4)
1574        elif fmt == "qe":
1575            from pymatgen.io.pwscf import PWInput
1576            return str(PWInput(self, pseudo={s: s + ".pseudo" for s in self.symbol_set}))
1577        elif fmt == "siesta":
1578            return structure2siesta(self)
1579        elif fmt in ("wannier90", "w90"):
1580            from abipy.wannier90.win import structure2wannier90
1581            return structure2wannier90(self)
1582        elif fmt.lower() == "poscar":
1583            # Don't call super for poscar because we need more significant_figures to
1584            # avoid problems with abinit space group routines where the default numerical tolerance is tight.
1585            from pymatgen.io.vasp import Poscar
1586            return Poscar(self).get_string(significant_figures=12)
1587        else:
1588            return super().to(fmt=fmt, **kwargs)
1589
1590    def displace(self, displ, eta, frac_coords=True, normalize=True):
1591        """
1592        Displace the sites of the structure along the displacement vector displ.
1593
1594        The displacement vector is first rescaled so that the maxium atomic displacement
1595        is one Angstrom, and then multiplied by eta. Hence passing eta=0.001, will move
1596        all the atoms so that the maximum atomic displacement is 0.001 Angstrom.
1597
1598        Args:
1599            displ: Displacement vector with 3*len(self) entries (fractional coordinates).
1600            eta: Scaling factor.
1601            frac_coords: Boolean stating whether the vector corresponds to fractional or cartesian coordinates.
1602        """
1603        # Get a copy since we are going to modify displ.
1604        displ = np.reshape(displ, (-1, 3)).copy()
1605
1606        if len(displ) != len(self):
1607            raise ValueError("Displ must contains 3 * natom entries")
1608        if np.iscomplexobj(displ):
1609            raise TypeError("Displacement cannot be complex")
1610
1611        if not frac_coords:
1612            # Convert to fractional coordinates.
1613            displ = np.reshape([self.lattice.get_fractional_coords(vec) for vec in displ], (-1,3))
1614
1615        # Normalize the displacement so that the maximum atomic displacement is 1 Angstrom.
1616        if normalize:
1617            dnorm = self.norm(displ, space="r")
1618            displ /= np.max(np.abs(dnorm))
1619
1620        # Displace the sites.
1621        for i in range(len(self)):
1622            self.translate_sites(indices=i, vector=eta * displ[i, :], frac_coords=True)
1623
1624    def get_smallest_supercell(self, qpoint, max_supercell):
1625        """
1626        Args:
1627            qpoint: q vector in reduced coordinates in reciprocal space
1628            max_supercell: vector with the maximum supercell size
1629
1630        Returns: the scaling matrix of the supercell
1631        """
1632        if np.allclose(qpoint, 0.0):
1633            scale_matrix = np.eye(3, 3)
1634            return scale_matrix
1635
1636        l = max_supercell
1637
1638        # Inspired from Exciting Fortran code phcell.F90
1639        # It should be possible to improve this coding.
1640        scale_matrix = np.zeros((3, 3), dtype=int)
1641        dmin = np.inf
1642        found = False
1643
1644        # Try to reduce the matrix
1645        rprimd = self.lattice.matrix
1646        for l1 in np.arange(-l[0], l[0] + 1):
1647            for l2 in np.arange(-l[1], l[1] + 1):
1648                for l3 in np.arange(-l[2], l[2] + 1):
1649                    lnew = np.array([l1, l2, l3])
1650                    ql = np.dot(lnew, qpoint)
1651                    # Check if integer and non zero !
1652                    if np.abs(ql - np.round(ql)) < 1e-6:
1653                        Rl = np.dot(lnew, rprimd)
1654                        # Normalize the displacement so that the maximum atomic displacement is 1 Angstrom.
1655                        dnorm = np.sqrt(np.dot(Rl, Rl))
1656                        if dnorm < (dmin-1e-6) and dnorm > 1e-6:
1657                            found = True
1658                            scale_matrix[:, 0] = lnew
1659                            dmin = dnorm
1660        if not found:
1661            raise ValueError('max_supercell is not large enough for this q-point')
1662
1663        found = False
1664        dmin = np.inf
1665        for l1 in np.arange(-l[0], l[0] + 1):
1666            for l2 in np.arange(-l[1], l[1] + 1):
1667                for l3 in np.arange(-l[2], l[2] + 1):
1668                    lnew = np.array([l1, l2, l3])
1669                    # Check if not parallel !
1670                    cp = np.cross(lnew, scale_matrix[:,0])
1671                    if np.dot(cp, cp) > 1e-6:
1672                        ql = np.dot(lnew, qpoint)
1673                        # Check if integer and non zero !
1674                        if np.abs(ql - np.round(ql)) < 1e-6:
1675                            Rl = np.dot(lnew, rprimd)
1676                            dnorm = np.sqrt(np.dot(Rl, Rl))
1677                            if dnorm < (dmin-1e-6) and dnorm > 1e-6:
1678                                found = True
1679                                scale_matrix[:, 1] = lnew
1680                                dmin = dnorm
1681        if not found:
1682            raise ValueError('max_supercell is not large enough for this q-point')
1683
1684        dmin = np.inf
1685        found = False
1686        for l1 in np.arange(-l[0], l[0] + 1):
1687            for l2 in np.arange(-l[1], l[1] + 1):
1688                for l3 in np.arange(-l[2], l[2] + 1):
1689                    lnew = np.array([l1, l2, l3])
1690                    # Check if not parallel!
1691                    cp = np.dot(np.cross(lnew, scale_matrix[:, 0]), scale_matrix[:, 1])
1692                    if cp > 1e-6:
1693                        # Should be positive as (R3 X R1).R2 > 0 for abinit!
1694                        ql = np.dot(lnew, qpoint)
1695                        # Check if integer and non zero!
1696                        if np.abs(ql - np.round(ql)) < 1e-6:
1697                            Rl = np.dot(lnew, rprimd)
1698                            dnorm = np.sqrt(np.dot(Rl,Rl))
1699                            if dnorm < (dmin-1e-6) and dnorm > 1e-6:
1700                                found = True
1701                                scale_matrix[:, 2] = lnew
1702                                dmin = dnorm
1703        if not found:
1704            raise ValueError('max_supercell is not large enough for this q-point')
1705
1706        # Fortran 2 python!!!
1707        return scale_matrix.T
1708
1709    def get_trans_vect(self, scale_matrix):
1710        """
1711        Returns the translation vectors for a given scale matrix.
1712
1713        Args:
1714            scale_matrix: Scale matrix defining the new lattice vectors in term of the old ones
1715
1716        Return: the translation vectors
1717        """
1718        scale_matrix = np.array(scale_matrix, np.int16)
1719        if scale_matrix.shape != (3, 3):
1720            scale_matrix = np.array(scale_matrix * np.eye(3), np.int16)
1721
1722        def range_vec(i):
1723            low = 0
1724            high = 0
1725            for z in scale_matrix[:, i]:
1726                if z > 0:
1727                    high += z
1728                else:
1729                    low += z
1730            return np.arange(low, high+1)
1731
1732        arange = range_vec(0)[:, None] * np.array([1, 0, 0])[None, :]
1733        brange = range_vec(1)[:, None] * np.array([0, 1, 0])[None, :]
1734        crange = range_vec(2)[:, None] * np.array([0, 0, 1])[None, :]
1735        all_points = arange[:, None, None] + brange[None, :, None] + crange[None, None, :]
1736        all_points = all_points.reshape((-1, 3))
1737
1738        # find the translation vectors (in terms of the initial lattice vectors)
1739        # that are inside the unit cell defined by the scale matrix
1740        # we're using a slightly offset interval from 0 to 1 to avoid numerical
1741        # precision issues
1742        #print(scale_matrix)
1743        inv_matrix = np.linalg.inv(scale_matrix)
1744
1745        frac_points = np.dot(all_points, inv_matrix)
1746        tvects = all_points[np.where(np.all(frac_points < 1-1e-10, axis=1)
1747                                     & np.all(frac_points >= -1e-10, axis=1))]
1748        assert len(tvects) == np.round(abs(np.linalg.det(scale_matrix)))
1749
1750        return tvects
1751
1752    def write_vib_file(self, xyz_file, qpoint, displ, do_real=True, frac_coords=True,
1753                       scale_matrix=None, max_supercell=None):
1754        """
1755        Write into the file descriptor xyz_file the positions and displacements of the atoms
1756
1757        Args:
1758            xyz_file: file_descriptor
1759            qpoint: qpoint to be analyzed
1760            displ: eigendisplacements to be analyzed
1761            do_real: True if you want to get only real part, False means imaginary part
1762            frac_coords: True if the eigendisplacements are given in fractional coordinates
1763            scale_matrix: Scale matrix for supercell
1764            max_supercell: Maximum size of supercell vectors with respect to primitive cell
1765        """
1766        if scale_matrix is None:
1767            if max_supercell is None:
1768                raise ValueError("If scale_matrix is not provided, please provide max_supercell!")
1769
1770            scale_matrix = self.get_smallest_supercell(qpoint, max_supercell=max_supercell)
1771
1772        old_lattice = self._lattice
1773        new_lattice = Lattice(np.dot(scale_matrix, old_lattice.matrix))
1774
1775        tvects = self.get_trans_vect(scale_matrix)
1776
1777        new_displ = np.zeros(3, dtype=float)
1778
1779        fmtstr = "{{}} {{:.{0}f}} {{:.{0}f}} {{:.{0}f}} {{:.{0}f}} {{:.{0}f}} {{:.{0}f}}\n".format(6)
1780
1781        for at, site in enumerate(self):
1782            for t in tvects:
1783                if do_real:
1784                    new_displ[:] = np.real(np.exp(2*1j*np.pi*(np.dot(qpoint,t)))*displ[at,:])
1785                else:
1786                    new_displ[:] = np.imag(np.exp(2*1j*np.pi*(np.dot(qpoint,t)))*displ[at,:])
1787                if frac_coords:
1788                    # Convert to fractional coordinates.
1789                    new_displ = self.lattice.get_cartesian_coords(new_displ)
1790
1791                # We don't normalize here !!!
1792                fcoords = site.frac_coords + t
1793
1794                coords = old_lattice.get_cartesian_coords(fcoords)
1795
1796                new_fcoords = new_lattice.get_fractional_coords(coords)
1797
1798                # New_fcoords -> map into 0 - 1
1799                new_fcoords = np.mod(new_fcoords, 1)
1800                coords = new_lattice.get_cartesian_coords(new_fcoords)
1801
1802                xyz_file.write(fmtstr.format(site.specie, coords[0], coords[1], coords[2],
1803                               new_displ[0], new_displ[1], new_displ[2]))
1804
1805    def frozen_2phonon(self, qpoint, displ1, displ2, eta=1, frac_coords=False, scale_matrix=None, max_supercell=None):
1806        """
1807        Creates the supercell needed for a given qpoint and adds the displacements.
1808        The displacements are normalized so that the largest atomic displacement will correspond to the
1809        value of eta in Angstrom.
1810
1811        Args:
1812            qpoint: q vector in reduced coordinate in reciprocal space.
1813            displ1: first displacement in real space of the atoms.
1814            displ2: second displacement in real space of the atoms.
1815            eta: pre-factor multiplying the displacement. Gives the value in Angstrom of the
1816                largest displacement.
1817            frac_coords: whether the displacements are given in fractional or cartesian coordinates
1818            scale_matrix: the scaling matrix of the supercell. If None a scaling matrix suitable for
1819                the qpoint will be determined.
1820            max_supercell: mandatory if scale_matrix is None, ignored otherwise. Defines the largest
1821                supercell in the search for a scaling matrix suitable for the q point.
1822
1823        Returns:
1824            A namedtuple with a Structure with the displaced atoms, a numpy array containing the
1825            displacements applied to each atom and the scale matrix used to generate the supercell.
1826        """
1827
1828        if scale_matrix is None:
1829            if max_supercell is None:
1830                raise ValueError("scale_matrix is not provided, please provide max_supercell!")
1831
1832            scale_matrix = self.get_smallest_supercell(qpoint, max_supercell=max_supercell)
1833
1834        scale_matrix = np.array(scale_matrix, np.int16)
1835        if scale_matrix.shape != (3, 3):
1836            scale_matrix = np.array(scale_matrix * np.eye(3), np.int16)
1837
1838        old_lattice = self._lattice
1839        new_lattice = Lattice(np.dot(scale_matrix, old_lattice.matrix))
1840
1841        tvects = self.get_trans_vect(scale_matrix)
1842
1843        if frac_coords:
1844            displ1 = np.array((old_lattice.get_cartesian_coords(d) for d in displ1))
1845            displ2 = np.array((old_lattice.get_cartesian_coords(d) for d in displ2))
1846        else:
1847            displ1 = np.array(displ1)
1848            displ2 = np.array(displ2)
1849
1850        # from here on displ are in cartesian coordinates
1851        norm_factor = np.linalg.norm(displ1+displ2, axis=1).max()
1852
1853        displ1 = eta * displ1 / norm_factor
1854        displ2 = eta * displ2 / norm_factor
1855
1856        new_displ1 = np.zeros(3, dtype=float)
1857        new_displ2 = np.zeros(3, dtype=float)
1858        new_sites = []
1859        displ_list = []
1860        for at,site in enumerate(self):
1861            for t in tvects:
1862                new_displ1[:] = np.real(np.exp(2*1j * np.pi * (np.dot(qpoint, t))) * displ1[at,:])
1863                new_displ2[:] = np.real(np.exp(2*1j * np.pi * (np.dot(qpoint, t))) * displ2[at,:])
1864
1865                displ_list.append(new_displ1 + new_displ2)
1866                coords = site.coords + old_lattice.get_cartesian_coords(t) + new_displ1 + new_displ2
1867                new_site = PeriodicSite(
1868                    site.species, coords, new_lattice,
1869                    coords_are_cartesian=True, properties=site.properties,
1870                    to_unit_cell=True)
1871                new_sites.append(new_site)
1872
1873        new_structure = self.__class__.from_sites(new_sites)
1874
1875        return dict2namedtuple(structure=new_structure, displ=np.array(displ_list), scale_matrix=scale_matrix)
1876
1877    def frozen_phonon(self, qpoint, displ, eta=1, frac_coords=False, scale_matrix=None, max_supercell=None):
1878        """
1879        Creates a supercell with displaced atoms for the specified q-point.
1880        The displacements are normalized so that the largest atomic displacement will correspond to the
1881        value of eta in Angstrom.
1882
1883        Args:
1884            qpoint: q vector in reduced coordinate in reciprocal space.
1885            displ: displacement in real space of the atoms.
1886            eta: pre-factor multiplying the displacement. Gives the value in Angstrom of the
1887                largest displacement.
1888            frac_coords: whether the displacements are given in fractional or cartesian coordinates
1889            scale_matrix: the scaling matrix of the supercell. If None a scaling matrix suitable for
1890                the qpoint will be determined.
1891            max_supercell: mandatory if scale_matrix is None, ignored otherwise. Defines the largest
1892                supercell in the search for a scaling matrix suitable for the q point.
1893
1894        Returns:
1895            A namedtuple with a Structure with the displaced atoms, a numpy array containing the
1896            displacements applied to each atom and the scale matrix used to generate the supercell.
1897        """
1898
1899        if scale_matrix is None:
1900            if max_supercell is None:
1901                raise ValueError("If scale_matrix is not provided, please provide max_supercell!")
1902
1903            scale_matrix = self.get_smallest_supercell(qpoint, max_supercell=max_supercell)
1904
1905        scale_matrix = np.array(scale_matrix, np.int16)
1906        if scale_matrix.shape != (3, 3):
1907            scale_matrix = np.array(scale_matrix * np.eye(3), np.int16)
1908        print("scale_matrix:", scale_matrix)
1909
1910        old_lattice = self._lattice
1911        new_lattice = Lattice(np.dot(scale_matrix, old_lattice.matrix))
1912
1913        tvects = self.get_trans_vect(scale_matrix)
1914        print("tvects", tvects)
1915
1916        if frac_coords:
1917            displ = np.array((old_lattice.get_cartesian_coords(d) for d in displ))
1918        else:
1919            displ = np.array(displ)
1920        # from here displ are in cartesian coordinates
1921
1922        displ = eta * displ / np.linalg.norm(displ, axis=1).max()
1923
1924        new_displ = np.zeros(3, dtype=float)
1925        new_sites = []
1926        displ_list = []
1927        for at, site in enumerate(self):
1928            for t in tvects:
1929                new_displ[:] = np.real(np.exp(2*1j*np.pi*(np.dot(qpoint,t)))*displ[at,:])
1930                displ_list.append(list(new_displ))
1931
1932                coords = site.coords + old_lattice.get_cartesian_coords(t) + new_displ
1933                new_site = PeriodicSite(
1934                    site.species, coords, new_lattice,
1935                    coords_are_cartesian=True, properties=site.properties,
1936                    to_unit_cell=True)
1937                new_sites.append(new_site)
1938
1939        new_structure = self.__class__.from_sites(new_sites)
1940
1941        return dict2namedtuple(structure=new_structure, displ=np.array(displ_list), scale_matrix=scale_matrix)
1942
1943    def calc_kptbounds(self):
1944        """Returns the suggested value for the ABINIT variable ``kptbounds``."""
1945        kptbounds = [k.frac_coords for k in self.hsym_kpoints]
1946        return np.reshape(kptbounds, (-1, 3))
1947
1948    def get_kpath_input_string(self, fmt="abinit", line_density=10):
1949        """
1950        Return string with input variables for band-structure calculations
1951        in the format used by code `fmt`.
1952        Use `line_density` points for the smallest segment (if supported by code).
1953        """
1954        lines = []; app = lines.append
1955        if fmt in ("abinit", "abivars"):
1956            app("# Abinit Structure")
1957            app(self.convert(fmt=fmt))
1958            app("\n# tolwfr 1e-20 iscf -2 # NSCF run")
1959            app('# To read previous DEN file, use: getden -1 or specify filename via getden_path "out_DEN"')
1960            app("\n# K-path in reduced coordinates:")
1961            app(" ndivsm %d" % line_density)
1962            app(" kptopt %d" % -(len(self.hsym_kpoints) - 1))
1963            app(" kptbounds")
1964            for k in self.hsym_kpoints:
1965                app("    {:+.5f}  {:+.5f}  {:+.5f}  # {kname}".format(*k.frac_coords, kname=k.name))
1966
1967        elif fmt in ("wannier90", "w90"):
1968            app("# Wannier90 structure")
1969            from abipy.wannier90.win import Wannier90Input
1970            win = Wannier90Input(self)
1971            win.set_kpath()
1972            app(win.to_string())
1973
1974        elif fmt == "siesta":
1975            app("# Siesta structure")
1976            app(self.convert(fmt=fmt))
1977            # Build normalized k-path.
1978            from .kpoints import Kpath
1979            vertices_names = [(k.frac_coords, k.name) for k in self.hsym_kpoints]
1980            kpath = Kpath.from_vertices_and_names(self, vertices_names, line_density=line_density)
1981            app("%block BandLines")
1982            prev_ik = 0
1983            for ik, k in enumerate(kpath):
1984                if not k.name: continue
1985                n = ik - prev_ik
1986                app("{}  {:+.5f}  {:+.5f}  {:+.5f}  # {kname}".format(n if n else 1, *k.frac_coords, kname=k.name))
1987                prev_ik = ik
1988            app("%endblock BandLines")
1989
1990        else:
1991            raise ValueError("Don't know how to generate string for code: `%s`" % str(fmt))
1992
1993        return "\n".join(lines)
1994
1995    def calc_ksampling(self, nksmall, symprec=0.01, angle_tolerance=5):
1996        """
1997        Return the k-point sampling from the number of divisions ``nksmall`` to be used for
1998        the smallest reciprocal lattice vector.
1999        """
2000        ngkpt = self.calc_ngkpt(nksmall)
2001        shiftk = self.calc_shiftk(symprec=symprec, angle_tolerance=angle_tolerance)
2002
2003        return AttrDict(ngkpt=ngkpt, shiftk=shiftk)
2004
2005    def calc_ngkpt(self, nksmall):
2006        """
2007        Compute the ABINIT variable ``ngkpt`` from the number of divisions used
2008        for the smallest lattice vector.
2009
2010        Args:
2011            nksmall (int): Number of division for the smallest lattice vector.
2012        """
2013        lengths = self.lattice.reciprocal_lattice.abc
2014        lmin = np.min(lengths)
2015
2016        ngkpt = np.ones(3, dtype=int)
2017        for i in range(3):
2018            ngkpt[i] = int(round(nksmall * lengths[i] / lmin))
2019            if ngkpt[i] == 0:
2020                ngkpt[i] = 1
2021
2022        return ngkpt
2023
2024    def calc_shiftk(self, symprec=0.01, angle_tolerance=5):
2025        """
2026        Find the values of ``shiftk`` and ``nshiftk`` appropriated for the sampling of the Brillouin zone.
2027
2028        When the primitive vectors of the lattice do NOT form a FCC or a BCC lattice,
2029        the usual (shifted) Monkhorst-Pack grids are formed by using nshiftk=1 and shiftk 0.5 0.5 0.5 .
2030        This is often the preferred k point sampling. For a non-shifted Monkhorst-Pack grid,
2031        use `nshiftk=1` and `shiftk 0.0 0.0 0.0`, but there is little reason to do that.
2032
2033        When the primitive vectors of the lattice form a FCC lattice, with rprim::
2034
2035                0.0 0.5 0.5
2036                0.5 0.0 0.5
2037                0.5 0.5 0.0
2038
2039        the (very efficient) usual Monkhorst-Pack sampling will be generated by using nshiftk= 4 and shiftk::
2040
2041            0.5 0.5 0.5
2042            0.5 0.0 0.0
2043            0.0 0.5 0.0
2044            0.0 0.0 0.5
2045
2046        When the primitive vectors of the lattice form a BCC lattice, with rprim::
2047
2048               -0.5  0.5  0.5
2049                0.5 -0.5  0.5
2050                0.5  0.5 -0.5
2051
2052        the usual Monkhorst-Pack sampling will be generated by using nshiftk= 2 and shiftk::
2053
2054                0.25  0.25  0.25
2055               -0.25 -0.25 -0.25
2056
2057        However, the simple sampling nshiftk=1 and shiftk 0.5 0.5 0.5 is excellent.
2058
2059        For hexagonal lattices with hexagonal axes, e.g. rprim::
2060
2061                1.0  0.0       0.0
2062               -0.5  sqrt(3)/2 0.0
2063                0.0  0.0       1.0
2064
2065        one can use nshiftk= 1 and shiftk 0.0 0.0 0.5
2066        In rhombohedral axes, e.g. using angdeg 3*60., this corresponds to shiftk 0.5 0.5 0.5,
2067        to keep the shift along the symmetry axis.
2068
2069        Returns:
2070            Suggested value of shiftk.
2071        """
2072        # Find lattice type.
2073        sym = SpacegroupAnalyzer(self, symprec=symprec, angle_tolerance=angle_tolerance)
2074        lattice_type, spg_symbol = sym.get_lattice_type(), sym.get_space_group_symbol()
2075
2076        # Check if the cell is primitive
2077        is_primitive = len(sym.find_primitive()) == len(self)
2078
2079        # Generate the appropriate set of shifts.
2080        shiftk = None
2081
2082        if is_primitive:
2083            if lattice_type == "cubic":
2084                if "F" in spg_symbol:
2085                    # FCC
2086                    shiftk = [0.5, 0.5, 0.5,
2087                              0.5, 0.0, 0.0,
2088                              0.0, 0.5, 0.0,
2089                              0.0, 0.0, 0.5]
2090
2091                elif "I" in spg_symbol:
2092                    # BCC
2093                    shiftk = [0.25,  0.25,  0.25,
2094                             -0.25, -0.25, -0.25]
2095                    #shiftk = [0.5, 0.5, 05])
2096
2097            elif lattice_type == "hexagonal":
2098                # Find the hexagonal axis and set the shift along it.
2099                for i, angle in enumerate(self.lattice.angles):
2100                    if abs(angle - 120) < 1.0:
2101                        j = (i + 1) % 3
2102                        k = (i + 2) % 3
2103                        hex_ax = [ax for ax in range(3) if ax not in [j,k]][0]
2104                        break
2105                else:
2106                    raise ValueError("Cannot find hexagonal axis")
2107
2108                shiftk = [0.0, 0.0, 0.0]
2109                shiftk[hex_ax] = 0.5
2110
2111            elif lattice_type == "tetragonal":
2112                if "I" in spg_symbol:
2113                    # BCT
2114                    shiftk = [0.25,  0.25,  0.25,
2115                             -0.25, -0.25, -0.25]
2116
2117        if shiftk is None:
2118            # Use default value.
2119            shiftk = [0.5, 0.5, 0.5]
2120
2121        return np.reshape(shiftk, (-1, 3))
2122
2123    def num_valence_electrons(self, pseudos):
2124        """
2125        Returns the number of valence electrons.
2126
2127        Args:
2128            pseudos: List of |Pseudo| objects or list of filenames.
2129        """
2130        nval, table = 0, PseudoTable.as_table(pseudos)
2131        for site in self:
2132            pseudo = table.pseudo_with_symbol(site.specie.symbol)
2133            nval += pseudo.Z_val
2134
2135        return int(nval) if int(nval) == nval else nval
2136
2137    def valence_electrons_per_atom(self, pseudos):
2138        """
2139        Returns the number of valence electrons for each atom in the structure.
2140
2141        Args:
2142            pseudos: List of |Pseudo| objects or list of filenames.
2143        """
2144        table = PseudoTable.as_table(pseudos)
2145        psp_valences = []
2146        for site in self:
2147            pseudo = table.pseudo_with_symbol(site.specie.symbol)
2148            psp_valences.append(pseudo.Z_val)
2149
2150        return psp_valences
2151
2152    def write_notebook(self, nbpath=None):
2153        """
2154        Write a jupyter_ notebook to ``nbpath``. If nbpath is None, a temporay file in the current
2155        working directory is created. Return path to the notebook.
2156        """
2157        nbformat, nbv, nb = self.get_nbformat_nbv_nb(title=None)
2158
2159        # Use pickle files for data persistence.
2160        # The notebook will reconstruct the object from this file
2161        _, tmpfile = tempfile.mkstemp(suffix='.pickle')
2162        with open(tmpfile, "wb") as fh:
2163            pickle.dump(self, fh)
2164
2165        nb.cells.extend([
2166            #nbv.new_markdown_cell("# This is a markdown cell"),
2167            nbv.new_code_cell("structure = abilab.Structure.from_file('%s')" % tmpfile),
2168            nbv.new_code_cell("print(structure)"),
2169            nbv.new_code_cell("print(structure.abi_string)"),
2170            nbv.new_code_cell("structure"),
2171            nbv.new_code_cell("print(structure.spget_summary())"),
2172            nbv.new_code_cell("if structure.abi_spacegroup is not None: print(structure.abi_spacegroup)"),
2173            nbv.new_code_cell("print(structure.hsym_kpoints)"),
2174            nbv.new_code_cell("structure.plot_bz();"),
2175            nbv.new_code_cell("#import panel as pn; pn.extension()\n#structure.get_panel()"),
2176            nbv.new_code_cell("# sanitized = structure.abi_sanitize(); print(sanitized)"),
2177            nbv.new_code_cell("# ase_atoms = structure.to_ase_atoms()"),
2178            nbv.new_code_cell("# structure.plot_atoms();"),
2179            nbv.new_code_cell("# jsmol_view = structure.get_jsmol_view(); jsmol_view"),
2180            nbv.new_code_cell("# ngl_view = structure.get_ngl_view(); ngl_view"),
2181            nbv.new_code_cell("# ctk_view = structure.get_crystaltk_view(); ctk_view"),
2182        ])
2183
2184        return self._write_nb_nbpath(nb, nbpath)
2185
2186
2187def dataframes_from_structures(struct_objects, index=None, symprec=1e-2, angle_tolerance=5,
2188                               with_spglib=True, cart_coords=False):
2189    """
2190    Build two pandas Dataframes_ with the most important geometrical parameters associated to
2191    a list of structures or a list of objects that can be converted into structures.
2192
2193    Args:
2194        struct_objects: List of objects that can be converted to structure.
2195            Support filenames, structure objects, Abinit input files, dicts and many more types.
2196            See ``Structure.as_structure`` for the complete list.
2197        index: Index of the |pandas-DataFrame|.
2198        symprec (float): Symmetry precision used to refine the structure.
2199        angle_tolerance (float): Tolerance on angles.
2200        with_spglib (bool): If True, spglib_ is invoked to get the spacegroup symbol and number.
2201        cart_coords: True if the ``coords`` dataframe should contain Cartesian cordinates
2202            instead of Reduced coordinates.
2203
2204    Return:
2205        namedtuple with two |pandas-DataFrames| named ``lattice`` and ``coords``
2206        ``lattice`` contains the lattice parameters. ``coords`` the atomic positions..
2207        The list of structures is available in the ``structures`` entry.
2208
2209    .. code-block:: python
2210
2211        dfs = dataframes_from_structures(files)
2212        dfs.lattice
2213        dfs.coords
2214        for structure in dfs.structures:
2215            print(structure)
2216    """
2217    structures = [Structure.as_structure(obj) for obj in struct_objects]
2218    # Build Frame with lattice parameters.
2219    # Use OrderedDict to have columns ordered nicely.
2220    odict_list = [(structure.get_dict4pandas(with_spglib=with_spglib, symprec=symprec,
2221                                             angle_tolerance=angle_tolerance)) for structure in structures]
2222
2223    import pandas as pd
2224    lattice_frame = pd.DataFrame(odict_list, index=index,
2225                                 columns=list(odict_list[0].keys()) if odict_list else None)
2226
2227    # Build Frame with atomic positions.
2228    vtos = lambda v: "%+0.6f %+0.6f %+0.6f" % (v[0], v[1], v[2])
2229    max_numsite = max(len(s) for s in structures)
2230    odict_list = []
2231    for structure in structures:
2232        if cart_coords:
2233            odict_list.append({i: (site.species_string, vtos(site.coords)) for i, site in enumerate(structure)})
2234        else:
2235            odict_list.append({i: (site.species_string, vtos(site.frac_coords)) for i, site in enumerate(structure)})
2236
2237    coords_frame = pd.DataFrame(odict_list, index=index,
2238                                columns=list(range(max_numsite)) if odict_list else None)
2239
2240    return dict2namedtuple(lattice=lattice_frame, coords=coords_frame, structures=structures)
2241
2242
2243class StructureModifier(object):
2244    """
2245    This object provides an easy-to-use interface for
2246    generating new structures according to some algorithm.
2247
2248    The main advantages of this approach are:
2249
2250        *) Client code does not have to worry about the fact
2251           that many methods of Structure modify the object in place.
2252
2253        *) One can render the interface more user-friendly. For example
2254           some arguments might have a unit that can be specified in input.
2255           For example one can pass a length in Bohr that will be automatically
2256           converted into Angstrom before calling the pymatgen methods
2257    """
2258    def __init__(self, structure):
2259        """
2260        Args:
2261            structure: Structure object.
2262        """
2263        # Get a copy to avoid any modification of the input.
2264        self._original_structure = structure.copy()
2265
2266    def copy_structure(self):
2267        """Returns a copy of the original structure."""
2268        return self._original_structure.copy()
2269
2270    def scale_lattice(self, vol_ratios):
2271        """
2272        Scale the lattice vectors so that length proportions and angles are preserved.
2273
2274        Args:
2275            vol_ratios: List with the ratios v/v0 where v0 is the volume of the original structure.
2276
2277        Return: List of new structures with desired volume.
2278        """
2279        vol_ratios = np.array(vol_ratios)
2280        new_volumes = self._original_structure.volume * vol_ratios
2281
2282        news = []
2283        for vol in new_volumes:
2284            new_structure = self.copy_structure()
2285            new_structure.scale_lattice(vol)
2286            news.append(new_structure)
2287
2288        return news
2289
2290    def make_supercell(self, scaling_matrix):
2291        """
2292        Create a supercell.
2293
2294        Args:
2295            scaling_matrix: A scaling matrix for transforming the lattice vectors.
2296                Has to be all integers. Several options are possible:
2297
2298                a. A full 3x3 scaling matrix defining the linear combination of the old lattice vectors.
2299                    E.g., [[2,1,0],[0,3,0],[0,0,1]] generates a new structure with lattice vectors
2300                    a' = 2a + b, b' = 3b, c' = c
2301                    where a, b, and c are the lattice vectors of the original structure.
2302                b. A sequence of three scaling factors. e.g., [2, 1, 1]
2303                   specifies that the supercell should have dimensions 2a x b x c.
2304                c. A number, which simply scales all lattice vectors by the same factor.
2305
2306        Returns:
2307            New structure.
2308        """
2309        new_structure = self.copy_structure()
2310        new_structure.make_supercell(scaling_matrix)
2311        return new_structure
2312
2313    def displace(self, displ, etas, frac_coords=True):
2314        """
2315        Displace the sites of the structure along the displacement vector displ.
2316
2317        The displacement vector is first rescaled so that the maxium atomic displacement
2318        is one Angstrom, and then multiplied by eta. Hence passing eta=0.001, will move
2319        all the atoms so that the maximum atomic displacement is 0.001 Angstrom.
2320
2321        Args:
2322            displ: Displacement vector with 3*len(self) entries (fractional coordinates).
2323            eta: Scaling factor.
2324            frac_coords: Boolean stating whether the vector corresponds to fractional or cartesian coordinates.
2325
2326        Returns:
2327            List of new structures with displaced atoms.
2328        """
2329        if not isinstance(etas, collections.abc.Iterable):
2330            etas = [etas]
2331
2332        news = []
2333        for eta in etas:
2334            new_structure = self.copy_structure()
2335            new_structure.displace(displ, eta, frac_coords=frac_coords)
2336            news.append(new_structure)
2337
2338        return news
2339
2340    def frozen_phonon(self, qpoint, displ, eta=1, frac_coords=False, scale_matrix=None, max_supercell=None):
2341
2342        return self._original_structure.frozen_phonon(qpoint, displ, eta, frac_coords, scale_matrix, max_supercell)
2343
2344    def frozen_2phonon(self, qpoint, displ1, displ2, eta=1, frac_coords=False, scale_matrix=None, max_supercell=None):
2345
2346        return self._original_structure.frozen_2phonon(qpoint, displ1, displ2, eta, frac_coords, scale_matrix,
2347                                                       max_supercell)
2348
2349
2350def diff_structures(structures, fmt="cif", mode="table", headers=(), file=sys.stdout):
2351    """
2352    Convert list of structure to string using format `fmt`, print diff to file `file`.
2353
2354    Args:
2355        structures: List of structures or list of objects that can be converted into structure e.g. filepaths
2356        fmt: Any output format supported by `structure.to` method. Non-case sensitive.
2357        mode: `table` to show results in tabular form or `diff` to show differences with unified diff.
2358        headers: can be an explicit list of column headers Otherwise a headerless table is produced
2359        file: Output Stream
2360    """
2361    outs = [s.convert(fmt=fmt).splitlines() for s in map(Structure.as_structure, structures)]
2362
2363    if mode == "table":
2364        try:
2365            from itertools import izip_longest as zip_longest  # Py2
2366        except ImportError:
2367            from itertools import zip_longest  # Py3k
2368
2369        table = [r for r in zip_longest(*outs, fillvalue=" ")]
2370        from tabulate import tabulate
2371        print(tabulate(table, headers=headers), file=file)
2372
2373    elif mode == "diff":
2374        import difflib
2375        fromfile, tofile = "", ""
2376        for i in range(1, len(outs)):
2377            if headers: fromfile, tofile = headers[0], headers[i]
2378            diff = "\n".join(difflib.unified_diff(outs[0], outs[i], fromfile=fromfile, tofile=tofile))
2379            print(diff, file=file)
2380
2381    else:
2382        raise ValueError("Unsupported mode: `%s`" % str(mode))
2383
2384
2385def structure2siesta(structure, verbose=0):
2386    """
2387    Return string with structural information in Siesta format from pymatgen structure
2388
2389    Args:
2390        structure: AbiPy structure.
2391        verbose: Verbosity level.
2392    """
2393
2394    if not structure.is_ordered:
2395        raise NotImplementedError("""\
2396Received disordered structure with partial occupancies that cannot be converted into a Siesta input
2397Please use OrderDisorderedStructureTransformation or EnumerateStructureTransformation
2398to build an appropriate supercell from partial occupancies or alternatively use the Virtual Crystal Approximation.""")
2399
2400    species_by_znucl = structure.species_by_znucl
2401
2402    lines = []
2403    app = lines.append
2404    app("NumberOfAtoms %d" % len(structure))
2405    app("NumberOfSpecies %d" % structure.ntypesp)
2406
2407    if verbose:
2408        app("# The species number followed by the atomic number, and then by the desired label")
2409    app("%block ChemicalSpeciesLabel")
2410    for itype, specie in enumerate(species_by_znucl):
2411        app("    %d %d %s" % (itype + 1, specie.number, specie.symbol))
2412    app("%endblock ChemicalSpeciesLabel")
2413
2414    # Write lattice vectors.
2415    # Set small values to zero. This usually happens when the CIF file
2416    # does not give structure parameters with enough digits.
2417    lvectors = np.where(np.abs(structure.lattice.matrix) > 1e-8, structure.lattice.matrix, 0.0)
2418    app("LatticeConstant 1.0 Ang")
2419    app("%block LatticeVectors")
2420    for r in lvectors:
2421        app("    %.10f %.10f %.10f" % (r[0], r[1], r[2]))
2422    app("%endblock LatticeVectors")
2423
2424    # Write atomic coordinates
2425    #% block AtomicCoordinatesAndAtomicSpecies
2426    #4.5000 5.0000 5.0000 1
2427    #5.5000 5.0000 5.0000 1
2428    #% endblock AtomicCoordinatesAndAtomicSpecies
2429    app("AtomicCoordinatesFormat Fractional")
2430    app("%block AtomicCoordinatesAndAtomicSpecies")
2431    for i, site in enumerate(structure):
2432        itype = species_by_znucl.index(site.specie)
2433        fc = np.where(np.abs(site.frac_coords) > 1e-8, site.frac_coords, 0.0)
2434        app("    %.10f %.10f %.10f %d" % (fc[0], fc[1], fc[2], itype + 1))
2435    app("%endblock AtomicCoordinatesAndAtomicSpecies")
2436
2437    return "\n".join(lines)
2438