Amber parsing

Parsers for extracting alchemical data from AMBER output files.

Some of the file parsing parts are adapted from alchemical-analysis.

Changed in version 1.0.0: Now raises ValueError when an invalid file is given to the parser. Now raises ValueError when inconsistency in MBAR states/data is found.

alchemlyb.parsing.amber.convert_to_pandas(file_datum: Any) DataFrame[source]

Convert the data structure from numpy to pandas format

class alchemlyb.parsing.amber.SectionParser(filename: str)[source]

A simple parser to extract data values from sections.

__init__(filename: str) None[source]

Opens a file according to its file type.

skip_lines(nlines: int) None | str[source]

Skip a given number of lines.

skip_after(pattern: str) bool[source]

Skip until after a line that matches a regex pattern.

extract_section(start: str, end: str, fields: list[str], limit: None | str = None, extra: str = '') list[float][source]

Extract data values (int, float) in fields from a section marked with start and end regexes. Do not read further than limit regex.

__next__() str[source]

Read next line of the filehandle and check for EOF.

close() None[source]

Close the filehandle.

class alchemlyb.parsing.amber.FEData[source]

A simple struct container to collect data from individual files.

__init__() None[source]
clambda
t0
dt
T
ntpr
bar_intervall
gradients: list[float]
mbar_energies: list[list[float]]
have_mbar
mbar_lambdas: list[str]
mbar_lambda_idx
alchemlyb.parsing.amber.file_validation(outfile: str) FEData[source]

Function that validate and parse an AMBER output file. ValueError are risen if inconsinstencies in the input file are found.

Parameters:

outfile (str) – Path to AMBER .out file to validate and extract data from.

Returns:

FEData object populated with data from the parsed AMBER output file.

Return type:

FEData

alchemlyb.parsing.amber.extract(outfile: str, T: float) dict[str, None | DataFrame][source]

Return reduced potentials u_nk and gradients dH/dl from AMBER outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed; needed to generated the reduced potential (in units of kT)

Returns:

A dictionary with keys of ‘u_nk’, which is a pandas DataFrame of reduced potentials for each alchemical state (k) for each frame (n), and ‘dHdl’, which is a Series of dH/dl as a function of time for this lambda window.

Return type:

dict

Added in version 1.0.0.

alchemlyb.parsing.amber.extract_dHdl(outfile: str, T: float) None | DataFrame[source]

Return gradients dH/dl from AMBER TI outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed

Returns:

dH/dl – dH/dl as a function of time for this lambda window.

Return type:

pandas.Series

Changed in version 0.5.0: The scipy.constants is used for parsers instead of the constants used by the corresponding MD engine.

alchemlyb.parsing.amber.extract_u_nk(outfile: str, T: float) None | DataFrame[source]

Return reduced potentials u_nk from AMBER outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed; needed to generated the reduced potential (in units of kT)

Returns:

u_nk – Reduced potential for each alchemical state (k) for each frame (n).

Return type:

pandas.DataFrame

Changed in version 0.5.0: The scipy.constants is used for parsers instead of the constants used by the corresponding MD engine.

The parsers featured in this module are constructed to properly parse Amber MD output files containing derivatives of the Hamiltonian and FEP (BAR/MBAR) data.

API Reference

This submodule includes these parsing functions:

alchemlyb.parsing.amber.extract_dHdl(outfile: str, T: float) None | DataFrame[source]

Return gradients dH/dl from AMBER TI outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed

Returns:

dH/dl – dH/dl as a function of time for this lambda window.

Return type:

pandas.Series

Changed in version 0.5.0: The scipy.constants is used for parsers instead of the constants used by the corresponding MD engine.

alchemlyb.parsing.amber.extract_u_nk(outfile: str, T: float) None | DataFrame[source]

Return reduced potentials u_nk from AMBER outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed; needed to generated the reduced potential (in units of kT)

Returns:

u_nk – Reduced potential for each alchemical state (k) for each frame (n).

Return type:

pandas.DataFrame

Changed in version 0.5.0: The scipy.constants is used for parsers instead of the constants used by the corresponding MD engine.

alchemlyb.parsing.amber.extract(outfile: str, T: float) dict[str, None | DataFrame][source]

Return reduced potentials u_nk and gradients dH/dl from AMBER outputfile.

Parameters:
  • outfile (str) – Path to AMBER .out file to extract data from.

  • T (float) – Temperature in Kelvin at which the simulations were performed; needed to generated the reduced potential (in units of kT)

Returns:

A dictionary with keys of ‘u_nk’, which is a pandas DataFrame of reduced potentials for each alchemical state (k) for each frame (n), and ‘dHdl’, which is a Series of dH/dl as a function of time for this lambda window.

Return type:

dict

Added in version 1.0.0.