![]() |
Chemical Data Processing Library Python API - Version 1.4.0
|
Contains classes and functions related to molecular force fields. More...
Classes | |
| class | AtomProperty |
| Provides keys for built-in Chem.Atom properties. More... | |
| class | BondProperty |
| Provides keys for built-in Chem.Bond properties. More... | |
| class | ElasticPotential |
| Stores parameters for an elastic potential (harmonic distance restraint) between a pair of atoms. More... | |
| class | ElasticPotentialList |
| Data structure for the storage of elastic potential parameter set records. More... | |
| class | Error |
| Base class of the ForceField subsystem exception hierarchy. More... | |
| class | InteractionFilterFunction3 |
| Generic wrapper for storing user-defined three-atom interaction filtering functions (see [FUNWRP]). More... | |
| class | InteractionFilterFunction4 |
| Generic wrapper for storing user-defined four-atom interaction filtering functions (see [FUNWRP]). More... | |
| class | InteractionType |
| Provides flags for the specification of a set of force field interaction types. More... | |
| class | MMFF94AngleBendingInteraction |
| Stores parameters for a single MMFF94 angle-bending interaction defined over an atom triplet. More... | |
| class | MMFF94AngleBendingInteractionList |
| Data structure for the storage of MMFF94 angle-bending interaction parameter set records. More... | |
| class | MMFF94AngleBendingInteractionParameterizer |
| Detects and parameterizes the MMFF94 angle-bending interactions of a molecular graph. More... | |
| class | MMFF94AngleBendingParameterTable |
| Data structure for the storage and lookup of MMFF94 angle-bending interaction parameters. More... | |
| class | MMFF94AromaticAtomTypeDefinitionTable |
| Data structure for the storage and lookup of MMFF94 aromatic atom type definitions. More... | |
| class | MMFF94AromaticSSSRSubset |
| Implements the extraction of all rings in the SSSR of a molecular graph that are considered aromatic according to the MMFF94 aromaticity model. More... | |
| class | MMFF94AtomChargeFunction |
| Generic wrapper class used to store a user-defined MMFF94 partial atom charge retrieval function. More... | |
| class | MMFF94AtomTypePropertyTable |
| Data structure for the storage and lookup of structural and chemical property data associated with numeric MMFF94 atom types. More... | |
| class | MMFF94AtomTyper |
| Assigns MMFF94 symbolic and numeric atom types to the atoms of a molecular graph. More... | |
| class | MMFF94BondChargeIncrementTable |
| Data structure for the storage and lookup of MMFF94 bond charge increments. More... | |
| class | MMFF94BondStretchingInteraction |
| Stores parameters for a single MMFF94 bond-stretching interaction between two bonded atoms. More... | |
| class | MMFF94BondStretchingInteractionList |
| Data structure for the storage of MMFF94 bond-stretching interaction parameter set records. More... | |
| class | MMFF94BondStretchingInteractionParameterizer |
| Detects and parameterizes the MMFF94 bond-stretching interactions of a molecular graph. More... | |
| class | MMFF94BondStretchingParameterTable |
| Data structure for the storage and lookup of MMFF94 bond-stretching interaction parameters. More... | |
| class | MMFF94BondStretchingRuleParameterTable |
| Data structure for the storage and lookup of MMFF94 bond-stretching interaction fallback parameters. More... | |
| class | MMFF94BondTypeIndexFunction |
| Generic wrapper class used to store a user-defined MMFF94 bond type index retrieval function. More... | |
| class | MMFF94BondTyper |
| Assigns MMFF94 bond type indices to the bonds of a molecular graph. More... | |
| class | MMFF94ChargeCalculator |
| Calculator for the MMFF94 partial atomic charges of a molecular graph. More... | |
| class | MMFF94DefaultStretchBendParameterTable |
| Data structure for the storage and lookup of MMFF94 stretch-bend interaction fallback parameters. More... | |
| class | MMFF94ElectrostaticInteraction |
| Stores parameters for a single MMFF94 electrostatic interaction between two non-bonded atoms. More... | |
| class | MMFF94ElectrostaticInteractionList |
| Data structure for the storage of MMFF94 electrostatic interaction parameter set records. More... | |
| class | MMFF94ElectrostaticInteractionParameterizer |
| Detects and parameterizes the MMFF94 electrostatic interactions of a molecular graph. More... | |
| class | MMFF94EnergyCalculator |
| Calculates the total MMFF94 force field energy for a set of atom 3D coordinates. More... | |
| class | MMFF94FormalAtomChargeDefinitionTable |
| Data structure for the storage and lookup of formal charge definitions used by the MMFF94 charge model. More... | |
| class | MMFF94GradientCalculator |
| Calculates the total MMFF94 force field energy and its gradient for a set of atom 3D coordinates. More... | |
| class | MMFF94HeavyToHydrogenAtomTypeMap |
| Data structure for the storage and lookup of heavy-to-hydrogen symbolic MMFF94 atom type mappings. More... | |
| class | MMFF94InteractionData |
| Container holding the full set of MMFF94 interaction parameters for a molecular graph. More... | |
| class | MMFF94InteractionParameterizer |
| Highly configurable implementation of the complete workflow for molecular graph MMFF94 interaction perception and parameterization. More... | |
| class | MMFF94NumericAtomTypeFunction |
| Generic wrapper class used to store a user-defined numeric MMFF94 atom type retrieval function. More... | |
| class | MMFF94OutOfPlaneBendingInteraction |
| Stores parameters for a single MMFF94 out-of-plane bending interaction at a trigonal center. More... | |
| class | MMFF94OutOfPlaneBendingInteractionList |
| Data structure for the storage of MMFF94 out-of-plane bending interaction parameter set records. More... | |
| class | MMFF94OutOfPlaneBendingInteractionParameterizer |
| Detects and parameterizes the MMFF94 out-of-plane bending interactions of a molecular graph. More... | |
| class | MMFF94OutOfPlaneBendingParameterTable |
| Data structure for the storage and lookup of MMFF94 out-of-plane bending interaction force constants. More... | |
| class | MMFF94ParameterSet |
| Provides constants for MMFF94 parameter set variant specification. More... | |
| class | MMFF94PartialBondChargeIncrementTable |
| Data structure for the storage and lookup of MMFF94 per-atom partial bond charge increment and formal charge adjustment factors. More... | |
| class | MMFF94PrimaryToParameterAtomTypeMap |
| Data structure for the storage and lookup of primary numeric MMFF94 atom type to fallback parameter atom type list mappings. More... | |
| class | MMFF94RingSetFunction |
| Generic wrapper class used to store a user-defined MMFF94 ring set retrieval function. More... | |
| class | MMFF94StretchBendInteraction |
| Stores paramters for a single MMFF94 stretch-bend coupling interaction. More... | |
| class | MMFF94StretchBendInteractionList |
| Data structure for the storage of MMFF94 stretch-bend coupling interaction parameter set records. More... | |
| class | MMFF94StretchBendInteractionParameterizer |
| class | MMFF94StretchBendParameterTable |
| Data structure for the storage and lookup of MMFF94 stretch-bend interaction parameters. More... | |
| class | MMFF94SymbolicAtomTypeFunction |
| Generic wrapper class used to store a user-defined symbolic MMFF94 atom type retrieval function. More... | |
| class | MMFF94SymbolicAtomTypePatternTable |
| Data structure storing SMARTS substructure patterns used to assign symbolic MMFF94 atom types during atom typing. More... | |
| class | MMFF94SymbolicToNumericAtomTypeMap |
| Lookup table mapping each symbolic MMFF94 atom type to its corresponding numeric MMFF94 atom type. More... | |
| class | MMFF94TorsionInteraction |
| Stores parameters for a single MMFF94 torsion interaction over an atom quadruplet i-j-k-l. More... | |
| class | MMFF94TorsionInteractionList |
| Data structure for the storage of MMFF94 torsion interaction parameter set records. More... | |
| class | MMFF94TorsionInteractionParameterizer |
| Detects and parameterizes the MMFF94 torsion interactions of a molecular graph. More... | |
| class | MMFF94TorsionParameterTable |
| Data structure for the storage and lookup of MMFF94 torsion interaction parameters. More... | |
| class | MMFF94VanDerWaalsAtomParameters |
| class | MMFF94VanDerWaalsInteraction |
| Stores parameters for a single MMFF94 Van der Waals interaction between two non-bonded atoms. More... | |
| class | MMFF94VanDerWaalsInteractionList |
| Data structure for the storage of MMFF94 Van der Waals interaction parameter set records. More... | |
| class | MMFF94VanDerWaalsInteractionParameterizer |
| Detects and parameterizes the MMFF94 Van der Waals interactions of a molecular graph. More... | |
| class | MMFF94VanDerWaalsParameterTable |
| Data structure for the storage and lookup of MMFF94 interaction parameters. More... | |
| class | MolecularGraphProperty |
| Provides keys for built-in Chem.MolecularGraph properties. More... | |
| class | ParameterizationFailed |
| Thrown when force field parameterization has failed. More... | |
| class | TopologicalAtomDistanceFunction |
| Generic wrapper class used to store a user-defined topological atom-pair distance function. More... | |
| class | UFFAtomType |
| Provides constants for the specification of numeric Universal Force Field (UFF) atom types. More... | |
| class | UFFAtomTypePropertyTable |
| Data structure for the storage and lookup of various properties associated with numeric UFF atom types. More... | |
Functions | |
| None | setMMFF94AromaticRings (Chem.MolecularGraph molgraph, Chem.FragmentList rings) |
| Sets the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph to rings. More... | |
| bool | hasMMFF94AromaticRings (Chem.MolecularGraph molgraph) |
| Tells whether the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph is set. More... | |
| Chem.FragmentList | getMMFF94AromaticRings (Chem.MolecularGraph molgraph) |
| Returns the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph. More... | |
| None | clearMMFF94AromaticRings (Chem.MolecularGraph molgraph) |
| Clears the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph. More... | |
| None | setMMFF94Charge (Chem.Atom atom, float charge) |
| Sets the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom to charge. More... | |
| bool | hasMMFF94Charge (Chem.Atom atom) |
| Tells whether the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom is set. More... | |
| float | getMMFF94Charge (Chem.Atom atom) |
| Returns the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom. More... | |
| None | clearMMFF94Charge (Chem.Atom atom) |
| Clears the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom. More... | |
| None | setMMFF94NumericType (Chem.Atom atom, int type) |
| Sets the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom to type. More... | |
| bool | hasMMFF94NumericType (Chem.Atom atom) |
| Tells whether the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom is set. More... | |
| int | getMMFF94NumericType (Chem.Atom atom) |
| Returns the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom. More... | |
| None | clearMMFF94NumericType (Chem.Atom atom) |
| Clears the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom. More... | |
| None | setMMFF94SymbolicType (Chem.Atom atom, str type) |
| Sets the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom to type. More... | |
| bool | hasMMFF94SymbolicType (Chem.Atom atom) |
| Tells whether the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom is set. More... | |
| str | getMMFF94SymbolicType (Chem.Atom atom) |
| Returns the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom. More... | |
| None | clearMMFF94SymbolicType (Chem.Atom atom) |
| Clears the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom. More... | |
| None | setMMFF94TypeIndex (Chem.Bond bond, int type_idx) |
| Sets the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond to type_idx. More... | |
| bool | hasMMFF94TypeIndex (Chem.Bond bond) |
| Tells whether the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond is set. More... | |
| int | getMMFF94TypeIndex (Chem.Bond bond) |
| Returns the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond. More... | |
| None | clearMMFF94TypeIndex (Chem.Bond bond) |
| Clears the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond. More... | |
| None | setUFFType (Chem.Atom atom, int type) |
| Sets the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom to type. More... | |
| bool | hasUFFType (Chem.Atom atom) |
| Tells whether the ForceField.AtomProperty.UFF_TYPE property of the atom atom is set. More... | |
| int | getUFFType (Chem.Atom atom) |
| Returns the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom. More... | |
| None | clearUFFType (Chem.Atom atom) |
| Clears the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom. More... | |
| None | assignMMFF94AtomTypes (Chem.MolecularGraph molgraph, bool strict, bool overwrite) |
| Assigns MMFF94 atom types to the atoms of the molecular graph molgraph. More... | |
| None | assignMMFF94BondTypeIndices (Chem.MolecularGraph molgraph, bool strict, bool overwrite) |
| Assigns MMFF94 bond type indices to the bonds of the molecular graph molgraph. More... | |
| None | assignUFFAtomTypes (Chem.MolecularGraph molgraph, bool overwrite) |
| Assigns UFF atom types to the atoms of the molecular graph molgraph. More... | |
| float | calcBondAngleCosDerivatives (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D term_atom1_deriv, Math.Vector3D ctr_atom_deriv, Math.Vector3D term_atom2_deriv) |
| Calculates the partial derivatives \( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_x}} \) of the of the cosine of the angle \( \vartheta_{ijk} \) between the bonds i-j and j-k. More... | |
| float | calcBondAngleCos (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, float r_ij, float r_jk) |
| Calculates the cosine of the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k. More... | |
| float | calcBondAngleCos (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos) |
| Calculates the cosine of the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k. More... | |
| float | calcBondAngle (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, float r_ij, float r_jk) |
| Calculates the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k. More... | |
| float | calcBondAngle (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos) |
| Calculates the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k. More... | |
| tuple | calcBondLengthsAndAngleCos (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos) |
| tuple | calcBondLengthsAndAngle (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos) |
| float | calcDihedralAngleCosDerivatives (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom1_pos, Math.Vector3D ctr_atom2_pos, Math.Vector3D term_atom2_pos, Math.Vector3D term_atom1_deriv, Math.Vector3D ctr_atom1_deriv, Math.Vector3D ctr_atom2_deriv, Math.Vector3D term_atom2_deriv) |
| Calculates the partial derivatives \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_x}} \) of the cosine of the angle \( \Phi_{ijkl} \) between the planes defined by the atom triplets i-j-k and j-k-l. More... | |
| float | calcDihedralAngleCos (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom1_pos, Math.Vector3D ctr_atom2_pos, Math.Vector3D term_atom2_pos) |
| Calculates the cosine of the dihedral angle \( \Phi_{ijkl} \) between the planes defined by the atom triplets i-j-k and j-k-l. More... | |
| float | calcDistanceDerivatives (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, Math.Vector3D atom1_deriv, Math.Vector3D atom2_deriv) |
| Calculates the partial derivatives \( \frac{\partial r_{ij}}{\partial \vec{p_x}} \) of the distance \( r_{ij} \) between two atoms i and j. More... | |
| float | calcDistance (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos) |
| Calculates the distance \( r_{ij} \) between two atoms i and j. More... | |
| float | calcElasticPotentialEnergy (ElasticPotential pot, Math.Vector3DArray coords) |
| Calculates the energy of a single elastic potential pot for the supplied atom coordinates coords. More... | |
| float | calcElasticPotentialEnergy (ElasticPotentialList list, Math.Vector3DArray coords) |
| float | calcElasticPotentialEnergy (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, float force_const, float ref_length) |
| Calculates the energy \( E_{ij} \) of an elastic potential applied on a pair of atoms i-j. More... | |
| float | calcElasticPotentialGradient (ElasticPotential pot, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the energy and gradient contribution of a single elastic potential pot for the supplied atom coordinates coords. More... | |
| float | calcElasticPotentialGradient (ElasticPotentialList list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcElasticPotentialGradient (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, Math.Vector3D atom1_grad, Math.Vector3D atom2_grad, float force_const, float ref_length) |
| Calculates the elastic potential energy gradient \( \nabla E_{ij} \) for a pair of atoms i-j. More... | |
| float | calcMMFF94AngleBendingEnergy (MMFF94AngleBendingInteraction iaction, Math.Vector3DArray coords) |
| Calculates the angle-bending interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94AngleBendingEnergy (MMFF94AngleBendingInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94AngleBendingEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, bool linear, float force_const, float ref_angle) |
| Calculates the angle bending interaction energy \( EA_{ijk} \) for two bonds i-j and j-k. More... | |
| float | calcMMFF94AngleBendingEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, float r_ij, float r_jk, bool linear, float force_const, float ref_angle) |
| Calculates the angle bending interaction energy \( EA_{ijk} \) for two bonds i-j and j-k. More... | |
| float | calcMMFF94AngleBendingGradient (MMFF94AngleBendingInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the angle-bending interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94AngleBendingGradient (MMFF94AngleBendingInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94AngleBendingGradient (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D term_atom1_grad, Math.Vector3D ctr_atom_grad, Math.Vector3D term_atom2_grad, bool linear, float force_const, float ref_angle) |
| Calculates the angle bending interaction energy gradient \( \nabla EA_{ijk} \) for two bonds i-j and j-k. More... | |
| None | calcMMFF94AtomCharges (Chem.MolecularGraph molgraph, bool strict, bool overwrite) |
| Calculates and assigns MMFF94 partial atomic charges to the atoms of the molecular graph molgraph. More... | |
| float | calcMMFF94BondStretchingEnergy (MMFF94BondStretchingInteraction iaction, Math.Vector3DArray coords) |
| Calculates the bond-stretching interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94BondStretchingEnergy (MMFF94BondStretchingInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94BondStretchingEnergy (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, float force_const, float ref_length) |
| Calculates the bond stretching interaction energy \( EB_{ij} \) for the bond i-j. More... | |
| float | calcMMFF94BondStretchingEnergy (float r_ij, float force_const, float ref_length) |
| Calculates the bond stretching interaction energy \( EB_{ij} \) for the bond i-j. More... | |
| float | calcMMFF94BondStretchingGradient (MMFF94BondStretchingInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the bond-stretching interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94BondStretchingGradient (MMFF94BondStretchingInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94BondStretchingGradient (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, Math.Vector3D atom1_grad, Math.Vector3D atom2_grad, float force_const, float ref_length) |
| Calculates the bond stretching interaction energy gradient \( \nabla EB_{ij} \) for the bond i-j. More... | |
| float | calcMMFF94ElectrostaticEnergy (MMFF94ElectrostaticInteraction iaction, Math.Vector3DArray coords) |
| Calculates the electrostatic interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94ElectrostaticEnergy (MMFF94ElectrostaticInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94ElectrostaticEnergy (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, float atom1_chg, float atom2_chg, float scale_fact, float de_const, float dist_expo) |
| Calculates the electrostatic interaction energy \( EQ_{ij} \) for the atom pair i-j. More... | |
| float | calcMMFF94ElectrostaticEnergy (float r_ij, float atom1_chg, float atom2_chg, float scale_fact, float de_const, float dist_expo) |
| Calculates the electrostatic interaction energy \( EQ_{ij} \) for the atom pair i-j. More... | |
| float | calcMMFF94ElectrostaticGradient (MMFF94ElectrostaticInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the electrostatic interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94ElectrostaticGradient (MMFF94ElectrostaticInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94ElectrostaticGradient (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, Math.Vector3D atom1_grad, Math.Vector3D atom2_grad, float atom1_chg, float atom2_chg, float scale_fact, float de_const, float dist_expo) |
| Calculates the electrostatic interaction energy gradient \( \nabla EQ_{ij} \) for the atom pair i-j. More... | |
| float | calcMMFF94OutOfPlaneBendingEnergy (MMFF94OutOfPlaneBendingInteraction iaction, Math.Vector3DArray coords) |
| Calculates the out-of-plane bending interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94OutOfPlaneBendingEnergy (MMFF94OutOfPlaneBendingInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94OutOfPlaneBendingEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D oop_atom_pos, float force_const) |
| Calculates the out-of-plane bending interaction energy \( EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k. More... | |
| float | calcMMFF94OutOfPlaneBendingEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D oop_atom_pos, float r_jl, float force_const) |
| Calculates the out-of-plane bending interaction energy \( EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k. More... | |
| float | calcMMFF94OutOfPlaneBendingGradient (MMFF94OutOfPlaneBendingInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the out-of-plane bending interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94OutOfPlaneBendingGradient (MMFF94OutOfPlaneBendingInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94OutOfPlaneBendingGradient (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D oop_atom_pos, Math.Vector3D term_atom1_grad, Math.Vector3D ctr_atom_grad, Math.Vector3D term_atom2_grad, Math.Vector3D oop_atom_grad, float force_const) |
| Calculates the out-of-plane bending interaction energy gradient \( \nabla EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k. More... | |
| float | calcMMFF94StretchBendEnergy (MMFF94StretchBendInteraction iaction, Math.Vector3DArray coords) |
| Calculates the stretch-bend coupling interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94StretchBendEnergy (MMFF94StretchBendInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94StretchBendEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, float ijk_force_const, float kji_force_const, float ref_angle, float ref_length1, float ref_length2) |
| Calculates the stretch-bend interaction energy \( EBA_{ijk} \) for two bonds i-j and j-k. More... | |
| float | calcMMFF94StretchBendEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, float r_ij, float r_jk, float ijk_force_const, float kji_force_const, float ref_angle, float ref_length1, float ref_length2) |
| Calculates the stretch-bend interaction energy \( EBA_{ijk} \) for two bonds i-j and j-k. More... | |
| float | calcMMFF94StretchBendGradient (MMFF94StretchBendInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the stretch-bend coupling interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94StretchBendGradient (MMFF94StretchBendInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94StretchBendGradient (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D term_atom1_grad, Math.Vector3D ctr_atom_grad, Math.Vector3D term_atom2_grad, float ijk_force_const, float kji_force_const, float ref_angle, float ref_length1, float ref_length2) |
| Calculates the stretch-bend interaction energy gradient \( \nabla EBA_{ijk} \) for two bonds i-j and j-k. More... | |
| float | calcMMFF94TorsionEnergy (MMFF94TorsionInteraction iaction, Math.Vector3DArray coords) |
| Calculates the torsion interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94TorsionEnergy (MMFF94TorsionInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94TorsionEnergy (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom1_pos, Math.Vector3D ctr_atom2_pos, Math.Vector3D term_atom2_pos, float tor_param1, float tor_param2, float tor_param3) |
| Calculates the torsion interaction energy \( ET_{ijkl} \) for the central bond j-k and the connected bonds i-j and k-l. More... | |
| float | calcMMFF94TorsionGradient (MMFF94TorsionInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the torsion interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94TorsionGradient (MMFF94TorsionInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94TorsionGradient (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom1_pos, Math.Vector3D ctr_atom2_pos, Math.Vector3D term_atom2_pos, Math.Vector3D term_atom1_grad, Math.Vector3D ctr_atom1_grad, Math.Vector3D ctr_atom2_grad, Math.Vector3D term_atom2_grad, float tor_param1, float tor_param2, float tor_param3) |
| Calculates the torsion interaction energy gradient \( \nabla ET_{ijkl} \) for the central bond j-k and the connected bonds i-j and k-l. More... | |
| float | calcMMFF94VanDerWaalsEnergy (MMFF94VanDerWaalsInteraction iaction, Math.Vector3DArray coords) |
| Calculates the Van der Waals interaction energy of iaction using the geometry from coords. More... | |
| float | calcMMFF94VanDerWaalsEnergy (MMFF94VanDerWaalsInteractionList ia_list, Math.Vector3DArray coords) |
| float | calcMMFF94VanDerWaalsEnergy (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, float e_IJ, float r_IJ, float r_IJ_7) |
| Calculates the Van der Waals interaction energy \( E_{vdW_{ij}} \) for the atom pair i-j. More... | |
| float | calcMMFF94VanDerWaalsEnergy (float r_ij, float e_IJ, float r_IJ, float r_IJ_7) |
| Calculates the Van der Waals interaction energy \( E_{vdW_{ij}} \) for the atom pair i-j. More... | |
| float | calcMMFF94VanDerWaalsGradient (MMFF94VanDerWaalsInteraction iaction, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| Calculates the Van der Waals interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad. More... | |
| float | calcMMFF94VanDerWaalsGradient (MMFF94VanDerWaalsInteractionList ia_list, Math.Vector3DArray coords, Math.Vector3DArray grad) |
| float | calcMMFF94VanDerWaalsGradient (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos, Math.Vector3D atom1_grad, Math.Vector3D atom2_grad, float e_IJ, float r_IJ, float r_IJ_7) |
| Calculates the Van der Waals interaction energy gradient \( \nabla E_{vdW_{ij}} \) for the atom pair i-j. More... | |
| float | calcOutOfPlaneAngleCosDerivatives (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D oop_atom_pos, Math.Vector3D term_atom1_deriv, Math.Vector3D ctr_atom_deriv, Math.Vector3D term_atom2_deriv, Math.Vector3D oop_atom_deriv) |
| Calculates the partial derivatives \( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_x}} \) of the cosine of the angle \( \omega_{ijk;l} \) between the bond j-l and the normal of the plane defined by the atoms i-j-k. More... | |
| float | calcOutOfPlaneAngle (Math.Vector3D term_atom1_pos, Math.Vector3D ctr_atom_pos, Math.Vector3D term_atom2_pos, Math.Vector3D oop_atom_pos, float r_jl) |
| Calculates the out-of-plane angle \( \chi_{ijk;l} \) between the bond j-l and the plane defined by the atoms i-j-k. More... | |
| float | calcSquaredDistance (Math.Vector3D atom1_pos, Math.Vector3D atom2_pos) |
| Calculates the squared distance \( r_{ij}^2 \) between two atoms i and j. More... | |
| None | filterInteractions (MMFF94InteractionData ia_data, MMFF94InteractionData filtered_ia_data, Util.BitSet inc_atom_mask) |
| Filters an MMFF94 interaction data set, retaining only those interactions that exclusively reference atoms in inc_atom_mask. More... | |
| Chem.FragmentList | perceiveMMFF94AromaticRings (Chem.MolecularGraph molgraph) |
| Perceives the list of MMFF94 aromatic rings of the molecular graph molgraph. More... | |
| Chem.FragmentList | perceiveMMFF94AromaticRings (Chem.MolecularGraph molgraph, bool overwrite) |
| Perceives and (optionally) stores the list of MMFF94 aromatic rings of the molecular graph molgraph. More... | |
| int | perceiveUFFType (Chem.Atom atom, Chem.MolecularGraph molgraph) |
| Perceives the UFF atom type of the atom atom in the context of the molecular graph molgraph. More... | |
Contains classes and functions related to molecular force fields.
| None CDPL.ForceField.setMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph, |
| Chem.FragmentList | rings | ||
| ) |
Sets the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph to rings.
| molgraph | The molecular graph for which to set the property value. |
| rings | The new list of MMFF94 aromatic rings. |
| bool CDPL.ForceField.hasMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph | ) |
Tells whether the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph is set.
| molgraph | The molecular graph for which to query the property value. |
True if the property is set, and False otherwise. | Chem.FragmentList CDPL.ForceField.getMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph | ) |
Returns the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph.
| molgraph | The molecular graph for which to return the property value. |
| None CDPL.ForceField.clearMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph | ) |
Clears the value of the ForceField.MolecularGraphProperty.MMFF94_AROMATIC_RINGS property of the molecular graph molgraph.
| molgraph | The molecular graph for which to clear the property value. |
| None CDPL.ForceField.setMMFF94Charge | ( | Chem.Atom | atom, |
| float | charge | ||
| ) |
Sets the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom to charge.
| atom | The atom for which to set the property value. |
| charge | The new MMFF94 partial atomic charge. |
| bool CDPL.ForceField.hasMMFF94Charge | ( | Chem.Atom | atom | ) |
Tells whether the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom is set.
| atom | The atom for which to query the property value. |
True if the MMFF94 partial atomic charge is set, and False otherwise. | float CDPL.ForceField.getMMFF94Charge | ( | Chem.Atom | atom | ) |
Returns the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom.
| atom | The atom for which to return the property value. |
| None CDPL.ForceField.clearMMFF94Charge | ( | Chem.Atom | atom | ) |
Clears the value of the ForceField.AtomProperty.MMFF94_CHARGE property of the atom atom.
| atom | The atom for which to clear the property value. |
| None CDPL.ForceField.setMMFF94NumericType | ( | Chem.Atom | atom, |
| int | type | ||
| ) |
Sets the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom to type.
| atom | The atom for which to set the property value. |
| type | The new numeric MMFF94 atom type. |
| bool CDPL.ForceField.hasMMFF94NumericType | ( | Chem.Atom | atom | ) |
Tells whether the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom is set.
| atom | The atom for which to query the property value. |
True if the numeric MMFF94 atom type is set, and False otherwise. | int CDPL.ForceField.getMMFF94NumericType | ( | Chem.Atom | atom | ) |
Returns the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom.
| atom | The atom for which to return the property value. |
| None CDPL.ForceField.clearMMFF94NumericType | ( | Chem.Atom | atom | ) |
Clears the value of the ForceField.AtomProperty.MMFF94_NUMERIC_TYPE property of the atom atom.
| atom | The atom for which to clear the property value. |
| None CDPL.ForceField.setMMFF94SymbolicType | ( | Chem.Atom | atom, |
| str | type | ||
| ) |
Sets the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom to type.
| atom | The atom for which to set the property value. |
| type | The new symbolic MMFF94 atom type. |
| bool CDPL.ForceField.hasMMFF94SymbolicType | ( | Chem.Atom | atom | ) |
Tells whether the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom is set.
| atom | The atom for which to query the property value. |
True if the symbolic MMFF94 atom type is set, and False otherwise. | str CDPL.ForceField.getMMFF94SymbolicType | ( | Chem.Atom | atom | ) |
Returns the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom.
| atom | The atom for which to return the property value. |
| None CDPL.ForceField.clearMMFF94SymbolicType | ( | Chem.Atom | atom | ) |
Clears the value of the ForceField.AtomProperty.MMFF94_SYMBOLIC_TYPE property of the atom atom.
| atom | The atom for which to clear the property value. |
| None CDPL.ForceField.setMMFF94TypeIndex | ( | Chem.Bond | bond, |
| int | type_idx | ||
| ) |
Sets the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond to type_idx.
| bond | The bond for which to set the property value. |
| type_idx | The new MMFF94 bond type index. |
| bool CDPL.ForceField.hasMMFF94TypeIndex | ( | Chem.Bond | bond | ) |
Tells whether the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond is set.
| bond | The bond for which to query the property value. |
True if the MMFF94 bond type index is set, and False otherwise. | int CDPL.ForceField.getMMFF94TypeIndex | ( | Chem.Bond | bond | ) |
Returns the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond.
| bond | The bond for which to return the property value. |
| None CDPL.ForceField.clearMMFF94TypeIndex | ( | Chem.Bond | bond | ) |
Clears the value of the ForceField.BondProperty.MMFF94_TYPE_INDEX property of the bond bond.
| bond | The bond for which to clear the property value. |
| None CDPL.ForceField.setUFFType | ( | Chem.Atom | atom, |
| int | type | ||
| ) |
Sets the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom to type.
| atom | The atom for which to set the property value. |
| type | The new numeric UFF atom type. |
| bool CDPL.ForceField.hasUFFType | ( | Chem.Atom | atom | ) |
Tells whether the ForceField.AtomProperty.UFF_TYPE property of the atom atom is set.
| atom | The atom for which to query the property value. |
True if the numeric UFF atom type is set, and False otherwise. | int CDPL.ForceField.getUFFType | ( | Chem.Atom | atom | ) |
Returns the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom.
| atom | The atom for which to return the property value. |
| None CDPL.ForceField.clearUFFType | ( | Chem.Atom | atom | ) |
Clears the value of the ForceField.AtomProperty.UFF_TYPE property of the atom atom.
| atom | The atom for which to clear the property value. |
| None CDPL.ForceField.assignMMFF94AtomTypes | ( | Chem.MolecularGraph | molgraph, |
| bool | strict, | ||
| bool | overwrite | ||
| ) |
Assigns MMFF94 atom types to the atoms of the molecular graph molgraph.
| molgraph | The molecular graph whose atoms will be typed. |
| strict | If True, atoms for which no MMFF94 type could be perceived cause an error to be reported. Otherwise, perception never fails and problematic atoms are assigned a suitable fallback type. |
| overwrite | Specifies whether already assigned atom type properties should be replaced. |
| None CDPL.ForceField.assignMMFF94BondTypeIndices | ( | Chem.MolecularGraph | molgraph, |
| bool | strict, | ||
| bool | overwrite | ||
| ) |
Assigns MMFF94 bond type indices to the bonds of the molecular graph molgraph.
| molgraph | The molecular graph whose bonds will be typed. |
| strict | If True, strict parameterization will be peformed that might fail. Otherwise, bonds with parameterization problems receive the type index 0. |
| overwrite | Specifies whether already assigned bond type index properties should be replaced. |
| None CDPL.ForceField.assignUFFAtomTypes | ( | Chem.MolecularGraph | molgraph, |
| bool | overwrite | ||
| ) |
Assigns UFF atom types to the atoms of the molecular graph molgraph.
| molgraph | The molecular graph whose atoms will be typed. |
| overwrite | Specifies whether already assigned atom type properties should be replaced. |
| float CDPL.ForceField.calcBondAngleCosDerivatives | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | term_atom1_deriv, | ||
| Math.Vector3D | ctr_atom_deriv, | ||
| Math.Vector3D | term_atom2_deriv | ||
| ) |
Calculates the partial derivatives \( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_x}} \) of the of the cosine of the angle \( \vartheta_{ijk} \) between the bonds i-j and j-k.
\( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_i}} = \frac{\vec{v_{jk}}}{r_{ji} \: r_{jk}} - \frac{\vec{v_{ji}} \: (\vec{v_{ji}} \cdot \vec{v_{jk}})}{r_{ji}^3 \: r_{jk}} \)
\( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_k}} = \frac{\vec{v_{ji}}}{r_{ji} \: r_{jk}} - \frac{\vec{v_{jk}} \: (\vec{v_{ji}} \cdot \vec{v_{jk}})}{r_{ji} \: r_{jk}^3} \)
\( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_j}} = -(\frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_i}} + \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_k}}) \)
\( \cos(\vartheta_{ijk}) = \frac{\vec{v_{ij}} \cdot \vec{v_{jk}}}{r_{ij} \: r_{jk}} \)
where
\( \vec{v_{ji}} = \vec{p_i} - \vec{p_j} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( r_{ji} = |\vec{v_{ji}}| \)
\( r_{jk} = |\vec{v_{jk}}| \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
| term_atom1_pos | The position \( \vec{p_i} \) of atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of atom k. |
| term_atom1_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_i}} \) at the given atom positions. |
| ctr_atom_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_j}} \) at the given atom positions. |
| term_atom2_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\vartheta_{ijk})}{\partial \vec{p_k}} \) at the given atom positions. |
| float CDPL.ForceField.calcBondAngleCos | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | r_ij, | ||
| float | r_jk | ||
| ) |
Calculates the cosine of the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k.
\( \cos(\vartheta_{ijk}) = \frac{\vec{v_{ij}} \cdot \vec{v_{jk}}}{|\vec{v_{ij}}| \: |\vec{v_{jk}}|} \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| r_ij | The length of the bond between atom i and j. |
| r_jk | The length of the bond between atom j and k. |
| float CDPL.ForceField.calcBondAngleCos | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos | ||
| ) |
Calculates the cosine of the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k.
\( \cos(\vartheta_{ijk}) = \frac{\vec{v_{ij}} \cdot \vec{v_{jk}}}{|\vec{v_{ij}}| \: |\vec{v_{jk}}|} \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| float CDPL.ForceField.calcBondAngle | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | r_ij, | ||
| float | r_jk | ||
| ) |
Calculates the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k.
\( \vartheta_{ijk} = \arccos(\frac{\vec{v_{ij}} \cdot \vec{v_{jk}}}{|\vec{v_{ij}}| \: |\vec{v_{jk}}|}) \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| r_ij | The length of the bond between atom i and j. |
| r_jk | The length of the bond between atom j and k. |
| float CDPL.ForceField.calcBondAngle | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos | ||
| ) |
Calculates the bond angle \( \vartheta_{ijk} \) between the two bonds i-j and j-k.
\( \vartheta_{ijk} = \arccos(\frac{\vec{v_{ij}} \cdot \vec{v_{jk}}}{|\vec{v_{ij}}| \: |\vec{v_{jk}}|}) \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| tuple CDPL.ForceField.calcBondLengthsAndAngleCos | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos | ||
| ) |
| term_atom1_pos | |
| ctr_atom_pos | |
| term_atom2_pos |
| tuple CDPL.ForceField.calcBondLengthsAndAngle | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos | ||
| ) |
| term_atom1_pos | |
| ctr_atom_pos | |
| term_atom2_pos |
| float CDPL.ForceField.calcDihedralAngleCosDerivatives | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom1_pos, | ||
| Math.Vector3D | ctr_atom2_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | term_atom1_deriv, | ||
| Math.Vector3D | ctr_atom1_deriv, | ||
| Math.Vector3D | ctr_atom2_deriv, | ||
| Math.Vector3D | term_atom2_deriv | ||
| ) |
Calculates the partial derivatives \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_x}} \) of the cosine of the angle \( \Phi_{ijkl} \) between the planes defined by the atom triplets i-j-k and j-k-l.
\( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_i}} = \vec{v_{jk}} \times \vec{a} \)
\( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_j}} = \vec{r_{ki}} \times \vec{a} - \vec{v_{lk}} \times \vec{b} \)
\( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_l}} = \vec{v_{jk}} \times \vec{b} \)
\( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_k}} = -(\frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_i}} + \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_j}} + \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_l}}) \)
\( \cos(\Phi_{ijkl}) = \frac{\vec{n_{ijk}} \cdot \vec{n_{jkl}}}{|\vec{n_{ijk}}| \: |\vec{n_{jkl}}|} \)
where
\( \vec{v_{ji}} = \vec{p_i} - \vec{p_j} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{v_{lk}} = \vec{p_k} - \vec{p_l} \)
\( \vec{r_{ki}} = \vec{p_i} - \vec{p_k} \)
\( \vec{n_{ijk}} = \vec{v_{ji}} \times \vec{v_{jk}} \)
\( \vec{n_{jkl}} = \vec{v_{jk}} \times \vec{v_{lk}} \)
\( \vec{a} = \frac{\frac{\vec{n_{jkl}}}{|\vec{n_{jkl}}|} - \cos(\Phi_{ijkl}) \: \frac{\vec{n_{ijk}}}{|\vec{n_{ijk}}|}}{|\vec{n_{ijk}}|} \)
\( \vec{b} = \frac{\frac{\vec{n_{ijk}}}{|\vec{n_{ijk}}|} - \cos(\Phi_{ijkl}) \: \frac{\vec{n_{jkl}}}{|\vec{n_{jkl}}|}}{|\vec{n_{jkl}}|} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
\( \vec{p_l} \) = coordinates of atom l.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom1_pos | The position \( \vec{p_j} \) of the central atom j. |
| ctr_atom2_pos | The position \( \vec{p_k} \) of the central atom k. |
| term_atom2_pos | The position \( \vec{p_l} \) of the terminal atom l. |
| term_atom1_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_i}} \) at the given atom positions. |
| ctr_atom1_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_j}} \) at the given atom positions. |
| ctr_atom2_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_k}} \) at the given atom positions. |
| term_atom2_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_l}} \) at the given atom positions. |
| float CDPL.ForceField.calcDihedralAngleCos | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom1_pos, | ||
| Math.Vector3D | ctr_atom2_pos, | ||
| Math.Vector3D | term_atom2_pos | ||
| ) |
Calculates the cosine of the dihedral angle \( \Phi_{ijkl} \) between the planes defined by the atom triplets i-j-k and j-k-l.
\( \cos(\Phi_{ijkl}) = \frac{\vec{n_{ijk}} \cdot \vec{n_{jkl}}}{|\vec{n_{ijk}}| \: |\vec{n_{jkl}}|} \)
where
\( \vec{v_{ji}} = \vec{p_i} - \vec{p_j} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{v_{lk}} = \vec{p_k} - \vec{p_l} \)
\( \vec{n_{ijk}} = \vec{v_{ji}} \times \vec{v_{jk}} \)
\( \vec{n_{jkl}} = \vec{v_{jk}} \times \vec{v_{lk}} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
\( \vec{p_l} \) = coordinates of atom l.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom1_pos | The position \( \vec{p_j} \) of the central atom j. |
| ctr_atom2_pos | The position \( \vec{p_k} \) of the central atom k. |
| term_atom2_pos | The position \( \vec{p_l} \) of the terminal atom l. |
| float CDPL.ForceField.calcDistanceDerivatives | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| Math.Vector3D | atom1_deriv, | ||
| Math.Vector3D | atom2_deriv | ||
| ) |
Calculates the partial derivatives \( \frac{\partial r_{ij}}{\partial \vec{p_x}} \) of the distance \( r_{ij} \) between two atoms i and j.
\( \frac{\partial r_{ij}}{\partial \vec{p_i}} = \frac{-\vec{v_{ij}}}{r_{ij}} \)
\( \frac{\partial r_{ij}}{\partial \vec{p_j}} = \frac{\vec{v_{ij}}}{r_{ij}} \)
\( r_{ij} = |\vec{v_{ij}}| \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| atom1_deriv | Output variable for the calculated partial derivative \( \frac{\partial r_{ij}}{\partial \vec{p_i}} \) at the given atom positions. |
| atom2_deriv | Output variable for the calculated partial derivative \( \frac{\partial r_{ij}}{\partial \vec{p_j}} \) at the given atom positions. |
| float CDPL.ForceField.calcDistance | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos | ||
| ) |
Calculates the distance \( r_{ij} \) between two atoms i and j.
\( r_{ij} = |\vec{v_{ij}}| \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| float CDPL.ForceField.calcElasticPotentialEnergy | ( | ElasticPotential | pot, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the energy of a single elastic potential pot for the supplied atom coordinates coords.
| pot | The parameters of the elastic potential to evaluate. |
| coords | The atom 3D coordinates for which the energy shall be calculated. |
| float CDPL.ForceField.calcElasticPotentialEnergy | ( | ElasticPotentialList | list, |
| Math.Vector3DArray | coords | ||
| ) |
| list | |
| coords |
| float CDPL.ForceField.calcElasticPotentialEnergy | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| float | force_const, | ||
| float | ref_length | ||
| ) |
Calculates the energy \( E_{ij} \) of an elastic potential applied on a pair of atoms i-j.
\( E_{ij} = k_{ij} \: \Delta r_{ij}^2 \)
where
\( k_{ij} \) = the force constant of the elastic potential.
\( \Delta r_{ij} \) = \( r_{ij} - r_{ij}^0 \), the difference between actual and reference distance of the atoms i and j.
| atom1_pos | The position of atom i. |
| atom2_pos | The position of atom j. |
| force_const | The force constant \( k_{ij} \). |
| ref_length | The reference distance \( r_{ij}^0 \). |
| float CDPL.ForceField.calcElasticPotentialGradient | ( | ElasticPotential | pot, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the energy and gradient contribution of a single elastic potential pot for the supplied atom coordinates coords.
| pot | The parameters of the elastic potential to evaluate. |
| coords | The atom 3D coordinates for which the energies/gradients shall be calculated. |
| grad | The atom gradient vector array receiving the accumulated contributions. |
| float CDPL.ForceField.calcElasticPotentialGradient | ( | ElasticPotentialList | list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| list | |
| coords | |
| grad |
| float CDPL.ForceField.calcElasticPotentialGradient | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| Math.Vector3D | atom1_grad, | ||
| Math.Vector3D | atom2_grad, | ||
| float | force_const, | ||
| float | ref_length | ||
| ) |
Calculates the elastic potential energy gradient \( \nabla E_{ij} \) for a pair of atoms i-j.
Energy function:
\( E_{ij} = k_{ij} \: \Delta r_{ij}^2 \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial E_{ij}}{\partial \vec{p_x}} = \frac{\partial E_{ij}}{\partial \Delta r_{ij}} \: \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} \)
\( \frac{\partial E_{ij}}{\partial \Delta r_{ij}} = 2 \: \Delta r_{ij} \: k_{ij} \)
for the calculation of the partial derivatives \( \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} \) see calcDistanceDerivatives().
where
\( k_{ij} \) = the force constant of the elastic potential.
\( \Delta r_{ij} \) = \( r_{ij} - r_{ij}^0 \), the difference between actual and reference distance of the atoms i and j.
\( \vec{p_x} \) = coordinates of the atoms i and j.
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| atom2_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| force_const | The force constant \( k_{ij} \). |
| ref_length | The reference distance \( r_{ij}^0 \). |
| float CDPL.ForceField.calcMMFF94AngleBendingEnergy | ( | MMFF94AngleBendingInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the angle-bending interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 angle-bending interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94AngleBendingEnergy | ( | MMFF94AngleBendingInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94AngleBendingEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| bool | linear, | ||
| float | force_const, | ||
| float | ref_angle | ||
| ) |
Calculates the angle bending interaction energy \( EA_{ijk} \) for two bonds i-j and j-k.
\( EA_{ijk} = 0.043844 \: \frac{ka_{IJK}}{2} \: \Delta \vartheta_{ijk}^2 \: (1 + cb \: \Delta \vartheta_{ijk}) \)
where
\( ka_{IJK} \) = angle bending force constant in \( \frac{md Ang}{rad^2} \) for the angle between atoms i, j and k of atom types I, J and K.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees (see calcBondAngle()).
\( cb \) = \( -0.007 \: deg^{-1} \), the "cubic-bend" constant.
For linear or near-linear bond angles such as those which occur in alkynes, nitriles, isonitriles, azides, and diazo compounds, the form used in DREIDING and UFF is employed:
\( EA_{ijk} = 143.9325 \: ka_{IJK} \:(1 + \cos(\vartheta_{ijk})) \)
where \( ka_{IJK} \) and \( \vartheta_{ijk} \) are defined as above.
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| linear | If True, the bond angle is linear. |
| force_const | The angle bending force constant \( ka_{IJK} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| float CDPL.ForceField.calcMMFF94AngleBendingEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | r_ij, | ||
| float | r_jk, | ||
| bool | linear, | ||
| float | force_const, | ||
| float | ref_angle | ||
| ) |
Calculates the angle bending interaction energy \( EA_{ijk} \) for two bonds i-j and j-k.
\( EA_{ijk} = 0.043844 \: \frac{ka_{IJK}}{2} \: \Delta \vartheta_{ijk}^2 \: (1 + cb \: \Delta \vartheta_{ijk}) \)
where
\( ka_{IJK} \) = angle bending force constant in \( \frac{md Ang}{rad^2} \) for the angle between atoms i, j and k of atom types I, J and K.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees (see calcBondAngle()).
\( cb \) = \( -0.007 \: deg^{-1} \), the "cubic-bend" constant.
For linear or near-linear bond angles such as those which occur in alkynes, nitriles, isonitriles, azides, and diazo compounds, the form used in DREIDING and UFF is employed:
\( EA_{ijk} = 143.9325 \: ka_{IJK} \:(1 + \cos(\vartheta_{ijk})) \)
where \( ka_{IJK} \) and \( \vartheta_{ijk} \) are defined as above.
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| r_ij | The length of the bond between atom i and j. |
| r_jk | The length of the bond between atom j and k. |
| linear | If True, the bond angle is linear. |
| force_const | The angle bending force constant \( ka_{IJK} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| float CDPL.ForceField.calcMMFF94AngleBendingGradient | ( | MMFF94AngleBendingInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the angle-bending interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 angle-bending interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94AngleBendingGradient | ( | MMFF94AngleBendingInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94AngleBendingGradient | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | term_atom1_grad, | ||
| Math.Vector3D | ctr_atom_grad, | ||
| Math.Vector3D | term_atom2_grad, | ||
| bool | linear, | ||
| float | force_const, | ||
| float | ref_angle | ||
| ) |
Calculates the angle bending interaction energy gradient \( \nabla EA_{ijk} \) for two bonds i-j and j-k.
Energy function employed for the non-linear case:
\( EA_{ijk} = 0.043844 \: \frac{ka_{IJK}}{2} \: \Delta \vartheta_{ijk}^2 \: (1 + cb \: \Delta \vartheta_{ijk}) \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EA_{ijk}}{\partial \vec{p_x}} = \frac{\partial EA_{ijk}}{\partial \vartheta_{ijk}} \: \frac{\partial \vartheta_{ijk}}{\partial \cos(\vartheta_{ijk})} \: \frac{\partial \cos(\vartheta_{ijk})}{\vec{p_x}} \)
\( \frac{\partial EA_{ijk}}{\partial \vartheta_{ijk}} = -ka_{IJK} \: (86.58992538 \: \vartheta_{ijk}^2 - 3.022558594 \: \vartheta_{ijk} \: \vartheta_{IJK}^0 - 143.9313616 \: \vartheta_{ijk} + 0.02637679965 \: \vartheta_{IJK}^{0^2} + 2.512076157 \: \vartheta_{IJK}^0) \)
\( \frac{\partial \vartheta_{ijk}}{\partial \cos(\vartheta_{ijk})} = \frac{-1}{\sqrt{1 - \cos(\vartheta_{ijk})^2}} \)
for the calculation of the partial derivatives \( \frac{\partial \cos(\vartheta_{ijk})}{\vec{p_x}} \) see calcBondAngleCosDerivatives().
For linear or near-linear bond angles such as those which occur in alkynes, nitriles, isonitriles, azides, and diazo compounds, the energy function form used in DREIDING and UFF is employed:
\( EA_{ijk} = 143.9325 \: ka_{IJK} \:(1 + \cos(\vartheta_{ijk})) \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EA_{ijk}}{\partial \vec{p_x}} = 143.9325 \: ka_{IJK} \: \frac{\partial \cos(\vartheta_{ijk})}{\vec{p_x}} \)
where
\( ka_{IJK} \) = angle bending force constant in \( \frac{md Ang}{rad^2} \) for the angle between atoms i, j and k of atom types I, J and K.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees.
\( cb \) = \( -0.007 \: deg^{-1} \), the "cubic-bend" constant.
\( \vec{p_x} \) = coordinates of the involved atoms i, j and k.
| term_atom1_pos | The position \( \vec{p_i} \) of atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of atom k. |
| term_atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| ctr_atom_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| term_atom2_grad | The output variable storing the accumulated energy gradient contributions for atom k. |
| linear | If True, the bond angle is linear. |
| force_const | The angle bending force constant \( ka_{IJK} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| None CDPL.ForceField.calcMMFF94AtomCharges | ( | Chem.MolecularGraph | molgraph, |
| bool | strict, | ||
| bool | overwrite | ||
| ) |
Calculates and assigns MMFF94 partial atomic charges to the atoms of the molecular graph molgraph.
| molgraph | The molecular graph for which to assign atom charges. |
| strict | If True, strict parameterization is performed (and may fail). Otherwise, in case of parameterization problems, fallback strategies take effect. |
| overwrite | Specifies whether already assigned charge properties should be replaced. |
| float CDPL.ForceField.calcMMFF94BondStretchingEnergy | ( | MMFF94BondStretchingInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the bond-stretching interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 bond-stretching interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94BondStretchingEnergy | ( | MMFF94BondStretchingInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94BondStretchingEnergy | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| float | force_const, | ||
| float | ref_length | ||
| ) |
Calculates the bond stretching interaction energy \( EB_{ij} \) for the bond i-j.
\( EB_{ij} = 143.9325 \: \frac{kb_{IJ}}{2} \: \Delta r_{ij}^2 \: (1 + cs \: \Delta r_{ij} + \frac{7}{12} \: cs^2 \: \Delta r_{ij}^2) \)
where
\( kb_{IJ} \) = the bond stretching force constant in \( \frac{md}{Ang} \) for bonded atoms i and j of types I and J.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J (see calcDistance()).
\( cs \) = \( -2 \: Ang^{-1} \), the "cubic stretch" constant.
Note: throughout this description, the indices i, j, k, ... represent atoms; I, J, K, ... denote the corresponding numerical MMFF atom types (or, occasionally, the atomic species).
| atom1_pos | The position of atom i. |
| atom2_pos | The position of atom j. |
| force_const | The bond stretching force constant \( kb_{IJ} \). |
| ref_length | The reference bond length \( r_{IJ}^0 \). |
| float CDPL.ForceField.calcMMFF94BondStretchingEnergy | ( | float | r_ij, |
| float | force_const, | ||
| float | ref_length | ||
| ) |
Calculates the bond stretching interaction energy \( EB_{ij} \) for the bond i-j.
\( EB_{ij} = 143.9325 \: \frac{kb_{IJ}}{2} \: \Delta r_{ij}^2 \: (1 + cs \: \Delta r_{ij} + \frac{7}{12} \: cs^2 \: \Delta r_{ij}^2) \)
where
\( kb_{IJ} \) = the bond stretching force constant in \( \frac{md}{Ang} \) for bonded atoms i and j of types I and J.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J (see calcDistance()).
\( cs \) = \( -2 \: Ang^{-1} \), the "cubic stretch" constant.
Note: throughout this description, the indices i, j, k, ... represent atoms; I, J, K, ... denote the corresponding numerical MMFF atom types (or, occasionally, the atomic species).
| r_ij | The length of the bond between atom i and j. |
| force_const | The bond stretching force constant \( kb_{IJ} \). |
| ref_length | The reference bond length \( r_{IJ}^0 \). |
| float CDPL.ForceField.calcMMFF94BondStretchingGradient | ( | MMFF94BondStretchingInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the bond-stretching interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 bond-stretching interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94BondStretchingGradient | ( | MMFF94BondStretchingInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94BondStretchingGradient | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| Math.Vector3D | atom1_grad, | ||
| Math.Vector3D | atom2_grad, | ||
| float | force_const, | ||
| float | ref_length | ||
| ) |
Calculates the bond stretching interaction energy gradient \( \nabla EB_{ij} \) for the bond i-j.
Energy function:
\( EB_{ij} = 143.9325 \: \frac{kb_{IJ}}{2} \: \Delta r_{ij}^2 \: (1 + cs \: \Delta r_{ij} + \frac{7}{12} \: cs^2 \: \Delta r_{ij}^2) \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EB_{ij}}{\partial \vec{p_x}} = \frac{\partial EB_{ij}}{\partial \Delta r_{ij}} \: \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} \)
\( \frac{\partial EB_{ij}}{\partial \Delta r_{ij}} = (167.92125 \: \Delta r_{ij}^3 \: cs^2 + 215.89875 \: \Delta r_{ij}^2 \: cs + 143.9325 \: \Delta r_{ij}) \: kb_{IJ} \)
for the calculation of the partial derivatives \( \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} \) see calcDistanceDerivatives().
where
\( kb_{IJ} \) = the bond stretching force constant in \( \frac{md}{Ang} \) for bonded atoms i and j of types I and J.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J.
\( cs \) = \( -2 \: Ang^{-1} \), the "cubic stretch" constant.
\( \vec{p_x} \) = coordinates of the involved atoms i and j.
Note: throughout this description, the indices i, j, k, ... represent atoms; I, J, K, ... denote the corresponding numerical MMFF atom types (or, occasionally, the atomic species).
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| atom2_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| force_const | The bond stretching force constant \( kb_{IJ} \). |
| ref_length | The reference bond length \( r_{IJ}^0 \). |
| float CDPL.ForceField.calcMMFF94ElectrostaticEnergy | ( | MMFF94ElectrostaticInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the electrostatic interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 electrostatic interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94ElectrostaticEnergy | ( | MMFF94ElectrostaticInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94ElectrostaticEnergy | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| float | atom1_chg, | ||
| float | atom2_chg, | ||
| float | scale_fact, | ||
| float | de_const, | ||
| float | dist_expo | ||
| ) |
Calculates the electrostatic interaction energy \( EQ_{ij} \) for the atom pair i-j.
\( EQ_{ij} = S \: 332.0716 \: \frac{q_i \: q_j}{D \: (R_{ij} + \delta)^n} \)
where
\( S \) = a scaling factor depending on the topological distance of i-j.
\( q_i \) and \( q_j \) = Partial atomic charges.
\( D \) = Dielectric constant.
\( R_{ij} \) = Interatomic distance (Ã…) (see calcDistance()).
\( \delta \) = Electrostatic buffering constant (0.05 Ã…).
\( n \) = Exponent (normally 1, but can be 2 for distance-dependent dielectric constant).
Note: 1-4 electrostatic interactions are scaled by 0.75 (thus, the electrostatic energy term becomes \( EQ_{14} \: 0.75 \)).
| atom1_pos | The position of atom i. |
| atom2_pos | The position of atom j. |
| atom1_chg | The partial atom charge \( q_i \) of atom i. |
| atom2_chg | The partial atom charge \( q_j \) of atom j. |
| scale_fact | The scaling factor for \( S \) depending on the topological i-j distance. |
| de_const | The dielectric constant \( D \). |
| dist_expo | The exponent \( n \). |
| float CDPL.ForceField.calcMMFF94ElectrostaticEnergy | ( | float | r_ij, |
| float | atom1_chg, | ||
| float | atom2_chg, | ||
| float | scale_fact, | ||
| float | de_const, | ||
| float | dist_expo | ||
| ) |
Calculates the electrostatic interaction energy \( EQ_{ij} \) for the atom pair i-j.
\( EQ_{ij} = S \: 332.0716 \: \frac{q_i \: q_j}{D \: (R_{ij} + \delta)^n} \)
where
\( S \) = a scaling factor depending on the topological distance of i-j.
\( q_i \) and \( q_j \) = Partial atomic charges.
\( D \) = Dielectric constant.
\( R_{ij} \) = Interatomic distance (Ã…) (see calcDistance()).
\( \delta \) = Electrostatic buffering constant (0.05 Ã…).
\( n \) = Exponent (normally 1, but can be 2 for distance-dependent dielectric constant).
Note: 1-4 electrostatic interactions are scaled by 0.75 (thus, the electrostatic energy term becomes \( EQ_{14} \: 0.75 \)).
| r_ij | The interatomic distance \( R_{ij} \) of atom i and atom j. |
| atom1_chg | The partial atom charge \( q_i \) of atom i. |
| atom2_chg | The partial atom charge \( q_j \) of atom j. |
| scale_fact | The scaling factor for \( S \) depending on the topological i-j distance. |
| de_const | The dielectric constant \( D \). |
| dist_expo | The exponent \( n \). |
| float CDPL.ForceField.calcMMFF94ElectrostaticGradient | ( | MMFF94ElectrostaticInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the electrostatic interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 electrostatic interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94ElectrostaticGradient | ( | MMFF94ElectrostaticInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94ElectrostaticGradient | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| Math.Vector3D | atom1_grad, | ||
| Math.Vector3D | atom2_grad, | ||
| float | atom1_chg, | ||
| float | atom2_chg, | ||
| float | scale_fact, | ||
| float | de_const, | ||
| float | dist_expo | ||
| ) |
Calculates the electrostatic interaction energy gradient \( \nabla EQ_{ij} \) for the atom pair i-j.
Energy function:
\( EQ_{ij} = S \: 332.0716 \: \frac{q_i \: q_j}{D \: (R_{ij} + \delta)^n} \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EQ_{ij}}{\partial \vec{p_x}} = \frac{\partial EQ_{ij}}{\partial R_{ij}} \: \frac{\partial R_{ij}}{\partial \vec{p_x}} \)
\( \frac{\partial EQ_{ij}}{\partial R_{ij}} = -S \: 332.0716 \: n \: \frac{q_i \: q_j}{D \: (R_{ij} + \delta)^{n + 1}} \)
for the calculation of the partial derivatives \( \frac{\partial R_{ij}}{\partial \vec{p_x}} \) see calcDistanceDerivatives().
where
\( S \) = a scaling factor depending on the topological distance of i-j.
\( q_i \) and \( q_j \) = partial atomic charges.
\( D \) = dielectric constant.
\( R_{ij} \) = interatomic distance (Ã…).
\( \delta \) = electrostatic buffering constant (0.05 Ã…).
\( n \) = exponent (normally 1, but can be 2 for distance-dependent dielectric constant).
\( \vec{p_x} \) = coordinates of the involved atoms i and j.
Note: 1-4 electrostatic interactions are scaled by 0.75 (thus, the electrostatic gradient term becomes \( EQ_{14} \: 0.75 \)).
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| atom2_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| atom1_chg | The partial atom charge \( q_i \) of atom i. |
| atom2_chg | The partial atom charge \( q_j \) of atom j. |
| scale_fact | The scaling factor for \( S \) depending on the topological i-j distance. |
| de_const | The dielectric constant \( D \). |
| dist_expo | The exponent \( n \). |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingEnergy | ( | MMFF94OutOfPlaneBendingInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the out-of-plane bending interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 out-of-plane bending interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingEnergy | ( | MMFF94OutOfPlaneBendingInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | oop_atom_pos, | ||
| float | force_const | ||
| ) |
Calculates the out-of-plane bending interaction energy \( EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k.
\( EOOP_{ijk;l} = 0.043844 \: \frac{koop_{IJK \colon L}}{2} \: \chi_{ijk;l}^2 \)
where
\( koop_{IJK \colon L} \) = out-of-plane bending force constant in \( \frac{md Ang}{rad^2} \).
\( \chi_{ijk;l} \) = angle in degrees between the bond j-l and the plane i-j-k, where j is the central atom (see calcOutOfPlaneAngle()).
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| oop_atom_pos | The position of the out-of-plane atom l. |
| force_const | The out-of-plane bending force constant \( koop_{IJK \colon L} \). |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | oop_atom_pos, | ||
| float | r_jl, | ||
| float | force_const | ||
| ) |
Calculates the out-of-plane bending interaction energy \( EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k.
\( EOOP_{ijk;l} = 0.043844 \: \frac{koop_{IJK \colon L}}{2} \: \chi_{ijk;l}^2 \)
where
\( koop_{IJK \colon L} \) = out-of-plane bending force constant in \( \frac{md Ang}{rad^2} \).
\( \chi_{ijk;l} \) = angle in degrees between the bond j-l and the plane i-j-k, where j is the central atom (see calcOutOfPlaneAngle()).
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| oop_atom_pos | The position of the out-of-plane atom l. |
| r_jl | The length of the bond between atom j and atom l. |
| force_const | The out-of-plane bending force constant \( koop_{IJK \colon L} \). |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingGradient | ( | MMFF94OutOfPlaneBendingInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the out-of-plane bending interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 out-of-plane bending interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingGradient | ( | MMFF94OutOfPlaneBendingInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94OutOfPlaneBendingGradient | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | oop_atom_pos, | ||
| Math.Vector3D | term_atom1_grad, | ||
| Math.Vector3D | ctr_atom_grad, | ||
| Math.Vector3D | term_atom2_grad, | ||
| Math.Vector3D | oop_atom_grad, | ||
| float | force_const | ||
| ) |
Calculates the out-of-plane bending interaction energy gradient \( \nabla EOOP_{ijk;l} \) for the bond j-l and the plane i-j-k.
Energy function:
\( EOOP_{ijk;l} = 0.043844 \: \frac{koop_{IJK \colon L}}{2} \: (\chi_{ijk;l} \: \frac{180}{\pi})^2 \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EOOP_{ijk;l}}{\partial \vec{p_x}} = \frac{\partial EOOP_{ijk;l}}{\partial \chi_{ijk;l}} \: \frac{\partial \chi_{ijk;l}}{\partial \cos(\alpha_{ijk;l})} \: \frac{\partial \cos(\alpha_{ijk;l})}{\partial \vec{p_x}} \)
\( \frac{\partial EOOP_{ijk;l}}{\partial \chi_{ijk;l}} = 0.043844 \: (\frac{180}{\pi})^2 \: \chi_{ijk;l} \: koop_{IJK \colon L} \)
\( \chi_{ijk;l} = \frac{\pi}{2} - \alpha_{ijk;l} \)
\( \frac{\partial \chi_{ijk;l}}{\partial \cos(\alpha_{ijk;l})} = \frac{-1}{\sqrt{1 - \cos(\alpha_{ijk;l})^2}} \)
for the calculation of the partial derivatives \( \frac{\partial \cos(\alpha_{ijk;l})}{\partial \vec{p_x}} \) see calcOutOfPlaneAngleCosDerivatives().
where
\( koop_{IJK \colon L} \) = out-of-plane bending force constant in \( \frac{md Ang}{rad^2} \).
\( \chi_{ijk;l} \) = angle in radians between the bond j-l and the plane i-j-k, where j is the central atom.
\( \vec{p_x} \) = coordinates of the involved atoms i, j, k and l.
| term_atom1_pos | The position \( \vec{p_i} \) of atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of atom k. |
| oop_atom_pos | The position \( \vec{p_l} \) of the out-of-plane atom l. |
| term_atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| ctr_atom_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| term_atom2_grad | The output variable storing the accumulated energy gradient contributions for atom k. |
| oop_atom_grad | The output variable storing the accumulated energy gradient contributions for atom l. |
| force_const | The out-of-plane bending force constant \( koop_{IJK \colon L} \). |
| float CDPL.ForceField.calcMMFF94StretchBendEnergy | ( | MMFF94StretchBendInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the stretch-bend coupling interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 stretch-bend interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94StretchBendEnergy | ( | MMFF94StretchBendInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94StretchBendEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | ijk_force_const, | ||
| float | kji_force_const, | ||
| float | ref_angle, | ||
| float | ref_length1, | ||
| float | ref_length2 | ||
| ) |
Calculates the stretch-bend interaction energy \( EBA_{ijk} \) for two bonds i-j and j-k.
\( EBA_{ijk} = 2.51210 \: (kba_{IJK} \: \Delta r_{ij} + kba_{KJI} \: \Delta r_{kj}) \: \Delta \vartheta_{ijk} \)
where
\( kba_{IJK} \) = force constant in \( \frac{md}{rad} \) for i-j stretch coupled to i-j-k bend.
\( kba_{KJI} \) = force constant in \( \frac{md}{rad} \) for k-j stretch coupled to i-j-k bend.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J.
\( \Delta r_{kj} \) = \( r_{kj} - r_{KJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms k and j of types K and J.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees.
For the calculation of \( r_{ij} \), \( r_{kj} \), and \( \vartheta_{ijk} \) see calcBondLengthsAndAngle().
Currently, stretch-bend interactions are omitted when the i-j-k interaction corresponds to a linear bond angle.
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| ijk_force_const | The stretch-bend force constant \( kba_{IJK} \). |
| kji_force_const | The stretch-bend force constant \( kba_{KJI} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| ref_length1 | The reference bond length \( r_{IJ}^0 \). |
| ref_length2 | The reference bond length \( r_{KJ}^0 \). |
| float CDPL.ForceField.calcMMFF94StretchBendEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | r_ij, | ||
| float | r_jk, | ||
| float | ijk_force_const, | ||
| float | kji_force_const, | ||
| float | ref_angle, | ||
| float | ref_length1, | ||
| float | ref_length2 | ||
| ) |
Calculates the stretch-bend interaction energy \( EBA_{ijk} \) for two bonds i-j and j-k.
\( EBA_{ijk} = 2.51210 \: (kba_{IJK} \: \Delta r_{ij} + kba_{KJI} \: \Delta r_{kj}) \: \Delta \vartheta_{ijk} \)
where
\( kba_{IJK} \) = force constant in \( \frac{md}{rad} \) for i-j stretch coupled to i-j-k bend.
\( kba_{KJI} \) = force constant in \( \frac{md}{rad} \) for k-j stretch coupled to i-j-k bend.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J.
\( \Delta r_{kj} \) = \( r_{kj} - r_{KJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms k and j of types K and J.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees.
For the calculation of \( r_{ij} \), \( r_{kj} \), and \( \vartheta_{ijk} \) see calcBondLengthsAndAngle().
Currently, stretch-bend interactions are omitted when the i-j-k interaction corresponds to a linear bond angle.
| term_atom1_pos | The position of atom i. |
| ctr_atom_pos | The position of the central atom j. |
| term_atom2_pos | The position of atom k. |
| r_ij | The length of the bond between atom i and j. |
| r_jk | The length of the bond between atom j and k. |
| ijk_force_const | The stretch-bend force constant \( kba_{IJK} \). |
| kji_force_const | The stretch-bend force constant \( kba_{KJI} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| ref_length1 | The reference bond length \( r_{IJ}^0 \). |
| ref_length2 | The reference bond length \( r_{KJ}^0 \). |
| float CDPL.ForceField.calcMMFF94StretchBendGradient | ( | MMFF94StretchBendInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the stretch-bend coupling interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 stretch-bend interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94StretchBendGradient | ( | MMFF94StretchBendInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94StretchBendGradient | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | term_atom1_grad, | ||
| Math.Vector3D | ctr_atom_grad, | ||
| Math.Vector3D | term_atom2_grad, | ||
| float | ijk_force_const, | ||
| float | kji_force_const, | ||
| float | ref_angle, | ||
| float | ref_length1, | ||
| float | ref_length2 | ||
| ) |
Calculates the stretch-bend interaction energy gradient \( \nabla EBA_{ijk} \) for two bonds i-j and j-k.
Energy function:
\( EBA_{ijk} = 2.51210 \: (kba_{IJK} \: \Delta r_{ij} + kba_{KJI} \: \Delta r_{kj}) \: \Delta \vartheta_{ijk} \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial EBA_{ijk}}{\partial \vec{p_x}} = 2.5121 \: \Delta \vartheta_{ijk} \: (kba_{IJK} \: \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} + kba_{KJI} \: \frac{\partial \Delta r_{kj}}{\partial \vec{p_x}}) + 2.5121 \: \frac{\partial \Delta \vartheta_{ijk}}{\partial \vec{p_x}} \: (kba_{IJK} \: \Delta r_{ij} + kba_{KJI} \: \Delta r_{kj}) \)
\( \frac{\partial \Delta \vartheta_{ijk}}{\partial \vec{p_x}} = \frac{\partial \Delta \vartheta_{ijk}}{\partial \vartheta_{ijk}} \: \frac{\partial \vartheta_{ijk}}{\partial \cos(\vartheta_{ijk})} \: \frac{\partial \cos(\vartheta_{ijk})}{\vec{p_x}} \)
\( \frac{\partial \Delta \vartheta_{ijk}}{\partial \vartheta_{ijk}} = \frac{180}{\pi} \)
\( \frac{\partial \vartheta_{ijk}}{\partial \cos(\vartheta_{ijk})} = \frac{-1}{\sqrt{1 - \cos(\vartheta_{ijk})^2}} \)
for the calculation of the partial derivatives \( \frac{\partial \cos(\vartheta_{ijk})}{\vec{p_x}} \) see calcBondAngleCosDerivatives() and for the calculation of \( \frac{\partial \Delta r_{ij}}{\partial \vec{p_x}} \) see calcDistanceDerivatives().
where
\( kba_{IJK} \) = force constant in \( \frac{md}{rad} \) for i-j stretch coupled to i-j-k bend.
\( kba_{KJI} \) = force constant in \( \frac{md}{rad} \) for k-j stretch coupled to i-j-k bend.
\( \Delta r_{ij} \) = \( r_{ij} - r_{IJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms i and j of types I and J.
\( \Delta r_{kj} \) = \( r_{kj} - r_{KJ}^0 \), the difference in angstroms between actual and reference bond lengths between bonded atoms k and j of types K and J.
\( \Delta \vartheta_{ijk} \) = \( \vartheta_{ijk} \: \frac{180}{\pi} - \vartheta_{IJK}^0 \), the difference between actual and reference i-j-k bond angles in degrees.
\( \vec{p_x} \) = coordinates of the involved atoms i, j and k.
Currently, stretch-bend interactions are omitted when the i-j-k interaction corresponds to a linear bond angle.
| term_atom1_pos | The position \( \vec{p_i} \) of atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of atom k. |
| term_atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| ctr_atom_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| term_atom2_grad | The output variable storing the accumulated energy gradient contributions for atom k. |
| ijk_force_const | The stretch-bend force constant \( kba_{IJK} \). |
| kji_force_const | The stretch-bend force constant \( kba_{KJI} \). |
| ref_angle | The reference bond angle \( \vartheta_{IJK}^0 \). |
| ref_length1 | The reference bond length \( r_{IJ}^0 \). |
| ref_length2 | The reference bond length \( r_{KJ}^0 \). |
| float CDPL.ForceField.calcMMFF94TorsionEnergy | ( | MMFF94TorsionInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the torsion interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 torsion interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94TorsionEnergy | ( | MMFF94TorsionInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94TorsionEnergy | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom1_pos, | ||
| Math.Vector3D | ctr_atom2_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| float | tor_param1, | ||
| float | tor_param2, | ||
| float | tor_param3 | ||
| ) |
Calculates the torsion interaction energy \( ET_{ijkl} \) for the central bond j-k and the connected bonds i-j and k-l.
\( ET_{ijkl} = 0.5 \: (V_1 \: (1 + \cos(\Phi_{ijkl})) + V_2 \: (1 - \cos(2 \: \Phi_{ijkl})) + V_3 \: (1 + \cos(3 \: \Phi_{ijkl}))) \)
where \( \Phi_{ijkl} \) is the i-j-k-l dihedral angle. The constants \( V_1 \), \( V_2 \) and \( V_3 \) depend on the atom types I, J, K and L for atoms i, j, k and l, where i-j, j-k and k-l are bonded pairs and i is not equal to l.
For the calculation of \( \cos(\Phi_{ijkl}) \) see calcDihedralAngleCos().
| term_atom1_pos | The position of the terminal atom i. |
| ctr_atom1_pos | The position of the central atom j. |
| ctr_atom2_pos | The position of the central atom k. |
| term_atom2_pos | The position of the terminal atom l. |
| tor_param1 | The torsion parameter \( V_1 \). |
| tor_param2 | The torsion parameter \( V_2 \). |
| tor_param3 | The torsion parameter \( V_3 \). |
| float CDPL.ForceField.calcMMFF94TorsionGradient | ( | MMFF94TorsionInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the torsion interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 torsion interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94TorsionGradient | ( | MMFF94TorsionInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94TorsionGradient | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom1_pos, | ||
| Math.Vector3D | ctr_atom2_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | term_atom1_grad, | ||
| Math.Vector3D | ctr_atom1_grad, | ||
| Math.Vector3D | ctr_atom2_grad, | ||
| Math.Vector3D | term_atom2_grad, | ||
| float | tor_param1, | ||
| float | tor_param2, | ||
| float | tor_param3 | ||
| ) |
Calculates the torsion interaction energy gradient \( \nabla ET_{ijkl} \) for the central bond j-k and the connected bonds i-j and k-l.
Energy function:
\( ET_{ijkl} = 0.5 \: (V_1 \: (1 + \cos(\Phi_{ijkl})) + V_2 \: (1 - \cos(2 \: \Phi_{ijkl})) + V_3 \: (1 + \cos(3 \: \Phi_{ijkl}))) \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial ET_{ijkl}}{\partial \vec{p_x}} = \frac{\partial ET_{ijkl}}{\partial \Phi_{ijkl}} \: \frac{\partial \Phi_{ijkl}}{\partial \cos(\Phi_{ijkl})} \: \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_x}} \)
\( \frac{\partial ET_{ijkl}}{\partial \Phi_{ijkl}} = V_2 \: \sin(2 \: \Phi_{ijkl}) - 0.5 \: V_1 \: \sin(\Phi_{ijkl}) - 1.5 \: V_3 \: \sin(3 \: \Phi_{ijkl}) \)
\( \frac{\partial \Phi_{ijkl}}{\partial \cos(\Phi_{ijkl})} = \frac{-1}{\sqrt{1 - \cos(\Phi_{ijkl})^2}} \)
for the calculation of the partial derivatives \( \frac{\partial \cos(\Phi_{ijkl})}{\partial \vec{p_x}} \) see calcDihedralAngleCosDerivatives().
where
\( \Phi_{ijkl} \) is the i-j-k-l dihedral angle. The constants \( V_1 \), \( V_2 \) and \( V_3 \) depend on the atom types I, J, K and L for atoms i, j, k and l, where i-j, j-k and k-l are bonded pairs and i is not equal to l.
\( \vec{p_x} \) = coordinates of the involved atoms i, j, k and l.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom1_pos | The position \( \vec{p_j} \) of the central atom j. |
| ctr_atom2_pos | The position \( \vec{p_k} \) of the central atom k. |
| term_atom2_pos | The position \( \vec{p_l} \) of the terminal atom l. |
| term_atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| ctr_atom1_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| ctr_atom2_grad | The output variable storing the accumulated energy gradient contributions for atom k. |
| term_atom2_grad | The output variable storing the accumulated energy gradient contributions for atom l. |
| tor_param1 | The torsion parameter \( V_1 \). |
| tor_param2 | The torsion parameter \( V_2 \). |
| tor_param3 | The torsion parameter \( V_3 \). |
| float CDPL.ForceField.calcMMFF94VanDerWaalsEnergy | ( | MMFF94VanDerWaalsInteraction | iaction, |
| Math.Vector3DArray | coords | ||
| ) |
Calculates the Van der Waals interaction energy of iaction using the geometry from coords.
| iaction | The MMFF94 Van der Waals interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| float CDPL.ForceField.calcMMFF94VanDerWaalsEnergy | ( | MMFF94VanDerWaalsInteractionList | ia_list, |
| Math.Vector3DArray | coords | ||
| ) |
| ia_list | |
| coords |
| float CDPL.ForceField.calcMMFF94VanDerWaalsEnergy | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| float | e_IJ, | ||
| float | r_IJ, | ||
| float | r_IJ_7 | ||
| ) |
Calculates the Van der Waals interaction energy \( E_{vdW_{ij}} \) for the atom pair i-j.
\( E_{vdW_{ij}} = \varepsilon_{IJ} \: (\frac{1.07 \: R_{IJ}^*}{(R_{ij} + 0.07 \: R_{IJ}^*)})^7 \: (\frac{1.12 \: R_{IJ}^{*^7}}{(R_{ij}^7 + 0.12 \: R_{IJ}^{*^7})} - 2) \;\;\;\; (1) \)
where
\( R_{ij} \) = the interatomic distance (see calcDistance()).
\( R_{II}^* = A_I \: \alpha_I^{PEXP} \;\;\;\; (2) \)
\( R_{IJ}^* = 0.5 \: (R_{II}^* + R_{JJ}^*) \: (1 + AFACT(1 - \exp(-BFACT \: \gamma_{IJ}^2))) \;\;\;\; (3) \)
\( \gamma_{IJ} = \frac{(R_{II}^* - R_{JJ}^*)}{(R_{II}^* + R_{JJ}^*)} \;\;\;\; (4) \)
\( \varepsilon_{IJ} = \frac{181.16 \: G_I \: GJ \: \alpha_I \: \alpha_J}{((\alpha_I / N_I)^{1/2} + (\alpha_J / N_J)^{1/2})} \: \frac{1}{R_{IJ}^{*^6}} \;\;\;\; (5) \)
MMFF employs a "Buffered 14-7" form (eq 1) together with an expression which relates the minimum-energy separation \( R_{II}^* \) to the atomic polarizability \( \alpha_I \) (eq 2), a specially formulated combination rule (eqs 3, 4), and a Slater-Kirkwood expression for the well depth \( \varepsilon_{IJ} \) (eq 5): The first non-comment line in the parameter file "MMFFVDW.PAR" contains five floating point numbers which define the variables PEXP, AFACT, BFACT, DARAD, and DAEPS, respectively. PEXP (currently 0.25), AFACT (currently 0.2) and BFACT (currently 12.0) are used in the equations shown above; DARAD and DAEPS are used as follows:
When either i or j is a hydrogen-bond donor, MMFF uses the simple arithmetic mean \( 0.5 \: (R_{II}^* + R_{JJ}^*) \) instead of eq 3 to obtain \( R_{IJ}^* \). If the i-j interaction is a donor-acceptor interaction, MMFF also scales \( R_{IJ}^* \) as given by eq 3 by DARAD (currently 0.8) and \( \varepsilon_{IJ} \) as given by eq 5 by DAEPS (currently 0.5).
| atom1_pos | The position of atom i. |
| atom2_pos | The position of atom j. |
| e_IJ | The precalculated value \( \varepsilon_{IJ} \). |
| r_IJ | The precalculated value \( R_{IJ}^* \). |
| r_IJ_7 | The precalculated value \( R_{IJ}^{*^7} \). |
| float CDPL.ForceField.calcMMFF94VanDerWaalsEnergy | ( | float | r_ij, |
| float | e_IJ, | ||
| float | r_IJ, | ||
| float | r_IJ_7 | ||
| ) |
Calculates the Van der Waals interaction energy \( E_{vdW_{ij}} \) for the atom pair i-j.
\( E_{vdW_{ij}} = \varepsilon_{IJ} \: (\frac{1.07 \: R_{IJ}^*}{(R_{ij} + 0.07 \: R_{IJ}^*)})^7 \: (\frac{1.12 \: R_{IJ}^{*^7}}{(R_{ij}^7 + 0.12 \: R_{IJ}^{*^7})} - 2) \;\;\;\; (1) \)
where
\( R_{ij} \) = the interatomic distance (see calcDistance()).
\( R_{II}^* = A_I \: \alpha_I^{PEXP} \;\;\;\; (2) \)
\( R_{IJ}^* = 0.5 \: (R_{II}^* + R_{JJ}^*) \: (1 + AFACT(1 - \exp(-BFACT \: \gamma_{IJ}^2))) \;\;\;\; (3) \)
\( \gamma_{IJ} = \frac{(R_{II}^* - R_{JJ}^*)}{(R_{II}^* + R_{JJ}^*)} \;\;\;\; (4) \)
\( \varepsilon_{IJ} = \frac{181.16 \: G_I \: GJ \: \alpha_I \: \alpha_J}{((\alpha_I / N_I)^{1/2} + (\alpha_J / N_J)^{1/2})} \: \frac{1}{R_{IJ}^{*^6}} \;\;\;\; (5) \)
MMFF employs a "Buffered 14-7" form (eq 1) together with an expression which relates the minimum-energy separation \( R_{II}^* \) to the atomic polarizability \( \alpha_I \) (eq 2), a specially formulated combination rule (eqs 3, 4), and a Slater-Kirkwood expression for the well depth \( \varepsilon_{IJ} \) (eq 5): The first non-comment line in the parameter file "MMFFVDW.PAR" contains five floating point numbers which define the variables PEXP, AFACT, BFACT, DARAD, and DAEPS, respectively. PEXP (currently 0.25), AFACT (currently 0.2) and BFACT (currently 12.0) are used in the equations shown above; DARAD and DAEPS are used as follows:
When either i or j is a hydrogen-bond donor, MMFF uses the simple arithmetic mean \( 0.5 \: (R_{II}^* + R_{JJ}^*) \) instead of eq 3 to obtain \( R_{IJ}^* \). If the i-j interaction is a donor-acceptor interaction, MMFF also scales \( R_{IJ}^* \) as given by eq 3 by DARAD (currently 0.8) and \( \varepsilon_{IJ} \) as given by eq 5 by DAEPS (currently 0.5).
| r_ij | The interatomic distance \( R_{ij} \) of atom i and atom j. |
| e_IJ | The precalculated value \( \varepsilon_{IJ} \). |
| r_IJ | The precalculated value \( R_{IJ}^* \). |
| r_IJ_7 | The precalculated value \( R_{IJ}^{*^7} \). |
| float CDPL.ForceField.calcMMFF94VanDerWaalsGradient | ( | MMFF94VanDerWaalsInteraction | iaction, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
Calculates the Van der Waals interaction energy of iaction and accumulates the corresponding atom-position gradient contribution into grad.
| iaction | The MMFF94 Van der Waals interaction record. |
| coords | The atom-coordinates array providing the geometry. |
| grad | The gradient vector to accumulate into. |
| float CDPL.ForceField.calcMMFF94VanDerWaalsGradient | ( | MMFF94VanDerWaalsInteractionList | ia_list, |
| Math.Vector3DArray | coords, | ||
| Math.Vector3DArray | grad | ||
| ) |
| ia_list | |
| coords | |
| grad |
| float CDPL.ForceField.calcMMFF94VanDerWaalsGradient | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos, | ||
| Math.Vector3D | atom1_grad, | ||
| Math.Vector3D | atom2_grad, | ||
| float | e_IJ, | ||
| float | r_IJ, | ||
| float | r_IJ_7 | ||
| ) |
Calculates the Van der Waals interaction energy gradient \( \nabla E_{vdW_{ij}} \) for the atom pair i-j.
Energy function:
\( E_{vdW_{ij}} = \varepsilon_{IJ} \: (\frac{1.07 \: R_{IJ}^*}{(R_{ij} + 0.07 \: R_{IJ}^*)})^7 \: (\frac{1.12 \: R_{IJ}^{*^7}}{(R_{ij}^7 + 0.12 \: R_{IJ}^{*^7})} - 2) \;\;\;\; (1) \)
The partial derivatives with respect to the atom coordinates \( \vec{p_x} \) are calculated by:
\( \frac{\partial E_{vdW_{ij}}}{\partial \vec{p_x}} = \frac{\partial E_{vdW_{ij}}}{\partial R_{ij}} \: \frac{\partial R_{ij}}{\partial \vec{p_x}} \)
\( \frac{\partial E_{vdW_{ij}}}{\partial R_{ij}} = \frac{-R_{IJ}^{*^7} \: \varepsilon_{IJ}}{(R_{ij} + 0.07 \: R_{IJ}^*)^8 \: (R_{ij}^7 + 0.12 \: R_{IJ}^{*^7})^2} \: (-22.48094067 \: R_{ij}^{14} + 19.78322779 \: R_{ij}^7 \: R_{IJ}^{*^7} + 0.8812528743 \: R_{ij}^6 \: R_{IJ}^{*^8} + 1.186993667 \: R_{IJ}^{*^{14}}) \)
for the calculation of the partial derivatives \( \frac{\partial R_{ij}}{\partial \vec{p_x}} \) see calcDistanceDerivatives().
where
\( R_{ij} \) = the interatomic distance.
\( R_{II}^* = A_I \: \alpha_I^{PEXP} \;\;\;\; (2) \)
\( R_{IJ}^* = 0.5 \: (R_{II}^* + R_{JJ}^*) \: (1 + AFACT(1 - \exp(-BFACT \: \gamma_{IJ}^2))) \;\;\;\; (3) \)
\( \gamma_{IJ} = \frac{(R_{II}^* - R_{JJ}^*)}{(R_{II}^* + R_{JJ}^*)} \;\;\;\; (4) \)
\( \varepsilon_{IJ} = \frac{181.16 \: G_I \: GJ \: \alpha_I \: \alpha_J}{((\alpha_I / N_I)^{1/2} + (\alpha_J / N_J)^{1/2})} \: \frac{1}{R_{IJ}^{*6}} \;\;\;\; (5) \)
\( \vec{p_x} \) = coordinates of the involved atoms i and j.
MMFF employs a "Buffered 14-7" form (eq 1) together with an expression which relates the minimum-gradient separation \( R_{II}^* \) to the atomic polarizability \( \alpha_I \) (eq 2), a specially formulated combination rule (eqs 3, 4), and a Slater-Kirkwood expression for the well depth \( \varepsilon_{IJ} \) (eq 5): The first non-comment line in the parameter file "MMFFVDW.PAR" contains five floating point numbers which define the variables PEXP, AFACT, BFACT, DARAD, and DAEPS, respectively. PEXP (currently 0.25), AFACT (currently 0.2) and BFACT (currently 12.0) are used in the equations shown above; DARAD and DAEPS are used as follows:
When either i or j is a hydrogen-bond donor, MMFF uses the simple arithmetic mean \( 0.5 \: (R_{II}^* + R_{JJ}^*) \) instead of eq 3 to obtain \( R_{IJ}^* \). If the i-j interaction is a donor-acceptor interaction, MMFF also scales \( R_{IJ}^* \) as given by eq 3 by DARAD (currently 0.8) and \( \varepsilon_{IJ} \) as given by eq 5 by DAEPS (currently 0.5).
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| atom1_grad | The output variable storing the accumulated energy gradient contributions for atom i. |
| atom2_grad | The output variable storing the accumulated energy gradient contributions for atom j. |
| e_IJ | The precalculated value \( \varepsilon_{IJ} \). |
| r_IJ | The precalculated value \( R_{IJ}^* \). |
| r_IJ_7 | The precalculated value \( R_{IJ}^{*^7} \). |
| float CDPL.ForceField.calcOutOfPlaneAngleCosDerivatives | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | oop_atom_pos, | ||
| Math.Vector3D | term_atom1_deriv, | ||
| Math.Vector3D | ctr_atom_deriv, | ||
| Math.Vector3D | term_atom2_deriv, | ||
| Math.Vector3D | oop_atom_deriv | ||
| ) |
Calculates the partial derivatives \( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_x}} \) of the cosine of the angle \( \omega_{ijk;l} \) between the bond j-l and the normal of the plane defined by the atoms i-j-k.
\( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_i}} = \frac{\vec{v_{jk}} \times \vec{v_{jl}}}{|\vec{n_{ijk}}| \: r_{jl}} - \cos(\omega_{ijk;l}) \: \frac{M_1 \cdot \vec{n_{ijk}}}{|\vec{n_{ijk}}|^2} \)
\( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_k}} = \frac{\vec{v_{jl}} \times \vec{v_{ji}}}{|\vec{n_{ijk}}| \: r_{jl}} - \cos(\omega_{ijk;l}) \: \frac{M_2 \cdot \vec{n_{ijk}}}{|\vec{n_{ijk}}|^2} \)
\( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_l}} = \frac{-1}{|\vec{n_{ijk}}| \: r_{jl}} \: (\frac{\vec{v_{jl}} (\vec{n_{ijk}} \cdot \vec{v_{jl}})}{r_{jl}^2} + \vec{r_{kl}} \times \vec{v_{il}} + \vec{v_{jl}} \times \vec{v_{ji}} + \vec{v_{jk}} \times \vec{v_{jl}}) \)
\( \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_j}} = -(\frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_i}} + \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_k}} + \frac{\partial \cos(\omega_{ijk;l})}{\partial \vec{p_l}}) \)
\( \cos(\omega_{ijk;l}) = \frac{\vec{n_{ijk}} \cdot \vec{v_{jl}}}{|\vec{n_{ijk}}| \: r_{jl}} \)
where
\( M_1 = \begin{vmatrix} 0 & -\vec{v_{jk}}.z & \vec{v_{jk}}.y \\ \vec{v_{jk}}.z & 0 & -\vec{v_{jk}}.x \\ -\vec{v_{jk}}.y & \vec{v_{jk}}.x & 0 \end{vmatrix} \)
\( M_2 = \begin{vmatrix} 0 & \vec{v_{ji}}.z & -\vec{v_{ji}}.y \\ -\vec{v_{ji}}.z & 0 & \vec{v_{ji}}.x \\ \vec{v_{ji}}.y & -\vec{v_{ji}}.x & 0 \end{vmatrix} \)
\( \vec{v_{ji}} = \vec{p_i} - \vec{p_j} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{v_{jl}} = \vec{p_l} - \vec{p_j} \)
\( \vec{r_{kl}} = \vec{p_l} - \vec{p_k} \)
\( \vec{v_{il}} = \vec{p_l} - \vec{p_i} \)
\( \vec{n_{ijk}} = \vec{v_{ji}} \times \vec{v_{jk}} \)
\( r_{jl} = |\vec{v_{jl}}| \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
\( \vec{p_l} \) = coordinates of atom l.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| oop_atom_pos | The position \( \vec{p_l} \) of the out-of-plane atom l. |
| term_atom1_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\omega_{ijkl})}{\partial \vec{p_i}} \) at the given atom positions. |
| ctr_atom_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\omega_{ijkl})}{\partial \vec{p_j}} \) at the given atom positions. |
| term_atom2_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\omega_{ijkl})}{\partial \vec{p_k}} \) at the given atom positions. |
| oop_atom_deriv | Output variable for the calculated partial derivative \( \frac{\partial \cos(\omega_{ijkl})}{\partial \vec{p_l}} \) at the given atom positions. |
| float CDPL.ForceField.calcOutOfPlaneAngle | ( | Math.Vector3D | term_atom1_pos, |
| Math.Vector3D | ctr_atom_pos, | ||
| Math.Vector3D | term_atom2_pos, | ||
| Math.Vector3D | oop_atom_pos, | ||
| float | r_jl | ||
| ) |
Calculates the out-of-plane angle \( \chi_{ijk;l} \) between the bond j-l and the plane defined by the atoms i-j-k.
\( \chi_{ijk;l} = \frac{\pi}{2} - \arccos(\frac{\vec{n_{ijk}} \cdot \vec{v_{jl}}}{|\vec{n_{ijk}}| \: |\vec{v_{jl}}|}) \)
where
\( \vec{v_{ji}} = \vec{p_i} - \vec{p_j} \)
\( \vec{v_{jk}} = \vec{p_k} - \vec{p_j} \)
\( \vec{v_{jl}} = \vec{p_l} - \vec{p_j} \)
\( \vec{n_{ijk}} = \vec{v_{ji}} \times \vec{v_{jk}} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
\( \vec{p_k} \) = coordinates of atom k.
\( \vec{p_l} \) = coordinates of atom l.
| term_atom1_pos | The position \( \vec{p_i} \) of the terminal atom i. |
| ctr_atom_pos | The position \( \vec{p_j} \) of the central atom j. |
| term_atom2_pos | The position \( \vec{p_k} \) of the terminal atom k. |
| oop_atom_pos | The position \( \vec{p_l} \) of the out-of-plane atom l. |
| r_jl | The length of the bond between atom j and atom l. |
| float CDPL.ForceField.calcSquaredDistance | ( | Math.Vector3D | atom1_pos, |
| Math.Vector3D | atom2_pos | ||
| ) |
Calculates the squared distance \( r_{ij}^2 \) between two atoms i and j.
\( r_{ij}^2 = |\vec{v_{ij}}|^2 \)
where
\( \vec{v_{ij}} = \vec{p_j} - \vec{p_i} \)
\( \vec{p_i} \) = coordinates of atom i.
\( \vec{p_j} \) = coordinates of atom j.
| atom1_pos | The position \( \vec{p_i} \) of atom i. |
| atom2_pos | The position \( \vec{p_j} \) of atom j. |
| None CDPL.ForceField.filterInteractions | ( | MMFF94InteractionData | ia_data, |
| MMFF94InteractionData | filtered_ia_data, | ||
| Util.BitSet | inc_atom_mask | ||
| ) |
Filters an MMFF94 interaction data set, retaining only those interactions that exclusively reference atoms in inc_atom_mask.
| ia_data | The input interaction data set. |
| filtered_ia_data | The output interaction data set receiving the retained interactions. |
| inc_atom_mask | A bit set whose set bits correspond to the indices of atoms that may appear in a retained interaction. |
| Chem.FragmentList CDPL.ForceField.perceiveMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph | ) |
Perceives the list of MMFF94 aromatic rings of the molecular graph molgraph.
| molgraph | The molecular graph. |
| Chem.FragmentList CDPL.ForceField.perceiveMMFF94AromaticRings | ( | Chem.MolecularGraph | molgraph, |
| bool | overwrite | ||
| ) |
Perceives and (optionally) stores the list of MMFF94 aromatic rings of the molecular graph molgraph.
| molgraph | The molecular graph to inspect/modify. |
| overwrite | Specifies whether an already existing value of the property should be replaced. |
| int CDPL.ForceField.perceiveUFFType | ( | Chem.Atom | atom, |
| Chem.MolecularGraph | molgraph | ||
| ) |
Perceives the UFF atom type of the atom atom in the context of the molecular graph molgraph.
| atom | The atom whose UFF type is to be perceived. |
| molgraph | The molecular graph providing the atom environment for type perception. |