1#!/usr/local/bin/python3.8
2"""
3Script to analyze/export/visualize the crystal structure saved in the netcdf files produced by ABINIT.
4"""
5import sys
6import os
7import argparse
8import numpy as np
9
10from pprint import pprint
11from tabulate import tabulate
12from monty.string import marquee
13from monty.functools import prof_main
14from monty.termcolor import cprint
15from pymatgen.io.vasp.outputs import Xdatcar
16from abipy import abilab
17from abipy.core.symmetries import AbinitSpaceGroup
18from abipy.core.kpoints import Ktables, Kpoint, IrredZone
19from abipy.core.structure import diff_structures
20from abipy.iotools.visualizer import Visualizer
21from abipy.iotools.xsf import xsf_write_structure
22from abipy.abio import factories
23
24
25def save_structure(structure, options):
26    """Save structure to file."""
27    if not options.savefile: return
28    print("Saving structure to %s" % options.savefile)
29    if os.path.exists(options.savefile):
30        backup = options.savefile + ".bkp"
31        print("%s already exists. Saving backup copy to: %s" % (options.savefile, backup))
32        os.rename(options.savefile, backup)
33
34    structure.to(filename=options.savefile)
35
36
37def check_ordered_structure(structure):
38    """Print a warning and sys.exit 1 if structure is disordered."""
39    if not structure.is_ordered:
40        cprint("""
41Cannot handle disordered structure with fractional site occupancies.
42Use OrderDisorderedStructureTransformation or EnumerateStructureTransformation
43to build an appropriate supercell from partial occupancies.""", color="magenta")
44        sys.exit(1)
45
46
47def get_epilog():
48    return """\
49Usage example:
50
51###################
52# Space group tools
53###################
54
55  abistruct.py spglib FILE                 => Read structure from FILE and analyze it with spglib.
56  abistruct.py abispg FILE                 => Read structure from FILE, and compute ABINIT space group.
57  abistruct.py primitive FILE              => Read structure from FILE, use pymatgen and spglib to find primitive structure.
58  abistruct.py abisanitize FILE            => Read structure from FILE, call abisanitize, compare structures
59                                              and save "abisanitized" structure to file.
60  abistruct.py conventional FILE           => Read structure from FILE, generate conventional structure
61                                              following doi:10.1016/j.commatsci.2010.05.010
62  abistruct.py proto FILE                 => Read structure from FILE, find possible crystallographic prototypes:
63                                             in the AFLOW library of crystallographic prototypes.
64                                             http://doi.org/10.1016/j.commatsci.2017.01.017
65
66##################
67# Conversion tools
68##################
69
70  abistruct.py convert FILE                => Print the ABINIT variables defining the structure.
71  abistruct.py convert FILE -f cif         => Read structure from FILE and output CIF file
72                                              (Use convert --help to get list of formats supported)
73  abistruct.py convert out_HIST.nc         => Read FINAL structure from the HIST file and
74                                              print the corresponding ABINIT variables.
75  abistruct.py supercell FILE -s 2 2 1     => Read structure from FILE and build [2, 2, 1] supercell,
76                                              print new structure using --format (default abivars).
77################
78# K-points tools
79################
80
81  abistruct.py kpath FILE                  => Read structure from FILE and print ABINIT variables for k-path.
82  abistruct.py bz FILE                     => Read structure from FILE, plot BZ with matplotlib.
83  abistruct.py ngkpt FILE -n 4             => Compute `ngkpt` and `shiftk` from the number of divisions used to sample
84                                              the smallest reciprocal lattice vector.
85  abistruct.py abikmesh FILE --ngkpt 2 2 2 => Read structure from FILE, call Abinit to sample the BZ
86                                              with a 2, 2, 2 k-mesh, print points in IBZ with weights.
87  abistruct.py ktables FILE -m 2 2 2       => Read structure from FILE, call spglib to sample the BZ
88                                              with a 2, 2, 2 k-mesh, print points in IBZ with weights.
89  abistruct.py lgk FILE -k 0.25 0 0        => Read structure from FILE, find little group of k-point,
90                                              print Bilbao character table.
91  abistruct.py kstar FILE -k 0.25 0 0      => Read structure from FILE, print star of k-point.
92  abistruct.py keq FILE -k 0.5 0 0 0 0.5 0  => Read structure from FILE, test whether k1 and k2 are
93                                               symmetry-equivalent k-points.
94
95###############
96# Miscelleanous
97###############
98
99  abistruct.py neighbors FILE              => Get neighbors for each atom in the unit cell, out to a distance radius.
100  abistruct.py interpolate FILE1 FILE2     => Interpolate between two structures. Useful for the construction of NEB inputs.
101  abistruct.py xrd FILE                    => X-ray diffraction plot.
102  abistruct.py oxistate FILE               => Estimate oxidation states with pymatgen bond valence analysis.
103
104###############
105# Visualization
106###############
107
108  abistruct.py visualize FILE -a vesta     => Visualize the structure with vesta (default if -a is not given)
109                                              Supports also ovito, xcrysden, vtk, mayavi, matplotlib See --help
110  abistruct.py ipython FILE                => Read structure from FILE and open it in the Ipython terminal.
111  abistruct.py notebook FILE               => Read structure from FILE and generate jupyter notebook.
112  abistruct.py panel FILE                  => Generate GUI in web browser to interact with the structure
113                                              Requires panel package.
114
115###########
116# Databases
117###########
118
119  abistruct.py cod_id 1526507              => Get structure from COD database (http://www.crystallography.net/cod).
120  abistruct.py cod_search MgB2             => Search for structures in the COD database.
121  abistruct.py mp_id mp-149                => Get structure from materials project database and print
122                                              JSON representation. Use e.g. `-f abivars` to change format.
123  abistruct.py mp_match FILE               => Read structure from FILE and find matching structures on the
124                                              Materials Project site. Use e.g. `-f cif` to change output format.
125  abistruct.py mp_search LiF               => Connect to the materials project database. Get structures corresponding
126                                              to a chemical system or formula e.g. `Fe2O3` or `Li-Fe-O` or
127                                              `Ir-O-*` for wildcard pattern matching.
128                                              Print info and Abinit input files. Use e.g. `-f POSCAR`
129                                              to change output format. `-f None` to disable structure output.
130  abistruct.py mp_pd FILE-or-elements      => Generate phase diagram with entries from the Materials Project.
131  abistruct.py mp_ebands FILE             => Fetch electron band structure from MP database. Print gaps.
132                                              Accept FILE with structure if ebands fro structure is wanted
133                                              or mp id e.g. "mp-149 or list of elements e.g `Li-Fe-O` or chemical formula.
134
135`FILE` is any file supported by abipy/pymatgen e.g Netcdf files, Abinit input/output, POSCAR, xsf ...
136Use `abistruct.py --help` for help and `abistruct.py COMMAND --help` to get the documentation for `COMMAND`.
137Use `-v` to increase verbosity level (can be supplied multiple times e.g -vv).
138"""
139
140
141def get_parser(with_epilog=False):
142
143    # Parent parser for commands that need to know the filepath
144    path_selector = argparse.ArgumentParser(add_help=False)
145    path_selector.add_argument('filepath', nargs="?",
146        help="File with the crystalline structure (Abinit Netcdf files, CIF, Abinit input/output files, POSCAR ...)")
147
148    # Parent parser for commands supporting (jupyter notebooks)
149    nb_parser = argparse.ArgumentParser(add_help=False)
150    nb_parser.add_argument('-nb', '--notebook', default=False, action="store_true", help='Generate jupyter notebook.')
151    nb_parser.add_argument('--classic-notebook', action='store_true', default=False,
152                           help="Use classic notebook instead of jupyterlab.")
153    nb_parser.add_argument('--no-browser', action='store_true', default=False,
154                           help=("Start the jupyter server to serve the notebook "
155                                 "but don't open the notebook in the browser.\n"
156                                 "Use this option to connect remotely from localhost to the machine running the kernel"))
157    nb_parser.add_argument('--foreground', action='store_true', default=False,
158        help="Run jupyter notebook in the foreground.")
159
160    parser = argparse.ArgumentParser(epilog=get_epilog() if with_epilog else "",
161                                     formatter_class=argparse.RawDescriptionHelpFormatter)
162    parser.add_argument('-V', '--version', action='version', version=abilab.__version__)
163
164    # Parser for commands that need to call spglib.
165    spgopt_parser = argparse.ArgumentParser(add_help=False)
166    spgopt_parser.add_argument('--symprec', default=1e-3, type=float,
167        help="""\
168symprec (float): Tolerance for symmetry finding. Defaults to 1e-3,
169which is fairly strict and works well for properly refined structures with atoms in the proper symmetry coordinates.
170For structures with slight deviations from their proper atomic positions (e.g., structures relaxed with electronic structure
171codes), a looser tolerance of 0.1 (the value used in Materials Project) is often needed.""")
172    spgopt_parser.add_argument('--angle-tolerance', default=5.0, type=float,
173        help="angle_tolerance (float): Angle tolerance for symmetry finding. Default: 5.0")
174    spgopt_parser.add_argument("--no-time-reversal", default=False, action="store_true", help="Don't use time-reversal.")
175    spgopt_parser.add_argument("--site-symmetry", default=False, action="store_true",
176                               help="Show site symmetries i.e. the point group operations that leave the site invariant.")
177
178    # Parent parser for common options.
179    copts_parser = argparse.ArgumentParser(add_help=False)
180    copts_parser.add_argument('-v', '--verbose', default=0, action='count', # -vv --> verbose=2
181        help='verbose, can be supplied multiple times to increase verbosity')
182    copts_parser.add_argument('--loglevel', default="ERROR", type=str,
183        help="Set the loglevel. Possible values: CRITICAL, ERROR (default), WARNING, INFO, DEBUG")
184
185    # Parent parser for commands that need to save structure to file.
186    savefile_parser = argparse.ArgumentParser(add_help=False)
187    savefile_parser.add_argument("-s", "--savefile", default="", type=str,
188        help="Save final structure to file. Format is detected from file extensions "
189             "e.g. out.abi for Abinit input, out.cif for CIF format.")
190
191    # Helper functions to construct sub-parsers.
192    def add_primitive_options(parser):
193        """Add --no-primitive and --primitive-standard options to a parser."""
194        group = parser.add_mutually_exclusive_group()
195        group.add_argument('--no-primitive', default=False, action='store_true', help="Do not enforce primitive cell.")
196        group.add_argument('--primitive-standard', default=False, action='store_true',
197            help="Enforce primitive standard cell.")
198
199    supported_formats = "(abivars, cif, xsf, poscar, qe, siesta, wannier90, cssr, json, None)"
200
201    def add_format_arg(parser, default, option=True, formats=None):
202        """Add --format option to a parser with default value `default`."""
203        formats = supported_formats if formats is None else formats
204        if option:
205            parser.add_argument("-f", "--format", default=default, type=str,
206                help="Output format. Default: %s. Accept: %s" % (default, formats))
207        else:
208            parser.add_argument('format', nargs="?", default=default, type=str,
209                help="Output format. Default: %s. Accept: %s" % (default, formats))
210
211    # Create the parsers for the sub-commands
212    subparsers = parser.add_subparsers(dest='command', help='sub-command help',
213        description="Valid subcommands, use command --help for help")
214
215    # Subparser for spglib command.
216    p_spglib = subparsers.add_parser('spglib', parents=[copts_parser, path_selector, spgopt_parser],
217        help="Analyze structure with spglib.")
218
219    # Subparser for abispg command.
220    p_abispg = subparsers.add_parser('abispg', parents=[copts_parser, path_selector, savefile_parser],
221        help="Extract/Compute Abinit space group from file with structure.")
222    p_abispg.add_argument("-t", "--tolsym", type=float, default=None, help="""\
223Gives the tolerance on the atomic positions (reduced coordinates), primitive vectors, or magnetization,
224to be considered equivalent, thanks to symmetry operations. This value is used by ABINIT in the recognition of the set
225of symmetries of the system, or the application of the symmetry operations to generate from a reduced set of atoms,
226The internal default is 1e-8. Setting tolsym to a value larger than 1e-8 will make Abinit detect the spacegroup within
227this tolerance and re-symmetrize the input structure. This option is useful if the structure has been taken from a CIF
228file that does not have enough significant digits.""")
229    p_abispg.add_argument("-d", "--diff-mode", type=str, default="table", choices=["table", "diff"],
230        help="Select diff output format.")
231
232    # Subparser for convert command.
233    p_convert = subparsers.add_parser('convert', parents=[copts_parser, path_selector],
234        help="Convert structure to the specified format.")
235    add_format_arg(p_convert, default="cif")
236
237    p_print = subparsers.add_parser('print', parents=[copts_parser, path_selector],
238                                    help="Print Structure to terminal.")
239
240    # Subparser for supercell command.
241    p_supercell = subparsers.add_parser('supercell', parents=[copts_parser, path_selector],
242        help="Generate supercell.")
243    p_supercell.add_argument("-s", "--scaling_matrix", nargs="+", required=True, type=int,
244        help="""\
245scaling_matrix: A scaling matrix for transforming the lattice vectors.
246Has to be all integers. Several options are possible:
247
248    a. A full 3x3 scaling matrix defining the linear combination
249       the old lattice vectors. E.g., -s 2,1,0 0,3,0, 0,0,1 generates a new structure with lattice vectors
250       a' = 2a + b, b' = 3b, c' = c where a, b, and c are the lattice vectors of the original structure.
251    b. An sequence of three scaling factors. E.g., -s 2, 1, 1 specifies that the supercell should
252       have dimensions 2a x b x c.
253    c. A number, which simply scales all lattice vectors by the same factor.
254    """)
255    add_format_arg(p_supercell, default="abivars")
256
257    # Subparser for abisanitize
258    p_abisanitize = subparsers.add_parser('abisanitize', parents=[copts_parser, path_selector, spgopt_parser, savefile_parser],
259        help="Sanitize structure with abi_sanitize, compare structures and save result to file.")
260    add_primitive_options(p_abisanitize)
261
262    p_primitive = subparsers.add_parser('primitive', parents=[copts_parser, path_selector, spgopt_parser, savefile_parser],
263        help="Use spglib to find a smaller unit cell than the input")
264
265    # Subparser for irefine
266    p_irefine = subparsers.add_parser('irefine', parents=[copts_parser, path_selector, spgopt_parser],
267        help="Refine structure with abi_sanitize iteratively, stop if target space group is obtained.")
268    p_irefine.add_argument("--target-spgnum", required=True, type=int, help="Target space group number.")
269    p_irefine.add_argument("--symprec-step", default=0.05, type=float, help='Increment for symprec.')
270    p_irefine.add_argument("--angle-tolerance-step", default=0.0, type=float, help='Increment for angle_tolerance.')
271    p_irefine.add_argument("--ntrial", default=50, type=int, help='Number of trials. Default 50.')
272    add_primitive_options(p_irefine)
273
274    # Subparser for conventional.
275    p_conventional = subparsers.add_parser('conventional', parents=[copts_parser, path_selector, spgopt_parser, savefile_parser],
276        help="Gives a structure with a conventional cell according to certain standards. "
277             "The standards are defined in doi:10.1016/j.commatsci.2010.05.010")
278
279    # Subparser for proto.
280    p_proto = subparsers.add_parser('proto', parents=[copts_parser, path_selector],
281        help=("Find prototype in the AFLOW LIBRARY OF CRYSTALLOGRAPHIC PROTOTYPES. "
282              "http://doi.org/10.1016/j.commatsci.2017.01.017"))
283    p_proto.add_argument("--ltol", default=0.2, type=float, help="fractional length tolerance.")
284    p_proto.add_argument("--stol", default=0.3, type=float, help="site tolerance.")
285    p_proto.add_argument("--angle-tol", default=5, type=float, help="angle tolerance.")
286
287    # Subparser for wyckoff.
288    p_wyckoff = subparsers.add_parser('wyckoff', parents=[copts_parser, spgopt_parser, path_selector],
289            help="Print wyckoff positions. WARNING: still under development!")
290    p_wyckoff.add_argument("--refine", default=False, action="store_true",
291                           help="Use spglib to refine structure before computation")
292
293    # Subparser for tensor_site.
294    p_tensor_site = subparsers.add_parser('tensor_site', parents=[copts_parser, spgopt_parser, path_selector],
295            help="Print symmetry properties of tensors due to site-symmetries. WARNING: still under development!")
296    p_tensor_site.add_argument("--refine", default=False, action="store_true",
297                               help="Use spglib to refine structure before computation")
298
299    # Subparser for neighbors.
300    p_neighbors = subparsers.add_parser('neighbors', parents=[copts_parser, path_selector],
301                                        help="Get neighbors for each atom in the unit cell, out to a distance radius.")
302    p_neighbors.add_argument("-r", "--radius", default=2, type=float, help="Radius of the sphere in Angstrom.")
303
304    # Subparser for interpolate.
305    p_interpolate = subparsers.add_parser('interpolate', parents=[copts_parser],
306        help=("Interpolate between two structures. Useful for the construction of NEB inputs."))
307    p_interpolate.add_argument("filepaths", nargs=2, help="Files with initial and final structures.")
308    p_interpolate.add_argument("-n", "--nimages", default=10, type=int, help="No. of interpolation images. Defaults to 10.")
309    p_interpolate.add_argument("--autosort_tol", default=0.5, type=float, help="""\
310A distance tolerance in Angstrom in which to automatically sort end_structure to match to the
311closest points in this particular structure. This is usually what you want in a NEB calculation.
3120 implies no sorting. Otherwise, a 0.5 value (default) usually works pretty well.""")
313    add_format_arg(p_interpolate, default="abivars")
314
315    # Subparser for xrd.
316    p_xrd = subparsers.add_parser('xrd', parents=[copts_parser, path_selector], help="X-ray diffraction plot.")
317    p_xrd.add_argument("-w", "--wavelength", default="CuKa", type=str, help=(
318        "The wavelength can be specified as a string. It must be one of the "
319        "supported definitions in the WAVELENGTHS dict declared in pymatgen/analysis/diffraction/xrd.py."
320        "Defaults to 'CuKa', i.e, Cu K_alpha radiation."))
321    p_xrd.add_argument("-s", "--symprec", default=0, type=float, help=(
322        "Symmetry precision for structure refinement. "
323        "If set to 0, no refinement is done. Otherwise, refinement is performed using spglib with provided precision."))
324    p_xrd.add_argument("-t", "--two-theta-range", default=(0, 90), nargs=2, help=(
325        "Tuple for range of two_thetas to calculate in degrees. Defaults to (0, 90)."))
326    p_xrd.add_argument("-nap", "--no-annotate-peaks", default=False, action="store_true",
327        help="Whether to annotate the peaks with plane information.")
328
329    # Subparser for oxistate.
330    p_oxistate = subparsers.add_parser('oxistate', parents=[copts_parser, path_selector],
331        help="Estimate oxidation states with pymatgen bond valence analysis.")
332    # Subparser for ipython.
333    p_ipython = subparsers.add_parser('ipython', parents=[copts_parser, path_selector],
334        help="Open IPython shell for advanced operations on structure object.")
335    # Subparser for notebook.
336    p_notebook = subparsers.add_parser('notebook', parents=[copts_parser, path_selector],
337        help="Read structure from file and generate jupyter notebook.")
338    p_notebook.add_argument('--classic-notebook', action='store_true', default=False,
339                            help="Use classic notebook instead of jupyterlab.")
340    p_notebook.add_argument('--no-browser', action='store_true', default=False,
341                            help=("Start the jupyter server to serve the notebook "
342                                  "but don't open the notebook in the browser.\n"
343                                  "Use this option to connect remotely from localhost to the machine running the kernel"))
344    p_notebook.add_argument('--foreground', action='store_true', default=False,
345        help="Run jupyter notebook in the foreground.")
346
347    p_panel = subparsers.add_parser('panel', parents=[copts_parser, path_selector],
348        help="Open GUI in web browser, requires panel package.")
349
350    # Subparser for kpath.
351    p_kpath = subparsers.add_parser('kpath', parents=[copts_parser, path_selector],
352        help="Read structure from file, generate k-path for band-structure calculations.")
353    add_format_arg(p_kpath, default="abinit", formats=["abinit", "wannier90", "siesta"])
354    # Subparser for bz.
355    p_bz = subparsers.add_parser('bz', parents=[copts_parser, path_selector],
356        help="Read structure from file, plot Brillouin zone with matplotlib.")
357    # Subparser for ngkpt.
358    p_ngkpt = subparsers.add_parser('ngkpt', parents=[copts_parser, path_selector],
359        help="Return the Abinit k-point sampling variables "
360             "from the number of divisions used to sample the smallest "
361             "lattice vector of the reciprocal lattice.")
362    p_ngkpt.add_argument("-n", "--nksmall", required=True, type=int,
363        help="Number of divisions used to sample the smallest reciprocal lattice vector.")
364    # Subparser for ktables.
365    p_ktables = subparsers.add_parser('ktables', parents=[copts_parser, path_selector],
366        help=("Read structure from filepath, call spglib to sample the BZ, "
367              "print k-points in the IBZ with weights."))
368    p_ktables.add_argument("-m", "--mesh", nargs=3, required=True, type=int, help="Mesh divisions e.g. 2 3 4")
369    p_ktables.add_argument("-s", "--is_shift", nargs="+", default=None, type=int,
370        help=("Three integers (spglib API). The kmesh is shifted along " +
371              "the axis in half of adjacent mesh points irrespective of the mesh numbers e.g. 1 1 1 " +
372              "Default: Unshifted mesh."))
373    p_ktables.add_argument("--no-time-reversal", default=False, action="store_true", help="Don't use time-reversal.")
374
375    # Subparser for abikmesh.
376    p_abikmesh = subparsers.add_parser('abikmesh', parents=[copts_parser, path_selector],
377        help=("Read structure from file, call Abinit to sample the BZ with ngkpt, shiftk, and kptopt. "
378              "Print k-points in the IBZ with weights."))
379    p_abikmesh.add_argument("--kppa", default=None, type=int,
380            help="Number of k-points per reciprocal atom. Mutually exclusive with ngkpt.")
381    p_abikmesh.add_argument("--ngkpt", nargs=3, default=None, type=int,
382            help="Mesh divisions e.g. 2 3 4")
383    p_abikmesh.add_argument("--shiftk", nargs="+", default=(0.5, 0.5, 0.5), type=float,
384       help="Kmesh shifts. Default: 0.5 0.5 0.5")
385    p_abikmesh.add_argument("--kptopt", default=1, type=int,
386            help="Kptopt input variable. Default: 1")
387
388    # Subparser for kmesh_jhu.
389    #p_kmesh_jhu = subparsers.add_parser('kmesh_jhu', parents=[copts_parser, path_selector], help="Foobar ")
390
391    # Subparser for lgk.
392    p_lgk = subparsers.add_parser('lgk', parents=[copts_parser, path_selector, spgopt_parser],
393        help="Read structure from file, find little group of k-point, print Bilbao character table.")
394    p_lgk.add_argument("-k", "--kpoint", nargs=3, required=True, type=float,
395        help="K-point in reduced coordinates e.g. 0.25 0 0")
396
397    # Subparser for kstar.
398    p_kstar = subparsers.add_parser('kstar', parents=[copts_parser, path_selector, spgopt_parser],
399        help="Read structure from file, print star of k-point.")
400    p_kstar.add_argument("-k", "--kpoint", nargs=3, required=True, type=float,
401        help="K-point in reduced coordinates e.g. 0.25 0 0")
402
403    # Subparser for keq.
404    p_keq = subparsers.add_parser('keq', parents=[copts_parser, path_selector, spgopt_parser],
405        help="Read structure from file, check whether two k-points are equivalent by symmetry.")
406    p_keq.add_argument("-k", "--kpoints", nargs=6, required=True, type=float,
407        help="K-points in reduced coordinates e.g. 0.25 0 0 0 0.25 0")
408
409    # Subparser for visualize command.
410    p_visualize = subparsers.add_parser('visualize', parents=[copts_parser, path_selector],
411        help=("Visualize the structure with the specified application. "
412              "Requires external app or optional python modules (mayavi, vtk)."))
413    p_visualize.add_argument("-a", "--appname", type=str, default="vesta",
414        help=("Application name. Possible options: %s, mpl (matplotlib), mayavi, vtk" % ", ".join(Visualizer.all_visunames())))
415
416    # Options for commands accessing the materials project database.
417    mp_rest_parser = argparse.ArgumentParser(add_help=False)
418    mp_rest_parser.add_argument("--mapi-key", default=None,
419        help="Pymatgen PMG_MAPI_KEY. Use value in .pmgrc.yaml if not specified.")
420    mp_rest_parser.add_argument("--endpoint", help="Pymatgen database.", default="https://www.materialsproject.org/rest/v2")
421    mp_rest_parser.add_argument("-b", "--browser", default=False, action='store_true',
422        help="Open materials-project webpages in browser")
423
424    # Subparser for mp_id command.
425    p_mpid = subparsers.add_parser('mp_id', parents=[copts_parser, mp_rest_parser],
426        help="Get structure from the pymatgen database. Export to format. Requires internet connection and PMG_MAPI_KEY.")
427    p_mpid.add_argument("mpid", type=str, default=None, help="Pymatgen identifier.")
428    add_format_arg(p_mpid, default="cif")
429
430    # Subparser for mp_match command.
431    p_mpmatch = subparsers.add_parser('mp_match', parents=[path_selector, mp_rest_parser, copts_parser, nb_parser],
432        help="Get structure from the pymatgen database. Requires internet connection and PMG_MAPI_KEY.")
433    add_format_arg(p_mpmatch, default="abivars")
434
435    # Subparser for mp_search command.
436    p_mpsearch = subparsers.add_parser('mp_search', parents=[mp_rest_parser, copts_parser, nb_parser],
437        help="Get structure from the pymatgen database. Requires internet connection and PMG_MAPI_KEY")
438    p_mpsearch.add_argument("chemsys_formula_id", type=str, default=None,
439        help="A chemical system (e.g., Li-Fe-O), or formula (e.g., Fe2O3) or materials_id (e.g., mp-149).")
440    p_mpsearch.add_argument("-s", "--select-spgnum", type=int, default=None,
441        help="Select structures with this space group number.")
442    add_format_arg(p_mpsearch, default="abivars")
443
444    # Subparser for mp_pd command.
445    p_mp_pda = subparsers.add_parser('mp_pd', parents=[mp_rest_parser, copts_parser],
446        help=("Generate phase diagram with entries from the Materials Project. "
447              "Requires internet connection and PMG_MAPI_KEY"))
448    p_mp_pda.add_argument("file_or_elements", type=str, default=None,
449        help="FILE with structure or elements e.g., Li-Fe-O).")
450    p_mp_pda.add_argument("-u", "--show-unstable", type=int, default=0,
451        help="""Whether unstable phases will be plotted as
452well as red crosses. If a number > 0 is entered, all phases with
453ehull < show_unstable will be shown.""")
454
455    # Subparser for mp_ebands command.
456    p_mp_ebands = subparsers.add_parser('mp_ebands', parents=[copts_parser, mp_rest_parser],
457        help="Get structure from the pymatgen database. Export to format. Requires internet connection and PMG_MAPI_KEY.")
458    p_mp_ebands.add_argument("chemsys_formula_id", type=str, default=None,
459        help="A chemical system (e.g., Li-Fe-O), or formula (e.g., Fe2O3) or materials_id (e.g., mp-149).")
460    #add_format_arg(p_mp_ebands, default="cif")
461
462    # Subparser for cod_search command.
463    p_codsearch = subparsers.add_parser('cod_search', parents=[copts_parser, nb_parser],
464        help="Get structure from COD database. Requires internet connection and mysql")
465    p_codsearch.add_argument("formula", type=str, default=None, help="formula (e.g., Fe2O3).")
466    p_codsearch.add_argument("-s", "--select-spgnum", type=int, default=None,
467        help="Select structures with this space group number.")
468    p_codsearch.add_argument('--primitive', default=False, action='store_true',
469        help="Convert COD cells into primitive cells.")
470    add_format_arg(p_codsearch, default="abivars")
471
472    # Subparser for cod_id command.
473    p_codid = subparsers.add_parser('cod_id', parents=[copts_parser],
474        help="Get structure from COD database. Requires internet connection and mysql")
475    p_codid.add_argument("cod_identifier", type=int, default=None, help="COD identifier.")
476    p_codid.add_argument('--primitive', default=False, action='store_true', help="Convert COD cell into primitive cell.")
477    add_format_arg(p_codid, default="abivars")
478
479    # Subparser for animate command.
480    p_animate = subparsers.add_parser('animate', parents=[copts_parser, path_selector],
481        help="Read structures from HIST.nc or XDATCAR. Print structures in Xcrysden AXSF format to stdout.")
482
483    return parser
484
485
486@prof_main
487def main():
488
489    def show_examples_and_exit(err_msg=None, error_code=1):
490        """Display the usage of the script."""
491        sys.stderr.write(get_epilog())
492        if err_msg: sys.stderr.write("Fatal Error\n" + err_msg + "\n")
493        sys.exit(error_code)
494
495    parser = get_parser(with_epilog=True)
496
497    # Parse command line.
498    try:
499        options = parser.parse_args()
500    except Exception as exc:
501        show_examples_and_exit(error_code=1)
502
503    if not options.command:
504        show_examples_and_exit(error_code=1)
505
506    if hasattr(options, "format"):
507        options.format = options.format.strip()
508
509    # loglevel is bound to the string value obtained from the command line argument.
510    # Convert to upper case to allow the user to specify --loglevel=DEBUG or --loglevel=debug
511    import logging
512    numeric_level = getattr(logging, options.loglevel.upper(), None)
513    if not isinstance(numeric_level, int):
514        raise ValueError('Invalid log level: %s' % options.loglevel)
515    logging.basicConfig(level=numeric_level)
516
517    if options.verbose > 2:
518        print(options)
519
520    if options.command == "spglib":
521        structure = abilab.Structure.from_file(options.filepath)
522        print(structure.spget_summary(symprec=options.symprec, angle_tolerance=options.angle_tolerance,
523                                      site_symmetry=options.site_symmetry, verbose=options.verbose))
524
525    elif options.command == "abispg":
526        structure = abilab.Structure.from_file(options.filepath)
527        check_ordered_structure(structure)
528        abi_spg = structure.abi_spacegroup
529
530        if abi_spg is not None and options.tolsym is None:
531            print(structure.spget_summary(verbose=options.verbose))
532        else:
533            # Here we compare Abinit wrt spglib. If abi_spg is None, we create a temporary
534            # task to run the code in dry-run mode.
535            if abi_spg is None:
536                print("FILE does not contain Abinit symmetry operations.")
537            cprint("Calling Abinit in --dry-run mode with chkprim = 0 to get space group.")
538            if options.tolsym is not None and options.tolsym > 1e-8:
539                cprint("Crystal structure will be re-symmetrized by Abinit with tolsym: %s" % options.tolsym, "yellow")
540
541            from abipy.data.hgh_pseudos import HGH_TABLE
542            gsinp = factories.gs_input(structure, HGH_TABLE, spin_mode="unpolarized")
543            gsinp["chkprim"] = 0
544            abistructure = gsinp.abiget_spacegroup(tolsym=options.tolsym)
545            print(abistructure.spget_summary(verbose=options.verbose))
546            print("")
547
548            diff_structures([structure, abistructure], mode=options.diff_mode,
549                            headers=["Input structure", "After Abinit symmetrization"], fmt="abivars")
550
551            # Save file.
552            save_structure(abistructure, options)
553
554    elif options.command == "convert":
555        fmt = options.format
556        if fmt == "cif" and options.filepath.endswith(".cif"): fmt = "abivars"
557        print(abilab.Structure.from_file(options.filepath).convert(fmt=fmt))
558
559    elif options.command == "print":
560        print(abilab.Structure.from_file(options.filepath).to_string(verbose=options.verbose))
561
562    elif options.command == "supercell":
563        structure = abilab.Structure.from_file(options.filepath)
564
565        options.scaling_matrix = np.array(options.scaling_matrix)
566        if len(options.scaling_matrix) == 9:
567            options.scaling_matrix.shape = (3, 3)
568        if options.verbose:
569            print("scaling matrix: ", options.scaling_matrix)
570
571        supcell = structure * options.scaling_matrix
572        #supcell = structure.make_supercell(scaling_matrix, to_unit_cell=True)
573        print(supcell.convert(fmt=options.format))
574
575    elif options.command == "abisanitize":
576        print("\nCalling abi_sanitize to get a new structure in which:")
577        print("    * Structure is refined.")
578        print("    * Reduced to primitive settings.")
579        print("    * Lattice vectors are exchanged if the triple product is negative\n")
580
581        structure = abilab.Structure.from_file(options.filepath)
582        sanitized = structure.abi_sanitize(symprec=options.symprec, angle_tolerance=options.angle_tolerance,
583                                           primitive=not options.no_primitive, primitive_standard=options.primitive_standard)
584        index = [options.filepath, "abisanitized"]
585        dfs = abilab.dataframes_from_structures([structure, sanitized], index=index, with_spglib=True)
586
587        abilab.print_dataframe(dfs.lattice, title="Lattice parameters:")
588        abilab.print_dataframe(dfs.coords, title="Atomic positions (columns give the site index):")
589
590        if not options.verbose:
591            print("\nUse -v for more info")
592        else:
593            #print("\nDifference between structures:")
594            if len(structure) == len(sanitized):
595                table = []
596                for line1, line2 in zip(str(structure).splitlines(), str(sanitized).splitlines()):
597                    table.append([line1, line2])
598                print(str(tabulate(table, headers=["Initial structure", "Abisanitized"])))
599
600            else:
601                print("\nInitial structure:")
602                print(structure)
603                print("\nabisanitized structure:")
604                print(sanitized)
605
606        # Save file.
607        save_structure(sanitized, options)
608
609    elif options.command == "primitive":
610        structure = abilab.Structure.from_file(options.filepath)
611        primitive = structure.get_primitive_structure(tolerance=0.25, use_site_props=False, constrain_latt=None)
612        separator = "\n" + 90 * "="
613        print("\nInitial structure:\n", structure, separator)
614        print("\nPrimitive structure returned by pymatgen:\n", primitive, separator)
615        print("\nAbinit input for primitive structure:\n", primitive.abi_string)
616        if structure != primitive:
617            print("\nInput structure is not primitive (according to pymatgen + spglib)")
618        else:
619            print("\nInput structure seems to be primitive (according to pymatgen + spglib)")
620
621    elif options.command == "irefine":
622        structure = abilab.Structure.from_file(options.filepath)
623        sanitized = structure.copy()
624        symprec, angle_tolerance = options.symprec, options.angle_tolerance
625        print("Calling abi_sanitize with increasing tolerances to reach target space group:", options.target_spgnum)
626        print("Using symprec_step: ", options.symprec_step, ", angle_tolerance_step:", options.angle_tolerance_step,
627              "ntrial", options.ntrial)
628        itrial = 0
629        while itrial < options.ntrial:
630            print(">>> Trying with symprec: %s, angle_tolerance: %s" % (symprec, angle_tolerance))
631            sanitized = sanitized.abi_sanitize(symprec=symprec, angle_tolerance=angle_tolerance,
632                primitive=not options.no_primitive, primitive_standard=options.primitive_standard)
633            spg_symb, spg_num = sanitized.get_space_group_info(symprec=symprec, angle_tolerance=angle_tolerance)
634            print(">>> Space-group number:", spg_symb, ", symbol:", spg_num, "for trial:", itrial)
635            if spg_num == options.target_spgnum:
636                print(2 * "\n", "# Final structure with space group number:", spg_symb, ", symbol:", spg_num, 2 * "\n")
637                print(sanitized.convert(fmt="cif"))
638                break
639
640            # Increment counter and tolerances.
641            itrial += 1
642            symprec += options.symprec_step
643            angle_tolerance += options.angle_tolerance_step
644        else:
645            print("Cannot find space group number:", options.target_spgnum, "after", options.ntrial, "iterations")
646            return 1
647
648        # Save file.
649        #save_structure(sanitized, options)
650
651    elif options.command == "conventional":
652        print("\nCalling get_conventional_standard_structure to get conventional structure:")
653        print("The standards are defined in Setyawan, W., & Curtarolo, S. (2010).")
654        print("High-throughput electronic band structure calculations: Challenges and tools. ")
655        print("Computational Materials Science, 49(2), 299-312. doi:10.1016/j.commatsci.2010.05.010\n")
656
657        structure = abilab.Structure.from_file(options.filepath)
658        conv = structure.get_conventional_standard_structure(international_monoclinic=True,
659                                           symprec=options.symprec, angle_tolerance=options.angle_tolerance)
660        index = [options.filepath, "conventional"]
661        dfs = abilab.dataframes_from_structures([structure, conv], index=index, with_spglib=True)
662
663        abilab.print_dataframe(dfs.lattice, title="Lattice parameters:")
664        if options.verbose:
665            abilab.print_dataframe(dfs.coords, title="Atomic positions (columns give the site index):")
666
667        if not options.verbose:
668            print("\nUse -v for more info")
669        else:
670            #print("\nDifference between structures:")
671            if len(structure) == len(conv):
672                table = []
673                for line1, line2 in zip(str(structure).splitlines(), str(conv).splitlines()):
674                    table.append([line1, line2])
675                print(str(tabulate(table, headers=["Initial structure", "Conventional"])))
676
677            else:
678                print("\nInitial structure:\n", structure)
679                print("\nConventional structure:\n", conv)
680
681        # Save file.
682        save_structure(conv, options)
683
684    elif options.command == "proto":
685        structure = abilab.Structure.from_file(options.filepath)
686        from pymatgen.analysis.aflow_prototypes import AflowPrototypeMatcher
687        m = AflowPrototypeMatcher(initial_ltol=options.ltol, initial_stol=options.stol,
688                                  initial_angle_tol=options.angle_tol)
689        dlist = m.get_prototypes(structure)
690        if not dlist:
691            cprint("Cannot find AFLOW prototype for structure.")
692            print(structure.to_string(verbose=options.verbose))
693            return 1
694        else:
695            cprint("Found %d matches" % len(dlist), "green")
696            for d in dlist:
697                if "snl" in d:
698                    snl = d.pop("snl")
699                    if options.verbose: pprint(snl.as_dict())
700                pprint(d)
701                url = "http://aflow.org/CrystalDatabase/%s.html" % d["tags"]["aflow"]
702                print("AFLOW url: %s\n" % url)
703            if not options.verbose:
704                print("Use --verbose to increase output level")
705
706    elif options.command == "wyckoff":
707        structure = abilab.Structure.from_file(options.filepath)
708        if options.refine:
709            print("Refining structure with symprec: %s, angle_tolerance: %s" % (options.symprec, options.angle_tolerance))
710            structure = structure.refine(symprec=options.symprec, angle_tolerance=options.angle_tolerance)
711        print(structure.spget_summary(verbose=options.verbose))
712        ss = structure.site_symmetries
713        df = ss.get_wyckoff_dataframe(verbose=options.verbose)
714        abilab.print_dataframe(df, title="\nWyckoff positions in reduced coordinates.")
715
716    elif options.command == "tensor_site":
717        structure = abilab.Structure.from_file(options.filepath)
718        if options.refine:
719            print("Refining structure with symprec: %s, angle_tolerance: %s" % (options.symprec, options.angle_tolerance))
720            structure = structure.refine(symprec=options.symprec, angle_tolerance=options.angle_tolerance)
721        print(structure.spget_summary(verbose=options.verbose))
722        ss = structure.site_symmetries
723        df = ss.get_tensor_rank2_dataframe(verbose=options.verbose)
724        abilab.print_dataframe(df, title="\nTensor components in reduced coordinates (rank 2, symmetric)")
725
726    elif options.command == "neighbors":
727        abilab.Structure.from_file(options.filepath).print_neighbors(radius=options.radius)
728
729    elif options.command == "interpolate":
730        initial_structure = abilab.Structure.from_file(options.filepaths[0])
731        end_structure = abilab.Structure.from_file(options.filepaths[1])
732        structures = initial_structure.interpolate(end_structure, nimages=options.nimages,
733                                                   interpolate_lattices=False, pbc=True,
734                                                   autosort_tol=options.autosort_tol)
735        structures = list(map(abilab.Structure.as_structure, structures))
736        for i, s in enumerate(structures):
737            print(marquee("Structure #%d" % i, mark="="))
738            print(s.convert(fmt=options.format))
739            print(" ")
740
741    elif options.command == "xrd":
742        structure = abilab.Structure.from_file(options.filepath)
743        two_theta_range = tuple(float(t) for t in options.two_theta_range)
744        structure.plot_xrd(wavelength=options.wavelength, two_theta_range=two_theta_range,
745                           symprec=options.symprec, annotate_peaks=not options.no_annotate_peaks)
746
747    elif options.command == "oxistate":
748        print(abilab.Structure.from_file(options.filepath).get_oxi_state_decorated())
749
750    elif options.command == "ipython":
751        structure = abilab.Structure.from_file(options.filepath)
752        print("Invoking Ipython, `structure` object will be available in the Ipython terminal")
753        import IPython
754        IPython.start_ipython(argv=[], user_ns={"structure": structure})
755
756    elif options.command == "notebook":
757        structure = abilab.Structure.from_file(options.filepath)
758        structure.make_and_open_notebook(nbpath=None, foreground=options.foreground,
759                                         classic_notebook=options.classic_notebook,
760                                         no_browser=options.no_browser)
761
762    elif options.command == "panel":
763        structure = abilab.Structure.from_file(options.filepath)
764        try:
765            import panel  # noqa: F401
766        except ImportError as exc:
767            cprint("Use `conda install panel` or `pip install panel` to install the python package.", "red")
768            raise exc
769
770        structure.get_panel().show()  #threaded=True)
771        return 0
772
773    elif options.command == "visualize":
774        structure = abilab.Structure.from_file(options.filepath)
775        print(structure)
776        print("Visualizing structure with:", options.appname)
777        structure.visualize(appname=options.appname)
778
779    elif options.command == "kpath":
780        structure = abilab.Structure.from_file(options.filepath)
781        print(structure.get_kpath_input_string(fmt=options.format, line_density=10))
782
783    elif options.command == "bz":
784        abilab.Structure.from_file(options.filepath).plot_bz()
785
786    elif options.command == "ngkpt":
787        d = abilab.Structure.from_file(options.filepath).calc_ksampling(options.nksmall)
788        print("ngkpt %d %d %d" % (d.ngkpt[0], d.ngkpt[1], d.ngkpt[2]))
789        print("nshiftk ", len(d.shiftk), "\nshiftk")
790        for s in d.shiftk:
791            print("  %s %s %s" % (s[0], s[1], s[2]))
792
793    elif options.command == "ktables":
794        structure = abilab.Structure.from_file(options.filepath)
795        k = Ktables(structure, options.mesh, options.is_shift, not options.no_time_reversal)
796        print(k)
797        print("")
798        print("NB: These results are obtained by calling spglib with the structure read from file.")
799        print("The k-points might differ from the ones expected by Abinit, especially if the space groups differ.")
800
801        if not options.verbose:
802            print("\nUse -v to obtain the BZ --> IBZ mapping.")
803        else:
804            print()
805            k.print_bz2ibz()
806
807    elif options.command == "abikmesh":
808        structure = abilab.Structure.from_file(options.filepath)
809        if options.kppa is None and options.ngkpt is None:
810            raise ValueError("Either ngkpt or kppa must be provided")
811
812        if options.kppa is not None:
813            print("Calling Abinit to compute the IBZ with kppa:", options.kppa, "and shiftk:", options.shiftk)
814            ibz = IrredZone.from_kppa(structure, options.kppa, options.shiftk,
815                                      kptopt=options.kptopt, verbose=options.verbose)
816        else:
817            print("Calling Abinit to compute the IBZ with ngkpt:", options.ngkpt, "and shiftk:", options.shiftk)
818            ibz = IrredZone.from_ngkpt(structure, options.ngkpt, options.shiftk,
819                                       kptopt=options.kptopt, verbose=options.verbose)
820
821        print(ibz.to_string(verbose=options.verbose))
822
823    #elif options.command == "kmesh_jhu":
824    #    structure = abilab.Structure.from_file(options.filepath)
825    #    from pymatgen.ext.jhu import get_kpoints
826    #    kpoints = get_kpoints(structure, min_distance=0, min_total_kpoints=1,
827    #                           kppra=None, gap_distance=7, remove_symmetry=None,
828    #                           include_gamma="auto", header="simple", incar=None)
829    #    #print(kpoints)
830
831    elif options.command == "lgk":
832        structure = abilab.Structure.from_file(options.filepath)
833        spgrp = structure.abi_spacegroup
834        if spgrp is None:
835            cprint("Your file does not contain Abinit symmetry operations.", "yellow")
836            cprint("Will call spglib to obtain the space group (assuming time-reversal: %s)" %
837                   (not options.no_time_reversal), "yellow")
838            spgrp = AbinitSpaceGroup.from_structure(structure, has_timerev=not options.no_time_reversal,
839                        symprec=options.symprec, angle_tolerance=options.angle_tolerance)
840        print()
841        print(marquee("Structure", mark="="))
842        print(structure.spget_summary(verbose=options.verbose))
843        print("\n")
844
845        print(marquee("Little Group", mark="="))
846        ltk = spgrp.find_little_group(kpoint=options.kpoint)
847        print(ltk.to_string(verbose=options.verbose))
848
849    elif options.command == "kstar":
850        structure = abilab.Structure.from_file(options.filepath)
851
852        # Call spglib to get spacegroup if Abinit spacegroup is not available.
853        if structure.abi_spacegroup is None:
854            structure.spgset_abi_spacegroup(has_timerev=not options.no_time_reversal)
855
856        kpoint = Kpoint(options.kpoint, structure.reciprocal_lattice)
857        kstar = kpoint.compute_star(structure.abi_spacegroup, wrap_tows=True)
858        print("Found %s points in the star of %s\n" % (len(kstar), repr(kpoint)))
859        for k in kstar:
860            print(4 * " ", repr(k))
861
862    elif options.command == "keq":
863        structure = abilab.Structure.from_file(options.filepath)
864
865        # Call spglib to get spacegroup if Abinit spacegroup is not available.
866        if structure.abi_spacegroup is None:
867            structure.spgset_abi_spacegroup(has_timerev=not options.no_time_reversal)
868
869        k1, k2 = options.kpoints[:3], options.kpoints[3:6]
870        k1tab = structure.abi_spacegroup.symeq(k1, k2)
871
872        if k1tab.isym != -1:
873            print("\nk1:", k1, "and k2:", k2, "are symmetry equivalent k-points\n")
874            print("Related by the symmetry operation (reduced coords):\n", k1tab.op)
875            print("With umklapp vector Go = TO(k1) - k2 =", k1tab.g0)
876        else:
877            print(k1, "and", k2, "are NOT symmetry equivalent")
878
879    elif options.command == "mp_id":
880        # Get the Structure corresponding to material_id.
881        structure = abilab.Structure.from_mpid(options.mpid, final=True,
882                                               api_key=options.mapi_key, endpoint=options.endpoint)
883        # Convert to format and print it.
884        print(structure.convert(fmt=options.format))
885
886    elif options.command == "mp_match":
887        mp = abilab.mp_match_structure(options.filepath)
888        if not mp.structures:
889            cprint("No structure found in MP database", "yellow")
890            return 1
891
892        if options.notebook:
893            return mp.make_and_open_notebook(foreground=options.foreground,
894                                             classic_notebook=options.classic_notebook,
895                                             no_browser=options.no_browser)
896        else:
897            mp.print_results(fmt=options.format, verbose=options.verbose)
898
899        if options.browser:
900            mp.open_browser(limit=None if options.verbose == 2 else 10)
901
902    elif options.command == "mp_search":
903        mp = abilab.mp_search(options.chemsys_formula_id)
904        if not mp.structures:
905            cprint("No structure found in Materials Project database", "yellow")
906            return 1
907        if options.select_spgnum: mp = mp.filter_by_spgnum(options.select_spgnum)
908
909        if options.notebook:
910            return mp.make_and_open_notebook(foreground=options.foreground,
911                                             classic_notebook=options.classic_notebook,
912                                             no_browser=options.no_browser)
913        else:
914            mp.print_results(fmt=options.format, verbose=options.verbose)
915
916        if options.browser:
917            mp.open_browser(limit=None if options.verbose == 2 else 10)
918
919    elif options.command == "mp_pd":
920        if os.path.exists(options.file_or_elements):
921            structure = abilab.Structure.from_file(options.file_or_elements)
922            elements = structure.symbol_set
923        else:
924            elements = options.file_or_elements.split("-")
925
926        if options.verbose > 1: print("Building phase-diagram for elements:", elements)
927        with abilab.restapi.get_mprester(api_key=options.mapi_key, endpoint=options.endpoint) as rest:
928            pdr = rest.get_phasediagram_results(elements)
929            pdr.print_dataframes(verbose=options.verbose)
930            pdr.plot(show_unstable=options.show_unstable)
931
932    elif options.command == "mp_ebands":
933        if os.path.exists(options.chemsys_formula_id):
934            mp = abilab.mp_match_structure(options.chemsys_formula_id)
935            for mpid in mp.ids:
936                ebands = abilab.ElectronBands.from_mpid(mpid, api_key=options.mapi_key, endpoint=options.endpoint)
937                print(ebands)
938        else:
939            if options.chemsys_formula_id.startswith("mp-"):
940                # Assume valid mp identifier.
941                mpid = options.chemsys_formula_id
942                ebands = abilab.ElectronBands.from_mpid(mpid, api_key=options.mapi_key, endpoint=options.endpoint)
943                print(ebands)
944            else:
945                mp = abilab.mp_search(options.chemsys_formula_id)
946                if not mp.structures:
947                    cprint("No structure found in Materials Project database", "yellow")
948                    return 1
949
950                for mpid, structure in zip(mp.ids, mp.structures):
951                    if structure is None:
952                        cprint("ignoring mpid %s because cannot find structure" % mpid, "red")
953                        continue
954                    ebands = abilab.ElectronBands.from_mpid(mpid, api_key=options.mapi_key, endpoint=options.endpoint)
955                    if ebands is None:
956                        cprint("Cannot get ebands for structure:\n%s" % str(structure), "red")
957                    else:
958                        print(ebands)
959
960    elif options.command == "cod_search":
961        cod = abilab.cod_search(options.formula, primitive=options.primitive)
962        if not cod.structures:
963            cprint("No structure found in COD database", "yellow")
964            return 1
965        if options.select_spgnum: cod = cod.filter_by_spgnum(options.select_spgnum)
966
967        if options.notebook:
968            return cod.make_and_open_notebook(foreground=options.foreground,
969                                              classic_notebook=options.classic_notebook,
970                                              no_browser=options.no_browser)
971        else:
972            cod.print_results(fmt=options.format, verbose=options.verbose)
973
974    elif options.command == "cod_id":
975        # Get the Structure from COD
976        structure = abilab.Structure.from_cod_id(options.cod_identifier, primitive=options.primitive)
977        # Convert to format and print it.
978        print(structure.convert(fmt=options.format))
979
980    elif options.command == "animate":
981        filepath = options.filepath
982        if any(filepath.endswith(ext) for ext in ("HIST", "HIST.nc")):
983            with abilab.abiopen(filepath) as hist:
984                structures = hist.structures
985
986        elif "XDATCAR" in filepath:
987            structures = Xdatcar(filepath).structures
988            if not structures:
989                raise RuntimeError("Your Xdatcar contains only one structure. Due to a bug "
990                    "in the pymatgen routine, your structures won't be parsed correctly"
991                    "Solution: Add another structure at the end of the file.")
992        else:
993            raise ValueError("Don't know how to handle file %s" % filepath)
994
995        xsf_write_structure(sys.stdout, structures)
996
997    else:
998        raise ValueError("Unsupported command: %s" % options.command)
999
1000    return 0
1001
1002
1003if __name__ == "__main__":
1004    sys.exit(main())
1005