EmbedMolecules#
- nvmolkit.embedMolecules.EmbedMolecules(
- molecules: list[Mol],
- params: EmbedParameters,
- confsPerMolecule: int = 1,
- maxIterations: int = -1,
- hardwareOptions: HardwareOptions | None = None,
Embed multiple molecules with multiple conformers on GPUs.
This function performs GPU-accelerated ETKDG conformer generation on multiple molecules. It uses CUDA for GPU acceleration and OpenMP for CPU parallelization to achieve high performance embedding of large molecule sets.
nvMolKit implements a subset of features specified in the EmbedParameters class. The following features are restricted:
useRandomCoords must be True
Bounds matrices are not supported (setBoundsMat)
Custom Coulomb potentials are not supported (SetCPCI)
Coordinate constraints are not supported (SetCoordMap)
embedFragmentsSeparately is not supported. All fragments will be embedded together.
- Parameters:
molecules – List of RDKit molecules to embed. Molecules should be prepared (sanitized, explicit hydrogens added if needed).
params – RDKit EmbedParameters object with embedding settings. Must have useRandomCoords=True for ETKDG.
confsPerMolecule – Number of conformers to generate per molecule (default: 1)
maxIterations – Maximum ETKDG iterations, -1 for automatic calculation (default: -1)
hardwareOptions – HardwareOptions with hardware settings. If None, uses defaults.
- Returns:
None. Input molecules are modified in-place with generated conformers.
- Raises:
ValueError – If any molecule in the input list is invalid, or if hardware configuration parameters are invalid
RuntimeError – If CUDA operations fail or embedding encounters errors
Example
>>> from rdkit import Chem >>> from rdkit.Chem.rdDistGeom import ETKDGv3 >>> from nvmolkit.types import HardwareOptions >>> import nvmolkit.embedMolecules as embed >>> >>> # Load molecules >>> mol1 = Chem.AddHs(Chem.MolFromSmiles('CCO')) >>> mol2 = Chem.AddHs(Chem.MolFromSmiles('CCC')) >>> >>> # Set up embedding parameters >>> params = ETKDGv3() >>> params.useRandomCoords = True # Required for nvMolKit ETKDG >>> >>> # Configure hardware options >>> hardware_opts = HardwareOptions( ... preprocessingThreads=8, ... batchSize=500, ... batchesPerGpu=4, ... gpuIds=[0, 1], ... ) >>> embed.EmbedMolecules([mol1, mol2], params, confsPerMolecule=5, hardwareOptions=hardware_opts) >>> >>> # Check conformers were generated >>> mol1.GetNumConformers() # Should be 5 >>> mol2.GetNumConformers() # Should be 5
Note
Input molecules are modified in-place with generated conformers
params.useRandomCoords must be True for ETKDG algorithm
If gpuIds is empty, all available GPUs (0 to N-1) will be used automatically