mpinterfaces package

Subpackages

Submodules

mpinterfaces.calibrate module

class mpinterfaces.calibrate.Calibrate(incar, poscar, potcar, kpoints, system=None, is_matrix=False, Grid_type='A', parent_job_dir='.', job_dir='Job', qadapter=None, job_cmd='qsub', wait=True, mappings_override=None, functional='PBE', database=None, magnetism=None, mag_init=None, reuse=None, reuse_override=None, reuse_incar=None, solvation=None, turn_knobs=OrderedDict([('ENCUT', []), ('KPOINTS', [])]), checkpoint_file=None, cal_logger=None)[source]

Bases: monty.json.MSONable

The base class for creating vasp work flows for calibrating the input parameters for different systems

A wrapper around Custodian

LOG_FILE = 'calibrate.json'
add_job(name='noname', job_dir='.')[source]

add a single job using the current incar, poscar, potcar and kpoints

as_dict()[source]
classmethod from_dict(d)[source]
key_to_name(key)[source]

convenient string mapping for the keys in the turn_knobs dict

Parameters:key – key to the knob dict
Returns:an appropriate string representation of the key so that the name doesnt clash with the filenames
kpoint_to_name(kpoint, grid_type)[source]

get a string representation for the given kpoint

Parameters:
  • kpoint – an iterable
  • grid_type – grid_type used for the KPOINTS
Returns:

Monkhorst Pack 2 2 2 will be named 2x2x2

Return type:

string representation for kpoint eg

potcar_to_name(mapping=None, functional=None)[source]

convert a symbol mapping and functional to a name that can be used for setting up the potcar jobs

Args:
mapping: example:- if mapping = {‘Pt’:’Pt_pv’,
‘Si’:’Si_GW’} then the name will be PBE_Pt_pv_Si_GW
with self.functional=”PBE”
Returns:string
recursive_jobs(n, keys, i)[source]

recursively setup the jobs: used by setup_matrix_job

Parameters:
  • n – total number of knobs aka parameters to be tuned
  • keys – list of knobs i.e parameter names
  • i – ith knob
run(job_cmd=None)[source]

run the vasp jobs through custodian if the job list is empty, run a single job with the initial input set

set_incar(param, val)[source]

set the incar paramter, param = val

set_kpoints(kpoint=None, poscar=None, ibzkpth=None)[source]

set the kpoint

set_poscar(scale=None, poscar=None)[source]

perturbs given structure by volume scaling factor or takes user defined variants of Poscar

Parameters:
  • scale – Volume Scaling parameter
  • poscar – Poscar object of user defined structure
  • the poscar (set) – volume scaled by the scale factor
set_potcar(mapping=None, functional='PBE')[source]

set the potcar: symbol to potcar type mapping

setup()[source]

set up the jobs for the given turn_knobs dict is_matrix = True implies that the params in the dict are interrelated. Otherwise calcs corresponding to each dict key is independent

setup_incar_jobs(param, val_list)[source]

set up incar jobs,, calls set_incar to set the value to param

Parameters:
  • param – Name of INCAR parameter
  • val_list – List of values to vary for the param
setup_kpoints_jobs(kpoints_list=None)[source]

setup the kpoint jobs

setup_matrix_job()[source]

set up jobs where the dict keys are interrelated mind: its an ordered dict, the order in which the keys are specified determines the nested directory structure

setup_poscar_jobs(scale_list=None, poscar_list=None)[source]

for scaling the latice vectors of the original structure, scale_list is volume scaling factor list

setup_potcar_jobs(mappings, functional_list)[source]

take a list of symbol mappings and setup the potcar jobs

val_to_name(val)[source]

convert a value to a string so that it can be used for naming the job directory the decimal points in floats are replaced with underscore character if the value is of type list, kpoint_to_name method is used since only kpoint values are expected to be of type list if the values is of type dict then potcar_to_name method is invoked

Parameters:val – knob value to be converted into an appropriate string representation
Returns:a string filename for the value
class mpinterfaces.calibrate.CalibrateBulk(incar, poscar, potcar, kpoints, system=None, is_matrix=False, Grid_type='A', parent_job_dir='.', job_dir='./Bulk', qadapter=None, job_cmd='qsub', wait=True, mappings_override=None, functional='PBE', turn_knobs={'ENCUT': [], 'KPOINTS': []}, checkpoint_file=None, cal_logger=None)[source]

Bases: mpinterfaces.calibrate.Calibrate

Calibrate parameters for Bulk calculations

class mpinterfaces.calibrate.CalibrateInterface(incar, poscar, potcar, kpoints, system=None, is_matrix=False, Grid_type='A', parent_job_dir='.', job_dir='./Interface', qadapter=None, job_cmd='qsub', wait=True, mappings_override=None, functional='PBE', turn_knobs={'VACUUM': [], 'THICKNESS': []}, from_ase=False, checkpoint_file=None, cal_logger=None)[source]

Bases: mpinterfaces.calibrate.CalibrateSlab

Calibrate paramters for interface calculations

create_interface()[source]

add params that you want to vary

