Manual

The complete WESTpy reference is reported below.

Atom

class westpy.atom.Atom(cell=None, ase_atom=None, symbol=None, cry_coord=None, abs_coord=None, **kwargs)

Bases: object

Class for representing a single atom.

An atom can be initialized by and exported to an ASE Atom object.

All internal quantities below are in atomic unit.

cell

the cell where the atom lives.

Type:

Cell

symbol

chemical symbol of the atom.

Type:

str

abs_coord

absolute coordinate.

Type:

float array

cry_coord

crystal coordinate.

Type:

float array

Extra attributes can be attached to an atom.

property cry_coord

BSE

class westpy.bse.bseresult.BSEResult(filename: str)

Bases: object

Parses BSE/TDDFT results.

Parameters:

filename (string) – Wbse or Westpp output file (JSON)

Example:

>>> from westpy.bse import *
>>> wbse = BSEResult("wbse.json")
plotSpectrum(ipol: str = None, energyRange: List[float] = [0.0, 10.0, 0.01], sigma: float = 0.1, n_extra: int = 0, fname: str = None)

Plots BSE/TDDFT absorption spectrum.

Parameters:
  • ipol (string) – which component to compute (“XX”, “XY”, “XZ”, “YX”, “YY”, “YZ”, “ZX”, “ZY”, “ZZ”, or “XYZ”)

  • energyRange (3-dim float) – energy range = min, max, step (eV)

  • sigma (float) – Broadening width (eV)

  • n_extra (int) – Number of extrapolation steps (Lanczos only)

  • fname (string) – Output file name

Example:

>>> from westpy.bse import *
>>> wbse = BSEResult("wbse.json")
>>> wbse.plotSpectrum(ipol="XYZ",energyRange=[0.0,10.0,0.01],sigma=0.1,n_extra=100000)
class westpy.bse.qbox2bse.Qbox2BSE(filename: str)

Bases: object

Parses Qbox output and generates files needed by WEST BSE.

Parameters:

filename (string) – Qbox output file (XML)

Example:

>>> from westpy.bse import *
>>> qb = Qbox2BSE("qb.out")
write_localization(filename: str = 'bis_info')

Reads localization from XML file then writes to text file.

Parameters:

filename (string) – name of Qbox bisection information file

Example:

>>> from westpy.bse import *
>>> qb = Qbox2BSE("qb.out")
>>> qb.write_localization()
write_wavefunction(filename: str = 'qb_wfc')

Reads wavefunctions from XML file then writes to HDF5 file.

Parameters:

filename (string) – name of Qbox wavefunction file

Example:

>>> from westpy.bse import *
>>> qb = Qbox2BSE("qb.out")
>>> qb.write_wavefunction()

Cell

class westpy.cell.Cell(ase_cell=None, R=None)

Bases: object

A cell that is defined by lattice constants and contains a list of atoms.

A cell can be initialized by and exported to an ASE Atoms object. A cell can also be exported to several useful formats such as Qbox and Quantum Espresso formats.

All internal quantities below are in atomic unit.

R

primitive lattice vectors (each row is one vector).

Type:

float 3x3 array

G

reciprocal lattice vectors (each row is one vector).

Type:

float 3x3 array

omega

cell volume.

Type:

float

atoms

list of atoms.

Type:

list of Atoms

property G
property R
property atoms
property isperiodic
property natoms
property nspecies
property omega
property species
update_lattice(R)

Data Container

class westpy.dataContainer.DataContainer

Bases: object

Class for representing an in-memory data container.

Example:

>>> from westpy import *
>>> dc = DataContainer()
checkKeys(printSummary=True)

Checks that all keys are described.

Parameters:

printSummary (boolean) – if True prints summary

Returns:

True if all keys are described, False otherwise.

Return type:

boolean

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertPoint({"a":1, "b":2},{"energy":-4.5})
>>> dc.upsertKey("a","the first letter")
>>> dc.upsertKey("b","another letter")
>>> dc.upsertKey("energy","a quantity")
>>> dc.checkKeys()
removeKey(key)

Removes the description of a key

Parameters:

key (string) – key

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertKey("a","the first letter")
>>> dc.removeKey("a")
removePoint(identifier)

Removes point with given identifier from the data container.

Parameters:

identifier (hashable object) – identifier

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertPoint({"a":1, "b":2},{"energy":-4.5})
>>> dc.removePoint({"a":1, "b":2})
showPoints()

