Skip to content

septop #

Automated BFE calculations using the seperated topology method

Classes:

Functions:

  • load_config

    Load a configuration from a YAML file.

  • create_state_dicts

    Map the lambda states specified in the configuration to a dictionary.

  • compute_ddg

    Computes the binding free energy from the complex and solution phase samples.

  • equilibrate_states

    Equilibrate the system at each lambda window.

  • run_hremd

    Perform replica exchange sampling for a system prepared for SepTop calculations.

  • setup_complex

    Prepares a system ready for running the SepTop method.

  • setup_solution

    Prepares a system ready for running the SepTop method.

  • run_solution_phase

    Run the solution phase of the SepTop calculation.

  • run_complex_phase

    Run the complex phase of the SepTop calculation.

  • submit_network

    Submits a set of SepTop calculations to an HPC queueing manager.

Attributes:

DEFAULT_LAMBDA_VDW_1_COMPLEX module-attribute #

DEFAULT_LAMBDA_VDW_1_COMPLEX = (
    [0.0] * 8 + [0.0, 0.0, 0.0] + tolist()
)

The default vdW lambda schedule of the first ligand in the complex phase.

DEFAULT_LAMBDA_CHARGES_1_COMPLEX module-attribute #

DEFAULT_LAMBDA_CHARGES_1_COMPLEX = (
    [0.0] * 8 + [0.25, 0.5, 0.75] + [1.0] * 8
)

The default charge lambda schedule of the first ligand in the complex phase.

DEFAULT_LAMBDA_VDW_2_COMPLEX module-attribute #

DEFAULT_LAMBDA_VDW_2_COMPLEX = (
    tolist() + [0.0, 0.0, 0.0] + [0.0] * 8
)

The default vdW lambda schedule of the second ligand in the complex phase.

DEFAULT_LAMBDA_CHARGES_2_COMPLEX module-attribute #

DEFAULT_LAMBDA_CHARGES_2_COMPLEX = (
    [1.0] * 8 + [0.75, 0.5, 0.25] + [0.0] * 8
)

The default charge lambda schedule of the second ligand in the complex phase.

DEFAULT_LAMBDA_BORESCH_LIGAND_1 module-attribute #

DEFAULT_LAMBDA_BORESCH_LIGAND_1 = (
    [0.0, 0.05, 0.1, 0.3, 0.5, 0.75, 1.0, 1.0]
    + [1.0] * 3
    + [1.0] * 8
)

The default lambda schedule of the Boresch restraint on the first ligand in the complex phase.

DEFAULT_LAMBDA_BORESCH_LIGAND_2 module-attribute #

DEFAULT_LAMBDA_BORESCH_LIGAND_2 = (
    [1.0] * 8
    + [1.0] * 3
    + [1.0, 0.95, 0.9, 0.7, 0.5, 0.25, 0.0, 0.0]
)

The default lambda schedule of the Boresch restraint on the second ligand in the complex phase.

DEFAULT_LAMBDA_VDW_1_SOLUTION module-attribute #

DEFAULT_LAMBDA_VDW_1_SOLUTION = [
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.12,
    0.24,
    0.36,
    0.48,
    0.6,
    0.7,
    0.77,
    0.85,
    1.0,
]

The default vdW lambda schedule of the first ligand in the solution phase.

DEFAULT_LAMBDA_CHARGES_1_SOLUTION module-attribute #

DEFAULT_LAMBDA_CHARGES_1_SOLUTION = [
    0.0,
    0.125,
    0.25,
    0.375,
    0.5,
    0.625,
    0.75,
    0.875,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
]

The default charge lambda schedule of the first ligand in the solution phase.

DEFAULT_LAMBDA_VDW_2_SOLUTION module-attribute #

DEFAULT_LAMBDA_VDW_2_SOLUTION = [
    1.0,
    0.85,
    0.77,
    0.7,
    0.6,
    0.48,
    0.36,
    0.24,
    0.12,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
    0.0,
]

The default vdW lambda schedule of the second ligand in the solution phase.

DEFAULT_LAMBDA_CHARGES_2_SOLUTION module-attribute #

DEFAULT_LAMBDA_CHARGES_2_SOLUTION = [
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    1.0,
    0.875,
    0.75,
    0.625,
    0.5,
    0.375,
    0.25,
    0.125,
    0.0,
]

The default charge lambda schedule of the second ligand in the solution phase.

DEFAULT_BORESCH_K_DISTANCE module-attribute #

DEFAULT_BORESCH_K_DISTANCE = 20.0 * _KCAL_PER_ANG_SQR

The default force constant of the Boresch distance restraint.

DEFAULT_BORESCH_K_THETA module-attribute #

