2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,2014,2015,2016,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_TOPOLOGY_ATOMS_H
39 #define GMX_TOPOLOGY_ATOMS_H
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
47 #include "gromacs/utility/unique_cptr.h"
51 /* The particle type */
62 /* The particle type names */
63 extern const char* ptype_str[eptNR + 1];
65 /* Enumerated type for pdb records. The other entries are ignored
66 * when reading a pdb file
87 real m, q; /* Mass and charge */
88 real mB, qB; /* Mass and charge for Free Energy calc */
89 unsigned short type; /* Atom type */
90 unsigned short typeB; /* Atom type for Free Energy calc */
91 int ptype; /* Particle type */
92 int resind; /* Index into resinfo (in t_atoms) */
93 int atomnumber; /* Atomic Number or 0 */
94 char elem[4]; /* Element name */
97 typedef struct t_resinfo
99 char** name; /* Pointer to the residue name */
100 int nr; /* Residue number */
101 unsigned char ic; /* Code for insertion of residues */
102 int chainnum; /* Iincremented at TER or new chain id */
103 char chainid; /* Chain identifier written/read to pdb */
104 char** rtp; /* rtp building block name (optional) */
107 typedef struct t_pdbinfo
109 int type; /* PDB record name */
110 int atomnr; /* PDB atom number */
111 char altloc; /* Alternate location indicator */
112 char atomnm[6]; /* True atom name including leading spaces */
113 real occup; /* Occupancy */
114 real bfac; /* B-factor */
115 gmx_bool bAnisotropic; /* (an)isotropic switch */
116 int uij[6]; /* Anisotropic B-factor */
119 //! Contains indices into group names for different groups.
120 using AtomGroupIndices = std::vector<int>;
122 typedef struct t_atoms
124 int nr; /* Nr of atoms */
125 t_atom* atom; /* Array of atoms (dim: nr) */
126 /* The following entries will not */
127 /* always be used (nres==0) */
128 char*** atomname; /* Array of pointers to atom name */
129 /* use: (*(atomname[i])) */
130 char*** atomtype; /* Array of pointers to atom types */
131 /* use: (*(atomtype[i])) */
132 char*** atomtypeB; /* Array of pointers to B atom types */
133 /* use: (*(atomtypeB[i])) */
134 int nres; /* The number of resinfo entries */
135 t_resinfo* resinfo; /* Array of residue names and numbers */
136 t_pdbinfo* pdbinfo; /* PDB Information, such as aniso. Bfac */
138 /* Flags that tell if properties are set for all nr atoms.
139 * For B-state parameters, both haveBState and the mass/charge/type
140 * flag should be TRUE.
142 gmx_bool haveMass; /* Mass available */
143 gmx_bool haveCharge; /* Charge available */
144 gmx_bool haveType; /* Atom type available */
145 gmx_bool haveBState; /* B-state parameters available */
146 gmx_bool havePdbInfo; /* pdbinfo available */
149 typedef struct t_atomtypes
151 int nr; /* number of atomtypes */
152 int* atomnumber; /* Atomic number, used for QM/MM */
155 #define PERTURBED(a) (((a).mB != (a).m) || ((a).qB != (a).q) || ((a).typeB != (a).type))
157 void init_atom(t_atoms* at);
158 void init_atomtypes(t_atomtypes* at);
159 void done_atom(t_atoms* at);
160 void done_and_delete_atoms(t_atoms* atoms);
161 void done_atomtypes(t_atomtypes* at);
163 void init_t_atoms(t_atoms* atoms, int natoms, gmx_bool bPdbinfo);
164 /* allocate memory for the arrays, set nr to natoms and nres to 0
165 * set pdbinfo to NULL or allocate memory for it */
167 void gmx_pdbinfo_init_default(t_pdbinfo* pdbinfo);
169 t_atoms* copy_t_atoms(const t_atoms* src);
170 /* copy an atoms struct from src to a new one */
172 void add_t_atoms(t_atoms* atoms, int natom_extra, int nres_extra);
173 /* allocate extra space for more atoms and or residues */
175 void t_atoms_set_resinfo(t_atoms* atoms,
177 struct t_symtab* symtab,
183 /* Set the residue name, number, insertion code and chain identifier
184 * of atom index atom_ind.
187 void pr_atoms(FILE* fp, int indent, const char* title, const t_atoms* atoms, gmx_bool bShownumbers);
188 void pr_atomtypes(FILE* fp, int indent, const char* title, const t_atomtypes* atomtypes, gmx_bool bShowNumbers);
190 /*! \brief Compare information in the t_atoms data structure.
192 * \param[in] fp Pointer to file to write to.
193 * \param[in] a1 Pointer to first data structure to compare.
194 * \param[in] a2 Pointer to second data structure or nullptr.
195 * \param[in] relativeTolerance Relative floating point comparison tolerance.
196 * \param[in] absoluteTolerance Absolute floating point comparison tolerance.
198 void compareAtoms(FILE* fp, const t_atoms* a1, const t_atoms* a2, real relativeTolerance, real absoluteTolerance);
200 /*! \brief Set mass for each atom using the atom and residue names using a database
202 * If atoms->haveMass = TRUE does nothing.
203 * If printMissingMasss = TRUE, prints details for first 10 missing masses
206 void atomsSetMassesBasedOnNames(t_atoms* atoms, gmx_bool printMissingMasses);
208 //! Deleter for t_atoms, needed until it has a proper destructor.
209 using AtomsDataPtr = gmx::unique_cptr<t_atoms, done_and_delete_atoms>;