interface_setup(turn_knobs=None)[source]
class mpinterfaces.calibrate.CalibrateMolecule(incar, poscar, potcar, kpoints, system=None, is_matrix=False, Grid_type='A', parent_job_dir='.', job_dir='./Molecule', qadapter=None, job_cmd='qsub', wait=True, mappings_override=None, functional='PBE', turn_knobs={'ENCUT': [], 'KPOINTS': []}, checkpoint_file=None, cal_logger=None)[source]

Bases: mpinterfaces.calibrate.Calibrate

Calibrate paramters for Molecule calculations

setup_kpoints_jobs(Grid_type='M', kpoints_list=None, conv_step=1)[source]
class mpinterfaces.calibrate.CalibrateSlab(incar, poscar, potcar, kpoints, system=None, is_matrix=False, Grid_type='A', parent_job_dir='.', job_dir='./Slab', qadapter=None, job_cmd='qsub', wait=True, mappings_override=None, functional='PBE', turn_knobs={'VACUUM': [], 'THICKNESS': []}, from_ase=False, checkpoint_file=None, cal_logger=None)[source]

Bases: mpinterfaces.calibrate.Calibrate

Calibrate paramters for Slab calculations

create_slab(vacuum=12, thickness=10)[source]

set the vacuum spacing, slab thickness and call sd_flags for top 2 layers

returns the poscar corresponding to the modified structure

set_reconstructed_surface(sites_to_add)[source]

Append sites as needed for reconstruction TODO

static set_sd_flags(interface=None, n_layers=2, top=True, bottom=True)[source]

set the relaxation flags for top and bottom layers of interface.

The upper and lower bounds of the z coordinate are determined based on the slab. All layers above and below the bounds will be relaxed. This means if there is a ligand on top of the slab, all of its atoms will also be relaxed.

setup_thickness_jobs(thickness_list)[source]

create slabs with the provided thickness settings

returns list of poscars

setup_vacuum_jobs(vacuum_list)[source]

create slabs with the provided vacuum settings

returns list of poscars

slab_setup(turn_knobs=None)[source]

invoke the set up methods corresponding to the dict keys: VACUUM and THICKNESS sets the POSCAR key in the turn_knobs

mpinterfaces.data_processor module

class mpinterfaces.data_processor.MPINTComputedEntry(structure, kpoints, incar, energy, correction=0.0, parameters=None, data=None, entry_id=None)[source]

Bases: pymatgen.entries.computed_entries.ComputedEntry

extend ComputedEntry to include structure as well as kpoints

as_dict()[source]
classmethod from_dict(d)[source]
class mpinterfaces.data_processor.MPINTVaspDrone(inc_structure=False, inc_incar_n_kpoints=False, parameters=None, data=None)[source]

Bases: pymatgen.apps.borg.hive.VaspToComputedEntryDrone

extend VaspToComputedEntryDrone to use the custom Vasprun: MPINTVasprun

as_dict()[source]
assimilate(path)[source]
classmethod from_dict(d)[source]
class mpinterfaces.data_processor.MPINTVasprun(filename, ionic_step_skip=None, ionic_step_offset=0, parse_dos=True, parse_eigen=True, parse_projected_eigen=False, parse_potcar_file=True)[source]

Bases: pymatgen.io.vasp.outputs.Vasprun

Extend Vasprun to use custom ComputedEntry: MPINTComputedEntry

get_computed_entry(inc_structure=False, inc_incar_n_kpoints=False, parameters=None, data=None)[source]

Returns a ComputedEntry from the vasprun.

Parameters:
  • inc_structure (bool) – Set to True if you want ComputedStructureEntries to be returned instead of ComputedEntries.
  • inc_incar_n_kpoints (bool) – along with inc_structure set to True if you want MPINTComputedEntries to be returned
  • parameters (list) – Input parameters to include. It has to be one of the properties supported by the Vasprun object. If parameters == None, a default set of parameters that are necessary for typical post-processing will be set.
  • data (list) – Output data to include. Has to be one of the properties supported by the Vasprun object.
Returns:

ComputedStructureEntry/ComputedEntry

mpinterfaces.database module

class mpinterfaces.database.MPINTVaspToDbTaskDrone(host='127.0.0.1', port=27017, database='vasp', user=None, password=None, collection='nanoparticles', parse_dos=False, compress_dos=False, simulate_mode=False, additional_fields=None, update_duplicates=True, mapi_key=None, use_full_uri=True, runs=None)[source]

Bases: matgendb.creator.VaspToDbTaskDrone

subclassing VaspToDbTaskDrone

generate_doc(dir_name, vasprun_files)[source]

Overridden

post_process(dir_name, d)[source]
customization:
adds system.json to the dictionary
mpinterfaces.database.analysis_and_error_checks(d, max_force_threshold=0.5, volume_change_threshold=0.2)[source]
mpinterfaces.database.get_uri(dir_name)[source]

Customized version of the original pymatgen-db version. Customization required because same job folder on hipergator gets different uri for different login nodes .

Returns the URI path for a directory. This allows files hosted on different file servers to have distinct locations. :param dir_name: A directory name.

Returns:/full/path/of/dir_name.
Return type:Full URI path, e.g., fileserver.host.com