Shows all points of the data container.

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertPoint({"a":1, "b":2},{"energy":-4.5})
>>> dc.showPoints()
upsertKey(key, description)

Updates or inserts a new key and its description.

Parameters:
  • key (string) – key

  • description (hashable object) – description

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertKey("a","the first letter")
upsertPoint(identifier, document, incremental_update=True)

Update or inserts an entry to the data container.

If identifier exists, update the document associated to the identifier, otherwise insert the document with the identifier.

Parameters:
  • identifier (hashable object) – identifier

  • document (hashable object) – document

  • incremental_update (boolean) – if the document exists, only update it, do not remove its other keys.

Example:

>>> from westpy import *
>>> dc = DataContainer()
>>> dc.upsertPoint({"a":1, "b":2},{"energy":-4.5})

Electronic Structure

class westpy.electronicStructure.ElectronicStructure

Bases: object

Class for representing an electronic structure calculation.

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
addDataPoint(ksb, key, what)

Adds datapoint to data.

Parameters:
  • ksb (3-dim tuple of int) – triplet of integers: k-point, spin, band (integer labels)

  • key (string) – key

  • what (hashable object) – content attached to key

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.addKey("eks","Kohn-Sham energy in eV")
>>> es.addDataPoint([1,1,1],"eks",-4.6789)
addKey(key, description)

Describes metadata key.

Parameters:
  • key (string) – key

  • description (hashable object) – description

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.addKey("eks","Kohn-Sham")
plotDOS(kk=[1], ss=[1], energyKeys=['eks'], sigma=0.1, weight=1.0, energyRange=[-20.0, 0.0, 0.01], fname='dos.png')

Plots desnity of states (DOS).

Parameters:
  • kk (list of int) – list of k-points

  • ss (list of int) – list of spin channels (must be [1], [2], or [1,2])

  • energyKeys (list of string) – energy keys (mush match the available keys)

  • sigma (float or string) – standard deviation of gaussian (eV), optional (mush match the available keys)

  • weight (float or string) – weight, optional (mush match the available keys)

  • energyRange (3-dim tuple of float) – energy range = min, max, step (eV), optional

  • fname (string) – output file name, optional

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.addKey("eks","Kohn-Sham energy in eV")
>>> es.addDataPoint([1,1,1],"eks",-4.6789)
>>> es.addDataPoint([1,1,2],"eks",-4.3456)
>>> es.addDataPoint([1,2,1],"eks",-4.4567)
>>> es.addDataPoint([1,2,2],"eks",-4.0123)
>>> es.plotDOS(kk=[1],ss=[1,2],energyKeys=["eks"],energyRange=[-5.,-3,0.01])
plotLDOS(kk=[1], ss=[1], energyKeys=['eks'], sigma=0.1, weight=1.0, energyRange=[-20.0, 0.0, 0.01], wfcKey='wfcFile', fname='ldos.png')

Plots desnity of states (DOS).

Parameters:
  • kk (list of int) – list of k-points

  • ss (list of int) – list of spin channels (must be [1], [2], or [1,2])

  • energyKeys (list of string) – energy keys (mush match the available keys)

  • sigma (float or string) – standard deviation of gaussian (eV), optional (mush match the available keys)

  • weight (float or string) – weight, optional (mush match the available keys)

  • energyRange (3-dim tuple of float) – energy range = min, max, step (eV), optional

  • wfcKey (string) – wavefunction file (mush match the available keys)

  • fname (string) – output file name, optional

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.addKey("eks","Kohn-Sham energy in eV")
>>> es.addDataPoint([1,1,1],"eks",-4.6789)
>>> es.addDataPoint([1,1,2],"eks",-4.3456)
>>> es.addDataPoint([1,2,1],"eks",-4.4567)
>>> es.addDataPoint([1,2,2],"eks",-4.0123)
>>> es.plotLDOS(kk=[1],ss=[1,2],energyKeys=["eks"],energyRange=[-5.,-3,0.01])
removeKey(key)

Removes key from metadata.

Parameters:

key (string) – key

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.addKey("eks","Kohn-Sham")
>>> es.removeKey("eks")
showKeys()

Shows keys in metadata.

Example:

>>> from westpy import *
>>> es = ElectronicStructure()
>>> es.showKeys()

Geometry

class westpy.geometry.Geometry(cell=None)

