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) 2013,2014,2015,2016, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/topology/atomprop.h"
47 #include "gromacs/topology/symtab.h"
48 #include "gromacs/utility/compare.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/utility/txtdump.h"
53 const char *ptype_str[eptNR+1] = {
54 "Atom", "Nucleus", "Shell", "Bond", "VSite", NULL
57 void init_atom(t_atoms *at)
68 at->haveCharge = FALSE;
70 at->haveBState = FALSE;
71 at->havePdbInfo = FALSE;
74 void init_atomtypes(t_atomtypes *at)
79 at->atomnumber = NULL;
84 void done_atom(t_atoms *at)
95 void done_atomtypes(t_atomtypes *atype)
100 sfree(atype->surftens);
101 sfree(atype->atomnumber);
102 sfree(atype->gb_radius);
106 void add_t_atoms(t_atoms *atoms, int natom_extra, int nres_extra)
112 srenew(atoms->atomname, atoms->nr+natom_extra);
113 srenew(atoms->atom, atoms->nr+natom_extra);
114 if (NULL != atoms->pdbinfo)
116 srenew(atoms->pdbinfo, atoms->nr+natom_extra);
118 if (NULL != atoms->atomtype)
120 srenew(atoms->atomtype, atoms->nr+natom_extra);
122 if (NULL != atoms->atomtypeB)
124 srenew(atoms->atomtypeB, atoms->nr+natom_extra);
126 for (i = atoms->nr; (i < atoms->nr+natom_extra); i++)
128 atoms->atomname[i] = NULL;
129 memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
130 if (NULL != atoms->pdbinfo)
132 std::memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
134 if (NULL != atoms->atomtype)
136 atoms->atomtype[i] = NULL;
138 if (NULL != atoms->atomtypeB)
140 atoms->atomtypeB[i] = NULL;
143 atoms->nr += natom_extra;
147 srenew(atoms->resinfo, atoms->nres+nres_extra);
148 for (i = atoms->nres; (i < atoms->nres+nres_extra); i++)
150 std::memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
152 atoms->nres += nres_extra;
156 void init_t_atoms(t_atoms *atoms, int natoms, gmx_bool bPdbinfo)
160 snew(atoms->atomname, natoms);
161 atoms->atomtype = NULL;
162 atoms->atomtypeB = NULL;
163 snew(atoms->resinfo, natoms);
164 snew(atoms->atom, natoms);
165 atoms->haveMass = FALSE;
166 atoms->haveCharge = FALSE;
167 atoms->haveType = FALSE;
168 atoms->haveBState = FALSE;
169 atoms->havePdbInfo = bPdbinfo;
170 if (atoms->havePdbInfo)
172 snew(atoms->pdbinfo, natoms);
176 atoms->pdbinfo = NULL;
180 void gmx_pdbinfo_init_default(t_pdbinfo *pdbinfo)
182 pdbinfo->type = epdbATOM;
184 pdbinfo->altloc = ' ';
185 pdbinfo->atomnm[0] = '\0';
186 pdbinfo->occup = 1.0;
188 pdbinfo->bAnisotropic = FALSE;
189 std::fill(pdbinfo->uij, pdbinfo->uij+6, 0.0);
192 t_atoms *copy_t_atoms(const t_atoms *src)
198 init_t_atoms(dst, src->nr, (NULL != src->pdbinfo));
200 if (NULL != src->atomname)
202 snew(dst->atomname, src->nr);
204 if (NULL != src->atomtype)
206 snew(dst->atomtype, src->nr);
208 if (NULL != src->atomtypeB)
210 snew(dst->atomtypeB, src->nr);
212 for (i = 0; (i < src->nr); i++)
214 dst->atom[i] = src->atom[i];
215 if (NULL != src->pdbinfo)
217 dst->pdbinfo[i] = src->pdbinfo[i];
219 if (NULL != src->atomname)
221 dst->atomname[i] = src->atomname[i];
223 if (NULL != src->atomtype)
225 dst->atomtype[i] = src->atomtype[i];
227 if (NULL != src->atomtypeB)
229 dst->atomtypeB[i] = src->atomtypeB[i];
232 dst->nres = src->nres;
233 for (i = 0; (i < src->nres); i++)
235 dst->resinfo[i] = src->resinfo[i];
240 void t_atoms_set_resinfo(t_atoms *atoms, int atom_ind, t_symtab *symtab,
241 const char *resname, int resnr, unsigned char ic,
242 int chainnum, char chainid)
246 ri = &atoms->resinfo[atoms->atom[atom_ind].resind];
247 ri->name = put_symtab(symtab, resname);
251 ri->chainnum = chainnum;
252 ri->chainid = chainid;
255 static void pr_atom(FILE *fp, int indent, const char *title, const t_atom *atom, int n)
259 if (available(fp, atom, indent, title))
261 indent = pr_title_n(fp, indent, title, n);
262 for (i = 0; i < n; i++)
264 pr_indent(fp, indent);
265 fprintf(fp, "%s[%6d]={type=%3hu, typeB=%3hu, ptype=%8s, m=%12.5e, "
266 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
267 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype],
268 atom[i].m, atom[i].q, atom[i].mB, atom[i].qB,
269 atom[i].resind, atom[i].atomnumber);
274 static void pr_strings2(FILE *fp, int indent, const char *title,
275 char ***nm, char ***nmB, int n, gmx_bool bShowNumbers)
279 if (available(fp, nm, indent, title))
281 indent = pr_title_n(fp, indent, title, n);
282 for (i = 0; i < n; i++)
284 pr_indent(fp, indent);
285 fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n",
286 title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
291 static void pr_resinfo(FILE *fp, int indent, const char *title, const t_resinfo *resinfo, int n,
292 gmx_bool bShowNumbers)
296 if (available(fp, resinfo, indent, title))
298 indent = pr_title_n(fp, indent, title, n);
299 for (i = 0; i < n; i++)
301 pr_indent(fp, indent);
302 fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
303 title, bShowNumbers ? i : -1,
304 *(resinfo[i].name), resinfo[i].nr,
305 (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
310 void pr_atoms(FILE *fp, int indent, const char *title, const t_atoms *atoms,
311 gmx_bool bShownumbers)
313 if (available(fp, atoms, indent, title))
315 indent = pr_title(fp, indent, title);
316 pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
317 pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
318 pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
319 pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
324 void pr_atomtypes(FILE *fp, int indent, const char *title, const t_atomtypes *atomtypes,
325 gmx_bool bShowNumbers)
328 if (available(fp, atomtypes, indent, title))
330 indent = pr_title(fp, indent, title);
331 for (i = 0; i < atomtypes->nr; i++)
333 pr_indent(fp, indent);
335 "atomtype[%3d]={radius=%12.5e, volume=%12.5e, gb_radius=%12.5e, surftens=%12.5e, atomnumber=%4d, S_hct=%12.5e)}\n",
336 bShowNumbers ? i : -1, atomtypes->radius[i], atomtypes->vol[i],
337 atomtypes->gb_radius[i],
338 atomtypes->surftens[i], atomtypes->atomnumber[i], atomtypes->S_hct[i]);
343 static void cmp_atom(FILE *fp, int index, const t_atom *a1, const t_atom *a2, real ftol, real abstol)
347 cmp_us(fp, "atom.type", index, a1->type, a2->type);
348 cmp_us(fp, "atom.ptype", index, a1->ptype, a2->ptype);
349 cmp_int(fp, "atom.resind", index, a1->resind, a2->resind);
350 cmp_int(fp, "atom.atomnumber", index, a1->atomnumber, a2->atomnumber);
351 cmp_real(fp, "atom.m", index, a1->m, a2->m, ftol, abstol);
352 cmp_real(fp, "atom.q", index, a1->q, a2->q, ftol, abstol);
353 cmp_us(fp, "atom.typeB", index, a1->typeB, a2->typeB);
354 cmp_real(fp, "atom.mB", index, a1->mB, a2->mB, ftol, abstol);
355 cmp_real(fp, "atom.qB", index, a1->qB, a2->qB, ftol, abstol);
359 cmp_us(fp, "atom.type", index, a1->type, a1->typeB);
360 cmp_real(fp, "atom.m", index, a1->m, a1->mB, ftol, abstol);
361 cmp_real(fp, "atom.q", index, a1->q, a1->qB, ftol, abstol);
365 void cmp_atoms(FILE *fp, const t_atoms *a1, const t_atoms *a2, real ftol, real abstol)
369 fprintf(fp, "comparing atoms\n");
373 cmp_int(fp, "atoms->nr", -1, a1->nr, a2->nr);
374 for (i = 0; (i < a1->nr); i++)
376 cmp_atom(fp, i, &(a1->atom[i]), &(a2->atom[i]), ftol, abstol);
381 for (i = 0; (i < a1->nr); i++)
383 cmp_atom(fp, i, &(a1->atom[i]), NULL, ftol, abstol);
388 void atomsSetMassesBasedOnNames(t_atoms *atoms, gmx_bool printMissingMasses)
392 /* We could decide to anyhow assign then or generate a fatal error,
393 * but it's probably most useful to keep the masses we have.
398 int maxWarn = (printMissingMasses ? 10 : 0);
401 gmx_atomprop_t aps = gmx_atomprop_init();
403 gmx_bool haveMass = TRUE;
404 for (int i = 0; i < atoms->nr; i++)
406 if (!gmx_atomprop_query(aps, epropMass,
407 *atoms->resinfo[atoms->atom[i].resind].name,
413 if (numWarn < maxWarn)
415 fprintf(stderr, "Can not find mass in database for atom %s in residue %d %s\n",
417 atoms->resinfo[atoms->atom[i].resind].nr,
418 *atoms->resinfo[atoms->atom[i].resind].name);
427 atoms->haveMass = haveMass;
429 gmx_atomprop_destroy(aps);