DEFAULT_BORESCH_K_THETA = 20.0 * _KCAL_PER_RAD_SQR

The default force constant of the Boresch angle restraint.

DEFAULT_RESTRAINT_MASK module-attribute #

DEFAULT_RESTRAINT_MASK = '!(:WAT,CL,NA,K) & !@/H'

The default Amber style selection mask to apply position restraints to.

DEFAULT_EQUILIBRATE_INTEGRATOR module-attribute #

DEFAULT_EQUILIBRATE_INTEGRATOR = LangevinIntegrator(
    timestep=2.0 * femtosecond, friction=1.0 / picosecond
)

The default integrator to use during equilibration.

DEFAULT_EQUILIBRATE_RESTRAINTS module-attribute #

DEFAULT_EQUILIBRATE_RESTRAINTS = {
    DEFAULT_RESTRAINT_MASK: FlatBottomRestraint(
        k=25.0 * _KCAL_PER_ANG_SQR, radius=1.5 * _ANGSTROM
    )
}

The default position restraints to apply during equilibration.

LAMBDA_BORESCH_LIGAND_1 module-attribute #

LAMBDA_BORESCH_LIGAND_1 = 'lambda_boresch_lig_1'

The name of the context variable used to control the Boresch-style restraints on the first ligand.

LAMBDA_BORESCH_LIGAND_2 module-attribute #

LAMBDA_BORESCH_LIGAND_2 = 'lambda_boresch_lig_2'

The name of the context variable used to control the Boresch-style restraints on the second ligand.

SepTopComplexRestraints pydantic-model #

Bases: BoreschRestraint

Configure the restraints to apply in the complex phase.

k_distance pydantic-field #

k_distance: OpenMMQuantity[_KCAL_PER_ANG_SQR]

Force constant [kcal/mol/Å^2] of the harmonic distance restraint between r3 and l1.

k_angle_a pydantic-field #

k_angle_a: OpenMMQuantity[_KCAL_PER_RAD_SQR]

Force constant [kcal/mol/rad^2] of the harmonic angle restraint on the angle formed by r2, r3, and l1.

k_angle_b pydantic-field #

k_angle_b: OpenMMQuantity[_KCAL_PER_RAD_SQR]

Force constant [kcal/mol/rad^2] of the harmonic angle restraint on the angle formed by r3, l1, and l2.

k_dihedral_a pydantic-field #

k_dihedral_a: OpenMMQuantity[_KCAL_PER_RAD_SQR]

Force constant [kcal/mol/rad^2] of the harmonic dihedral restraint on the dihedral angle formed by r1, r2, r3, and l1.

k_dihedral_b pydantic-field #

k_dihedral_b: OpenMMQuantity[_KCAL_PER_RAD_SQR]

Force constant [kcal/mol/rad^2] of the harmonic dihedral restraint on the dihedral angle formed by r2, r3, l1, and l2.

k_dihedral_c pydantic-field #

k_dihedral_c: OpenMMQuantity[_KCAL_PER_RAD_SQR]

Force constant [kcal/mol/rad^2] of the harmonic dihedral restraint on the dihedral angle formed by r3, l1, l2, and l3.

SepTopSolutionRestraints pydantic-model #

Bases: BaseModel

Configure the restraints to apply in the solution phase.

SepTopSetupStage pydantic-model #

Bases: BaseModel

Configure how the complex will be solvated and restrained prior to equilibration

SepTopStates pydantic-model #

Bases: BaseModel

Configure the lambda schedules.

SepTopEquilibrateStage pydantic-model #

Bases: BaseModel

Configure how the system will be equilibrated prior to replica exchange.

SepTopSamplingStage pydantic-model #

Bases: HREMD

Configure how the system will be sampled using Hamiltonian replica exchange.

temperature pydantic-field #

temperature: OpenMMQuantity[kelvin] = DEFAULT_TEMPERATURE

The temperature to sample at.

n_warmup_steps pydantic-field #

n_warmup_steps: int = 150000

The number of steps to run each replica for before starting hremd trials. All energies gathered during this period will be discarded.

n_steps_per_cycle pydantic-field #

n_steps_per_cycle: int = 1000

The number of steps to propagate the system by before attempting an exchange.

n_cycles pydantic-field #

n_cycles: int = 2500

The number of cycles of 'propagate the system' -> 'exchange replicas' to run.

max_step_retries pydantic-field #

max_step_retries: int = 5

The maximum number of times to attempt to step if a NaN is encountered before raising an exception

swap_mode pydantic-field #

swap_mode: HREMDSwapModeLiteral | None = value

The mode in which to propose state swaps between replicas. This can either be: 'neighbours', only try and swap adjacent states or ii. 'all', try and swap all states stochastically. If None, no replica exchanges will be attempted.