Bases: object

Class for representing a set of atoms in a periodic cell.

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.setCell( (1,0,0), (0,1,0), (0,0,1) )
>>> geom.addAtom( symbol="Si", abs_coord=(0,0,0) )
>>> geom.addSpecies( "Si", "http://www.quantum-simulation.org/potentials/sg15_oncv/upf/Si_ONCV_PBE-1.1.upf" )

Note

Vectors are set in a.u. by default. If you set units=Angstrom a coversion to a.u. will be made.

Bohr = 1.0
addAtom(symbol, abs_coord, units=1.0)

Adds a single atom.

Parameters:
  • symbol (string) – chemical symbol

  • abs_coord (3-dim tuple of float) – absolute coordinates

  • units ("Bohr" or "Angstrom") – Units, optional

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtom( symbol="Si", abs_coord=(0,0,0) )
addAtomsFromOnlineXYZ(url)

Adds atoms from XYZ file (only one image) located at url.

Parameters:

url (string) – url

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomsFromOnlineXYZ( "https://west-code.org/database/gw100/xyz/CH4.xyz" )
addAtomsFromXYZFile(fname)

Adds atoms from XYZ file (only one image).

Parameters:

fname (string) – file name

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomFromXYZFile( "CH4.xyz" )
addFracCoordAtom(symbol, frac_coord)

adds a single atom by fractional coords

Parameters:
  • symbol (string) – chemical symbol

  • frac_coord (3-dim tuple of float) – fractional coordinates

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addFracCoordAtom( "Si", (0,1/3.0,2/3.0) )
addSpecies(symbol, url)

Adds a species.

Parameters:
  • symbol (string) – chemical symbol

  • url (string) – url

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addSpecies( "Si", "http://www.quantum-simulation.org/potentials/sg15_oncv/upf/Si_ONCV_PBE-1.1.upf" )

Note

You can use this method to add either upf or xml pseudopotentials. However it is forbidden to mix them.

downloadPseudopotentials()

Download Pseudopotentials.

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomsFromOnlineXYZ( "https://west-code.org/database/gw100/xyz/CH4.xyz" )
>>> geom.addSpecies( "C", "http://www.quantum-simulation.org/potentials/sg15_oncv/xml/C_ONCV_PBE-1.0.xml")
>>> geom.addSpecies( "H", "http://www.quantum-simulation.org/potentials/sg15_oncv/xml/H_ONCV_PBE-1.0.xml")
>>> geom.downloadPseudopotentials()

Note

Pseudopotential files will be downloaded in the current directory.

getNumberOfAtoms()

Returns number of atoms.

Returns:

number of atoms

Return type:

int

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomsFromOnlineXYZ( "https://west-code.org/database/gw100/xyz/CH4.xyz" )
>>> nat = geom.getNumberOfAtoms()
>>> print( nat )
5
getNumberOfElectrons()

Returns number of electrons.

Returns:

number of electrons

Return type:

int

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomsFromOnlineXYZ( "https://west-code.org/database/gw100/xyz/CH4.xyz" )
>>> geom.addSpecies( "C", "http://www.quantum-simulation.org/potentials/sg15_oncv/xml/C_ONCV_PBE-1.0.xml")
>>> geom.addSpecies( "H", "http://www.quantum-simulation.org/potentials/sg15_oncv/xml/H_ONCV_PBE-1.0.xml")
>>> nelec = geom.getNumberOfElectrons()
>>> print( nelec )
8
getNumberOfSpecies()

Returns number of species.

Returns:

number of species

Return type:

int

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.addAtomsFromOnlineXYZ( "https://west-code.org/database/gw100/xyz/CH4.xyz" )
>>> ntyp = geom.getNumberOfSpecies()
>>> print( ntyp )
2
isValid()

Checks if geometry is valid

The method checks that:
  • the cell is set

  • at least one atom has been added

  • the pseudopotentials of all species are defined

  • the pseudopotentials do not contain a mix of upf and xml formats

setCell(a1=(0, 0, 0), a2=(0, 0, 0), a3=(0, 0, 0), units=1.0)

Sets cell, given the three vectors \(a_1\), \(a_2\), \(a_3\).

Parameters:
  • a1 (3-dim tuple of float) – \(a_1\)

  • a2 (3-dim tuple of float) – \(a_2\)

  • a3 (3-dim tuple of float) – \(a_3\)

  • units ("Bohr" or "Angstrom") – Units, optional

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.setCell( (1,0,0), (0,1,0), (0,0,1) )
view(style='stick', width=800, height=800, ix=1, iy=1, iz=1, debug=False)