mpinterfaces.default_logger module

Creates a default logger for outputing text.

mpinterfaces.default_logger.get_default_logger(name, output_stream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)[source]

mpinterfaces.firetasks module

class mpinterfaces.firetasks.MPINTCalibrateTask(*args, **kwargs)[source]

Bases: fireworks.core.firework.FireTaskBase

Calibration Task

optional_params = ['que_params']
run_task(fw_spec)[source]

launch jobs to the queue

class mpinterfaces.firetasks.MPINTDatabaseTask(*args, **kwargs)[source]

Bases: fireworks.core.firework.FireTaskBase, fireworks.utilities.fw_serializers.FWSerializable

submit data to the database firetask

optional_params = ['dbase_params']
required_params = ['measure_dir']
run_task(fw_spec)[source]

go through the measurement job dirs and put the measurement jobs in the database

class mpinterfaces.firetasks.MPINTMeasurementTask(*args, **kwargs)[source]

Bases: fireworks.core.firework.FireTaskBase, fireworks.utilities.fw_serializers.FWSerializable

Measurement Task

optional_params = ['que_params', 'job_cmd', 'other_params', 'fw_id']
required_params = ['measurement']
run_task(fw_spec)[source]

setup up a measurement task using the prior calibration jobs and run

mpinterfaces.firetasks.get_cal_obj(d)[source]

construct a calibration object from the input dictionary, d

returns a calibration object

mpinterfaces.firetasks.load_class(mod, name)[source]

load class named name from module mod this function is adapted from materialsproject

mpinterfaces.instrument module

class mpinterfaces.instrument.MPINTJob(job_cmd, name='noname', output_file='job.out', parent_job_dir='.', job_dir='untitled', suffix='', final=True, gzipped=False, backup=False, vis=None, auto_npar=True, settings_override=None, wait=True, vjob_logger=None)[source]

Bases: custodian.custodian.Job

defines a job i.e setup the required input files and launch the job

Parameters:
  • job_cmd – a list, the command to be issued in each job_dir eg: [‘qsub’, ‘submit_job’]
  • job_dir – the directory from which the jobs will be launched
as_dict()[source]
classmethod from_dict(d)[source]
name()[source]
postprocess()[source]
run()[source]

move to the job_dir, launch the job and back to the parent job directory

setup()[source]

write the input files to the job_dir

class mpinterfaces.instrument.MPINTVaspErrors[source]

Bases: custodian.custodian.ErrorHandler

handles restarting of jobs that exceed the walltime employs the check + correct method of custodian ErrorHandler

class mpinterfaces.instrument.MPINTVaspInputSet(name, incar, poscar, potcar, kpoints, qadapter=None, script_name='submit_script', vis_logger=None, reuse_path=None, **kwargs)[source]

Bases: pymatgen.io.vasp.sets.DictSet

defines the set of input required for a vasp job i.e create INCAR, POSCAR, POTCAR & KPOINTS files

as_dict()[source]
classmethod from_dict(d)[source]
write_input(job_dir, make_dir_if_not_present=True, write_cif=False)[source]

the input files are written to the job_dir process(if needed) and write the input files in each directory structures read from the poscar files in the directory

class mpinterfaces.instrument.MPINTVaspJob(job_cmd, name='noname', output_file='job.out', parent_job_dir='.', job_dir='untitled', suffix='', final=True, gzipped=False, backup=False, vis=None, auto_npar=True, settings_override=None, wait=True, vjob_logger=None)[source]

Bases: mpinterfaces.instrument.MPINTJob

defines a vasp job i.e setup the required input files and launch the job

Parameters:
  • job_cmd – a list, the command to be issued in each job_dir eg: [‘qsub’, ‘submit_job’]
  • job_dir – the directory from which the jobs will be launched
get_final_energy()[source]

mpinterfaces.interface module

class mpinterfaces.interface.Interface(strt, hkl=[1, 1, 1], min_thick=10, min_vac=10, supercell=[1, 1, 1], name=None, adsorb_on_species=None, adatom_on_lig=None, ligand=None, displacement=1.0, surface_coverage=None, scell_nmax=10, coverage_tol=0.25, solvent=None, start_from_slab=False, validate_proximity=False, to_unit_cell=False, coords_are_cartesian=False, primitive=True, from_ase=False, lll_reduce=False, center_slab=True, max_normal_search=None, force_normalize=False, x_shift=0, y_shift=0, rot=[0, 0, 0])[source]

Bases: pymatgen.core.surface.Slab

Interface = slab + ligand + environment(solvent) Creates a Slab - Ligand Interface of given coverage and given slab-ligand displacement