max_swaps pydantic-field #

max_swaps: int | None = None

The maximum number of swap proposals to make if running in 'all' mode. This variable does nothing when running in 'neighbours' mode.

trajectory_interval pydantic-field #

trajectory_interval: int | None = None

The number of cycles to run before saving the current replica states to DCD trajectory files. If None, no trajectories will be saved.

trajectory_enforce_pbc pydantic-field #

trajectory_enforce_pbc: bool = False

Whether to apply periodic boundary conditions when retrieving coordinates for writing to trajectory files.

checkpoint_interval pydantic-field #

checkpoint_interval: int | None = None

The number of cycles to run before saving the current replica states to checkpoint files. If None, no checkpoints will be saved.

SepTopPhaseConfig pydantic-model #

Bases: BaseModel

Configure one phase (i.e. complex or solution) of a separated topology FE calculation.

SepTopConfig pydantic-model #

Bases: BaseModel

Configuration a separated topology FE calculation.

load_config #

load_config(path: Path) -> SepTopConfig

Load a configuration from a YAML file.

Parameters:

  • path (Path) –

    The path to the YAML configuration.

Returns:

Source code in femto/fe/septop/_config.py
def load_config(path: pathlib.Path) -> SepTopConfig:
    """Load a configuration from a YAML file.

    Args:
        path: The path to the YAML configuration.

    Returns:
        The loaded configuration.
    """
    return SepTopConfig(**yaml.safe_load(path.read_text()))

create_state_dicts #

create_state_dicts(
    config: SepTopStates, system: System
) -> list[dict[str, float]]

Map the lambda states specified in the configuration to a dictionary.

Parameters:

  • config (SepTopStates) –

    The configuration.

  • system (System) –

    The system being simulated.

Returns:

  • list[dict[str, float]]

    The dictionary of lambda states.

Source code in femto/fe/septop/_utils.py
def create_state_dicts(
    config: "femto.fe.septop.SepTopStates", system: openmm.System
) -> list[dict[str, float]]:
    """Map the lambda states specified in the configuration to a dictionary.

    Args:
        config: The configuration.
        system: The system being simulated.

    Returns:
        The dictionary of lambda states.
    """
    from femto.fe.septop import LAMBDA_BORESCH_LIGAND_1, LAMBDA_BORESCH_LIGAND_2

    states = [
        {
            LAMBDA_VDW_LIGAND_1: config.lambda_vdw_ligand_1[i],
            LAMBDA_CHARGES_LIGAND_1: config.lambda_charges_ligand_1[i],
            **(
                {LAMBDA_VDW_LIGAND_2: config.lambda_vdw_ligand_2[i]}
                if config.lambda_vdw_ligand_2 is not None
                else {}
            ),
            **(
                {LAMBDA_CHARGES_LIGAND_2: config.lambda_charges_ligand_2[i]}
                if config.lambda_charges_ligand_2 is not None
                else {}
            ),
            **(
                {LAMBDA_BORESCH_LIGAND_1: config.lambda_boresch_ligand_1[i]}
                if config.lambda_boresch_ligand_1 is not None
                else {}
            ),
            **(
                {LAMBDA_BORESCH_LIGAND_2: config.lambda_boresch_ligand_2[i]}
                if config.lambda_boresch_ligand_2 is not None
                else {}
            ),
            **(
                {femto.md.rest.REST_CTX_PARAM: config.bm_b0[i]}
                if config.bm_b0 is not None
                else {}
            ),
        }
        for i in range(len(config.lambda_vdw_ligand_1))
    ]

    return [
        femto.md.utils.openmm.evaluate_ctx_parameters(state, system) for state in states
    ]

compute_ddg #

compute_ddg(
    config: SepTopConfig,
    complex_u_kn: ndarray,
    complex_n_k: ndarray,
    complex_system: System,
    solution_u_kn: ndarray,
    solution_n_k: ndarray,
    solution_system: System,
) -> DataFrame

Computes the binding free energy from the complex and solution phase samples.

Parameters:

  • config (SepTopConfig) –

    The configuration.

  • complex_u_kn (ndarray) –

    The complex phase samples.

  • complex_n_k (ndarray) –

    The complex phase sample counts.

  • complex_system (System) –

    The complex phase system.

  • solution_u_kn (ndarray) –

    The solution phase samples.

  • solution_n_k (ndarray) –

    The solution phase sample counts.

  • solution_system (System) –

    The solution phase system.

Returns:

  • DataFrame

    A pandas DataFrame containing the total binding free energy and its components.