Display simulation box geom in Angstrom. ix, iy, iz is the perodic display. style can be line, stick, sphere.

Ground State

class westpy.groundState.GroundState(geom, xc, ecut)

Bases: object

Class for representing a ground state calculation with DFT.

Parameters:
  • geom (Geometry) – geometry (cell, atoms, species)

  • xc (string) – exchange-correlation functional

  • ecut (float) – energy cutoff for the wavefunction (Rydberg)

Example:

>>> from westpy import *
>>> geom = Geometry()
>>> geom.setCell( (1,0,0), (0,1,0), (0,0,1) )
>>> geom.addAtom( "Si", (0,0,0) )
>>> geom.addSpecies( "Si", "http://www.quantum-simulation.org/potentials/sg15_oncv/upf/Si_ONCV_PBE-1.1.upf" )
>>> gs = GroundState(geom,"PBE",30.0)

Note

Vectors are set in a.u. by default. If you set units=Angstrom a coversion to a.u. will be made.

downloadPseudopotentials()

Download Pseudopotentials.

Example:

>>> gs.downloadPseudopotentials()

Note

Pseudopotential files will be downloaded in the current directory.

generateInputPW(fname='pw.in')

Generates input file for pwscf. Valid only for QuantumEspresso calculations.

Parameters:

fname (string) – fname, optional

Example:

>>> gs.generateInputPW("pw.in")
generateInputQbox(fname='qbox.in')

Generates input file for qbox. Valid only for Qbox calculations.

Parameters:

fname (string) – fname, optional

Example:

>>> gs.generateInputQbox("qbox.in")
setCollinearSpin(tot_magnetization=0.0)

Sets collinear spin.

Parameters:

tot_magnetization (float) – Total majority spin charge - minority spin charge, optional

Example:

>>> gs.setCollinearSpin()
setIsolated()

Sets isolated system. Valid only for QuantumEspresso calculations.

Example:

>>> gs.setIsolated()
setKmesh(kmesh)

Sets the uniform grid for k-points sampling.

Parameters:

kmesh (3-dim tuple of int) – kmesh

Example:

>>> gs.setKmesh((2,2,2))
setNempty(nempty)

Sets the number of empty bands.

Parameters:

nempty (int) – number of empty bands

Example:

>>> gs.setNempty(10)
setNonCollinearSpin(lspinorb=False)

Sets non collinear spin. Requires fully relativistic pseudopotentials. Valid only for QuantumEspresso calculations.

Optionally spin-orbit can be turned on.

Parameters:

lspinorb (boolean) – spin-orbit, optional

Example:

>>> gs.setNonCollinearSpin()
updateSpecies(symbol, url)

Update a species.

Parameters:
  • symbol (string) – chemical symbol

  • url (string) – url

Example:

>>> geom.addSpecies( "Si", "http://www.quantum-simulation.org/potentials/sg15_oncv/upf/Si_ONCV_PBE-1.1.upf" )

Note

You can use this method to add either upf or xml pseudopotentials. However it is forbidden to mix them.

Life Time

westpy.lifetime.radiative_lifetime(westpp_file, ispin, band1, band2, n=None, e_zpl=None)

Computes radiative lifetime.

Parameters:
  • westpp_file (string) – The JSON output file of Westpp

  • ispin (int) – spin index

  • band1 (int) – band index (transition from band1 to band2 is computed)

  • band2 (int) – band index (transition from band1 to band2 is computed)

  • n (float or function of e_zpl) – refractive index

  • e_zpl (float) – zero-phonon line (ZPL) energy (Rydberg)

Example:

>>> from westpy import *
>>> tau = radiative_lifetime("westpp.json",2,101,102,2.0,1.25)

QDET

class westpy.qdet.qdetresult.QDETResult(filename: str, point_group: PointGroup | None = None, wfct_filenames: list | None = None, symmetrize: Dict[str, bool] = {})

Bases: object

solve(nelec: Tuple = None, nroots: int = 10, verbose: bool = True) Dict

Build and diagonalize effective Hamiltonians for given active space.

Parameters:
  • nelec (2-dim tuple of int) – Number of electrons in each spin channel

  • nroots (int) – Number of roots for FCI calculations

  • verbose (boolean) – If True, write detailed info for FCI calculations

