pymatgen.io.vasp.sets module

class DictSet(structure, config_dict, files_to_transfer=None, user_incar_settings=None, user_kpoints_settings=None, user_potcar_settings=None, constrain_total_magmom=False, sort_structure=True, potcar_functional='PBE', force_gamma=False, reduce_structure=None, vdw=None)[source]

Bases: pymatgen.io.vasp.sets.VaspInputSet

Concrete implementation of VaspInputSet that is initialized from a dict settings. This allows arbitrary settings to be input. In general, this is rarely used directly unless there is a source of settings in yaml format (e.g., from a REST interface). It is typically used by other VaspInputSets for initialization.

Special consideration should be paid to the way the MAGMOM initialization for the INCAR is done. The initialization differs depending on the type of structure and the configuration settings. The order in which the magmom is determined is as follows:

  1. If the site itself has a magmom setting, that is used.
  2. If the species on the site has a spin setting, that is used.
  3. If the species itself has a particular setting in the config file, that is used, e.g., Mn3+ may have a different magmom than Mn4+.
  4. Lastly, the element symbol itself is checked in the config file. If there are no settings, VASP’s default of 0.6 is used.
Parameters:
  • structure (Structure) – The Structure to create inputs for.
  • config_dict (dict) – The config dictionary to use.
  • files_to_transfer (dict) – A dictionary of {filename: filepath}. This allows the transfer of files from a previous calculation.
  • user_incar_settings (dict) – User INCAR settings. This allows a user to override INCAR settings, e.g., setting a different MAGMOM for various elements or species. Note that in the new scheme, ediff_per_atom and hubbard_u are no longer args. Instead, the config_dict supports EDIFF_PER_ATOM and EDIFF keys. The former scales with # of atoms, the latter does not. If both are present, EDIFF is preferred. To force such settings, just supply user_incar_settings={“EDIFF”: 1e-5, “LDAU”: False} for example. The keys ‘LDAUU’, ‘LDAUJ’, ‘LDAUL’ are special cases since pymatgen defines different values depending on what anions are present in the structure, so these keys can be defined in one of two ways, e.g. either {“LDAUU”:{“O”:{“Fe”:5}}} to set LDAUU for Fe to 5 in an oxide, or {“LDAUU”:{“Fe”:5}} to set LDAUU to 5 regardless of the input structure.
  • user_kpoints_settings (dict or Kpoints) – Allow user to override kpoints setting by supplying a dict E.g., {“reciprocal_density”: 1000}. User can also supply Kpoints object. Default is None.
  • (dict (user_potcar_settings) – Allow user to override POTCARs. E.g., {“Gd”: “Gd_3”}. This is generally not recommended. Default is None.
  • constrain_total_magmom (bool) – Whether to constrain the total magmom (NUPDOWN in INCAR) to be the sum of the expected MAGMOM for all species. Defaults to False.
  • sort_structure (bool) – Whether to sort the structure (using the default sort order of electronegativity) before generating input files. Defaults to True, the behavior you would want most of the time. This ensures that similar atomic species are grouped together.
  • potcar_functional (str) – Functional to use. Default (None) is to use the functional in Potcar.DEFAULT_FUNCTIONAL. Valid values: “PBE”, “LDA”, “PW91”, “LDA_US”
  • force_gamma (bool) – Force gamma centered kpoint generation. Default (False) is to use the Automatic Density kpoint scheme, which will use the Gamma centered generation scheme for hexagonal cells, and Monkhorst-Pack otherwise.
  • reduce_structure (None/str) – Before generating the input files, generate the reduced structure. Default (None), does not alter the structure. Valid values: None, “niggli”, “LLL”.
  • vdw – Adds default parameters for van-der-Waals functionals supported by VASP to INCAR. Supported functionals are: DFT-D2, undamped DFT-D3, DFT-D3 with Becke-Jonson damping, Tkatchenko-Scheffler, Tkatchenko-Scheffler with iterative Hirshfeld partitioning, MBD@rSC, dDsC, Dion’s vdW-DF and DF2, optPBE, optB88, and optB86b.
incar
kpoints

Writes out a KPOINTS file using the fully automated grid method. Uses Gamma centered meshes for hexagonal cells and Monk grids otherwise.

Algorithm:
Uses a simple approach scaling the number of divisions along each reciprocal lattice vector proportional to its length.
nelect

Gets the default number of electrons for a given structure.

poscar
write_input(output_dir, make_dir_if_not_present=True, include_cif=False)[source]
class MITMDSet(structure, start_temp, end_temp, nsteps, time_step=2, spin_polarized=False, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MITRelaxSet

Class for writing a vasp md run. This DOES NOT do multiple stage runs.

Parameters:
  • structure (Structure) – Input structure.
  • start_temp (int) – Starting temperature.
  • end_temp (int) – Final temperature.
  • nsteps (int) – Number of time steps for simulations. The NSW parameter.
  • time_step (int) – The time step for the simulation. The POTIM parameter. Defaults to 2fs.
  • spin_polarized (bool) – Whether to do spin polarized calculations. The ISPIN parameter. Defaults to False.
  • **kwargs – Other kwargs supported by DictSet.
kpoints
class MITNEBSet(structures, unset_encut=False, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MITRelaxSet

Class for writing NEB inputs. Note that EDIFF is not on a per atom basis for this input set.

Parameters:
  • unset_encut (bool) – Whether to unset ENCUT.
  • **kwargs – Other kwargs supported by DictSet.
poscar
poscars
write_input(output_dir, make_dir_if_not_present=True, write_cif=False, write_path_cif=False, write_endpoint_inputs=False)[source]

NEB inputs has a special directory structure where inputs are in 00, 01, 02, ….

Parameters:
  • output_dir (str) – Directory to output the VASP input files
  • make_dir_if_not_present (bool) – Set to True if you want the directory (and the whole path) to be created if it is not present.
  • write_cif (bool) – If true, writes a cif along with each POSCAR.
  • write_path_cif (bool) – If true, writes a cif for each image.
  • write_endpoint_inputs (bool) – If true, writes input files for running endpoint calculations.
class MITRelaxSet(structure, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.DictSet

Standard implementation of VaspInputSet utilizing parameters in the MIT High-throughput project. The parameters are chosen specifically for a high-throughput project, which means in general pseudopotentials with fewer electrons were chosen.

Please refer:

A Jain, G. Hautier, C. Moore, S. P. Ong, C. Fischer, T. Mueller,
K. A. Persson, G. Ceder. A high-throughput infrastructure for density
functional theory calculations. Computational Materials Science,
2011, 50(8), 2295-2310. doi:10.1016/j.commatsci.2011.02.023
CONFIG = {'INCAR': {'ALGO': 'FAST', 'EDIFF': 1e-05, 'ENCUT': 520, 'IBRION': 2, 'ICHARG': 1, 'ISIF': 3, 'ISMEAR': -5, 'ISPIN': 2, 'ISYM': 0, 'LDAU': True, 'LDAUJ': {'F': {'Ag': 0, 'Co': 0, 'Cr': 0, 'Cu': 0, 'Fe': 0, 'Mn': 0, 'Mo': 0, 'Nb': 0, 'Ni': 0, 'Re': 0, 'Ta': 0, 'V': 0, 'W': 0}, 'O': {'Ag': 0, 'Co': 0, 'Cr': 0, 'Cu': 0, 'Fe': 0, 'Mn': 0, 'Mo': 0, 'Nb': 0, 'Ni': 0, 'Re': 0, 'Ta': 0, 'V': 0, 'W': 0}, 'S': {'Fe': 0, 'Mn': 0}}, 'LDAUL': {'F': {'Ag': 2, 'Co': 2, 'Cr': 2, 'Cu': 2, 'Fe': 2, 'Mn': 2, 'Mo': 2, 'Nb': 2, 'Ni': 2, 'Re': 2, 'Ta': 2, 'V': 2, 'W': 2}, 'O': {'Ag': 2, 'Co': 2, 'Cr': 2, 'Cu': 2, 'Fe': 2, 'Mn': 2, 'Mo': 2, 'Nb': 2, 'Ni': 2, 'Re': 2, 'Ta': 2, 'V': 2, 'W': 2}, 'S': {'Fe': 2, 'Mn': 2.5}}, 'LDAUTYPE': 2, 'LDAUU': {'F': {'Ag': 1.5, 'Co': 3.4, 'Cr': 3.5, 'Cu': 4, 'Fe': 4.0, 'Mn': 3.9, 'Mo': 4.38, 'Nb': 1.5, 'Ni': 6, 'Re': 2, 'Ta': 2, 'V': 3.1, 'W': 4.0}, 'O': {'Ag': 1.5, 'Co': 3.4, 'Cr': 3.5, 'Cu': 4, 'Fe': 4.0, 'Mn': 3.9, 'Mo': 4.38, 'Nb': 1.5, 'Ni': 6, 'Re': 2, 'Ta': 2, 'V': 3.1, 'W': 4.0}, 'S': {'Fe': 1.9, 'Mn': 2.5}}, 'LDAUPRINT': 1, 'LORBIT': '11', 'LREAL': 'AUTO', 'LWAVE': False, 'NELM': 200, 'NELMIN': 6, 'NSW': 99, 'PREC': 'Accurate', 'SIGMA': 0.05, 'MAGMOM': {'Co': 5, 'Co3+': 0.6, 'Co4+': 1, 'Cr': 5, 'Fe': 5, 'Mn': 5, 'Mn3+': 4, 'Mn4+': 3, 'Mo': 5, 'Ni': 5, 'V': 5, 'W': 5, 'Ce': 5, 'Eu': 10, 'Eu2+': 7, 'Eu3+': 6}}, 'KPOINTS': {'length': 25}, 'POTCAR': {'Ac': {'symbol': 'Ac', 'hash': 'd6854224d20e3de6e6fd7399503791d1'}, 'Ag': {'symbol': 'Ag', 'hash': 'e8ffa02fe3f3a51338ac1ac91ae968b9'}, 'Al': {'symbol': 'Al', 'hash': 'a6fd9a46aec185f4ad2acd0cbe4ae2fa'}, 'Ar': {'symbol': 'Ar', 'hash': 'e782fc6292623b396091bf8b871c272f'}, 'As': {'symbol': 'As', 'hash': '8005364db225a254e52cba350bedd032'}, 'Au': {'symbol': 'Au', 'hash': 'a9182d436a13194b744640ac940ab9b0'}, 'B': {'symbol': 'B', 'hash': '18ed2875dfa6305324cec3d7d59273ae'}, 'Ba': {'symbol': 'Ba_sv', 'hash': 'c0477913afb63dfae3439f3534fbf0ed'}, 'Be': {'symbol': 'Be', 'hash': 'fb974e44d56a8c62c6bbd1a1eb70c3a7'}, 'Bi': {'symbol': 'Bi', 'hash': 'e29661c79d59abae3b3ba69eae24b1a5'}, 'Br': {'symbol': 'Br', 'hash': '40f9594b4506684a69158c8975cfb9d6'}, 'C': {'symbol': 'C', 'hash': 'c0a8167dbb174fe492a3db7f5006c0f8'}, 'Ca': {'symbol': 'Ca_sv', 'hash': 'eb006721e214c04b3c13146e81b3a27d'}, 'Cd': {'symbol': 'Cd', 'hash': '0506b2d0ac28d5fe2b5ced77a701aa86'}, 'Ce': {'symbol': 'Ce', 'hash': 'ff3a09f2ff91798e58eb4b9854e9be4a'}, 'Cl': {'symbol': 'Cl', 'hash': '779b9901046c78fe51c5d80224642aeb'}, 'Co': {'symbol': 'Co', 'hash': 'b169bca4e137294d2ab3df8cbdd09083'}, 'Cr': {'symbol': 'Cr', 'hash': '82c14307937c7509fda4e9bc023d243d'}, 'Cs': {'symbol': 'Cs_sv', 'hash': '096b53a7d80cc0086976bcda50d536e5'}, 'Cu': {'symbol': 'Cu', 'hash': '8ca4e43a30de0c397e51f16bbb20d678'}, 'Dy': {'symbol': 'Dy_3', 'hash': 'd4a05220ab0a2d4c03a76872ea724a1e'}, 'Er': {'symbol': 'Er_3', 'hash': 'daa65a04877317f8c3c593ddeaa8a132'}, 'Eu': {'symbol': 'Eu', 'hash': 'd466d046adf21f6146ee9644049ea268'}, 'F': {'symbol': 'F', 'hash': '180141c33d032bfbfff30b3bea9d23dd'}, 'Fe': {'symbol': 'Fe', 'hash': '9530da8244e4dac17580869b4adab115'}, 'Ga': {'symbol': 'Ga', 'hash': '6e0b9d58412b1bfcd7252aff13d476c2'}, 'Gd': {'symbol': 'Gd', 'hash': '1f0d42b1e5f6769d319d3f247992aeb9'}, 'Ge': {'symbol': 'Ge', 'hash': '79e788788c31e196a460553010512d3f'}, 'H': {'symbol': 'H', 'hash': 'bb43c666e3d36577264afe07669e9582'}, 'He': {'symbol': 'He', 'hash': '47f9434aa3db96c85d7c4b3e4c2df09b'}, 'Hf': {'symbol': 'Hf', 'hash': 'b113f150cbf9c736f8244a6c25b0482e'}, 'Hg': {'symbol': 'Hg', 'hash': 'c2f15dfb5fd53396c5427635e5019160'}, 'Ho': {'symbol': 'Ho_3', 'hash': '661891464a27e87cf7e1324dd1893b77'}, 'I': {'symbol': 'I', 'hash': 'f4ff16a495dd361ff5824ee61b418bb0'}, 'In': {'symbol': 'In', 'hash': '7df38c0cdb4e6d9a9b93f09d690bb3ae'}, 'Ir': {'symbol': 'Ir', 'hash': 'dbcf7dcc6f4fb40df7b3d26904f60a66'}, 'K': {'symbol': 'K_sv', 'hash': '3e84f86d37f203a4fb01de36af57e430'}, 'Kr': {'symbol': 'Kr', 'hash': '39b9b85ae3982e6c012fb549b2840ce5'}, 'La': {'symbol': 'La', 'hash': '9b3ce03d18f7c0b40471a817ff91b287'}, 'Li': {'symbol': 'Li', 'hash': '65e83282d1707ec078c1012afbd05be8'}, 'Lu': {'symbol': 'Lu_3', 'hash': 'd40a90babf1224b88ffb4c3273ac3848'}, 'Mg': {'symbol': 'Mg', 'hash': '1771eb72adbbfa6310d66e7517e49930'}, 'Mn': {'symbol': 'Mn', 'hash': 'd082dba29b57ab59b3165e605dbf71b8'}, 'Mo': {'symbol': 'Mo_pv', 'hash': '84e18fd84a98e3d7fa8f055952410df0'}, 'N': {'symbol': 'N', 'hash': 'b98fd027ddebc67da4063ff2cabbc04b'}, 'Na': {'symbol': 'Na', 'hash': '1a89e79f7e21d99e8cf5788979f6a987'}, 'Nb': {'symbol': 'Nb_pv', 'hash': '7bcee99a4dc3094be0f9fd7961c02966'}, 'Nd': {'symbol': 'Nd', 'hash': '0c64e63070cee837c967283fffa001df'}, 'Ne': {'symbol': 'Ne', 'hash': '52064eee378b9e37a295a674f1c278f0'}, 'Ni': {'symbol': 'Ni', 'hash': '653f5772e68b2c7fd87ffd1086c0d710'}, 'Np': {'symbol': 'Np', 'hash': '20cb30b714200c4db870550b288ac4cd'}, 'O': {'symbol': 'O', 'hash': '7a25bc5b9a5393f46600a4939d357982'}, 'Os': {'symbol': 'Os', 'hash': '35c2cb48d48a9c38c40fb82bbe70626d'}, 'P': {'symbol': 'P', 'hash': '7dc3393307131ae67785a0cdacb61d5f'}, 'Pa': {'symbol': 'Pa', 'hash': 'a1fdb1089d0727f415416ec8082246ba'}, 'Pb': {'symbol': 'Pb', 'hash': '704c2c967247d7f84090d2536c91877d'}, 'Pd': {'symbol': 'Pd', 'hash': 'a395eb3aaf2fcab12fac3030a1146f61'}, 'Pm': {'symbol': 'Pm', 'hash': 'a2c9485ea86b2a7cf175077e6e5c7b3e'}, 'Pr': {'symbol': 'Pr', 'hash': '92f191499bf5346ea652bb806350ad87'}, 'Pt': {'symbol': 'Pt', 'hash': 'a604ea3c6a9cc23c739b762f625cf449'}, 'Pu': {'symbol': 'Pu', 'hash': 'f1d01e845dccc52d448679911f301a73'}, 'Rb': {'symbol': 'Rb_pv', 'hash': 'e447c648d870b066b3514e6b800727ab'}, 'Re': {'symbol': 'Re', 'hash': '72385e193c92a8acfe17ea49004c2be1'}, 'Rh': {'symbol': 'Rh', 'hash': '2c3dba3fcc6058ca1b1cfa75e45084bc'}, 'Ru': {'symbol': 'Ru_pv', 'hash': '7925f4d4b68076d70af7cd86eef9ba8d'}, 'S': {'symbol': 'S', 'hash': 'd368db6899d8839859bbee4811a42a88'}, 'Sb': {'symbol': 'Sb', 'hash': 'd82c022b02fc5344e85bd1909f9ee3e7'}, 'Sc': {'symbol': 'Sc_sv', 'hash': 'dc386f505ad0c43385a7715b4111cb75'}, 'Se': {'symbol': 'Se', 'hash': '67a8804ede9f1112726e3d136978ef19'}, 'Si': {'symbol': 'Si', 'hash': 'b2b0ea6feb62e7cde209616683b8f7f5'}, 'Sm': {'symbol': 'Sm_3', 'hash': 'e5e274e7cd99602ca81d146155abdf88'}, 'Sn': {'symbol': 'Sn_d', 'hash': '849b0795e148f93113a06be8fd5f5001'}, 'Sr': {'symbol': 'Sr_sv', 'hash': 'ca6a5429c120a0ab705824386a76fe5b'}, 'Ta': {'symbol': 'Ta', 'hash': 'd4e2cfe9338ef80da592d5bb9dc782c7'}, 'Tb': {'symbol': 'Tb_3', 'hash': '0790955c547003956c0fd4f080f7f508'}, 'Tc': {'symbol': 'Tc', 'hash': '9592642886319309a39d55c5717c6f48'}, 'Te': {'symbol': 'Te', 'hash': '72719856e22fb1d3032df6f96d98a0f2'}, 'Th': {'symbol': 'Th', 'hash': 'aea79f322180fa6f0bfa74cb2a156dcf'}, 'Ti': {'symbol': 'Ti', 'hash': 'c617e8b539c3f44a0ab6e8da2a92d318'}, 'Tl': {'symbol': 'Tl', 'hash': '2aa0d5406aaab7ebfbc761da382f1352'}, 'Tm': {'symbol': 'Tm_3', 'hash': '94a07cb7949b01305cb161da0cbfb492'}, 'U': {'symbol': 'U', 'hash': '72702eabbb1bc02b4167590dc848ed5d'}, 'V': {'symbol': 'V_pv', 'hash': '7f1297a2e1d963e2a4d81b61f85e4ded'}, 'W': {'symbol': 'W_pv', 'hash': '2a33e0d5c700640535f60ac0a12177ab'}, 'Xe': {'symbol': 'Xe', 'hash': '338472e581f58b41d37c002a5e22353b'}, 'Y': {'symbol': 'Y_sv', 'hash': '4ed187e77cd54f198bb88020278b143d'}, 'Yb': {'symbol': 'Yb', 'hash': '9f472bd422f640710f7d93e2d9ce89f4'}, 'Zn': {'symbol': 'Zn', 'hash': 'e35ee27f8483a63bb68dbc236a343af3'}, 'Zr': {'symbol': 'Zr', 'hash': 'd221d2c0bac4f8e81af2f5c42a314274'}}}
class MPHSEBSSet(structure, user_incar_settings=None, added_kpoints=None, mode='Uniform', reciprocal_density=None, kpoints_line_density=20, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPHSERelaxSet

Implementation of a VaspInputSet for HSE band structure computations. Remember that HSE band structures must be self-consistent in VASP. A band structure along symmetry lines for instance needs BOTH a uniform grid with appropriate weights AND a path along the lines with weight 0.

Thus, the “Uniform” mode is just like regular static SCF but allows adding custom kpoints (e.g., corresponding to known VBM/CBM) to the uniform grid that have zero weight (e.g., for better gap estimate).

The “Line” mode is just like Uniform mode, but additionally adds k-points along symmetry lines with zero weight.

Parameters:
  • structure (Structure) – Structure to compute
  • user_incar_settings (dict) – A dict specifying additional incar settings
  • added_kpoints (list) – a list of kpoints (list of 3 number list) added to the run. The k-points are in fractional coordinates
  • mode (str) – “Line” - generate k-points along symmetry lines for bandstructure. “Uniform” - generate uniform k-points grid
  • reciprocal_density (int) – k-point density to use for uniform mesh
  • kpoints_line_density (int) – k-point density for high symmetry lines
  • **kwargs (dict) – Any other parameters to pass into DictVaspInputSet
classmethod from_prev_calc(prev_calc_dir, mode='gap', reciprocal_density=50, copy_chgcar=True, **kwargs)[source]

Generate a set of Vasp input files for HSE calculations from a directory of previous Vasp run. if mode==”gap”, it explicitly adds VBM and CBM of the prev run to the k-point list of this run.

Parameters:
  • prev_calc_dir (str) – Directory containing the outputs (vasprun.xml and OUTCAR) of previous vasp run.
  • mode (str) – Either “uniform”, “gap” or “line”
  • reciprocal_density (int) – density of k-mesh
  • copy_chgcar (bool) – whether to copy CHGCAR of previous run
  • **kwargs – All kwargs supported by MPHSEBSStaticSet, other than prev_structure which is determined from the previous calc dir.
kpoints
class MPHSERelaxSet(structure, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.DictSet

Same as the MPRelaxSet, but with HSE parameters.

CONFIG = {'INCAR': {'ALGO': 'All', 'EDIFF_PER_ATOM': 5e-05, 'ENCUT': 520, 'HFSCREEN': 0.2, 'IBRION': 2, 'ICHARG': 1, 'ISIF': 3, 'ISMEAR': 0, 'ISPIN': 2, 'LHFCALC': True, 'LORBIT': 11, 'LREAL': 'AUTO', 'LWAVE': False, 'NELM': 100, 'NSW': 99, 'PREC': 'Accurate', 'PRECFOCK': 'Fast', 'SIGMA': 0.05, 'MAGMOM': {'Co': 5, 'Co3+': 0.6, 'Co4+': 1, 'Cr': 5, 'Fe': 5, 'Mn': 5, 'Mn3+': 4, 'Mn4+': 3, 'Mo': 5, 'Ni': 5, 'V': 5, 'W': 5, 'Ce': 5, 'Eu': 10, 'Eu2+': 7, 'Eu3+': 6}}, 'KPOINTS': {'reciprocal_density': 50}, 'POTCAR': {'Ac': 'Ac', 'Ag': 'Ag', 'Al': 'Al', 'Ar': 'Ar', 'As': 'As', 'Au': 'Au', 'B': 'B', 'Ba': 'Ba_sv', 'Be': 'Be_sv', 'Bi': 'Bi', 'Br': 'Br', 'C': 'C', 'Ca': 'Ca_sv', 'Cd': 'Cd', 'Ce': 'Ce', 'Cl': 'Cl', 'Co': 'Co', 'Cr': 'Cr_pv', 'Cs': 'Cs_sv', 'Cu': 'Cu_pv', 'Dy': 'Dy_3', 'Er': 'Er_3', 'Eu': 'Eu', 'F': 'F', 'Fe': 'Fe_pv', 'Ga': 'Ga_d', 'Gd': 'Gd', 'Ge': 'Ge_d', 'H': 'H', 'He': 'He', 'Hf': 'Hf_pv', 'Hg': 'Hg', 'Ho': 'Ho_3', 'I': 'I', 'In': 'In_d', 'Ir': 'Ir', 'K': 'K_sv', 'Kr': 'Kr', 'La': 'La', 'Li': 'Li_sv', 'Lu': 'Lu_3', 'Mg': 'Mg_pv', 'Mn': 'Mn_pv', 'Mo': 'Mo_pv', 'N': 'N', 'Na': 'Na_pv', 'Nb': 'Nb_pv', 'Nd': 'Nd_3', 'Ne': 'Ne', 'Ni': 'Ni_pv', 'Np': 'Np', 'O': 'O', 'Os': 'Os_pv', 'P': 'P', 'Pa': 'Pa', 'Pb': 'Pb_d', 'Pd': 'Pd', 'Pm': 'Pm_3', 'Pr': 'Pr_3', 'Pt': 'Pt', 'Pu': 'Pu', 'Rb': 'Rb_sv', 'Re': 'Re_pv', 'Rh': 'Rh_pv', 'Ru': 'Ru_pv', 'S': 'S', 'Sb': 'Sb', 'Sc': 'Sc_sv', 'Se': 'Se', 'Si': 'Si', 'Sm': 'Sm_3', 'Sn': 'Sn_d', 'Sr': 'Sr_sv', 'Ta': 'Ta_pv', 'Tb': 'Tb_3', 'Tc': 'Tc_pv', 'Te': 'Te', 'Th': 'Th', 'Ti': 'Ti_pv', 'Tl': 'Tl_d', 'Tm': 'Tm_3', 'U': 'U', 'V': 'V_pv', 'W': 'W_pv', 'Xe': 'Xe', 'Y': 'Y_sv', 'Yb': 'Yb_2', 'Zn': 'Zn', 'Zr': 'Zr_sv'}}
class MPNonSCFSet(structure, prev_incar=None, mode='line', nedos=601, reciprocal_density=100, sym_prec=0.1, kpoints_line_density=20, optics=False, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPRelaxSet

Init a MPNonSCFSet. Typically, you would use the classmethod from_prev_calc to initialize from a previous SCF run.

Parameters:
  • structure (Structure) – Structure to compute
  • prev_incar (Incar/string) – Incar file from previous run.
  • mode (str) – Line or Uniform mode supported.
  • nedos (int) – nedos parameter. Default to 601.
  • reciprocal_density (int) – density of k-mesh by reciprocal volume (defaults to 100)
  • sym_prec (float) – Symmetry precision (for Uniform mode).
  • kpoints_line_density (int) – Line density for Line mode.
  • optics (bool) – whether to add dielectric function
  • **kwargs – kwargs supported by MPVaspInputSet.
classmethod from_prev_calc(prev_calc_dir, copy_chgcar=True, nbands_factor=1.2, standardize=False, sym_prec=0.1, international_monoclinic=True, reciprocal_density=100, kpoints_line_density=20, small_gap_multiply=None, **kwargs)[source]

Generate a set of Vasp input files for NonSCF calculations from a directory of previous static Vasp run.

Parameters:
  • prev_calc_dir (str) – The directory contains the outputs( vasprun.xml and OUTCAR) of previous vasp run.
  • copy_chgcar – Whether to copy the old CHGCAR. Defaults to True.
  • nbands_factor (float) – Multiplicative factor for NBANDS. Choose a higher number if you are doing an LOPTICS calculation.
  • standardize (float) – Whether to standardize to a primitive standard cell. Defaults to False.
  • sym_prec (float) – Tolerance for symmetry finding. If not 0, the final structure from the previous run will be symmetrized to get a primitive standard cell. Set to 0 if you don’t want that.
  • international_monoclinic (bool) – Whether to use international convention (vs Curtarolo) for monoclinic. Defaults True.
  • reciprocal_density (int) – density of k-mesh by reciprocal volume in uniform mode (defaults to 100)
  • kpoints_line_density (int) – density of k-mesh in line mode (defaults to 20)
  • small_gap_multiply ([float, float]) – If the gap is less than 1st index, multiply the default reciprocal_density by the 2nd index.
  • **kwargs – All kwargs supported by MPNonSCFSet, other than structure, prev_incar and prev_chgcar which are determined from the prev_calc_dir.
incar
kpoints
class MPRelaxSet(structure, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.DictSet

Implementation of VaspInputSet utilizing parameters in the public Materials Project. Typically, the pseudopotentials chosen contain more electrons than the MIT parameters, and the k-point grid is ~50% more dense. The LDAUU parameters are also different due to the different psps used, which result in different fitted values.

CONFIG = {'INCAR': {'ALGO': 'FAST', 'EDIFF_PER_ATOM': 5e-05, 'ENCUT': 520, 'IBRION': 2, 'ICHARG': 1, 'ISIF': 3, 'ISMEAR': -5, 'ISPIN': 2, 'LDAU': True, 'LDAUJ': {'F': {'Co': 0, 'Cr': 0, 'Fe': 0, 'Mn': 0, 'Mo': 0, 'Ni': 0, 'V': 0, 'W': 0}, 'O': {'Co': 0, 'Cr': 0, 'Fe': 0, 'Mn': 0, 'Mo': 0, 'Ni': 0, 'V': 0, 'W': 0}}, 'LDAUL': {'F': {'Co': 2, 'Cr': 2, 'Fe': 2, 'Mn': 2, 'Mo': 2, 'Ni': 2, 'V': 2, 'W': 2}, 'O': {'Co': 2, 'Cr': 2, 'Fe': 2, 'Mn': 2, 'Mo': 2, 'Ni': 2, 'V': 2, 'W': 2}}, 'LDAUTYPE': 2, 'LDAUU': {'F': {'Co': 3.32, 'Cr': 3.7, 'Fe': 5.3, 'Mn': 3.9, 'Mo': 4.38, 'Ni': 6.2, 'V': 3.25, 'W': 6.2}, 'O': {'Co': 3.32, 'Cr': 3.7, 'Fe': 5.3, 'Mn': 3.9, 'Mo': 4.38, 'Ni': 6.2, 'V': 3.25, 'W': 6.2}}, 'LDAUPRINT': 1, 'LORBIT': 11, 'LREAL': 'AUTO', 'LWAVE': False, 'NELM': 100, 'NSW': 99, 'PREC': 'Accurate', 'SIGMA': 0.05, 'MAGMOM': {'Co': 5, 'Co3+': 0.6, 'Co4+': 1, 'Cr': 5, 'Fe': 5, 'Mn': 5, 'Mn3+': 4, 'Mn4+': 3, 'Mo': 5, 'Ni': 5, 'V': 5, 'W': 5, 'Ce': 5, 'Eu': 10, 'Eu2+': 7, 'Eu3+': 6}}, 'KPOINTS': {'reciprocal_density': 64}, 'POTCAR': {'Ac': 'Ac', 'Ag': 'Ag', 'Al': 'Al', 'Ar': 'Ar', 'As': 'As', 'Au': 'Au', 'B': 'B', 'Ba': 'Ba_sv', 'Be': 'Be_sv', 'Bi': 'Bi', 'Br': 'Br', 'C': 'C', 'Ca': 'Ca_sv', 'Cd': 'Cd', 'Ce': 'Ce', 'Cl': 'Cl', 'Co': 'Co', 'Cr': 'Cr_pv', 'Cs': 'Cs_sv', 'Cu': 'Cu_pv', 'Dy': 'Dy_3', 'Er': 'Er_3', 'Eu': 'Eu', 'F': 'F', 'Fe': 'Fe_pv', 'Ga': 'Ga_d', 'Gd': 'Gd', 'Ge': 'Ge_d', 'H': 'H', 'He': 'He', 'Hf': 'Hf_pv', 'Hg': 'Hg', 'Ho': 'Ho_3', 'I': 'I', 'In': 'In_d', 'Ir': 'Ir', 'K': 'K_sv', 'Kr': 'Kr', 'La': 'La', 'Li': 'Li_sv', 'Lu': 'Lu_3', 'Mg': 'Mg_pv', 'Mn': 'Mn_pv', 'Mo': 'Mo_pv', 'N': 'N', 'Na': 'Na_pv', 'Nb': 'Nb_pv', 'Nd': 'Nd_3', 'Ne': 'Ne', 'Ni': 'Ni_pv', 'Np': 'Np', 'O': 'O', 'Os': 'Os_pv', 'P': 'P', 'Pa': 'Pa', 'Pb': 'Pb_d', 'Pd': 'Pd', 'Pm': 'Pm_3', 'Pr': 'Pr_3', 'Pt': 'Pt', 'Pu': 'Pu', 'Rb': 'Rb_sv', 'Re': 'Re_pv', 'Rh': 'Rh_pv', 'Ru': 'Ru_pv', 'S': 'S', 'Sb': 'Sb', 'Sc': 'Sc_sv', 'Se': 'Se', 'Si': 'Si', 'Sm': 'Sm_3', 'Sn': 'Sn_d', 'Sr': 'Sr_sv', 'Ta': 'Ta_pv', 'Tb': 'Tb_3', 'Tc': 'Tc_pv', 'Te': 'Te', 'Th': 'Th', 'Ti': 'Ti_pv', 'Tl': 'Tl_d', 'Tm': 'Tm_3', 'U': 'U', 'V': 'V_pv', 'W': 'W_pv', 'Xe': 'Xe', 'Y': 'Y_sv', 'Yb': 'Yb_2', 'Zn': 'Zn', 'Zr': 'Zr_sv'}}
class MPSOCSet(structure, saxis=(0, 0, 1), prev_incar=None, reciprocal_density=100, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPStaticSet

Init a MPSOCSet.

Parameters:
  • structure (Structure) – the structure must have the ‘magmom’ site property and each magnetic moment value must have 3 components. eg:- magmom = [[0,0,2], …]
  • saxis (tuple) – magnetic moment orientation
  • prev_incar (Incar) – Incar file from previous run.
  • reciprocal_density (int) – density of k-mesh by reciprocal volume (defaults to 100)
  • **kwargs – kwargs supported by MPVaspInputSet.
classmethod from_prev_calc(prev_calc_dir, copy_chgcar=True, nbands_factor=1.2, standardize=False, sym_prec=0.1, international_monoclinic=True, reciprocal_density=100, small_gap_multiply=None, **kwargs)[source]

Generate a set of Vasp input files for SOC calculations from a directory of previous static Vasp run. SOC calc requires all 3 components for MAGMOM for each atom in the structure.

Parameters:
  • prev_calc_dir (str) – The directory contains the outputs( vasprun.xml and OUTCAR) of previous vasp run.
  • copy_chgcar – Whether to copy the old CHGCAR. Defaults to True.
  • nbands_factor (float) – Multiplicative factor for NBANDS. Choose a higher number if you are doing an LOPTICS calculation.
  • standardize (float) – Whether to standardize to a primitive standard cell. Defaults to False.
  • sym_prec (float) – Tolerance for symmetry finding. If not 0, the final structure from the previous run will be symmetrized to get a primitive standard cell. Set to 0 if you don’t want that.
  • international_monoclinic (bool) – Whether to use international convention (vs Curtarolo) for monoclinic. Defaults True.
  • reciprocal_density (int) – density of k-mesh by reciprocal volume (defaults to 100)
  • small_gap_multiply ([float, float]) – If the gap is less than 1st index, multiply the default reciprocal_density by the 2nd index.
  • **kwargs – All kwargs supported by MPSOCSet, other than structure, prev_incar and prev_chgcar which are determined from the prev_calc_dir.
incar
class MPStaticSet(structure, prev_incar=None, prev_kpoints=None, lepsilon=False, lcalcpol=False, reciprocal_density=100, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPRelaxSet

Run a static calculation.

Parameters:
  • structure (Structure) – Structure from previous run.
  • prev_incar (Incar) – Incar file from previous run.
  • prev_kpoints (Kpoints) – Kpoints from previous run.
  • lepsilon (bool) – Whether to add static dielectric calculation
  • reciprocal_density (int) – For static calculations, we usually set the reciprocal density by volume. This is a convenience arg to change that, rather than using user_kpoints_settings. Defaults to 100, which is ~50% more than that of standard relaxation calculations.
  • **kwargs – kwargs supported by MPRelaxSet.
classmethod from_prev_calc(prev_calc_dir, standardize=False, sym_prec=0.1, international_monoclinic=True, reciprocal_density=100, small_gap_multiply=None, **kwargs)[source]

Generate a set of Vasp input files for static calculations from a directory of previous Vasp run.

Parameters:
  • prev_calc_dir (str) – Directory containing the outputs( vasprun.xml and OUTCAR) of previous vasp run.
  • standardize (float) – Whether to standardize to a primitive standard cell. Defaults to False.
  • sym_prec (float) – Tolerance for symmetry finding. If not 0, the final structure from the previous run will be symmetrized to get a primitive standard cell. Set to 0 if you don’t want that.
  • international_monoclinic (bool) – Whether to use international convention (vs Curtarolo) for monoclinic. Defaults True.
  • reciprocal_density (int) – density of k-mesh by reciprocal volume (defaults to 100)
  • small_gap_multiply ([float, float]) – If the gap is less than 1st index, multiply the default reciprocal_density by the 2nd index.
  • **kwargs – All kwargs supported by MPStaticSet, other than prev_incar and prev_structure and prev_kpoints which are determined from the prev_calc_dir.
incar
kpoints
class MVLElasticSet(structure, potim=0.015, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPRelaxSet

MVL denotes VASP input sets that are implemented by the Materials Virtual Lab (http://www.materialsvirtuallab.org) for various research.

This input set is used to calculate elastic constants in VASP. It is used in the following work:

Z. Deng, Z. Wang, I.-H. Chu, J. Luo, S. P. Ong.
“Elastic Properties of Alkali Superionic Conductor Electrolytes
from First Principles Calculations”, J. Electrochem. Soc.
2016, 163(2), A67-A74. doi: 10.1149/2.0061602jes

To read the elastic constants, you may use the Outcar class which parses the elastic constants.

Parameters:
  • scale (float) – POTIM parameter. The default of 0.015 is usually fine, but some structures may require a smaller step.
  • user_incar_settings (dict) – A dict specifying additional incar settings.
class MVLGBSet(structure, k_product=40, slab_mode=False, is_metal=True, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPRelaxSet

Class for writing a vasp input files for grain boundary calculations, slab or bulk.

Parameters:
  • structure (Structure) – provide the structure
  • k_product – Kpoint number * length for a & b directions, also for c direction in bulk calculations. Default to 40.
  • slab_mode (bool) – Defaults to False. Use default (False) for a bulk supercell. Use True if you are performing calculations on a slab-like (i.e., surface) of the GB, for example, when you are calculating the work of separation.
  • is_metal (bool) – Defaults to True. This determines whether an ISMEAR of 1 is used (for metals) or not (for insulators and semiconductors) by default. Note that it does not override user_incar_settings, which can be set by the user to be anything desired.
  • **kwargs – Other kwargs supported by DictVaspInputSet.
incar
kpoints

k_product, default to 40, is kpoint number * length for a & b directions, also for c direction in bulk calculations Automatic mesh & Gamma is the default setting.

class MVLGWSet(structure, prev_incar=None, nbands=None, potcar_functional='PBE_54', reciprocal_density=100, mode='STATIC', **kwargs)[source]

Bases: pymatgen.io.vasp.sets.DictSet

MVL denotes VASP input sets that are implemented by the Materials Virtual Lab (http://www.materialsvirtuallab.org) for various research. This is a flexible input set for GW calculations.

Note that unlike all other input sets in this module, the PBE_54 series of functional is set as the default. These have much improved performance for GW calculations.

A typical sequence is mode=”STATIC” -> mode=”DIAG” -> mode=”GW” -> mode=”BSE”. For all steps other than the first one (static), the recommendation is to use from_prev_calculation on the preceding run in the series.

Parameters:
  • structure (Structure) – Input structure.
  • prev_incar (Incar/string) – Incar file from previous run.
  • mode (str) – Supported modes are “STATIC” (default), “DIAG”, “GW”, and “BSE”.
  • nbands (int) – For subsequent calculations, it is generally recommended to perform NBANDS convergence starting from the NBANDS of the previous run for DIAG, and to use the exact same NBANDS for GW and BSE. This parameter is used by from_previous_calculation to set nband.
  • potcar_functional (str) – Defaults to “PBE_54”.
  • **kwargs – All kwargs supported by DictSet. Typically, user_incar_settings is a commonly used option.
CONFIG = {'INCAR': {'ALGO': 'Normal', 'EDIFF': 1e-08, 'IBRION': -1, 'ICHARG': 1, 'ISMEAR': 0, 'ISPIN': 2, 'LORBIT': 11, 'LREAL': 'AUTO', 'LWAVE': True, 'NELM': 100, 'PREC': 'Accurate', 'SIGMA': 0.01, 'MAGMOM': {'Co': 5, 'Co3+': 0.6, 'Co4+': 1, 'Cr': 5, 'Fe': 5, 'Mn': 5, 'Mn3+': 4, 'Mn4+': 3, 'Mo': 5, 'Ni': 5, 'V': 5, 'W': 5, 'Ce': 5, 'Eu': 10, 'Eu2+': 7, 'Eu3+': 6}}, 'KPOINTS': {'reciprocal_density': 100}, 'POTCAR': {'Ac': 'Ac', 'Ag': 'Ag_sv_GW', 'Al': 'Al_GW', 'Ar': 'Ar_GW', 'As': 'As_GW', 'At': 'At_d_GW', 'Au': 'Au_sv_GW', 'B': 'B_GW', 'Ba': 'Ba_sv_GW', 'Be': 'Be_sv_GW', 'Bi': 'Bi_d_GW', 'Br': 'Br_GW', 'C': 'C_GW', 'Ca': 'Ca_sv_GW', 'Cd': 'Cd_sv_GW', 'Ce': 'Ce_GW', 'Cl': 'Cl_GW', 'Co': 'Co_sv_GW', 'Cr': 'Cr_sv_GW', 'Cs': 'Cs_sv_GW', 'Cu': 'Cu_sv_GW', 'Dy': 'Dy_3', 'Er': 'Er_3', 'Eu': 'Eu', 'F': 'F_GW', 'Fe': 'Fe_sv_GW', 'Ga': 'Ga_d_GW', 'Gd': 'Gd', 'Ge': 'Ge_d_GW', 'H': 'H_GW', 'He': 'He_GW', 'Hf': 'Hf_sv_GW', 'Hg': 'Hg_sv_GW', 'Ho': 'Ho_3', 'I': 'I_GW', 'In': 'In_d_GW', 'Ir': 'Ir_sv_GW', 'K': 'K_sv_GW', 'Kr': 'Kr_GW', 'La': 'La_GW', 'Li': 'Li_sv_GW', 'Lu': 'Lu_3', 'Mg': 'Mg_sv_GW', 'Mn': 'Mn_sv_GW', 'Mo': 'Mo_sv_GW', 'N': 'N_GW', 'Na': 'Na_sv_GW', 'Nb': 'Nb_sv_GW', 'Nd': 'Nd_3', 'Ne': 'Ne_GW', 'Ni': 'Ni_sv_GW', 'Np': 'Np', 'O': 'O_GW', 'Os': 'Os_sv_GW', 'P': 'P_GW', 'Pa': 'Pa', 'Pb': 'Pb_d_GW', 'Pd': 'Pd_sv_GW', 'Pm': 'Pm_3', 'Po': 'Po_d_GW', 'Pr': 'Pr_3', 'Pt': 'Pt_sv_GW', 'Pu': 'Pu', 'Rb': 'Rb_sv_GW', 'Re': 'Re_sv_GW', 'Rh': 'Rh_sv_GW', 'Rn': 'Rn_d_GW', 'Ru': 'Ru_sv_GW', 'S': 'S_GW', 'Sb': 'Sb_d_GW', 'Sc': 'Sc_sv_GW', 'Se': 'Se_GW', 'Si': 'Si_GW', 'Sm': 'Sm_3', 'Sn': 'Sn_d_GW', 'Sr': 'Sr_sv_GW', 'Ta': 'Ta_sv_GW', 'Tb': 'Tb_3', 'Tc': 'Tc_sv_GW', 'Te': 'Te_GW', 'Th': 'Th', 'Ti': 'Ti_sv_GW', 'Tl': 'Tl_d_GW', 'Tm': 'Tm_3', 'U': 'U', 'V': 'V_sv_GW', 'W': 'W_sv_GW', 'Xe': 'Xe_GW', 'Y': 'Y_sv_GW', 'Yb': 'Yb_2', 'Zn': 'Zn_sv_GW', 'Zr': 'Zr_sv_GW'}}
SUPPORTED_MODES = ('DIAG', 'GW', 'STATIC', 'BSE')
classmethod from_prev_calc(prev_calc_dir, copy_wavecar=True, mode='DIAG', nbands_factor=5, ncores=16, **kwargs)[source]

Generate a set of Vasp input files for GW or BSE calculations from a directory of previous Exact Diag Vasp run.

Parameters:
  • prev_calc_dir (str) – The directory contains the outputs( vasprun.xml of previous vasp run.
  • copy_wavecar – Whether to copy the old WAVECAR, WAVEDER and associated files. Defaults to True.
  • mode (str) – Supported modes are “STATIC”, “DIAG” (default), “GW”, and “BSE”.
  • nbands_factor (int) – Multiplicative factor for NBANDS. Only applies if mode==”DIAG”. Need to be tested for convergence.
  • ncores (int) – numbers of cores you do calculations. VASP will alter NBANDS if it was not dividable by ncores. Only applies if mode==”DIAG”.
  • **kwargs – All kwargs supported by MVLGWSet, other than structure, prev_incar and mode, which are determined from the prev_calc_dir.
incar
kpoints

Generate gamma center k-points mesh grid for GW calc, which is requested by GW calculation.

class MVLNPTMDSet(structure, start_temp, end_temp, nsteps, time_step=2, spin_polarized=False, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MITMDSet

Class for writing a vasp md run in NPT ensemble.

Notes

To eliminate Pulay stress, the default ENCUT is set to a rather large value of ENCUT, which is 1.5 * ENMAX.

Parameters:
  • structure (Structure) – input structure.
  • start_temp (int) – Starting temperature.
  • end_temp (int) – Final temperature.
  • nsteps (int) – Number of time steps for simulations. The NSW parameter.
  • time_step (int) – The time step for the simulation. The POTIM parameter. Defaults to 2fs.
  • spin_polarized (bool) – Whether to do spin polarized calculations. The ISPIN parameter. Defaults to False.
  • **kwargs – Other kwargs supported by DictSet.
class MVLSlabSet(structure, k_product=50, bulk=False, auto_dipole=False, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.MPRelaxSet

Class for writing a set of slab vasp runs, including both slabs (along the c direction) and orient unit cells (bulk), to ensure the same KPOINTS, POTCAR and INCAR criterion.

Parameters:
  • k_product – default to 50, kpoint number * length for a & b directions, also for c direction in bulk calculations
  • bulk (bool) – Set to True for bulk calculation. Defaults to False.
  • **kwargs – Other kwargs supported by DictSet.
kpoints
k_product, default to 50, is kpoint number * length for a & b
directions, also for c direction in bulk calculations

Automatic mesh & Gamma is the default setting.

class VaspInputSet[source]

Bases: monty.json.MSONable

Base class representing a set of Vasp input parameters with a structure supplied as init parameters. Typically, you should not inherit from this class. Start from DictSet or MPRelaxSet or MITRelaxSet.

all_input

Returns all input files as a dict of {filename – vasp object}

Returns:object}, e.g., {‘INCAR’: Incar object, …}
Return type:dict of {filename
as_dict(verbosity=2)[source]
incar

Incar object

kpoints

Kpoints object

poscar

Poscar object

potcar

Potcar object.

potcar_symbols

List of POTCAR symbols.

write_input(output_dir, make_dir_if_not_present=True, include_cif=False)[source]

Writes a set of VASP input to a directory.

Parameters:
  • output_dir (str) – Directory to output the VASP input files
  • make_dir_if_not_present (bool) – Set to True if you want the directory (and the whole path) to be created if it is not present.
  • include_cif (bool) – Whether to write a CIF file in the output directory for easier opening by VESTA.
batch_write_input(structures, vasp_input_set=<class 'pymatgen.io.vasp.sets.MPRelaxSet'>, output_dir='.', make_dir_if_not_present=True, subfolder=None, sanitize=False, include_cif=False, **kwargs)[source]

Batch write vasp input for a sequence of structures to output_dir, following the format output_dir/{group}/{formula}_{number}.

Parameters:
  • structures ([Structure]) – Sequence of Structures.
  • vasp_input_set (VaspInputSet) – VaspInputSet class that creates vasp input files from structures. Note that a class should be supplied. Defaults to MPRelaxSet.
  • output_dir (str) – Directory to output files. Defaults to current directory “.”.
  • make_dir_if_not_present (bool) – Create the directory if not present. Defaults to True.
  • subfolder (callable) – Function to create subdirectory name from structure. Defaults to simply “formula_count”.
  • sanitize (bool) – Boolean indicating whether to sanitize the structure before writing the VASP input files. Sanitized output are generally easier for viewing and certain forms of analysis. Defaults to False.
  • include_cif (bool) – Whether to output a CIF as well. CIF files are generally better supported in visualization programs.
  • **kwargs – Additional kwargs are passed to the vasp_input_set class in addition to structure.
get_structure_from_prev_run(vasprun, outcar=None, sym_prec=0.1, international_monoclinic=True)[source]

Process structure from previous run.

Parameters:
  • vasprun (Vasprun) – Vasprun that contains the final structure from previous run.
  • outcar (Outcar) – Outcar that contains the magnetization info from previous run.
  • sym_prec (float) – Tolerance for symmetry finding for standardization. If no standardization is desired, set to 0 or a False.
  • international_monoclinic (bool) – Whether to use international convention (vs Curtarolo) for monoclinic. Defaults True.
Returns:

Returns the magmom-decorated structure that can be passed to get Vasp input files, e.g. get_kpoints.

get_vasprun_outcar(path, parse_dos=True, parse_eigen=True)[source]