**** NAMD **** .. py:class:: NAMD This interface currently works as both an MD driver and an MM energy evaluator for MD simulations with the `NAMD `_ molecular dynamics simulation package (version 2.14 or later) A NAMD :ref:`Theory ` object can be specified as follows:: my_namd = NAMD(frag=..., par=..., pdb=..., psf=..., ) which would then be called as ``theory=namd`` and ``driver='namd'`` in a subsequent :py:class:`MD` task:: md = MD(theory=namd, freq_dcd=..., freq_chk=..., driver='namd', minimise=..., nsteps=..., nsnapshots=..., timestep=..., temperature=..., boundary=..., ) md.run(dryrun=False) This page explains the NAMD-specific keyword options, including those for NAMD as roles of both MD and MM drivers. For more generic MD options see :py:class:`MD`. .. note:: The NAMD interface is not currently compatible with other tasks such as ``SP()`` or ``Opt()``. This support will be added in a future release of ChemShell. General options =============== .. py:attribute:: binary (default: ``True``) NAMD outputs in binary format .. py:attribute:: coor (default: ``'_namd.coor'``) Filename of binary coordinates a restarted job runs from .. py:attribute:: input (default: ``'_namd.inp'``) Input script's filename .. py:attribute:: margin (default: ``4.72432`` Bohr) This parameters affects only the performance. Defaulted to ``0.0`` in NAMD, see: `NAMD Manual `__ .. py:attribute:: output (default: ``'_namd.out'``) Output filename .. py:attribute:: par (default: ``''``) Forcefield parameters filename(s). Allowed types of value are: * A Python :py:class:`str` of a single filename * A Python :py:class:`list` of multiple filenames. They will be concatenated by ChemShell to generate a single file .. py:attribute:: pdb (default: ``''``) PDB filename for input structure .. py:attribute:: prefix (default: ``'_namd'``) Prefix for the .coor (coordinates) and .vel (velocity) filenames .. py:attribute:: prefix_restart (default: ``'_namd'``) Prefix for files (.coor, .vel, and .xsc) from which a job is to restart .. py:attribute:: psf (default: ``''``) Name of PSF file to use. This does not need to be specified if the :py:meth:`NAMD.psfgen` method has been run beforehand .. py:attribute:: psfgen_options (default: ``{}``) Keyword option(s) for the PSFGEN interface (see :py:meth:`NAMD.psfgen` below). The PSFGEN keywords must be given as :py:class:`str`, for example, ``psfgen_options={'residue_aliases':{'CMO':'CO'},'topologies':['camphor.rtf']}`` .. py:attribute:: vel (default: ``'_namd.vel'``) Filename of the velocity file .. py:attribute:: xsc (default: ``'_namd.xsc'``) Filename of the extended system configuration (XSC) .. py:attribute:: xst (default: ``'_namd.xst'``) Filename of the extended system trajectory (XST) Forcefield options ================== .. py:attribute:: cutoff (default: ``22.6767`` Bohr) Cutoff distance for electrostatic and van der Waals interactions. Default equals to 12.0 Å. .. py:attribute:: exclude Choose a type to be excluded from non-bonded interactions. Accpeted values are: ``'none'``, ``'1-2'``, ``'1-3'``, ``'1-4'``, or ``'scaled1-4'`` .. py:attribute:: ff_dir (default: ``''``) Directory containing topology and parameter files. If left blank it will be chosen to be the built-in library of ChemShell's default :py:attr:`forcefield ` CHARMM. .. py:attribute:: freq_nonbonded (default: ``1``) Frequency to evaluate non-bonded interactions .. py:attribute:: freq_full_elect (default: ``1``) Frequence to evaluate electrostatic interactions .. py:attribute:: merge_cross (default: ``True``) Merge the crossterm (or CMAP) contribution into the dihedral energy .. py:attribute:: pairlist_dist (default: ``26.4562`` Bohr) Distance criterion to generate pairs of atoms with electrostatic and van der Waals interactions. Default equals to 14.0 Å. .. py:attribute:: scaling14 (default: ``1.0``) Scaling factor for 1-4 interactions .. py:attribute:: switching (default: ``True``) Apply smoothing functions to the electrostatic and van der Waals interactions .. py:attribute:: switch_dist (default: ``18.8973`` Bohr) Distance to start applying the :py:attr:`smooothing functions `. Default equals to 10.0 Å. MD options ========== .. py:attribute:: nsteps_per_cycle (default: ``1``) Number of timesteps between cycles of atom reassignments .. py:attribute:: seed (default: ``12345``) Seed number for the random number generator if :py:attr:`temperature` or :py:attr:`langevin` is selected Constrains ---------- .. py:attribute:: constraints (default: ``''``) Filename (in PDB format) of constraints definition .. py:attribute:: constraints_ref (default: ``''``) Filename (in PDB format) of a reference structure without :py:attr:`constraints` .. py:attribute:: fixed_atoms (default: ``''``) Filename (in PDB format) of a structure with fixed atoms (b-factors are 1.0) .. py:attribute:: fixed_atoms_forces (default: ``True``) Toggle computation of interactions between fixed atoms Pressure control ---------------- .. py:attribute:: constant_area (default: ``False``) Keep the ratio of the unit cell in the x-y plane constant while allowing fluctuations along all axes. See `NAMD Manual `__ .. py:attribute:: flexible_cell (default: ``False``) ``True`` allows the three orthogonal dimensions of the periodic cell to fluctuate independently. See `NAMD Manual `__ .. py:attribute:: group_pressure (default: ``False``) Calculate the pressure using hydrogen--group---based pseudo-molecular virial and kinetic energy. Only valied when the rigid bonds model (SHAKE) is used Periodic boundary conditions ---------------------------- .. py:attribute:: pme (default: ``True``) Toggles the particle mesh Ewald (PME) method .. py:attribute:: pme_grid_sizes (default: ``[0,0,0]``) Can be a :py:class:`list`, :py:class:`tuple`, or `NumPy array `_ of size 3. The values must be divisible by 2, 3, and 5 at the same time. If left undefined they will be determined by ChemShell. .. py:attribute:: wrap_all (default: ``True``) Output coordinates by wrapping all atoms within the periodic simulation box .. py:attribute:: wrap_water (default: ``True``) Output coordinates by wrapping water molecules within the periodic simulation box PSFGEN options ============== The PSFGEN plugin is available with NAMD. To use it in ChemShell, run the member function of a created ``NAMD`` instance:: my_namd.psfgen(my_biofrag, ...) .. py:method:: NAMD.psfgen(bioFrag, allow_undefined=False, atom_aliases={}, dryrun=False, ff='charmm', ff_dir='', guess_coord=True, haem='', ignore_undefined=True, input='_psfgen.inp', merge_solvents=False, no_auto=[], output='_psfgen.out', output_pdb='', output_psf='', output_dir='_psfgen', override_standard=False, patch=True, patch_disulphides=True, patch_titratables=True, path='', residue_aliases={}, solvents=[], topologies=[], use_newer_version=True) :param bioFrag: The biomolecular fragment (typically created from PDB/PQR/PSF) :type bioFrag: :py:class:`Fragment` :param allow_undefined: Tolerate with undefined residues :type allow_undefined: bool, optional :param atom_aliases: Atom names to be aliased :type atom_aliases: dict, optional :param dryrun: To generate the input script only without really executing PSFGEN :type dryrun: bool, optional :param ff: Forcefield scheme name (we currently only support CHARMM-type forcefield) :type ff: str, optional :param ff_dir: Path to forcefield data if the ChemShell built-in ones should be overridden :type ff_dir: str, optional :param guess_coord: Let PSFGEN guess positions of the missing atoms :type guess_coord: bool, optional :param haem: Residue name of the haem in case there is one and ChemShell cannot recognise it :type haem: str, optional :param input: Filename of the PSFGEN input script :type input: str, optional :param merge_solvents: Merge solvents to a single segment (Use with caution because it may cause the resIDs to exceed the 9999 limit!) :type merge_solvents: bool, optional :param no_auto: Name(s) of segments for which automatic angles/dihedrals assignment should not apply. They are solvents by default (``no_auto=[]``) :type no_auto: list, optional :param output: Name of the PSFGEN output file :type output: str, optional :param output_pdb: Name of the resulting PDB file. It is *_TITLE_psfgen.pdb* by default where *TITLE* is *bioFrag*'s title string :type output_pdb: str, optional :param output_psf: Name of the resulting PSF file. It is *_TITLE_psfgen.psf* by default where *TITLE* is *bioFrag*'s title string :type output_psf: str, optional :param output_dir: The subdirectory where the PSFGEN output files are generated :type output_dir: str, optional :param override_standard: Data from ChemShell's built-in library will override the standard forcefield data (currently CHARMM-type only) if ``override_standard=True`` :type override_standard: bool, optional :param patch: Patch residues :type patch: bool, optional :param patch_disulphides: Patch the disulphide bonds :type patch_disulphides: bool, optional :param patch_titratables: Patch the titratable amino acids :type patch_titratables: bool, optional :param patch: Apply residue patches via PSFGEN :type patch: bool, optional :param path: Path to the PSFGEN executable if the current ChemShell installation does not have NAMD :type path: str, optional :param solvent: Extra name(s) of solvent segments in addition to the ChemShell-recognisable names (currently water only such as 'CW', 'SOLA', 'SOLB',...) :type solvent: list, optional :param topologies: Extra topology file(s) (.rtf) :type topologies: str or list, optional :param use_newer_version: Use the one of the newer version if multiple :py:attr:`topologies` (.rtf) files containing definitions of same residues are chosen. For example, ChemShell may suggest *top_all22_prot.rtf* and *top_all36_prot.rtf* at the same time, but the latter will be used if ``use_newer_version=True`` :type use_newer_version: bool, optional