class westpy.qdet.symm.PointGroup(name: str, operations: Sequence[PointGroupOperation], ctable: Dict[str, Sequence])

Bases: object

compute_rep_on_orbitals(orbitals: Sequence[VData], orthogonalize: bool = False) Tuple[PointGroupRep, List[str]]
Compute representation matrix on the Hilbert space spanned by a set

of orbitals.

Parameters:
  • orbitals – a set of orbitals.

  • orthogonalize – if True, orthorgonalize representation matrix.

Returns:

(matrix representation, symmetries)

class westpy.qdet.symm.PointGroupInversion(origin: Sequence[float] | ndarray = (0.0, 0.0, 0.0), cell: ndarray = None)

Bases: PointGroupOperation

class westpy.qdet.symm.PointGroupOperation(T: ndarray, origin: Sequence[float] | ndarray = None, cell: ndarray = None)

Bases: object

An operation in the point group.

Parameters:
  • T – 4x4 affine transformation matrix for point group operation

  • origin – origin of operation, defined as follows: let [x,y,z] denote the coordinates of the origin in crystal (fractional) unit, and [nx,ny,nz] denote the shape of the volumetric data, i.e., the dimension of the real space grid, then origin = [x*nx,y*ny,z*nz]

  • cell – transformation matrix to the volumetric cell data, defined as follows let A = np.array([x1,y1,z1], [x2,y2,z2], [x3,y3,z3]).T denote the crystal vector matrix, where [xi,yi,zi] is the ith crystal basis vector in any unit, and B = 1 / np.array([nx,0,0], [0,ny,0], [0,0,nz]), where [nx,ny,nz] denote the shape of the volumetric data, i.e., the dimension of the real space grid, then cell = A @ B

property inv

Inverse operator

class westpy.qdet.symm.PointGroupReflection(normal: Sequence[float] | ndarray, origin: Sequence[float] | ndarray = (0.0, 0.0, 0.0), cell: ndarray = None)

Bases: PointGroupOperation

class westpy.qdet.symm.PointGroupRep(point_group: PointGroup, orbitals: Sequence[VData], orthogonalize: bool = False)

Bases: object

class westpy.qdet.symm.PointGroupRotateReflection(rotvec: Sequence[float] | ndarray, origin: Sequence[float] | ndarray = (0.0, 0.0, 0.0), multiple: int = 1, cell: ndarray = None)

Bases: PointGroupOperation

class westpy.qdet.symm.PointGroupRotation(rotvec: Sequence[float] | ndarray, origin: Sequence[float] | ndarray = (0.0, 0.0, 0.0), cell: ndarray = None)

Bases: PointGroupOperation

Relaxation

class westpy.relaxation.bfgs_iter(run_pw: str, run_wbse: str, run_nscf: str = None, run_wbse_init: str = None, pp: list = [], pw_input: str = 'pw.in', nscf_input: str = 'nscf.in', wbse_input: str = 'wbse.in', wbse_init_input: str = 'wbse_init.in', l_copy_save_dir: bool = True, l_restart: bool = False, energy_thr: float = 0.0001, grad_thr: float = 0.001, maxiter: int = 100, w1: float = 0.01, w2: float = 0.5, bfgs_ndim: int = 1, trust_radius_ini: float = 0.5, trust_radius_min: float = 0.0002, trust_radius_max: float = 0.8)

Bases: object

Class for carrying out BFGS geometry relaxation.

Parameters:
  • run_pw (string) – Full command to run pwscf, e.g., mpirun -n 2 /path/to/qe/bin/pw.x

  • run_wbse (string) – Full command to run wbse, e.g., mpirun -n 4 /path/to/qe/bin/wbse.x -nb 4

  • run_nscf (string) – Full command to run nscf, e.g., mpirun -n 2 /path/to/qe/bin/pw.x

  • run_wbse_init (string) – Full command to run wbse_init, e.g., mpirun -n 4 /path/to/qe/bin/wbse_init.x

  • pp (string) – List of pseudopotential files

  • pw_input (string) – pw.x input file name

  • nscf_input (string) – pw.x input file name for nscf calculation

  • wbse_input (string) – wbse.x input file name

  • wbse_init_input (string) – wbse_init.x input file name

  • l_copy_save_dir (boolean) – If False, does not copy .save dir (True if startingpot/wfc=’file’)

  • l_restart (boolean) – If True, restart an unfinished run

  • energy_thr (float) – Convergence threshold on total energy (Ry) for ionic minimization

  • grad_thr (float) – Convergence threshold on forces (Ry/Bohr) for ionic minimization

  • maxiter (int) – Maximum number of BFGS steps

  • w1 (float) – Parameters used in line search based on the Wolfe conditions

  • w2 (float) – Parameters used in line search based on the Wolfe conditions

  • bfgs_ndim (int) – Dimension of BFGS. Only bfgs_ndim == 1 implemented

  • trust_radius_ini (float) – Initial ionic displacement in the structural relaxation

  • trust_radius_min (float) – Minimum ionic displacement in the structural relaxation

  • trust_radius_max (float) – Maximum ionic displacement in the structural relaxation