Source code in femto/fe/septop/_analyze.py
def compute_ddg(
    config: "femto.fe.septop.SepTopConfig",
    complex_u_kn: numpy.ndarray,
    complex_n_k: numpy.ndarray,
    complex_system: openmm.System,
    solution_u_kn: numpy.ndarray,
    solution_n_k: numpy.ndarray,
    solution_system: openmm.System,
) -> "pandas.DataFrame":
    """Computes the binding free energy from the complex and solution phase samples.

    Args:
        config: The configuration.
        complex_u_kn: The complex phase samples.
        complex_n_k: The complex phase sample counts.
        complex_system: The complex phase system.
        solution_u_kn: The solution phase samples.
        solution_n_k: The solution phase sample counts.
        solution_system: The solution phase system.

    Returns:
        A pandas DataFrame containing the total binding free energy and its components.
    """
    import pandas
    import femto.fe.septop

    samples = {
        "complex": (complex_u_kn, complex_n_k),
        "solution": (solution_u_kn, solution_n_k),
    }
    results = {}

    for phase in "complex", "solution":
        phase_config: femto.fe.septop.SepTopPhaseConfig = getattr(config, phase)
        phase_u_kn, phase_n_k = samples[phase]

        estimated, _ = femto.fe.ddg.estimate_ddg(
            phase_u_kn, phase_n_k, phase_config.sample.temperature
        )
        del estimated["ddG_0_kcal_mol"]
        del estimated["ddG_0_error_kcal_mol"]

        results.update({f"{phase}_{k}": v for k, v in estimated.items()})

    results.update(compute_complex_correction(config.complex, complex_system))
    results.update(compute_solution_correction(config.solution, solution_system))

    results["ddG_kcal_mol"] = (
        results["complex_ddG_kcal_mol"] + results["complex_ddG_correction_kcal_mol"]
    ) - (results["solution_ddG_kcal_mol"] + results["solution_ddG_correction_kcal_mol"])

    results["ddG_error_kcal_mol"] = numpy.sqrt(
        results["complex_ddG_error_kcal_mol"] ** 2
        + results["solution_ddG_error_kcal_mol"] ** 2
    )

    return pandas.DataFrame([results])

equilibrate_states #

equilibrate_states(
    system: System,
    topology: Structure,
    states: SepTopStates,
    config: SepTopEquilibrateStage,
    platform: OpenMMPlatform,
    reporter: Reporter | None = None,
) -> list[State]

Equilibrate the system at each lambda window.

Parameters:

  • system (System) –

    The system to simulate.

  • topology (Structure) –

    The topology of the system to simulate.

  • states (SepTopStates) –

    The states of the system to simulate.

  • config (SepTopEquilibrateStage) –

    Configuration settings.

  • platform (OpenMMPlatform) –

    The accelerator to use.

  • reporter (Reporter | None, default: None ) –

    The (optional) reporter to use to record system statistics such as volume and energy.

Returns:

  • list[State]

    The final equilibrated state.

Source code in femto/fe/septop/_equilibrate.py
def equilibrate_states(
    system: openmm.System,
    topology: parmed.Structure,
    states: "femto.fe.septop.SepTopStates",
    config: "femto.fe.septop.SepTopEquilibrateStage",
    platform: femto.md.constants.OpenMMPlatform,
    reporter: femto.md.reporting.Reporter | None = None,
) -> list[openmm.State]:
    """Equilibrate the system at each lambda window.

    Args:
        system: The system to simulate.
        topology: The topology of the system to simulate.
        states: The states of the system to simulate.
        config: Configuration settings.
        platform: The accelerator to use.
        reporter: The (optional) reporter to use to record system statistics such as
            volume and energy.

    Returns:
        The final equilibrated state.
    """
    import femto.fe.septop

    reporter = femto.md.reporting.NullReporter() if reporter is None else reporter

    openmm_reporter = femto.md.reporting.openmm.OpenMMStateReporter(
        reporter, "equilibration", config.report_interval
    )

    state_dicts = femto.fe.septop.create_state_dicts(states, system)

    equilibrated_coords = femto.md.simulate.simulate_states(
        system, topology, state_dicts, config.stages, platform, openmm_reporter
    )
    return equilibrated_coords

run_hremd #

run_hremd(
    system: System,
    topology: Structure,
    coords: list[State],
    states: SepTopStates,
    config: SepTopSamplingStage,
    platform: OpenMMPlatform,
    output_dir: Path,
    reporter: Reporter | None = None,
)

Perform replica exchange sampling for a system prepared for SepTop calculations.

Parameters:

  • system (System) –

    The system.

  • topology (Structure) –

    The topology associated with the system.

  • coords (list[State]) –

    The starting coordinates for each state.

  • states (SepTopStates) –

    The lambda states to sample.

  • config (SepTopSamplingStage) –

    Configuration settings.

  • platform (OpenMMPlatform) –

    The platform to run on.

  • output_dir (Path) –

    The directory to store the sampled energies and statistics to, and any trajectory files if requested.

  • reporter (Reporter | None, default: None ) –

    The reporter to log statistics such as online estimates of the free energy to.

Source code in femto/fe/septop/_sample.py
def run_hremd(
    system: openmm.System,
    topology: parmed.Structure,
    coords: list[openmm.State],
    states: "femto.fe.septop.SepTopStates",
    config: "femto.fe.septop.SepTopSamplingStage",
    platform: femto.md.constants.OpenMMPlatform,
    output_dir: pathlib.Path,
    reporter: femto.md.reporting.Reporter | None = None,
):
    """Perform replica exchange sampling for a system prepared for SepTop calculations.

    Args:
        system: The system.
        topology: The topology associated with the system.
        coords: The starting coordinates for each state.
        states: The lambda states to sample.
        config: Configuration settings.
        platform: The platform to run on.
        output_dir: The directory to store the sampled energies and statistics to, and
            any trajectory files if requested.
        reporter: The reporter to log statistics such as online estimates of the
            free energy to.
    """
    import femto.fe.septop

    system = copy.deepcopy(system)

    n_barostats = sum(
        1
        for force in system.getForces()
        if isinstance(force, openmm.MonteCarloBarostat)
    )
    if n_barostats > 0:
        raise RuntimeError("the system should not contain a barostat already")

    if config.pressure is not None:
        barostat = openmm.MonteCarloBarostat(
            config.pressure, config.temperature, config.barostat_frequency
        )
        system.addForce(barostat)

    femto.md.utils.openmm.assign_force_groups(system)

    state_dicts = femto.fe.septop.create_state_dicts(states, system)

    integrator = femto.md.utils.openmm.create_integrator(
        config.integrator, config.temperature
    )
    simulation = femto.md.utils.openmm.create_simulation(
        system, topology, coords[0], integrator, state_dicts[0], platform
    )

    analysis_fn = None

    if reporter is not None and config.analysis_interval is not None:
        analysis_fn = functools.partial(_analyze, config=config, reporter=reporter)

    return femto.md.hremd.run_hremd(
        simulation,
        state_dicts,
        config,
        output_dir,
        initial_coords=coords,
        analysis_fn=analysis_fn,
        analysis_interval=config.analysis_interval,
    )

setup_complex #

setup_complex(
    config: SepTopSetupStage,
    receptor: AmberParm,
    ligand_1: AmberParm,
    ligand_2: AmberParm | None,
    receptor_ref_query: tuple[str, str, str] | None = None,
    ligand_1_ref_query: tuple[str, str, str] | None = None,
    ligand_2_ref_query: tuple[str, str, str] | None = None,
) -> tuple[Structure, System]

Prepares a system ready for running the SepTop method.

Returns:

  • tuple[Structure, System]

    The prepared topology and OpenMM system object.

Source code in femto/fe/septop/_setup.py
def setup_complex(
    config: "femto.fe.septop.SepTopSetupStage",
    receptor: parmed.amber.AmberParm,
    ligand_1: parmed.amber.AmberParm,
    ligand_2: parmed.amber.AmberParm | None,
    receptor_ref_query: tuple[str, str, str] | None = None,
    ligand_1_ref_query: tuple[str, str, str] | None = None,
    ligand_2_ref_query: tuple[str, str, str] | None = None,
) -> tuple[parmed.Structure, openmm.System]:
    """Prepares a system ready for running the SepTop method.

    Returns:
        The prepared topology and OpenMM system object.
    """
    import femto.fe.septop

    config = copy.deepcopy(config)

    restraint_config = config.restraints

    if not isinstance(restraint_config, femto.fe.septop.SepTopComplexRestraints):
        raise ValueError("invalid restraint config")

    system, topology, ligand_1_ref_idxs, ligand_2_ref_idxs = _setup_system(
        config, ligand_1, ligand_2, receptor, ligand_1_ref_query, ligand_2_ref_query
    )

    ligands = [ligand_1] + ([] if ligand_2 is None else [ligand_2])
    idx_offset = sum(len(ligand.atoms) for ligand in ligands)

    _LOGGER.info("applying restraints.")

    if receptor_ref_query is None:
        _LOGGER.info("selecting receptor reference atoms")

        receptor_ref_idxs_1 = femto.fe.reference.select_receptor_idxs(
            receptor, ligand_1, ligand_1_ref_idxs
        )
        receptor_ref_idxs_2 = receptor_ref_idxs_1

        if ligand_2 is not None and not femto.fe.reference.check_receptor_idxs(
            receptor, receptor_ref_idxs_1, ligand_2, ligand_1_ref_idxs
        ):
            _LOGGER.info("selecting alternate receptor reference atoms for ligand 2")
            receptor_ref_idxs_2 = femto.fe.reference.select_receptor_idxs(
                receptor, ligand_2, ligand_2_ref_idxs
            )

    else:
        receptor_ref_idxs_1 = femto.fe.reference.queries_to_idxs(
            receptor, receptor_ref_query
        )
        receptor_ref_idxs_2 = receptor_ref_idxs_1

    _LOGGER.info(f"receptor ref idxs for ligand 1={receptor_ref_idxs_1}")
    receptor_ref_idxs_1 = tuple(i + idx_offset for i in receptor_ref_idxs_1)

    _apply_complex_restraints(
        topology,
        receptor_ref_idxs_1,
        ligand_1_ref_idxs,
        restraint_config,
        system,
        LAMBDA_BORESCH_LIGAND_1,
    )

    if ligand_2 is not None:
        _LOGGER.info(f"receptor ref idxs for ligand 2={receptor_ref_idxs_2}")
        receptor_ref_idxs_2 = tuple(i + idx_offset for i in receptor_ref_idxs_2)

        _apply_complex_restraints(
            topology,
            receptor_ref_idxs_2,
            ligand_2_ref_idxs,
            restraint_config,
            system,
            LAMBDA_BORESCH_LIGAND_2,
        )

    femto.md.utils.openmm.assign_force_groups(system)

    return topology, system

setup_solution #

setup_solution(
    config: SepTopSetupStage,
    ligand_1: AmberParm,
    ligand_2: AmberParm | None,
    ligand_1_ref_query: tuple[str, str, str] | None = None,
    ligand_2_ref_query: tuple[str, str, str] | None = None,
) -> tuple[Structure, System]

Prepares a system ready for running the SepTop method.

Returns:

  • tuple[Structure, System]

    The prepared topology and OpenMM system object.

Source code in femto/fe/septop/_setup.py
def setup_solution(
    config: "femto.fe.septop.SepTopSetupStage",
    ligand_1: parmed.amber.AmberParm,
    ligand_2: parmed.amber.AmberParm | None,
    ligand_1_ref_query: tuple[str, str, str] | None = None,
    ligand_2_ref_query: tuple[str, str, str] | None = None,
) -> tuple[parmed.Structure, openmm.System]:
    """Prepares a system ready for running the SepTop method.

    Returns:
        The prepared topology and OpenMM system object.
    """
    import femto.fe.septop

    config = copy.deepcopy(config)

    if config.solvent.box_padding is None:
        raise NotImplementedError

    restraint_config = config.restraints

    if not isinstance(restraint_config, femto.fe.septop.SepTopSolutionRestraints):
        raise ValueError("invalid restraint config")

    ligand_offset = None

    if ligand_2 is not None:
        ligand_offset = _compute_ligand_offset(ligand_1, ligand_2)
        _offset_ligand(ligand_2, ligand_offset)

    system, topology, ligand_1_ref_idxs, ligand_2_ref_idxs = _setup_system(
        config,
        ligand_1,
        ligand_2,
        None,
        ligand_1_ref_query,
        ligand_2_ref_query,
        None if ligand_2 is None else -ligand_offset,
    )

    if ligand_2 is not None:
        _apply_solution_restraints(
            topology,
            ligand_1_ref_idxs[1],
            ligand_2_ref_idxs[1],
            config.restraints,
            system,
        )
    femto.md.utils.openmm.assign_force_groups(system)

    return topology, system

run_solution_phase #

run_solution_phase(
    config: SepTopConfig,
    ligand_1_coords: Path,
    ligand_1_params: Path,
    ligand_2_coords: Path | None,
    ligand_2_params: Path | None,
    output_dir: Path,
    report_dir: Path | None = None,
    ligand_1_ref_atoms: tuple[str, str, str] | None = None,
    ligand_2_ref_atoms: tuple[str, str, str] | None = None,
)

Run the solution phase of the SepTop calculation.

Parameters:

  • config (SepTopConfig) –

    The configuration.

  • ligand_1_coords (Path) –

    The coordinates of the first ligand.

  • ligand_1_params (Path) –

    The parameters of the first ligand.

  • ligand_2_coords (Path | None) –

    The coordinates of the second ligand.

  • ligand_2_params (Path | None) –

    The parameters of the second ligand.

  • output_dir (Path) –

    The directory to store all outputs in.

  • report_dir (Path | None, default: None ) –

    The directory to store the report in.

  • ligand_1_ref_atoms (tuple[str, str, str] | None, default: None ) –

    The AMBER style query masks that select the first ligands reference atoms.

  • ligand_2_ref_atoms (tuple[str, str, str] | None, default: None ) –

    The AMBER style query masks that select the second ligands reference atoms.

Source code in femto/fe/septop/_runner.py
def run_solution_phase(
    config: "femto.fe.septop.SepTopConfig",
    ligand_1_coords: pathlib.Path,
    ligand_1_params: pathlib.Path,
    ligand_2_coords: pathlib.Path | None,
    ligand_2_params: pathlib.Path | None,
    output_dir: pathlib.Path,
    report_dir: pathlib.Path | None = None,
    ligand_1_ref_atoms: tuple[str, str, str] | None = None,
    ligand_2_ref_atoms: tuple[str, str, str] | None = None,
):
    """Run the solution phase of the SepTop calculation.

    Args:
        config: The configuration.
        ligand_1_coords: The coordinates of the first ligand.
        ligand_1_params: The parameters of the first ligand.
        ligand_2_coords: The coordinates of the second ligand.
        ligand_2_params: The parameters of the second ligand.
        output_dir: The directory to store all outputs in.
        report_dir: The directory to store the report in.
        ligand_1_ref_atoms: The AMBER style query masks that select the first ligands
            reference atoms.
        ligand_2_ref_atoms: The AMBER style query masks that select the second ligands
            reference atoms.
    """

    prepare_fn = functools.partial(
        _prepare_solution_phase,
        config.solution,
        ligand_1_coords,
        ligand_1_params,
        ligand_2_coords,
        ligand_2_params,
        ligand_1_ref_atoms,
        ligand_2_ref_atoms,
    )
    _run_phase(config.solution, prepare_fn, output_dir, report_dir)

run_complex_phase #

run_complex_phase(
    config: SepTopConfig,
    ligand_1_coords: Path,
    ligand_1_params: Path,
    ligand_2_coords: Path | None,
    ligand_2_params: Path | None,
    receptor_coords: Path,
    receptor_params: Path | None,
    output_dir: Path,
    report_dir: Path | None = None,
    ligand_1_ref_atoms: tuple[str, str, str] | None = None,
    ligand_2_ref_atoms: tuple[str, str, str] | None = None,
    receptor_ref_atoms: tuple[str, str, str] | None = None,
)

Run the complex phase of the SepTop calculation.

Parameters:

  • config (SepTopConfig) –

    The configuration.

  • ligand_1_coords (Path) –

    The coordinates of the first ligand.

  • ligand_1_params (Path) –

    The parameters of the first ligand.

  • ligand_2_coords (Path | None) –

    The coordinates of the second ligand.

  • ligand_2_params (Path | None) –

    The parameters of the second ligand.

  • receptor_coords (Path) –

    The coordinates of the receptor.

  • receptor_params (Path | None) –

    The parameters of the receptor.

  • output_dir (Path) –

    The directory to store all outputs in.

  • report_dir (Path | None, default: None ) –

    The directory to store the logs / reports in.

  • ligand_1_ref_atoms (tuple[str, str, str] | None, default: None ) –

    The AMBER style query masks that select the first ligands reference atoms.

  • ligand_2_ref_atoms (tuple[str, str, str] | None, default: None ) –

    The AMBER style query masks that select the second ligands reference atoms.

  • receptor_ref_atoms (tuple[str, str, str] | None, default: None ) –

    The AMBER style query mask that selects the receptor atoms used to align the ligand.

Source code in femto/fe/septop/_runner.py
def run_complex_phase(
    config: "femto.fe.septop.SepTopConfig",
    ligand_1_coords: pathlib.Path,
    ligand_1_params: pathlib.Path,
    ligand_2_coords: pathlib.Path | None,
    ligand_2_params: pathlib.Path | None,
    receptor_coords: pathlib.Path,
    receptor_params: pathlib.Path | None,
    output_dir: pathlib.Path,
    report_dir: pathlib.Path | None = None,
    ligand_1_ref_atoms: tuple[str, str, str] | None = None,
    ligand_2_ref_atoms: tuple[str, str, str] | None = None,
    receptor_ref_atoms: tuple[str, str, str] | None = None,
):
    """Run the complex phase of the SepTop calculation.

    Args:
        config: The configuration.
        ligand_1_coords: The coordinates of the first ligand.
        ligand_1_params: The parameters of the first ligand.
        ligand_2_coords: The coordinates of the second ligand.
        ligand_2_params: The parameters of the second ligand.
        receptor_coords: The coordinates of the receptor.
        receptor_params: The parameters of the receptor.
        output_dir: The directory to store all outputs in.
        report_dir: The directory to store the logs / reports in.
        ligand_1_ref_atoms: The AMBER style query masks that select the first ligands
            reference atoms.
        ligand_2_ref_atoms: The AMBER style query masks that select the second ligands
            reference atoms.
        receptor_ref_atoms: The AMBER style query mask that selects the receptor atoms
            used to align the ligand.
    """

    prepare_fn = functools.partial(
        _prepare_complex_phase,
        config.complex,
        ligand_1_coords,
        ligand_1_params,
        ligand_2_coords,
        ligand_2_params,
        receptor_coords,
        receptor_params,
        ligand_1_ref_atoms,
        ligand_2_ref_atoms,
        receptor_ref_atoms,
    )
    _run_phase(config.complex, prepare_fn, output_dir, report_dir)

submit_network #

submit_network(
    config: SepTopConfig,
    network: Network,
    output_dir: Path,
    queue_options: SLURMOptions,
    mpi_command: list[str] | None = None,
) -> list[tuple[str, str, str]]

Submits a set of SepTop calculations to an HPC queueing manager.

Parameters:

  • config (SepTopConfig) –

    The configuration.

  • network (Network) –

    The network of edges to run.

  • output_dir (Path) –

    The directory to store any outputs in.

  • queue_options (SLURMOptions) –

    The options to use when submitting the jobs.

  • mpi_command (list[str] | None, default: None ) –

    The mpi runner command to use. The default is "srun --mpi=pmix".

Returns:

  • list[tuple[str, str, str]]

    The ids of the submitted jobs.

Source code in femto/fe/septop/_runner.py
def submit_network(
    config: "femto.fe.septop.SepTopConfig",
    network: femto.fe.inputs.Network,
    output_dir: pathlib.Path,
    queue_options: femto.fe.utils.queue.SLURMOptions,
    mpi_command: list[str] | None = None,
) -> list[tuple[str, str, str]]:
    """Submits a set of SepTop calculations to an HPC queueing manager.

    Args:
        config: The configuration.
        network: The network of edges to run.
        output_dir: The directory to store any outputs in.
        queue_options: The options to use when submitting the jobs.
        mpi_command: The mpi runner command to use. The default is
            ``"srun --mpi=pmix"``.

    Returns:
        The ids of the submitted jobs.
    """

    mpi_command = mpi_command if mpi_command is not None else ["srun", "--mpi=pmix"]

    output_dir.mkdir(exist_ok=True, parents=True)

    date_str = datetime.datetime.now().strftime("%Y-%m-%d--%H-%M-%S")
    config_path = output_dir / f"config-{date_str}.yaml"
    config_path.write_text(config.model_dump_yaml(sort_keys=False))

    femto_command = ["femto", "septop", "--config", config_path]

    slurm_job_ids = []

    for edge in network.edges:
        edge_dir = output_dir / f"{edge.ligand_1.name}~{edge.ligand_2.name}"

        complex_output_dir = edge_dir / "complex"
        solution_output_dir = edge_dir / "solution"

        ligand_args = [
            *_create_run_flags(edge.ligand_1, "ligand-1"),
            *_create_run_flags(edge.ligand_2, "ligand-2"),
        ]

        run_solution_id = femto.fe.utils.queue.submit_slurm_job(
            [
                *mpi_command,
                *femto_command,
                "run-solution",
                *ligand_args,
                f"--output-dir={solution_output_dir}",
                f"--report-dir={solution_output_dir}",
            ],
            queue_options,
            edge_dir / f"run-solution-{date_str}.out",
        )
        run_complex_id = femto.fe.utils.queue.submit_slurm_job(
            [
                *mpi_command,
                *femto_command,
                "run-complex",
                *_create_run_flags(network.receptor, "receptor"),
                *ligand_args,
                f"--output-dir={complex_output_dir}",
                f"--report-dir={complex_output_dir}",
            ],
            queue_options,
            edge_dir / f"run-complex-{date_str}.out",
        )

        analyze_id = femto.fe.utils.queue.submit_slurm_job(
            [
                *femto_command,
                "analyze",
                "--complex-samples",
                complex_output_dir / "_sample/samples.arrow",
                "--complex-system",
                complex_output_dir / "_setup/system.xml",
                "--solution-samples",
                solution_output_dir / "_sample/samples.arrow",
                "--solution-system",
                solution_output_dir / "_setup/system.xml",
                "--output",
                edge_dir / "ddg.csv",
            ],
            queue_options,
            edge_dir / f"analyze-{date_str}.out",
            [run_solution_id, run_complex_id],
        )

        slurm_job_ids.append((run_solution_id, run_complex_id, analyze_id))

    return slurm_job_ids