sisl.io.siesta.ncSileSiesta
- class sisl.io.siesta.ncSileSiesta(filename, mode='r', lvl=0, access=1, *args, **kwargs)
Bases:
SileCDFSiesta
Generic NetCDF output file containing a large variety of information
Plotting
Plotting functions for the
ncSileSiesta
class.plot.geometry
(*[, axes, atoms, ...])Calls
read_geometry
and creates aGeometryPlot
from its output.plot.grid
(name[, index, ...])Calls
read_grid
and creates aGridPlot
from its output.Methods
base_directory
([relative_to])Retrieve the base directory of the file, relative to the path relative_to
close
()dir_file
([filename, filename_base])File of the current Sile
iter
([group, dimension, variable, levels, root])Iterator on all groups, variables and dimensions.
read
(*args, **kwargs)Generic read method which should be overloaded in child-classes
Returns a set of atoms corresponding to the basis-sets in the nc file
read_density_matrix
(**kwargs)Returns a density matrix from the underlying NetCDF file
read_dynamical_matrix
(**kwargs)Returns a dynamical matrix from the underlying NetCDF file
read_energy_density_matrix
(**kwargs)Returns energy density matrix from the underlying NetCDF file
Returns the fermi-level
Returns a vector with final forces contained.
Reads the force-constant stored in the nc file
Returns Geometry object from a Siesta.nc file
read_grid
(name[, index])Reads a grid in the current Siesta.nc file
read_hamiltonian
(**kwargs)Returns a Hamiltonian from the underlying NetCDF file
Reads the force-constant stored in the nc file
Returns a Lattice object from a Siesta.nc file
Returns number of supercell connections
read_overlap
(**kwargs)Returns a overlap matrix from the underlying NetCDF file
write
(*args, **kwargs)Generic write method which should be overloaded in child-classes
write_basis
(atoms)Write the current atoms orbitals as the basis
write_density_matrix
(DM, **kwargs)Writes density matrix model to file
write_dynamical_matrix
(D, **kwargs)Writes dynamical matrix model to file
write_energy_density_matrix
(EDM, **kwargs)Writes energy density matrix model to file
write_geometry
(geometry)Creates the NetCDF file and writes the geometry information
write_hamiltonian
(H, **kwargs)Writes Hamiltonian model to file
write_overlap
(S, **kwargs)Write the overlap matrix to the NetCDF file
Attributes
File of the current Sile
File of the current Sile
Return a list of available grids in this file.
- base_directory(relative_to='.')
Retrieve the base directory of the file, relative to the path relative_to
- close()
- dir_file(filename=None, filename_base='')
File of the current Sile
- iter(group=True, dimension=True, variable=True, levels=-1, root=None)
Iterator on all groups, variables and dimensions.
This iterator iterates through all groups, variables and dimensions in the
Dataset
The generator sequence will _always_ be:
Group
Dimensions in group
Variables in group
As the dimensions are generated before the variables it is possible to copy groups, dimensions, and then variables such that one always ensures correct dependencies in the generation of a new
SileCDF
.- Parameters:
group (
bool
(True)) – whether the iterator yields Group instancesdimension (
bool
(True)) – whether the iterator yields Dimension instancesvariable (
bool
(True)) – whether the iterator yields Variable instanceslevels (
int
(-1)) – number of levels to traverse, with respect toroot
variable, i.e. number of sub-groups this iterator will return.root (
str
(None)) – the base root to start iterating from.
Examples
Script for looping and checking each instance.
>>> for gv in self.iter(): ... if self.isGroup(gv): ... # is group ... elif self.isDimension(gv): ... # is dimension ... elif self.isVariable(gv): ... # is variable
- plot.geometry(*, axes: Axes = ['x', 'y', 'z'], atoms: AtomsIndex = None, atoms_style: Sequence[AtomsStyleSpec] = [], atoms_scale: float = 1.0, atoms_colorscale: Colorscale | None = None, drawing_mode: Literal['scatter', 'balls', None] = None, bind_bonds_to_ats: bool = True, points_per_bond: int = 20, bonds_style: StyleSpec = {}, bonds_scale: float = 1.0, bonds_colorscale: Colorscale | None = None, show_atoms: bool = True, show_bonds: bool = True, show_cell: Literal['box', 'axes', False] = 'box', cell_style: StyleSpec = {}, nsc: tuple[int, int, int] = (1, 1, 1), atoms_ndim_scale: tuple[float, float, float] = (16, 16, 1), bonds_ndim_scale: tuple[float, float, float] = (1, 1, 10), dataaxis_1d: np.ndarray | Callable | None = None, arrows: Sequence[AtomArrowSpec] = (), backend='plotly') GeometryPlot
Calls
read_geometry
and creates aGeometryPlot
from its output.- Parameters:
axes – The axes to project the geometry to.
atoms – The atoms to plot. If None, all atoms are plotted.
atoms_style – List of style specifications for the atoms. See the showcase notebooks for examples.
atoms_scale – Scaling factor for the size of all atoms.
atoms_colorscale – Colorscale to use for the atoms in case the color attribute is an array of values. If None, the default colorscale is used for each backend.
drawing_mode – The method used to draw the atoms.
bind_bonds_to_ats – Whether to display only bonds between atoms that are being displayed.
points_per_bond – When the points are drawn using points instead of lines (e.g. in some frameworks to draw multicolor bonds), the number of points used per bond.
bonds_style – Style specification for the bonds. See the showcase notebooks for examples.
bonds_scale – Scaling factor for the width of all bonds.
bonds_colorscale – Colorscale to use for the bonds in case the color attribute is an array of values. If None, the default colorscale is used for each backend.
show_atoms – Whether to display the atoms.
show_bonds – Whether to display the bonds.
show_cell – Mode to display the cell. If False, the cell is not displayed.
cell_style – Style specification for the cell. See the showcase notebooks for examples.
nsc – Number of unit cells to display in each direction.
atoms_ndim_scale – Scaling factor for the size of the atoms for different dimensionalities (1D, 2D, 3D).
bonds_ndim_scale – Scaling factor for the width of the bonds for different dimensionalities (1D, 2D, 3D).
dataaxis_1d – Only meaningful for 1D plots. The data to plot on the Y axis.
arrows – List of arrow specifications to display. See the showcase notebooks for examples.
backend – The backend to use to generate the figure.
See also
GeometryPlot
The plot class used to generate the plot.
read_geometry
The method called to get the data.
- plot.grid(name, index=0, *, data_kwargs={}, axes: Axes = ['z'], represent: Literal['real', 'imag', 'mod', 'phase', 'deg_phase', 'rad_phase'] = 'real', transforms: Sequence[str | Callable] = (), reduce_method: Literal['average', 'sum'] = 'average', boundary_mode: str = 'grid-wrap', nsc: tuple[int, int, int] = (1, 1, 1), interp: tuple[int, int, int] = (1, 1, 1), isos: Sequence[dict] = [], smooth: bool = False, colorscale: Colorscale | None = None, crange: tuple[float, float] | None = None, cmid: float | None = None, show_cell: Literal['box', 'axes', False] = 'box', cell_style: dict = {}, x_range: Sequence[float] | None = None, y_range: Sequence[float] | None = None, z_range: Sequence[float] | None = None, plot_geom: bool = False, geom_kwargs: dict = {}, backend: str = 'plotly') GridPlot
Calls
read_grid
and creates aGridPlot
from its output.- Parameters:
name (
str
) – name of the grid variable to readindex (
int
orarray_like
, optional) – the spin-index for retrieving one of the components. If a vector is passed it refers to the fraction per indexed component. I.e.[0.5, 0.5]
will return sum of half the first two components. Default to the first component.spin (optional) – same as index argument. spin argument has precedence.
axes – The axes to project the grid to.
represent – The representation of the grid to plot.
transforms – List of transforms to apply to the grid before plotting.
reduce_method – The method used to reduce the grid axes that are not displayed.
boundary_mode – The method used to deal with the boundary conditions. Only used if the grid is to be orthogonalized. See scipy docs for more info on the possible values.
nsc – The number of unit cells to display in each direction.
interp – The interpolation factor to use for each axis to make the grid smoother.
isos – List of isosurfaces or isocontours to plot. See the showcase notebooks for examples.
smooth – Whether to ask the plotting backend to make an attempt at smoothing the grid display.
colorscale – Colorscale to use for the grid display in the 2D representation. If None, the default colorscale is used for each backend.
crange – Min and max values for the colorscale.
cmid – The value at which the colorscale is centered.
show_cell – Method used to display the unit cell. If False, the cell is not displayed.
cell_style – Style specification for the cell. See the showcase notebooks for examples.
x_range – The range of the x axis to take into account. Even if the X axis is not displayed! This is important because the reducing operation will only be applied on this range.
y_range – The range of the y axis to take into account. Even if the Y axis is not displayed! This is important because the reducing operation will only be applied on this range.
z_range – The range of the z axis to take into account. Even if the Z axis is not displayed! This is important because the reducing operation will only be applied on this range.
plot_geom – Whether to plot the associated geometry (if any).
geom_kwargs – Keyword arguments to pass to the geometry plot of the associated geometry.
backend – The backend to use to generate the figure.
See also
GridPlot
The plot class used to generate the plot.
read_grid
The method called to get the data.
- read(*args, **kwargs)
Generic read method which should be overloaded in child-classes
- Parameters:
kwargs – keyword arguments will try and search for the attribute
read_<>
and call it with the remaining**kwargs
as arguments.
- read_density_matrix(**kwargs) DensityMatrix [source]
Returns a density matrix from the underlying NetCDF file
- read_dynamical_matrix(**kwargs) DynamicalMatrix [source]
Returns a dynamical matrix from the underlying NetCDF file
This assumes that the dynamical matrix is stored in the field “H” as would the Hamiltonian. This is counter-intuitive but is required when using PHtrans.
- read_energy_density_matrix(**kwargs) EnergyDensityMatrix [source]
Returns energy density matrix from the underlying NetCDF file
- read_force_constant()
Reads the force-constant stored in the nc file
- Returns:
force constants (
numpy.ndarray with 5 dimensions containing all the forces. The 2nd dimensions contains
) – contains the directions, and 3rd dimensions contains -/+ displacements.
- read_grid(name, index=0, **kwargs) Grid [source]
Reads a grid in the current Siesta.nc file
Enables the reading and processing of the grids created by Siesta
- Parameters:
name (
str
) – name of the grid variable to readindex (
int
orarray_like
, optional) – the spin-index for retrieving one of the components. If a vector is passed it refers to the fraction per indexed component. I.e.[0.5, 0.5]
will return sum of half the first two components. Default to the first component.spin (optional) – same as index argument. spin argument has precedence.
- read_hamiltonian(**kwargs) Hamiltonian [source]
Returns a Hamiltonian from the underlying NetCDF file
- read_hessian()[source]
Reads the force-constant stored in the nc file
- Returns:
force constants (
numpy.ndarray with 5 dimensions containing all the forces. The 2nd dimensions contains
) – contains the directions, and 3rd dimensions contains -/+ displacements.
- write(*args, **kwargs)
Generic write method which should be overloaded in child-classes
- Parameters:
**kwargs – keyword arguments will try and search for the attribute write_ and call it with the remaining
**kwargs
as arguments.
- write_basis(atoms: Atoms)[source]
Write the current atoms orbitals as the basis
- Parameters:
atoms – atom specifications to write.
- write_density_matrix(DM, **kwargs)[source]
Writes density matrix model to file
- Parameters:
DM (
DensityMatrix
) – the model to be saved in the NC file
- write_dynamical_matrix(D, **kwargs)[source]
Writes dynamical matrix model to file
- Parameters:
D (
DynamicalMatrix
) – the model to be saved in the NC file
- write_energy_density_matrix(EDM, **kwargs)[source]
Writes energy density matrix model to file
- Parameters:
EDM (
EnergyDensityMatrix
) – the model to be saved in the NC file
- write_hamiltonian(H, **kwargs)[source]
Writes Hamiltonian model to file
- Parameters:
H (
Hamiltonian
) – the model to be saved in the NC fileEf (
float
, optional) – the Fermi level of the electronic structure (in eV), default to 0.
- property base_file
File of the current Sile
- property file
File of the current Sile
- property grids
Return a list of available grids in this file.
- plot
Plotting functions for the
ncSileSiesta
class.