abipy.dynamics package
Submodules
abipy.dynamics.analyzer module
Tools to analyze MD trajectories and compute diffusion coefficients.
- abipy.dynamics.analyzer.read_structure_postac_ucmats(traj_filepath, step_skip)[source]
Read all configurations from an ASE trajectory file.
- Parameters:
- Returns:
initial Structure, (nsteps, natom, 3) array with the Cartesian coords, (nsteps,3,3) array with cell vectors.
- Return type:
tuple with
- abipy.dynamics.analyzer.parse_lammps_input(filepath, verbose=0)[source]
Extract parameters, atoms and structure from a LAMMPS input file. Returns namedtuple with results.
- class abipy.dynamics.analyzer.MdAnalyzer(structure, temperature, times, cart_positions, ucmats, engine, pos_order='tac', evp_df=None | pandas.core.frame.DataFrame)[source]
Bases:
HasPickleIO
High-level interface to read MD trajectories and metadata from external files, compute the MSQD and plot the results.
- classmethod from_abiml_dir(directory, step_skip=1)[source]
Build an instance from a directory containing an ASE trajectory file and a JSON file with the MD parameters as produced by the abiml.py md script.
- Return type:
- classmethod from_hist_file(hist_filepath, step_skip=1)[source]
Build an instance from an ABINIT HIST.nc file.
- Return type:
- classmethod from_vaspruns(filepaths)[source]
Build an instance from a list of Vasprun files (must be ordered in sequence of MD simulation).
- Return type:
- classmethod from_qe_input(filepath, step_skip=1)[source]
Build an instance from a CP/PW input file.
- Conventions for Quantum ESPRESSO input files:
“.pwi” -> pw.x “.cpi” -> cp.x
- classmethod from_lammps_dir(directory, step_skip=1, basename='in.lammps')[source]
Build an instance from a directory containing a LAMMPS input file.
- Return type:
- classmethod from_lammps_input(input_filepath, step_skip=1)[source]
Build an instance from a LAMMPS input file.
- Parameters:
- Return type:
- classmethod from_lammpstrj(traj_filepath, input_filepath, step_skip=1)[source]
Build an instance from a LAMMPS trajectory file and a log file.
- Parameters:
- Return type:
- resample_step(start_at_step, take_every)[source]
Resample the trajectory. Start at iteration start_at_step and increase the timestep by taking every take_every iteration.
- Return type:
- resample_time(start_time, new_timestep)[source]
Resample the trajectory. Start at time start_time and use new timestep new_timestep.
- Return type:
- set_color_symbol(dict_or_string)[source]
Set the dictionary mapping chemical_symbol –> color used in the matplotlib plots.
- get_it_ts(t0)[source]
Return the index of time t0 in self.times and the array with the time values.
- iatoms_with_symbol(symbol, atom_inds=None)[source]
Array with the index of the atoms with the given chemical symbol. If atom_inds is not None, filter sites accordingly.
- Return type:
- get_sqdt_iatom(iatom, it0=0)[source]
Compute the square displacement vs time for a given atomic index starting from time index it0.
- Return type:
array
- get_sqdt_symbol(symbol, it0=0, atom_inds=None)[source]
Compute the square displacement vs time averaged over atoms with the same chemical symbol starting from time index it0. atoms_inds adds an additional filter on the site index.
- Return type:
array
- get_dw_symbol(symbol, t0=0.0, tmax=None, atom_inds=None)[source]
Compute diffusion coefficient by performing a naive linear regression of the raw MQST.
- get_msdtt0_symbol_tmax(symbol, tmax, atom_inds=None, nprocs=None)[source]
Calculates the MSD for every possible pair of time points using the formula:
$$MSD(t,t_0) = frac{1}{N} sum_{i=1}^{N} (vec{r}_i(t+t_0) - vec{r}_i(t_0))^2$$
where $N$ is the number of particles with the given symbol, and $vec{r}_i(t)$ is the position vector.
- plot_sqdt_atoms(symbols='all', t0=0.0, atom_inds=None, ax=None, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot the square displacement of atoms vs time.
- Parameters:
symbols – List of chemical symbols to consider.
t0 (
float
) – Initial time in ps.atom_inds – List of atom indices to include. None to disable filtering.
ax –
matplotlib.axes.Axes
or None if a new figure should be created.xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles.
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_sqdt_symbols(symbols, t0=0.0, atom_inds=None, with_dw=0, ax=None, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot the square displacement averaged over all atoms of the same specie vs time.
- Parameters:
symbols – List of chemical symbols to consider. “all” for all symbols in structure.
t0 (
float
) – Initial time in ps.atom_inds – List of atom indices to include. None to disable filtering.
with_dw – If != 0 compute diffusion coefficient via least-squares fit in the time-interval [t0, with_dw]. If with_dw < 0 e.g. -1 the time-interval is set to [t0, tmax].
ax –
matplotlib.axes.Axes
or None if a new figure should be created.xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles.
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_sqdt_symbols_tmax(symbols, tmax, atom_inds=None, nprocs=None, ax=None, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot the square displacement averaged over all atoms of the same specie vs time.
- Parameters:
symbols – List of chemical symbols to consider. “all” for all symbols in structure.
tmax (
float
) – Max time in ps.atom_inds – List of atom indices to include. None to disable filtering.
nprocs – Number of procs to use.
ax –
matplotlib.axes.Axes
or None if a new figure should be created.xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plot_lattices(what_list=('abc', 'angles', 'volume'), ax_list=None, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot lattice lengths/angles/volume as a function of time.
- Parameters:
what_list – List of strings specifying the quantities to plot. Default all
ax_list – List of axis or None if a new figure should be created.
xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- class abipy.dynamics.analyzer.Msdtt0(*, index_tmax, symbol, arr_tt0, mda)[source]
Bases:
object
This object stores:
$$MSD(t,t_0) = frac{1}{N} sum_{i=1}^{N} (vec{r}_i(t+t_0) - vec{r}_i(t_0))^2$$
where $N$ is the number of particles of a particular chemical symbol and $vec{r}_i(t)$ is the position vector.
-
mda:
MdAnalyzer
- plot(ax=None, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot <msd($t, t_0$)>$ averaged over the initial time t0.
- Parameters:
ax –
matplotlib.axes.Axes
or None if a new figure should be created.xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- get_diffusion_with_sigma(block_size1, block_size2, fit_time_start, fit_time_stop, sigma_berend)[source]
Compute diffusion coefficient with uncertainty.
- Return type:
- plot_mat(cmap='jet', fontsize=8, ax=None, **kwargs)[source]
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
-
mda:
- class abipy.dynamics.analyzer.Msdtt0List(iterable=(), /)[source]
Bases:
list
A list of Msdtt0 objects.
- plot(sharex=True, sharey=True, fontsize=8, **kwargs)[source]
Plot all Msdtt0 objects on a grid.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- class abipy.dynamics.analyzer.SigmaBerend(*, temperature, latex_formula, it1, time1, block_sizes1, sigmas1, delta_sigmas1, it2, time2, block_sizes2, sigmas2, delta_sigmas2)[source]
Bases:
object
Stores the variance of correlated data as function of block number.
- plot(fontsize=8, ax_list=None, **kwargs)[source]
Plot variance of correlated data as function of block number.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- class abipy.dynamics.analyzer.DiffusionData(*, diffusion_coeff, err_diffusion_coeff, conductivity, err_conductivity, temperature, symbol, latex_formula, avg_volume, ncarriers, block_size1, block_size2, fit_time_start, fit_time_stop, min_angcoeff, max_angcoeff, best_angcoeff, engine, msd_t, err_msd, sigma_berend, times, quote, timeArrScorrelated, msdSScorrelated, errMSDScorrelated)[source]
Bases:
HasPickleIO
Diffusion results for a given temperature.
-
sigma_berend:
SigmaBerend
- plot(ax=None, fontsize=8, **kwargs)[source]
Plot MDS(t) with errors.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
-
sigma_berend:
- class abipy.dynamics.analyzer.DiffusionDataList(iterable=(), /)[source]
Bases:
list
A list of DiffusionData objects.
- get_dataframe(add_keys=None)[source]
Dataframe with diffusion results.
- Parameters:
add_keys – optional list of attributes to add.
- Return type:
- plot_sigma_berend(**kwargs)[source]
Plot variance of correlated data as function of block number for all objects stored in DiffusionDataList.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- plot(**kwargs)[source]
Plot … for all objects stored DiffusionDataList.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- get_arrhenius_nmtuples(df, hue)[source]
Return list of namedtuple objects.
- Parameters:
df –
pandas.DataFrame
.hue – Variable defining how to group data. If None, no grouping is performed.
- class abipy.dynamics.analyzer.MultiMdAnalyzer(mdas, temp_colormap='jet')[source]
Bases:
HasPickleIO
High-level interface to analyze multiple MD trajectories
- classmethod from_abiml_dirs(directories, step_skip=1, pmode='processes')[source]
Build an instance from a list of directories produced by abiml.py md
- Return type:
- classmethod from_lammps_dirs(directories, step_skip=1, basename='in.lammps', pmode='processes')[source]
Build an instance from a list of directories containing LAMMPS results.
- Return type:
- classmethod from_vaspruns(vasprun_filepaths, step_skip=1, pmode='processes')[source]
Build an instance from a list of vasprun files.
- Return type:
- classmethod from_hist_files(hist_filepaths, step_skip=1, pmode='processes')[source]
Build an instance from a list of ABINIT HIST.nc files.
- Return type:
- get_params_dataframe()[source]
Dataframe with the parameters of the different MdAnalyzers.
- Return type:
- plot_sqdt_symbols(symbols, t0=0.0, xy_log=None, fontsize=8, xlims=None, **kwargs)[source]
Plot the square displacement averaged over all atoms of the same specie vs time for the different temperatures.
- Parameters:
symbols – List of chemical symbols to consider. “all” for all symbols in structure.
t0 (
float
) – Initial time in ps.ax –
matplotlib.axes.Axes
or None if a new figure should be created.xy_log – None or empty string for linear scale. “x” for log scale on x-axis. “xy” for log scale on x- and y-axis. “x:semilog” for semilog scale on x-axis.
fontsize – fontsize for legends and titles.
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- abipy.dynamics.analyzer.msd_tt0_from_tac(pos_tac, size_t, nprocs=None)[source]
Calculates the MSD for every possible pair of time points in the input array, using the formula:
$$MSD(t,t_0) = frac{1}{N} sum_{i=1}^{N} (vec{r}_i(t+t_0) - vec{r}_i(t_0))^2$$
where $N$ is the number of particles, and $vec{r}_i(t)$ is the position vector.
- abipy.dynamics.analyzer.block_mean_var(data, data_mean, n_block)[source]
Perform the block mean and the block variance of data.
- abipy.dynamics.analyzer.sigma_berend(nblock_step, tot_block, data)[source]
- Parameters:
- Return type:
Return: (data_in_block, sigma, delta_sigma)
- class abipy.dynamics.analyzer.ArrheniusEntry(*, key, symbol, composition, temps, diffusions, err_diffusions, volumes, mpl_style)[source]
Bases:
object
- classmethod from_file(filepath, key, mpl_style)[source]
Read data in CSV format. Assuming header with at least the following entries: temperature,diffusion,err_diffusion,volume,symbol,composition
- Return type:
- class abipy.dynamics.analyzer.ArrheniusPlotter(entries=None)[source]
Bases:
object
This object stores conductivities D(T) computed for different structures and/or different ML-potentials and allows one to produce Arrnehnius plots on the same figure. Internally, the results are indexed by a unique key that is be used as label in the matplotlib plot. The style for each key can be customized by setting the mpl_style dict with the options that will be passed to ax.plot.
In the simplest case, one reads the data from external files in CSV format.
Example
from abipy.dynamics.analyzer import ArrheniusPlotter
- key_path = {
“matgl-MD”: “diffusion_cLLZO-matgl.csv”, “m3gnet-MD”: “diffusion_cLLZO-m3gnet.csv”,
}
- mpl_style_key = {
“matgl-MD” : dict(c=’blue’), “m3gnet-MD”: dict(c=’purple’),
}
plotter = ArrheniusPlotter() for key, path in key_path.items():
plotter.add_entry_from_file(path, key, mpl_style=mpl_style_key[key])
# temperature grid refined thinvt_arange = (0.6, 2.5, 0.01) xlims, ylims = (0.5, 2.0), (-7, -4)
- plotter.plot(thinvt_arange=thinvt_arange,
xlims=xlims, ylims=ylims, text=’LLZO cubic’, savefig=None)
- plot(thinvt_arange=None, what='diffusion', ncar=None, colormap='jet', with_t=True, text=None, ax=None, fontsize=8, xlims=None, ylims=None, **kwargs)[source]
Arrhenius plot.
- Parameters:
thinvt_arange – start, stop, step for 1000/T mesh. If None, the mesh is automatically computed.
what – Selects the quantity to plot. Possibile values: “diffusion”, “sigma”, “tsigma”.
ncar – Number of carriers. Required if what is “sigma” or “tsigma”.
colormap – Colormap used to select the color if entry.mpl_style does not provide it.
with_t – True to add a twin axes with the value of T
text
ax –
matplotlib.axes.Axes
or None if a new figure should be created.fontsize – fontsize for legends and titles.
xlims – Set the data limits for the x-axis. Accept tuple e.g.
(left, right)
or scalar e.g.left
. If left (right) is None, default values are used.ylims – Similar to xlims but for the y-axis.
- Return type:
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
abipy.dynamics.cpx module
Parser for the output files produced by the CP code.
See: https://www.quantum-espresso.org/Doc/INPUT_CP.html
prefix.pos : atomic positions prefix.vel : atomic velocities prefix.for : atomic forces prefix.cel : cell parameters prefix.str : stress tensors prefix.evp : energies prefix.hrs : Hirshfeld effective volumes (ts-vdw) prefix.eig : eigen values prefix.nos : Nose-Hoover variables prefix.spr : spread of Wannier orbitals prefix.wfc : center of Wannier orbitals prefix.ncg : number of Poisson CG steps (PBE0) prefix_ndw.save/ : write restart folder prefix_ndr.save/ : read restart folder
- abipy.dynamics.cpx.parse_file_with_blocks(filepath, block_len)[source]
Parse a QE file whose format is:
20 0.00241888
0.31882936800000E+01 0.14832370390000E+02 0.12288296100000E+01 0.78323146900000E+01 0.67870403900000E+01 0.12288296100000E+01 0.20744346700000E+01 0.59953799200000E+01 0.47375825000000E+01
20 0.00241888
0.31889894893205E+01 0.14833320999242E+02 0.12271083047604E+01 0.78345550025202E+01 0.67881458009071E+01 0.12273445304341E+01 0.20762325213327E+01 0.59955571558384E+01 0.47335647385293E+01
i.e. a header followed by block_len lines
- Parameters:
- Return type:
Return: (number_of_steps, array_to_be_reshaped, list_of_headers)
- abipy.dynamics.cpx.parse_file_with_header(filepath)[source]
Parse a file with an header followed by csv-like columns e.g. an .evp file
- Return type:
- class abipy.dynamics.cpx.EvpFile(filepath)[source]
Bases:
TextFile
,NotebookWriter
The evp file is a file that contains the electronic and ionic energies and pressures for each time step of a CP simulation. The evp file has the following format:
- # nfi tps(ps) ekinc Tcell(K) Tion(K) etot enthal econs
0 0 300 -2099.0164 7.4066063 -2091.6098 19363.353 2188.2406 5.0978969 10 0.01 557.43004 -2105.3429 13.762216 -2091.5807 16917.626 2188.2406 5.0978969
…
NB: Energies are in Hartree and not in Ry.
Inheritance Diagram
- COLS_DICT = {'econs': Key(name='econs', info='etot+kinetic_energy_due_to_ions_moving', color='b'), 'econt': Key(name='econt', info='econs+ekinc+(thermostat_total_energy)', color='b'), 'ekinc': Key(name='ekinc', info='electron fake kinetic energy', color='b'), 'enthal': Key(name='enthal', info='etot+external_pressure*volume', color='b'), 'etot': Key(name='etot', info='dft total energy (without ionic and cell kinetic energy)', color='b'), 'nfi': Key(name='nfi', info='number of iterations', color='b'), 'temphc': Key(name='temphc', info='temperature due to “cell” kinetic energy in K', color='b'), 'tempp': Key(name='tempp', info='temperature due to ionic displacements within the cell', color='b'), 'time(ps)': Key(name='time(ps)', info='time', color='b'), 'tps(ps)': Key(name='tps(ps)', info='time', color='b')}
- plot(**kwargs)[source]
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- plot_hist(**kwargs)[source]
Plot one histogram per column.
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- Return type:
- abipy.dynamics.cpx.traj_to_qepos(traj_filepath, pos_filepath)[source]
Convert an ASE trajectory file to QE POS file.
- class abipy.dynamics.cpx.Qe2Extxyz(code, qe_input, pos_filepath, for_filepath, cell_filepath, evp_filepath, str_filepath=None)[source]
Bases:
object
Convert QE/CP output files into ASE extended xyz format.
Example
from abipy.dynamics.cpx import Qe2Extxyz converter = Qe2Extxyz.from_input(“cp.in”, code=”cp”) converter.write_xyz(“extended.xyz”, take_every=1)
- classmethod from_input(qe_input, code, prefix=None)[source]
Build object from QE/CP input assuming all the other output files are located in the same directory with the given prefix.
- abipy.dynamics.cpx.downsample_xyz(input_xyz, take_every, output_xyz, skip_head=None, verbose=1)[source]
Downsample an XYZ file. Return the number of configurations written in the new file.
- Parameters:
input_xyz – Input XYZ file.
take_every – Extract configurations from input_xyz every take_every step.
output_xyz – Output XYZ file.
skip_head – If not None, skip the first skip_head configurations.
verbose – Verbosity level.
- Return type:
abipy.dynamics.hist module
History file with structural relaxation results.
- class abipy.dynamics.hist.HistFile(filepath)[source]
Bases:
AbinitNcFile
,NotebookWriter
File with the history of a structural relaxation or molecular dynamics calculation.
Usage example:
with HistFile("foo_HIST") as hist: hist.plot()
Inheritance Diagram
- get_fstats_dict(step)[source]
Return
monty.collections.AttrDict
with stats on the forces at the givenstep
.- Return type:
AttrDict
- property initial_structure: Structure
The initial
abipy.core.structure.Structure
.
- property final_structure: Structure
The
abipy.core.structure.Structure
of the last iteration.
- structures()[source]
List of
abipy.core.structure.Structure
objects at the different steps.
- etotals()[source]
numpy.ndarray
with total energies in eV at the different steps.
- get_relaxation_analyzer()[source]
Return a pymatgen
RelaxationAnalyzer
object to analyze the relaxation in a calculation.- Return type:
RelaxationAnalyzer
- to_xdatcar(filepath=None, groupby_type=True, to_unit_cell=False, **kwargs)[source]
Return Xdatcar pymatgen object. See write_xdatcar for the meaning of arguments.
- Parameters:
to_unit_cell (bool) – Whether to translate sites into the unit cell.
kwargs – keywords arguments passed to Xdatcar constructor.
- write_xdatcar(filepath='XDATCAR', groupby_type=True, overwrite=False, to_unit_cell=False)[source]
Write Xdatcar file with unit cell and atomic positions to file
filepath
.- Parameters:
filepath – Xdatcar filename. If None, a temporary file is created.
groupby_type – If True, atoms are grouped by type. Note that this option may change the order of the atoms. This option is needed because there are post-processing tools (e.g. ovito) that do not work as expected if the atoms in the structure are not grouped by type.
overwrite – raise RuntimeError, if False and filepath exists.
to_unit_cell (bool) – Whether to translate sites into the unit cell.
- Return type:
- Returns:
path to Xdatcar file.
- visualize(appname='ovito', to_unit_cell=False)[source]
Visualize the crystalline structure with visualizer. See
Visualizer
for the list of applications and formats supported.- Parameters:
to_unit_cell (bool) – Whether to translate sites into the unit cell.
- plot_ax(ax, what, fontsize=8, **kwargs)[source]
Helper function to plot quantity
what
on axisax
with matplotlib.- Parameters:
fontsize – fontsize for legend.
method. (kwargs are passed to matplotlib plot)
- Return type:
- plotly_traces(fig, what, rcd=None, fontsize=8, showlegend=False, **kwargs)[source]
Helper function to plot quantity
what
on figurefig
with plotly.- Parameters:
rcd – If
fig
has subplots,rcd
is used to add traces on these subplots.fontsize – fontsize for legend.
method. (kwargs are passed to fig.add_scatter)
- plot(what_list=None, ax_list=None, fontsize=8, **kwargs)[source]
Plot the evolution of structural parameters (lattice lengths, angles and volume) as well as pressure, info on forces and total energy with matplotlib.
- Parameters:
what_list
ax_list – List of
matplotlib.axes.Axes
. If None, a new figure is created.fontsize – fontsize for legend
- Return type:
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- plotly(what_list=None, fig=None, fontsize=12, **kwargs)[source]
Plot the evolution of structural parameters (lattice lengths, angles and volume) as well as pressure, info on forces and total energy with plotly.
- Parameters:
what_list
fig – The fig for plot and the DOS plot. If None, a new figure is created.
fontsize – fontsize for legend
Returns:
plotly.graph_objects.Figure
Keyword arguments controlling the display of the figure: ================ ==================================================================== kwargs Meaning ================ ==================================================================== title Title of the plot (Default: None). show True to show the figure (default: True). hovermode True to show the hover info (default: False) savefig “abc.png” , “abc.jpeg” or “abc.webp” to save the figure to a file. write_json Write plotly figure to write_json JSON file.
Inside jupyter-lab, one can right-click the write_json file from the file menu and open with “Plotly Editor”. Make some changes to the figure, then use the file menu to save the customized plotly plot. Requires jupyter labextension install jupyterlab-chart-editor. See https://github.com/plotly/jupyterlab-chart-editor
- renderer (str or None (default None)) –
A string containing the names of one or more registered renderers (separated by ‘+’ characters) or None. If None, then the default renderers specified in plotly.io.renderers.default are used. See https://plotly.com/python-api-reference/generated/plotly.graph_objects.Figure.html
config (dict) A dict of parameters to configure the figure. The defaults are set in plotly.js. chart_studio True to push figure to chart_studio server. Requires authenticatios.
Default: False.
- template Plotly template. See https://plotly.com/python/templates/
- [“plotly”, “plotly_white”, “plotly_dark”, “ggplot2”,
“seaborn”, “simple_white”, “none”]
Default is None that is the default template is used.
- plotly_energies(fig=None, fontsize=12, **kwargs)[source]
Plot the total energies as function of the iteration step with plotly.
Keyword arguments controlling the display of the figure: ================ ==================================================================== kwargs Meaning ================ ==================================================================== title Title of the plot (Default: None). show True to show the figure (default: True). hovermode True to show the hover info (default: False) savefig “abc.png” , “abc.jpeg” or “abc.webp” to save the figure to a file. write_json Write plotly figure to write_json JSON file.
Inside jupyter-lab, one can right-click the write_json file from the file menu and open with “Plotly Editor”. Make some changes to the figure, then use the file menu to save the customized plotly plot. Requires jupyter labextension install jupyterlab-chart-editor. See https://github.com/plotly/jupyterlab-chart-editor
- renderer (str or None (default None)) –
A string containing the names of one or more registered renderers (separated by ‘+’ characters) or None. If None, then the default renderers specified in plotly.io.renderers.default are used. See https://plotly.com/python-api-reference/generated/plotly.graph_objects.Figure.html
config (dict) A dict of parameters to configure the figure. The defaults are set in plotly.js. chart_studio True to push figure to chart_studio server. Requires authenticatios.
Default: False.
- template Plotly template. See https://plotly.com/python/templates/
- [“plotly”, “plotly_white”, “plotly_dark”, “ggplot2”,
“seaborn”, “simple_white”, “none”]
Default is None that is the default template is used.
- plot_energies(ax=None, fontsize=8, **kwargs)[source]
Plot the total energies as function of the iteration step with matplotlib.
- Parameters:
ax –
matplotlib.axes.Axes
or None if a new figure should be created.fontsize – Legend and title fontsize.
- Return type:
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- yield_figs(**kwargs)[source]
This function generates a predefined list of matplotlib figures with minimal input from the user.
- yield_plotly_figs(**kwargs)[source]
This function generates a predefined list of matplotlib figures with minimal input from the user.
- mvplot_trajectories(colormap='hot', sampling=1, figure=None, show=True, with_forces=True, **kwargs)[source]
Call mayavi to plot atomic trajectories and the variation of the unit cell.
- get_panel(**kwargs)[source]
Build panel with widgets to interact with the
abipy.dynamics.hist.HistFile
either in a notebook or in panel app.
- class abipy.dynamics.hist.HistRobot(*args)[source]
Bases:
Robot
This robot analyzes the results contained in multiple HIST.nc files.
Inheritance Diagram
- EXT = 'HIST'
- get_dataframe(with_geo=True, index=None, abspath=False, with_spglib=True, funcs=None, **kwargs)[source]
Return a
pandas.DataFrame
with the most important final results and the filenames as index.- Parameters:
with_geo – True if structure info should be added to the dataframe
abspath – True if paths in index should be absolute. Default: Relative to getcwd().
index – Index of the dataframe, if None, robot labels are used
with_spglib – If True, spglib is invoked to get the space group symbol and number
- Return type:
- kwargs:
- attrs:
List of additional attributes of the
abipy.dynamics.hist.HistFile
to add to thepandas.DataFrame
.- funcs: Function or list of functions to execute to add more data to the DataFrame.
Each function receives a
abipy.dynamics.hist.HistFile
object and returns a tuple (key, value) where key is a string with the name of column and value is the value to be inserted.
- gridplot(what_list=None, sharex='row', sharey='row', fontsize=8, **kwargs)[source]
Plot the
what
value extracted from multiple HIST.nc files on a grid.- Parameters:
what_list – List of quantities to plot. Must be in [“energy”, “abc”, “angles”, “volume”, “pressure”, “forces”]
sharex – True if xaxis should be shared.
sharey – True if yaxis should be shared.
fontsize – fontsize for legend.
- Return type:
Returns:
matplotlib.figure.Figure
Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- combiplot(what_list=None, colormap='jet', fontsize=6, **kwargs)[source]
Plot multiple HIST.nc files on a grid. One plot for each
what
value.- Parameters:
what_list – List of strings with the quantities to plot. If None, all quanties are plotted.
colormap – matplotlib color map.
fontsize – fontisize for legend.
- Return type:
Returns:
matplotlib.figure.Figure
.Keyword arguments controlling the display of the figure:
kwargs
Meaning
title
Title of the plot (Default: None).
show
True to show the figure (default: True).
savefig
“abc.png” or “abc.eps” to save the figure to a file.
size_kwargs
Dictionary with options passed to fig.set_size_inches e.g. size_kwargs=dict(w=3, h=4)
tight_layout
True to call fig.tight_layout (default: False)
ax_grid
True (False) to add (remove) grid from all axes in fig. Default: None i.e. fig is left unchanged.
ax_annotate
Add labels to subplots e.g. (a), (b). Default: False
fig_close
Close figure. Default: False.
plotly
Try to convert mpl figure to plotly.
- class abipy.dynamics.hist.HistReader(path)[source]
Bases:
ETSF_Reader
This object reads data from the HIST file.
Inheritance Diagram
- read_eterms(unit='eV')[source]
monty.collections.AttrDict
with the decomposition of the total energy in unitsunit
- Return type:
AttrDict
- read_cart_forces(unit='eV ang^-1')[source]
Read and return a
numpy.ndarray
with the cartesian forces in unitunit
. Shape (num_steps, natom, 3)- Return type:
- read_reduced_forces()[source]
Read and return a
numpy.ndarray
with the forces in reduced coordinates Shape (num_steps, natom, 3)- Return type: