Loading [Contrib]/a11y/accessibility-menu.js
Skip to content
Scm logo
Scm logo
  • Amsterdam Modeling Suite
    Amsterdam Modeling Suite
    Atomistic Scale
    Electronic Structure
    ADF

    Understand and predict chemical properties with our fast and accurate molecular DFT code.

    Periodic DFT

    BAND & Quantum Espresso: Calculate reactivity, band gaps, optical response, and other properties for periodic systems.

    DFTB & MOPAC

    Model larger molecules and periodic systems, or prescreen many candidates, with the fast electronic structure methods DFTB and MOPAC.

    Interatomic Potentials
    ReaxFF

    Study large, chemically evolving systems with ReaxFF molecular dynamics.

    Machine Learning Potentials

    Use preparametrized ML potentials M3GNET, ANI-1ccx or your own models.

    Force Fields

    GFN-FF, Apple&P, UFF, and more- (polarizable) force fields.

    Meso- & Macroscale
    kMC and Microkinetics

    Predict catalytic turn-over frequencies with microkinetics and kinetic Monte Carlo.

    Bumblebee: OLED stacks

    3D kinetic Monte Carlo for simulating OLED device-level physics

    Fluid Thermodynamics
    COSMO-RS

    Quick physical property predictions, thermodynamic properties in solution, and solvent screening.

    Amsterdam Modeling Suite: computational chemistry with expert support to advance your chemistry & materials R&D

    Discover the Suite Pricing & licensing
  • Applications
  • Tools
    Tools
    Workflows and Utilities
    OLED workflows

    Automatic workflows to simulate physical vapor deposition and calculate properties for OLED device modeling.

    ChemTraYzer2

    Automatically extract reaction pathways and reaction rates from reactive MD trajectories.

    Conformers

    Easily generate, screen, refine, and select conformers. Pass on to other modules for conformational averaging.

    Reactions Discovery

    Predict chemical (side) reactions from nothing but constituent molecules.

    AMS Driver
    Properties

    Calculate frequencies, phonons, and more. Use forces and energies from AMS or external engines.

    PES Exploration

    Minimize structures, find transitions states, scan multiple coordinates.

    Molecular Dynamics

    Use advanced thermo- and barostats, non-equilibrium and accelerated MD, molecule gun.

    Monte Carlo

    Grand Canonical Monte Carlo to study absorption, (dis)charge processes.

    Interfaces
    ParAMS

    Versatile graphical and python scripting tools to create training sets and parametrize DFTB, ReaxFF, and machine learned potentials.

    PLAMS

    Versatile python scripting interface to create your own computational chemistry workflows

    GUI

    Powerful graphical interface to set up, run, and analyze calculations. Even across different platforms.

    VASP

    Interface to popular plane-wave code VASP. Easily set up PES Scans to create training data.

    The SCM team wants to make computational chemistry work for you!

    Check out the tutorials Questions? Contact us!
  • Docs & Support
    Docs & Support
    Downloads
    Windows

    ams2025.102

    Mac

    ams2025.102

    Linux

    ams2025.102

    See all
    Documentation
    Overview

    Documentation links for all our modules and tools

    Tutorials

    Get started quickly with our Tutorials!

    Installation Manual

    Quick-start guide and extensive installation manual

    Brochures

    Brochure and flyers for different applications

    Other Resources
    Changelog

    Latest changes to our binaries

    Webinars

    Workshops

    Knowledgebank

    Research highlights

    FAQ

    General FAQs on licensing.

    Pricing and licensing

    Price and licensing information.

    Consulting & Support
    How can we help?
  • Company
    • About us
    • Team
    • Projects & Collaborations
    • Partners & Contributors
    • News
    • Events
    • Careers
    • Pricing & Licensing
    • Contact us
Search
  • Free trial
  • Contact

Home > Documentation

Navigate to
  • Documentation
  • Tutorials
  • Installation

Table of contents

  • 1. Introduction
    • 1.1. What is PLAMS
    • 1.2. What can be done with PLAMS
    • 1.3. Simple example
    • 1.4. What PLAMS is not
    • 1.5. About this documentation
  • 2. Getting started
    • 2.1. Library contents
    • 2.2. Installing PLAMS
    • 2.3. Updating PLAMS
    • 2.4. Running PLAMS
    • 2.5. Defaults file
    • 2.6. The launch script
      • 2.6.1. Working folder location
      • 2.6.2. Passing variables
      • 2.6.3. Importing past jobs
      • 2.6.4. Restarting failed script
      • 2.6.5. Multiple input scripts
  • 3. Components overview
    • 3.1. Settings
      • 3.1.1. Tree-like structure
      • 3.1.2. Dot notation
      • 3.1.3. Case sensitivity
      • 3.1.4. Global settings
      • 3.1.5. API
    • 3.2. Jobs
      • 3.2.1. Preparing a job
        • 3.2.1.1. Contents of job settings
        • 3.2.1.2. Default settings
      • 3.2.2. Running a job
        • 3.2.2.1. Name conflicts
        • 3.2.2.2. Prerun and postrun methods
        • 3.2.2.3. Preview mode
      • 3.2.3. Job API
      • 3.2.4. Single jobs
        • 3.2.4.1. Subclassing SingleJob
      • 3.2.5. Multijobs
        • 3.2.5.1. Using MultiJob
    • 3.3. Results
      • 3.3.1. Files in the job folder
      • 3.3.2. Synchronization of parallel job executions
        • 3.3.2.1. Examples
      • 3.3.3. Cleaning job folder
        • 3.3.3.1. Cleaning for multijobs
      • 3.3.4. API
    • 3.4. Job runners
      • 3.4.1. Local job runner
      • 3.4.2. Remote job runner
    • 3.5. Job manager
      • 3.5.1. Rerun prevention
      • 3.5.2. Pickling
      • 3.5.3. Restarting scripts
      • 3.5.4. API
    • 3.6. Public functions
      • 3.6.1. Logging
      • 3.6.2. Binding decorators
    • 3.7. Molecule
      • 3.7.1. Molecule
        • 3.7.1.1. Atom labeling
      • 3.7.2. Atom
      • 3.7.3. Bond
      • 3.7.4. RDKit interface
      • 3.7.5. ASE interface
    • 3.8. Utilities
      • 3.8.1. Periodic Table
      • 3.8.2. Units
      • 3.8.3. Geometry tools
  • 4. Interfaces
    • 4.1. Amsterdam Modeling Suite
      • 4.1.1. AMS driver
        • 4.1.1.1. Preparing input
        • 4.1.1.2. Preparing runscript
        • 4.1.1.3. Molecule handling
        • 4.1.1.4. AMSJob API
        • 4.1.1.5. AMSResults API
      • 4.1.2. ADF
        • 4.1.2.1. Preparing input
        • 4.1.2.2. Preparing runscript
        • 4.1.2.3. API
      • 4.1.3. ReaxFF
      • 4.1.4. Analysis tools: Densf, FCF
      • 4.1.5. KF files
    • 4.2. Other programs
      • 4.2.1. Dirac
        • 4.2.1.1. Preparing a calculation
        • 4.2.1.2. Results extraction
        • 4.2.1.3. API
      • 4.2.2. Crystal
        • 4.2.2.1. Preparing a calculation
        • 4.2.2.2. Molecule parsing
        • 4.2.2.3. Results extraction
        • 4.2.2.4. Example
        • 4.2.2.5. API
      • 4.2.3. DFTB+
        • 4.2.3.1. Preparing a calculation
        • 4.2.3.2. Results extraction
        • 4.2.3.3. Example
        • 4.2.3.4. API
      • 4.2.4. CP2K
        • 4.2.4.1. Settings
        • 4.2.4.2. Molecule parsing
        • 4.2.4.3. Loading jobs
        • 4.2.4.4. Molecule loading
        • 4.2.4.5. API
      • 4.2.5. ORCA
        • 4.2.5.1. API
      • 4.2.6. MOPAC (standalone program)
        • 4.2.6.1. Preparing input
        • 4.2.6.2. API
  • 5. Examples
    • 5.1. Simple examples
      • 5.1.1. Charge transfer integrals with ADF
      • 5.1.2. Dictionary MultiJob
    • 5.2. Advanced examples
      • 5.2.1. Tuning the range separation
    • 5.3. Recipes
      • 5.3.1. ADF fragment job
      • 5.3.2. NBO with ADF
      • 5.3.3. Numerical gradients
      • 5.3.4. Numerical Hessian
      • 5.3.5. ReaxFF molecule gun
      • 5.3.6. Global Minimum Search
      • 5.3.7. Vibrational analysis with ASE
        • 5.3.7.1. Example
        • 5.3.7.2. API
PLAMS
  • Documentation/
  • PLAMS/
  • 5. Examples/
  • 5.3.6. Global Minimum Search

5.3.6. Global Minimum Search¶

(contributed by Bas van Beek)

This module implements a scheme for finding/approaching the conformational global minimum of a Molecule. The script accomplishes this by systematically varying all dihedral angles, going through the following steps in the process:

  1. A list of bonds is created, filtering out all bonds are either terminal, part of a ring or have a bond order larger than 1.
  2. The geometry of the molecule is optimized with the first dihedral angle set to 120, 0 and -120 degree; the lowest energy conformer is returned.
  3. This process is repeated in an incremental manner for all valid dihedral angles found in step 1., each step starting from the lowest energy conformer of the previous step.
  4. After all dihedral angles have been exhausted, the final geometry is returned.

Optimizations are possible at various levels of theory, RDKit UFF being a cheap default. Alternatively, the geometry can be optimized at an arbitrary level of theory with the help of the PLAMS JobRunner. Besides the input molecule an additional two arguments, at minimum, are required: A type object of a class derived from Job and a dictionary of all keyword arguments that should be passed to aforementioned job (e.g. the job Settings). See below for an exampling using ADFJob (DFT level):

s = Settings()
s.input.Geometry
s.input.basis.type = 'DZP'
s.input.XC.GGA = 'PBE'
s.input.NumericalQuality = 'Good'

mol_in = Molecule('my_mol.xyz')
job_kwarg = {'settings': s, 'name': 'my_first_DFTJob'}
mol_out = global_minimum(mol_in, job_type=ADFJob, **job_kwarg)

Depending on the Job class, it may be necessary to manually bind the get_energy() and get_main_molecule() functions to the jobs matching Results class, these two functions being used for reading energies and geometries, respectively. See below for an exampling using AMSJob (DFTB level):

s = Settings()
s.input.ams.Task = 'GeometryOptimization'
s.input.DFTB.Model = 'DFTB3'
s.input.DFTB.ResourcesDir = 'DFTB.org/3ob-3-1'

mol_in = Molecule('my_mol.xyz')
job_kwarg = {'settings': s, 'name': 'my_first_AMSJob'}
mol_out = global_minimum(mol_in, job_type=AMSJob, **job_kwarg)

Lastly, by tweaking the job settings including but not excluded to: single points, constrained geometry optimizations, linear transits or a transition state search. The only requirement is that the job yields both an energy and a geometry which can be read with the get_energy() and get_main_molecule() functions, respectively. See below for an exampling using ADFJob (PBE level). In this example the optimizer searches for a transition state along a reaction coordinate defined by atoms 1 & 2, using the TSRC method, while simultaneously varying all dihedral angles:

s = Settings()
s.input.Geometry = 'TransitionState'
s.input.TSRC.Dist = '1 2 1.0’
s.input.basis.type = 'DZP'
s.input.XC.GGA = 'PBE'
s.input.NumericalQuality = 'Good'

mol_in = Molecule('my_mol.xyz')
job_kwarg = {'settings': s, 'name': 'my_second_ADFJob'}
mol_out = global_minimum(mol_in, job_type=ADFJob, **job_kwarg)

The source code:


__all__ = ['global_minimum']

import sys

try:
    from rdkit.Chem import AllChem, rdForceFieldHelpers
    from ..interfaces.molecule.rdkit import to_rdmol, from_rdmol
except ImportError:
    pass

from ..core.functions import init, finish


def global_minimum(mol, n_scans=1, no_h=True, no_ring=True, bond_orders=[1.0], job_type=False, path='.', **kwarg):
    """
    Find the global minimum of the ligand (RDKit UFF or user-defined PLAMS |Job|) by systematically varying dihedral angles within the molecule.

    :param |Molecule| mol: The input molecule
    :param int n_scans:: The number of times the global minimum search should be repeated
    :param bool no_h: If hydrogen-containing bonds should ignored
    :param bool no_ring: If bonds in ring systems should ignored
    :param list bond_orders: A list of accepted bond orders (floats); a bond will be ignored if its bond order is not in *bond_orders*
    :param type or bool job_type: A type object of a class derived from |Job|. Set to ''False'' to use RDKit UFF
    :param str path: The path where the PLAMS working directory will be stored
    :param dict kwarg: Keyword arguments for *job_type*
    :return |Molecule|: A copy of *mol* with a newly optimized geometry
    """
    # Creates guess bonds if no bonds are present
    if len(mol.bonds) == 0:
        mol.guess_bonds()

    # Create a list of 2-tuples (i.e. atomic indices) representing (valid) bonds within the molecule
    bond_list = find_bond(mol, no_h=no_h, no_ring=no_ring, bond_orders=bond_orders)

    # Search for the global minimum with RDKit UFF or with PLAMS at an user-defined level of theory
    if not job_type:
        if 'rdkit.Chem.AllChem' not in sys.modules or 'rdkit.Chem.rdForceFieldHelpers' not in sys.modules:
            raise ImportError('rdkit.chem module not found, aborting RDKit UFF optimization')
        if not no_ring:
            raise TypeError('no_ring=False is not supported in combination with RDKit UFF')
        if not rdForceFieldHelpers.UFFHasAllMoleculeParams(to_rdmol(mol)):
            raise ValueError('No UFF parameters found for one or more atoms')

        for i in range(n_scans):
            for bond in bond_list:
                mol = global_minimum_scan_rdkit(mol, bond)

        # Optimize the molecule even if no valid bonds are found
        if not bond_list:
            rdmol = to_rdmol(mol)
            AllChem.UFFGetMoleculeForceField(rdmol).Minimize()
            mol = from_rdmol(rdmol)

    else:
        init(path=path)
        for i in range(n_scans):
            for bond in bond_list:
                mol = global_minimum_scan_plams(mol, bond, job_type, **kwarg)

        # Optimize the molecule even if no valid bonds are found
        if not bond_list:
            job = job_type(molecule=mol, **kwarg)
            results = job.run()
            mol = results.get_main_molecule()
        finish()

    return mol


def find_bond(mol, no_h=True, no_ring=True, bond_orders=[1.0]):
    """
    Create a list of bonds. Each entry is a tuple with indices of atoms forming a dihedral.
    Consider only diherdals with axis being a single bond, so that rotation is possible.

    :param |Molecule| mol: The input molecule
    :param bool no_h: If hydrogen-containing bonds should ignored
    :param bool no_ring: If bonds in ring systems should ignored
    :param list bond_orders: A list of accepted bond orders (floats); A bond will be ignored if its bond order is not in *bond_orders*
    :return list: A list of 2-tuples containing the atomic indices of valid bonds
    """
    mol.set_atoms_id()

    # Mark atoms that can form an "axis" of a diherdal, i.e atoms with more than one (non-hydrogen) neighbor
    for atom in mol:
        neighbors = mol.neighbors(atom)
        if no_h:
            neighbors = [at for at in mol.neighbors(atom) if at.atnum != 1]
        if no_ring:
            neighbors = [mol.in_ring(at) for at in neighbors]
        atom.mark = (len(neighbors) > 1)

    # For each bond with both ends marked add one bond to the list
    ret = []
    for i, bond in enumerate(mol.bonds):
        if bond.atom1.mark and bond.atom2.mark and bond.order in bond_orders:
            at1, at2 = bond.atom1, bond.atom2
            ret.append((at1.id-1 + 1, at2.id-1 + 1))

    # Clean up the molecule
    mol.unset_atoms_id()
    for atom in mol:
        del atom.mark

    return ret


def global_minimum_scan_plams(mol, bond_tuple, job_type, **kwarg):
    """
    Optimize the molecule (A PLAMS |Job|) with 3 different values for the given dihedral angle and find the lowest energy conformer.
    The matching PLAMS |Results| object must have access to the |get_energy()| and |get_main_molecule()| functions.
    If required, functions can be added manually to a class with the |add_to_class()| function.

    :param |Molecule| mol: The input molecule
    :param tuple bond_tuple: A 2-tuple containing the atomic indices of valid bonds
    :param type job_type: A type object of a class derived from |Job|
    :param dict kwarg: Keyword arguments for *job_type*
    :return |Molecule|: A copy of *mol* with a newly optimized geometry
    """
    # Define a number of variables and create 3 copies of the ligand
    angles = (-120, 0, 120)
    mol_list = [mol.copy() for i in range(3)]
    for angle, mol in zip(angles, mol_list):
        bond = mol[bond_tuple]
        atom = mol[bond_tuple[0]]
        mol.rotate_bond(bond, atom, angle, unit='degree')

    # Optimize the geometry for all dihedral angles in angle_list
    # The geometry that yields the minimum energy is returned
    energy_list = []
    for mol in mol_list:
        job = job_type(**kwarg)
        job.molecule = mol
        results = job.run()
        energy_list.append(results.get_energy())
        mol_new = results.get_main_molecule()
        for at, at_new in zip(mol, mol_new):
            at.coords = at_new.coords
    minimum = energy_list.index(min(energy_list))
    return mol_list[minimum]


def global_minimum_scan_rdkit(mol, bond_tuple):
    """
    Optimize the molecule (RDKit UFF) with 3 different values for the given dihedral angle and find the lowest energy conformer.

    :param |Molecule| mol: The input molecule
    :param tuple bond_tuple: A 2-tuples containing the atomic indices of valid bonds
    :return |Molecule|: A copy of *mol* with a newly optimized geometry
    """
    # Define a number of variables and create 3 copies of the ligand
    uff = AllChem.UFFGetMoleculeForceField
    angles = (-120, 0, 120)
    mol_list = [mol.copy() for i in range(3)]
    for angle, mol in zip(angles, mol_list):
        bond = mol[bond_tuple]
        atom = mol[bond_tuple[0]]
        mol.rotate_bond(bond, atom, angle, unit='degree')

    # Optimize the geometry for all dihedral angles in angle_list
    # The geometry that yields the minimum energy is returned
    mol_list = [to_rdmol(mol, properties=False) for mol in mol_list]
    for rdmol in mol_list:
        uff(rdmol).Minimize()
    energy_list = [uff(rdmol).CalcEnergy() for rdmol in mol_list]
    minimum = energy_list.index(min(energy_list))
    return from_rdmol(mol_list[minimum])
Next Previous
AMS Modules
Electronic Structure
ADF: molecular DFT Periodic DFT DFTB & MOPAC
Interatomic Potentials
ReaxFF ML Potentials Force Fields
Kinetics
kMC and Microkinetics Bumblebee: OLEDs
Macroscale
COSMO-RS
Application Areas
Research Topics
Batteries Biotechnology Bonding Analysis Catalysis Heavy Elements Inorganic Chemistry Materials Science Nanoscience Oil & Gas OLEDs Perovskites Polymers Semiconductors Spectroscopy
Where to use AMS?
Industry Government Lab National Supercomputer Academic Research Teaching
Tools
Workflows
Conformers OLED workflows Reaction analysis Reaction discovery
AMS Driver
Hybrid Engine Molecular Dynamics Monte Carlo PES Exploration Properties
Python Utilities
ACErxn ParAMS PLAMS pyZacros
Interfaces
GUI VASP Parametrization
Documentation & Support
Downloads Documentation Videos Release notes Changelog Previous releases Webinars Workshops AMS Literature Brochure
Company
About us Careers Contact us Events News Our team Partners & Contributors Projects & Collaborations
Pricing & Licensing
Get a price quote Pricing structure Ordering License terms Resellers FAQ

© SCM – Software Chemistry & Materials 2025 – All Rights Reserved

Privacy Disclaimer Cookies Terms of Use Copyright
Manage Cookie Consent
We use cookies to optimise site functionality and give you the best possible experience.
See our cookie statement for all information.
Functional cookies Always active
The technical storage or access is strictly necessary for the legitimate purpose of enabling the use of a specific service explicitly requested by the subscriber or user, or for the sole purpose of carrying out the transmission of a communication over an electronic communications network.
Preferences
The technical storage or access is necessary for the legitimate purpose of storing preferences that are not requested by the subscriber or user.
Statistics
The technical storage or access that is used exclusively for statistical purposes. The technical storage or access that is used exclusively for anonymous statistical purposes. Without a subpoena, voluntary compliance on the part of your Internet Service Provider, or additional records from a third party, information stored or retrieved for this purpose alone cannot usually be used to identify you.
Marketing
The technical storage or access is required to create user profiles to send advertising, or to track the user on a website or across several websites for similar marketing purposes.
Manage options Manage services Manage {vendor_count} vendors Read more about these purposes
View preferences
{title} {title} {title}