Parameters:
  • strt – Starting Structure Object for Slab of the Interface
  • hkl – Miller Index of Slab
  • min_thick – Minimum Slab Thickness in Angstroms desired
  • min_vac – Minimum Vacuum Spacing (Padding top and bottom, each) in Angstroms
  • supercell – Trial supercell to start with to enforce coverage, default 1x1x1
  • name – System name to specify database entry (can be a combination of miller indices of slab and ligand and solvent) eg: “PbS [1,1,1] + Hydrazine in DMF (epsilon = 37.5)”
  • adsorb_on_species – Reference atom on slab to adsorb on
  • adatom_on_lig – bonding atom on ligand
  • ligand – structure object for ligand
  • displacement – initial adsorption distance desired above the adsorb_on_species
  • surface_coverage – Number of ligands desired per surface area of slab, in ligands per square angstroms
  • scell_max – Maximum number of supercells to create (used for finding supercell for the given coverage requirement
  • coverage_tol – Tolerance for coverage calculation in Ligands per square Angstroms
  • solvent – Name of solvent to be added for the run
  • start_from_slab – Whether slab is given as input. Useful when custom reconstructed slabs are to be used
  • validate_proximity – Check whether any atoms are too close (using pymatgen default of 0.01 Angstroms)
  • to_unit_cell – Pymatgen Slab routine to find unit cell
  • coords_are_cartesian – Whether the input coordinates are in cartesian
  • from_ase – Whether to create Slab using python-ase for producing slabs that have orthogonal lattice vectors
  • NOTE
  • starting from the bulk structure, create slab (if) –
  • note – if the starting structure is a slab, the vaccum extension
  • not possible (is) –
as_dict()[source]
calc_energy(model='Coulomb')[source]

calculates energy of the interface according to a defined energy model, useful for pre-screening candidate structures, returns 5 lowest energy structures :param model: energy model to compute according to, default is

Coulomb

Returns: energy of structure according to model

cover_surface(site_indices)[source]

puts the ligand molecule on the given list of site indices

create_interface()[source]

creates the interface i.e creates a slab of given thicknes and vacuum space. It ensures that the cell is big enough and have enough ligands to satify the surface coverage criterion also sets the slab on which the ligand is adsorbed

enforce_coverage()[source]

adjusts the supercell size and the number of adsorbed ligands so as to meet the surface coverage criterion within the given tolerance limit(specified as fraction of the required surface coverage)

returns the number of ligands and the supercell size that satisfies the criterion

get_index(species_string)[source]

get the first site index of the atomic species

get_reduced_scell()[source]

enforces the surface coverage criterion and generates the list all reduced lattice vectors that correspond to the computed supercell size and returns the one with similar lattice vector norms

set_slab()[source]

set the slab on to which the ligand is adsorbed

set_top_atoms()[source]

set the list of top and bottom atoms indices

class mpinterfaces.interface.Ligand(mols, cm_dist=[], angle={}, link={}, remove=[], charge=0, spin_multiplicity=None, validate_proximity=False)[source]

Bases: pymatgen.core.structure.Molecule

Construct ligand from molecules

as_dict()[source]
copy()[source]
create_ligand()[source]

create the ligand by assembling the provided individual molecules and removeing the specified atoms from the molecules

get_perp_vec(vec1, vec2)[source]

returns the vector that is perpendicular to the vec1 and vec2 if the vectors are parllel, then perp_vec = (0, -z, y)

link the molecules together connect the specified atoms of mol to the atoms of other molecules in the list connection means putting the atomm of the mol at a position that is the average of the position of the atoms of the molecules given in the list

position_mols()[source]

position the center of masses of the molecules wrt each other first movement is in the x direction

rotate_mols()[source]

rotate the molecules wrt each other using the provided info

set_distance_matrix(mol)[source]

sets the distance matrix for the molecule

set_mol_vecs()[source]

get the start and end indices to define the vector that defines the molecule sets the vectors that point from the start index atom to the farthest atom for each molecule

mpinterfaces.lammps module

class mpinterfaces.lammps.CalibrateLammps(parameters, structure=None, parent_job_dir='.', job_dir='./cal_lammps', qadapter=None, job_cmd='qsub', wait=True, is_matrix=False, turn_knobs=OrderedDict([('STRUCTURES', []), ('PAIR_COEFF', [])]), checkpoint_file=None, cal_logger=None)[source]

Bases: mpinterfaces.calibrate.Calibrate

Defines LAMMPS workflow consisting of LAMMPS jobs

add_job(name='noname', job_dir='.')[source]

add lammps job given MPINTLammps object

as_dict()[source]
classmethod from_dict(d)[source]
key_to_name(key)[source]
set_paircoeff(structure, pcoeff)[source]
setup_genericparam_jobs(key, vals)[source]
setup_paircoeff_jobs(paircoeff_list)[source]
setup_params_jobs(params_list)[source]
setup_structure_jobs(structures, paircoeff)[source]
val_to_name(val)[source]
class mpinterfaces.lammps.MPINTLammps(structure, parameters={}, label='mpintlmp', specorder=None, always_triclinic=False, no_data_file=False)[source]

Bases: ase.calculators.lammpsrun.LAMMPS, monty.json.MSONable

setup LAMMPS for given structure and parameters extends ase.calculators.lammpsrun.LAMMPS

as_dict()[source]
classmethod from_dict(d)[source]
write_lammps_data(f)[source]

write lammps structure data from ase with custom modifications

write_lammps_in(lammps_in=None, lammps_trj=None, lammps_data=None)[source]

write lammps input file from ase with custom modifications

class mpinterfaces.lammps.MPINTLammpsInput(mplmp, qadapter=None, vis_logger=None)[source]

Bases: monty.json.MSONable

create inputs(data and input file) for LAMMPS

as_dict()[source]
classmethod from_dict(d)[source]
write_input(job_dir, make_dir_if_not_present=True)[source]

write LAMMPS input set to job_dir

class mpinterfaces.lammps.MPINTLammpsJob(job_cmd, name='mpintlmpjob', output_file='job.out', parent_job_dir='.', job_dir='untitled', final=True, gzipped=False, backup=False, vis=None, settings_override=None, wait=True, vjob_logger=None)[source]

Bases: mpinterfaces.instrument.MPINTJob

Define LAMMPS job

as_dict()[source]
classmethod from_dict(d)[source]
get_final_energy(lammps_log='log.lammps')[source]

return the final total energy

mpinterfaces.measurement module

class mpinterfaces.measurement.Measurement(cal_objs, parent_job_dir='.', job_dir='./Measurement')[source]

Bases: object

Takes in calibrate objects and use that to perform various measurement calculations such as solvation, ligand binding energy etc. The default behaviour is to setup and run static calculations for all the given calibrate jobs.

Serves as Base Class. Override this class for custom measurements.

Parameters:
  • cal_objs – List of Calibration Object Names
  • parent_job_dir – Directory in which measuremnt is setup
  • job_dir – Path name to directory for running the Measurement modules
get_energy(cal)[source]

measures the energy of a single cal object a single cal object can have multiple calculations returns energies lists

make_measurements()[source]

gets the energies and processes it override this for custom measurements

run(job_cmd=None)[source]

run jobs

setup()[source]

setup static jobs for all the calibrate objects copies CONTCAR to POSCAR sets NSW = 0

class mpinterfaces.measurement.MeasurementInterface(cal_objs, parent_job_dir='.', job_dir='./MeasurementInterface')[source]

Bases: mpinterfaces.measurement.Measurement

Interface

Takes list of Calibration Objects of Interface, Slab and Ligand and separates them

Parameters:cal_objs – List of Calibration Objects
make_measurements()[source]

get the slab, ligand and interface energies compute binding energies

setup()[source]

setup static jobs for the calibrate objects copies CONTCAR to POSCAR sets NSW = 0 write system.json file for database crawler

class mpinterfaces.measurement.MeasurementSolvation(cal_obj, parent_job_dir='.', job_dir='./MeasurementSolvation', sol_params=None)[source]

Bases: mpinterfaces.measurement.Measurement

Solvation with poisson-boltzmann(test verison)

make_measurements()[source]

get solvation energies

setup()[source]

setup solvation jobs for the calibrate objects copies WAVECAR and sets the solvation params in the incar file also dumps system.json file in each directory for the database crawler mind: works only for cal objects that does only single calculations

mpinterfaces.nanoparticle module

class mpinterfaces.nanoparticle.Nanoparticle(structure, rmax=15, hkl_family=((1, 0, 0), (1, 1, 1)), surface_energies=(28, 25))[source]

Bases: pymatgen.core.structure.Molecule

Construct nanoparticle using wulff construction

create()[source]

creates the nanoparticle by chopping of the corners normal to the specified surfaces. the distance to the surface from the center of the particel = normalized surface energy * max radius

get_centered_molecule()[source]
get_normals()[source]

get the normal to the plane (h,k,l)

set_miller_family()[source]

get all miller indices for the given maximum index get the list of indices that correspond to the given family of indices

mpinterfaces.rest module

exception mpinterfaces.rest.MWRestError[source]

Bases: Exception

Exception class for MWRestAdaptor. Raised when the query has problems, e.g., bad query format.

class mpinterfaces.rest.MWRester(api_key=None, endpoint='https://www.materialsweb.org/rest')[source]

Bases: object

A class to conveniently interface with the MaterialsWeb REST interface. The recommended way to use MWRester is with the “with” context manager to ensure that sessions are properly closed after usage:

with MWRester("API_KEY") as m:
    do_something

MWRester uses the “requests” package, which provides for HTTP connection pooling. All connections are made via https for security.

Parameters:
  • api_key (str) – A String API key for accessing the MaterialsWeb REST interface. Please obtain your API key at https://www.materialsweb.org.
  • endpoint (str) – Url of endpoint to access the MaterialsWeb REST interface. Defaults to the standard MaterialsWeb REST address, but can be changed to other urls implementing a similar interface.
get_data(chemsys_formula_id, data_type='vasp', prop='')[source]

Flexible method to get any data using the MaterialsWeb REST interface. Generally used by other methods for more specific queries. Format of REST return is always a list of dict (regardless of the number of pieces of data returned. The general format is as follows: [{“material_id”: material_id, “property_name” : value}, ...] :param chemsys_formula_id: A chemical system (e.g., Li-Fe-O),

or formula (e.g., Fe2O3) or materials_id (e.g., mp-1234).
Parameters:
  • data_type (str) – Type of data to return. Currently can either be “vasp” or “exp”.
  • prop (str) – Property to be obtained. Should be one of the MWRester.supported_task_properties. Leave as empty string for a general list of useful properties.
get_structure_by_material_id(material_id, final=True)[source]

Get a Structure corresponding to a material_id. :param material_id: MaterialsWeb material_id (a string,

e.g., mp-1234).
Parameters:final (bool) – Whether to get the final structure, or the initial (pre-relaxation) structure. Defaults to True.
Returns:Structure object.
supported_properties = ('energy', 'energy_per_atom', 'volume', 'formation_energy_per_atom', 'nsites', 'unit_cell_formula', 'pretty_formula', 'is_hubbard', 'elements', 'nelements', 'e_above_hull', 'hubbards', 'is_compatible', 'spacegroup', 'task_ids', 'band_gap', 'density', 'icsd_id', 'icsd_ids', 'cif', 'total_magnetization', 'material_id', 'oxide_type', 'tags', 'elasticity')
supported_task_properties = ('energy', 'energy_per_atom', 'volume', 'formation_energy_per_atom', 'nsites', 'unit_cell_formula', 'pretty_formula', 'is_hubbard', 'elements', 'nelements', 'e_above_hull', 'hubbards', 'is_compatible', 'spacegroup', 'band_gap', 'density', 'icsd_id', 'cif')

mpinterfaces.transformations module

mpinterfaces.transformations.generate_all_configs(mat2d, substrate, nlayers_2d=2, nlayers_substrate=2, seperation=5)[source]

For the given lattice matched 2D material and substrate structures, this functions computes all unique sites in the interface layers and subsequently generates all possible unique 2d/substrate interfaces and writes the corresponding poscar files

Parameters:
  • mat2d – Lattice and symmetry-matched 2D material structure
  • substrate – Lattice and symmetry-matched 2D substrate structure
  • nlayers_substrate – number of substrate layers
  • nlayers_2d – number of 2d material layers
  • seperation – seperation between the substrate and the 2d material
Returns:

None

TODO: give additional random placement of 2D material on substrate

mpinterfaces.transformations.get_aligned_lattices(slab_sub, slab_2d, max_area=200, max_mismatch=0.05, max_angle_diff=1, r1r2_tol=0.2)[source]

given the 2 slab structures and the alignment paramters, return slab structures with lattices that are aligned with respect to each other

mpinterfaces.transformations.get_angle(a, b)[source]

angle between lattice vectors a and b in degrees

mpinterfaces.transformations.get_area(uv)[source]

area of the parallelogram

mpinterfaces.transformations.get_matching_lattices(iface1, iface2, max_area=100, max_mismatch=0.01, max_angle_diff=1, r1r2_tol=0.02)[source]

computes a list of matching reduced lattice vectors that satify the max_area, max_mismatch and max_anglele_diff criteria

mpinterfaces.transformations.get_mismatch(a, b)[source]

percentage mistmatch between the lattice vectors a and b

mpinterfaces.transformations.get_r_list(area1, area2, max_area, tol=0.02)[source]

returns a list of r1 and r2 values that satisfies: r1/r2 = area2/area1 with the constraints: r1 <= Area_max/area1 and r2 <= Area_max/area2 r1 and r2 corresponds to the supercell sizes of the 2 interfaces that align them

mpinterfaces.transformations.get_reduced_uv(uv, tm)[source]

returns reduced lattice vectors

mpinterfaces.transformations.get_trans_matrices(n)[source]

yields a list of 2x2 transformation matrices for the given supercell n: size

mpinterfaces.transformations.get_uniq_layercoords(struct, nlayers, top=True)[source]

returns the coordinates of unique sites in the top or bottom nlayers of the given structure.

Parameters:
  • struct – input structure
  • nlayers – number of layers
  • top – top or bottom layers, default is top layer
Returns:

numpy array of unique coordinates

mpinterfaces.transformations.get_uv(ab, t_mat)[source]

return u and v, the supercell lattice vectors obtained through the transformation matrix

mpinterfaces.transformations.reduced_supercell_vectors(ab, n)[source]

returns all possible reduced in-plane lattice vectors and transition matrices for the given starting unit cell lattice vectors(ab) and the supercell size n

mpinterfaces.utils module

mpinterfaces.utils.add_vacuum(structure, vacuum)[source]

Adds padding to a slab or 2D material.

Parameters:
  • structure (Structure) – Structure to add vacuum to
  • vacuum (float) – Vacuum thickness to add in Angstroms
Returns:

Structure object with vacuum added.

mpinterfaces.utils.align_axis(structure, axis='c', direction=(0, 0, 1))[source]

Rotates a structure so that the specified axis is along the [001] direction. This is useful for adding vacuum, and in general for using vasp compiled with no z-axis relaxation.

Parameters:
  • structure (Structure) – Pymatgen Structure object to rotate.
  • axis – Axis to be rotated. Can be ‘a’, ‘b’, ‘c’, or a 1x3 vector.
  • direction (vector) – Final axis to be rotated to.
Returns:

structure. Rotated to align axis along direction.

mpinterfaces.utils.center_slab(structure)[source]

Centers the atoms in a slab structure around 0.5 fractional height.

Parameters:structure (Structure) – Structure to center
Returns:Centered Structure object.
mpinterfaces.utils.ensure_vacuum(structure, vacuum)[source]

Adds padding to a slab or 2D material until the desired amount of vacuum is reached.

Parameters:
  • structure (Structure) – Structure to add vacuum to
  • vacuum (float) – Final desired vacuum thickness in Angstroms
Returns:

Structure object with vacuum added.

mpinterfaces.utils.get_ase_slab(pmg_struct, hkl=(1, 1, 1), min_thick=10, min_vac=10)[source]

takes in the intial structure as pymatgen Structure object uses ase to generate the slab returns pymatgen Slab object

Parameters:
  • pmg_struct – pymatgen structure object
  • hkl – hkl index of surface of slab to be created
  • min_thick – minimum thickness of slab in Angstroms
  • min_vac – minimum vacuum spacing
mpinterfaces.utils.get_convergence_data(jfile, params=('ENCUT', 'KPOINTS'))[source]

returns data dict in the following format {‘Al’:

{‘ENCUT’: [ [500,1.232], [600,0.8798] ],
‘KPOINTS’:[ [], [] ]

},

‘W’: ...

}

Note: processes only INCAR parmaters and KPOINTS

mpinterfaces.utils.get_convergence_data_custom(jfile, params=('ENCUT', 'KPOINTS'))[source]

returns data dict in the following format {‘Al’:

{‘ENCUT’: [ [500,1.232], [600,0.8798] ],
‘KPOINTS’:[ [], [] ]

},

‘W’: ...

}

Note: processes only INCAR parmaters and KPOINTS Sufficient tagging of the data assumed from species,Poscar comment line and potcar functional

mpinterfaces.utils.get_job_state(job)[source]
Parameters:job – job
Returns:the job state and the job output file name
mpinterfaces.utils.get_logger(log_file_name)[source]

writes out logging file. Very useful project logging, recommended for use to monitor the start and completion of steps in the workflow Arg:

log_file_name: name of the log file, log_file_name.log
mpinterfaces.utils.get_magmom_afm(poscar, database=None)[source]

returns the magmom string which is an N length list

mpinterfaces.utils.get_magmom_mae(poscar, mag_init)[source]

mae

mpinterfaces.utils.get_magmom_string(structure)[source]

Based on a POSCAR, returns the string required for the MAGMOM setting in the INCAR. Initializes transition metals with 6.0 bohr magneton and all others with 0.5. TEST: integration of mat2d function with mpinterfaces calibrate.py Consider moving to mpinterfaces.utils

Parameters:structure (Structure) – Pymatgen Structure object
Returns:string with INCAR setting for MAGMOM according to mat2d database calculations
mpinterfaces.utils.get_markovian_path(points)[source]

Calculates the shortest path connecting an array of 2D points. Useful for sorting linemode k-points.

Parameters:points (list) – list/array of points of the format [[x_1, y_1, z_1], [x_2, y_2, z_2], ...]
Returns:A sorted list of the points in order on the markovian path.
Return type:list
mpinterfaces.utils.get_opt_params(data, species, param='ENCUT', ev_per_atom=0.001)[source]

return optimum parameter default: 1 meV/atom

mpinterfaces.utils.get_opt_params_custom(data, tag, param='ENCUT', ev_per_atom=1.0)[source]
Parameters:
  • data – dictionary of convergence data
  • tag – key to dictionary of convergence dara
  • param – parameter to be optimized
  • ev_per_atom – minimizing criterion in eV per unit
Returns
[list] optimum parameter set consisting of tag, potcar object, poscar object, list of convergence data energies sorted according to param

default criterion: 1 meV/atom

mpinterfaces.utils.get_rotation_matrix(axis, theta)[source]

Find the rotation matrix associated with counterclockwise rotation about the given axis by theta radians. Credit: http://stackoverflow.com/users/190597/unutbu

Parameters:
  • axis (list) – rotation axis of the form [x, y, z]
  • theta (float) – rotational angle in radians
Returns:

array. Rotation matrix.

mpinterfaces.utils.get_run_cmmnd(nnodes=1, ntasks=16, walltime='10:00:00', job_bin=None, job_name=None, mem=None)[source]

returns the fireworks CommonAdapter based on the queue system specified by mpint_config.yaml and the submit file template also specified in mpint_config.yaml NOTE: for the job_bin, please specify the mpi command as well:

Eg: mpiexec /path/to/binary
mpinterfaces.utils.get_spacing(structure)[source]

Returns the interlayer spacing for a 2D material or slab.

Parameters:
  • structure (Structure) – Structure to check spacing for.
  • cut (float) – a fractional z-coordinate that must be within the vacuum region.
Returns:

float. Spacing in Angstroms.

mpinterfaces.utils.get_structure_type(structure, write_poscar_from_cluster=False)[source]

This is a topology-scaling algorithm used to describe the periodicity of bonded clusters in a bulk structure.

Parameters:
  • structure (structure) – Pymatgen structure object to classify.
  • write_poscar_from_cluster (bool) – Set to True to write a POSCAR from the sites in the cluster.
Returns:

string. ‘molecular’ (0D), ‘chain’, ‘layered’, ‘heterogeneous’

(intercalated 3D), or ‘conventional’ (3D)

mpinterfaces.utils.is_converged(directory)[source]

Check if a relaxation has converged.

Parameters:directory (str) – path to directory to check.
Returns:boolean. Whether or not the job is converged.
mpinterfaces.utils.jobs_from_file(filename='calibrate.json')[source]

read in json file of format caibrate.json(the default logfile created when jobs are run through calibrate) and return the list of job objects.

Parameters:filename – checkpoint file name
Returns:list of all jobs
mpinterfaces.utils.launch_daemon(steps, interval, handlers=None, ld_logger=None)[source]

run all the ‘steps’ in daemon mode checks job status every ‘interval’ seconds also runs all the error handlers

mpinterfaces.utils.partition_jobs(turn_knobs, max_jobs)[source]

divide turn_knobs into smaller turn_knobs so that each one of them has smaller max_jobs jobs

mpinterfaces.utils.print_exception()[source]

Error exception catching function for debugging can be a very useful tool for a developer move to utils and activate when debug mode is on

mpinterfaces.utils.remove_z_kpoints()[source]

Strips all linemode k-points from the KPOINTS file that include a z-component, since these are not relevant for 2D materials and slabs.

mpinterfaces.utils.set_sd_flags(poscar_input=None, n_layers=2, top=True, bottom=True, poscar_output='POSCAR2')[source]

set the relaxation flags for top and bottom layers of interface. The upper and lower bounds of the z coordinate are determined based on the slab. :param poscar_input: input poscar file name :param n_layers: number of layers to be relaxed :param top: whether n_layers from top are be relaxed :param bottom: whether n_layers from bottom are be relaxed :param poscar_output: output poscar file name

Returns:None writes the modified poscar file
mpinterfaces.utils.slab_from_file(hkl, filename)[source]

reads in structure from the file and returns slab object. useful for reading in 2d/substrate structures from file. :param hkl: miller index of the slab in the input file. :param filename: structure file in any format

supported by pymatgen
Returns:Slab object
mpinterfaces.utils.update_checkpoint(job_ids=None, jfile=None, **kwargs)[source]

rerun the jobs with job ids in the job_ids list. The jobs are read from the json checkpoint file, jfile. If no job_ids are given then the checkpoint file will be updated with corresponding final energy

Parameters:
  • job_ids – list of job ids to update or q resolve
  • jfile – check point file
mpinterfaces.utils.update_submission_template(default_template, qtemplate)[source]

helper function for writing a CommonAdapter template fireworks submission file based on a provided default_template which contains hpc resource allocation information and the qtemplate which is a yaml of commonly modified user arguments

mpinterfaces.utils.write_circle_mesh_kpoints(center=[0, 0, 0], radius=0.1, resolution=20)[source]

Create a circular mesh of k-points centered around a specific k-point and write it to the KPOINTS file. Non-circular meshes are not supported, but would be easy to code. All k-point weights are set to 1.

Parameters:
  • center (list) – x, y, and z coordinates of mesh center. Defaults to Gamma.
  • radius (float) – Size of the mesh in inverse Angstroms.
  • resolution (int) – Number of mesh divisions along the radius in the 3 primary directions.
mpinterfaces.utils.write_pbs_runjob(name, nnodes, nprocessors, pmem, walltime, binary)[source]

writes a runjob based on a name, nnodes, nprocessors, walltime, and binary. Designed for runjobs on the Hennig group_list on HiperGator 1 (PBS).

Parameters:
  • name (str) – job name.
  • nnodes (int) – number of requested nodes.
  • nprocessors (int) – number of requested processors.
  • pmem (str) – requested memory including units, e.g. ‘1600mb’.
  • walltime (str) – requested wall time, hh:mm:ss e.g. ‘2:00:00’.
  • binary (str) – absolute path to binary to run.
mpinterfaces.utils.write_potcar(pot_path=None, types='None')[source]

Writes a POTCAR file based on a list of types.

Parameters:
  • pot_path (str) – can be changed to override default location of POTCAR files.
  • types (list) – list of same length as number of elements containing specifications for the kind of potential desired for each element, e.g. [‘Na_pv’, ‘O_s’]. If left as ‘None’, uses the defaults in the ‘potcar_symbols.yaml’ file in the package root.
mpinterfaces.utils.write_slurm_runjob(name, ntasks, pmem, walltime, binary)[source]

writes a runjob based on a name, nnodes, nprocessors, walltime, and binary. Designed for runjobs on the Hennig group_list on HiperGator 2 (SLURM).

Parameters:
  • name (str) – job name.
  • ntasks (int) – total number of requested processors.
  • pmem (str) – requested memory including units, e.g. ‘1600mb’.
  • walltime (str) – requested wall time, hh:mm:ss e.g. ‘2:00:00’.
  • binary (str) – absolute path to binary to run.

Module contents

mpinterfaces.get_struct_from_mp(formula, MAPI_KEY='', all_structs=False)[source]

fetches the structure corresponding to the given formula from the materialsproject database.

Note: Get the api key from materialsproject website. The one used here is nolonger valid.

Note: for the given formula there are many structures available, this function returns the one with the lowest energy above the hull unless all_structs is set to True