Example:

>>> from westpy import *
>>> run_pw = "mpirun -n 1 pw.x"
>>> run_wbse = "mpirun -n 4 wbse.x -nb 4"
>>> run_wbse_init = "mpirun -n 4 wbse_init.x -ni 4"
>>> bfgs = bfgs_iter(run_pw=run_pw, run_wbse=run_wbse, grad_thr=1e-4, maxiter=30)
>>> bfgs.solve()
solve()

Run BFGS geometry relaxation.

Session

class westpy.session.Session(emailId)

Bases: object

Class for setting up a session, connected to a remove server via rest APIs.

Example:

>>> from westpy import *
>>> session = Session("your.email@domain.edu")
getToken()

Returns the token of the session.

Example:

>>> from westpy import *
>>> session = Session("your.email@domain.edu")
>>> token = session.getToken()
>>> print(token)
run(executable=None, inputFile=None, outputFile=None, downloadUrl=[], number_of_cores=2)

Runs the executable on the remote server.

Parameters:
  • executable (string) – name of executable

  • inputFile (string) – name of input file

  • outputFile (string) – name of output file

  • downloadUrl (string) – URLs to be downloaded

  • number_of_cores (int) – number of cores

Example:

>>> from westpy import *
>>> session = Session("your.email@domain.edu")
>>> session.run( "pw", "pw.in", "pw.out", ["http://www.quantum-simulation.org/potentials/sg15_oncv/upf/C_ONCV_PBE-1.0.upf"] , 2 )
>>> session.stop()
status()

Returns whether the session is active and time left.

Example:

>>> from westpy import *
>>> session = Session("your.email@domain.edu")
>>> session.status()
stop()

Stops the session and clears the remote workspace.

Example:

>>> from westpy import *
>>> session = Session("your.email@domain.edu")
>>> session.stop()

Units

Westpy uses Hartree atomic units.

class westpy.units.Units(*args, **kwargs)

Bases: dict

Dictionary for units that supports .attribute access.

westpy.units.set_units()

Sets Rydberg atomic units.

Available units are:
  • Bohr

  • Angstrom, Ang

  • Rydberg, Ry

  • Hartree, Ha

  • eV

  • Joule

Note

westpy operates in Rydberg atomic units.

Utilities

Set of utilities.

westpy.utils.bool2str(logical)

Converts a boolean type into a string .TRUE. or .FALSE. .

Parameters:

logical (boolean) – logical

Returns:

.TRUE. or .FALSE.

Return type:

string

Example:

>>> from westpy import *
>>> t = bool2str(True)
>>> f = bool2str(False)
>>> print(t,f)
.TRUE. .FALSE.
westpy.utils.convertYaml2Json(fyml, fjson)

Converts the file from YAML to JSON.

Parameters:
  • fyml (string) – Name of YAML file

  • fjson (string) – Name of JSON file

Example:

>>> from westpy import *
>>> convertYaml2Json("file.yml","file.json")

Note

The file fjson will be created, fyml will not be overwritten.

westpy.utils.download(url, fname=None)

Downloads a file from url.

Parameters:
  • url (string) – url

  • fname (string) – file name

Example:

>>> from westpy import *
>>> download("https://west-code.org/database/gw100/xyz/CH4.xyz")

Note

The file will be downloaded in the current directory.

westpy.utils.extractFileNamefromUrl(url)

Extracts a file name from url.

Parameters:

url (string) – url

Returns:

file name

Return type:

string

Example:

>>> from westpy import *
>>> extractFileNamefromUrl("https://west-code.org/database/gw100/xyz/CH4.xyz")
westpy.utils.gaussian(x, mu, sig)

return normal distribution at point x.

\(f(x;\mu,\sigma) = \frac{1}{\sigma\sqrt{2\pi}}e^{-\frac{(x-\mu)^2}{2\sigma^2}}\)

Parameters:
Returns:

\(f(x;\mu,\sigma)\)

Return type:

float

Example:

>>> from westpy import *
>>> gaussian(1.0,2.0,3.0)
westpy.utils.listLinesWithKeyfromOnlineText(url, key)

List lines from text file located at url, with key.

Parameters:
  • url (string) – url

  • key (string) – keyword

Returns:

list of lines

Return type:

list

Example:

>>> from westpy import *
>>> url = "http://www.quantum-simulation.org/potentials/sg15_oncv/upf/Si_ONCV_PBE-1.1.upf"
>>> key = "z_valence"
>>> l = listLinesWithKeyfromOnlineText(url,key)
>>> print(l)
['       z_valence="    4.00"']

Note

Can be used to grep values from a UPF file.

westpy.utils.listValuesWithKeyFromOnlineXML(url, key)

List values from XML file located at url, with key.

Parameters:
  • url (string) – url

  • key (string) – keyword

Returns:

list of values

Return type:

list

Example:

>>> from westpy import *
>>> url = "http://www.quantum-simulation.org/potentials/sg15_oncv/xml/Si_ONCV_PBE-1.1.xml"
>>> key = "valence_charge"
>>> l = listLinesWithKeyfromOnlineXML(url,key)
>>> print(l)
['4']

Note

Can be used to grep values from an XML file.

westpy.utils.readJsonFile(fname)

Reads data from file using the JSON format.

Parameters:

fname (string) – file name

Returns:

data

Return type:

dict/list

Example:

>>> from westpy import *
>>> data = readJsonFile("mass.json")

Note

The file will be read from the current directory.

westpy.utils.read_cube(fname)

Read cube file into numpy array

Parameters:

fname (string) – file name of cube file

Returns:

(data, metadata)

Return type:

(np.array, dict)

westpy.utils.read_imcube(rfname, ifname='')

Convenience function to read in two cube files at once, where one contains the real part and the other contains the imag part. If only one file name given, other file name is inferred.

Parameters:
  • rfname (string) – file name of cube file of real part

  • ifname (string) – optional, file name of cube file of imag part

Returns:

(data, metadata), where data is (real part + j*imag part)

Return type:

(np.array, dict)

westpy.utils.wfreq2df(fname='wfreq.json', dfKeys=['eks', 'eqpLin', 'eqpSec', 'sigmax', 'sigmac_eks', 'sigmac_eqpLin', 'sigmac_eqpSec', 'vxcl', 'vxcnl', 'hf'])

Loads the wfreq JSON output into a pandas dataframe.

Parameters:
  • fname (string) – file name of JSON output file

  • dfKeys (list of string) – energy keys to be added to dataframe

Returns:

(dataframe, data)

Return type:

(pd.DataFrame, dict)

westpy.utils.writeJsonFile(fname, data)

Writes data to file using the JSON format.

Parameters:
  • fname (string) – file name

  • data (dict/list) – data

Example:

>>> from westpy import *
>>> data = {}
>>> data["mass"] = 1.0
>>> writeJsonFile("mass.json",data)

Note

The file will be generated in the current directory.

westpy.utils.write_cube(data, meta, fname)

Write volumetric data to cube file along

Parameters:
  • data (list of float) – volumetric data consisting of real values

  • meta (dict) – dict containing metadata with following keys: - atoms: list of atoms in the form (mass, [position]) - org: origin - xvec,yvec,zvec: lattice vector basis

  • fname (string) – file name of cubefile (existing files overwritten)

westpy.utils.write_imcube(data, meta, rfname, ifname='')

Convenience function to write two cube files from complex valued volumetric data, one for the real part and one for the imaginary part. Data about atoms, origin and lattice vectors are kept same for both. If only one file name given, other file name is inferred.

Parameters:
  • data (list of float) – volumetric data consisting of complex values

  • meta (dict) – dict containing metadata with following keys: - atoms: list of atoms in the form (mass, [position]) - org: origin - xvec,yvec,zvec: lattice vector basis

  • rfname (string) – file name of cubefile of real part

  • ifname (string) – file name of cubefile of imag part