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,2017,2018,2019, 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] = { "Atom", "Nucleus", "Shell", "Bond", "VSite", nullptr };
55 void init_atom(t_atoms* at)
60 at->resinfo = nullptr;
61 at->atomname = nullptr;
62 at->atomtype = nullptr;
63 at->atomtypeB = nullptr;
64 at->pdbinfo = nullptr;
66 at->haveCharge = FALSE;
68 at->haveBState = FALSE;
69 at->havePdbInfo = FALSE;
72 void init_atomtypes(t_atomtypes* at)
75 at->atomnumber = nullptr;
78 void done_atom(t_atoms* at)
89 void done_and_delete_atoms(t_atoms* atoms)
95 void done_atomtypes(t_atomtypes* atype)
98 sfree(atype->atomnumber);
101 void add_t_atoms(t_atoms* atoms, int natom_extra, int nres_extra)
107 srenew(atoms->atomname, atoms->nr + natom_extra);
108 srenew(atoms->atom, atoms->nr + natom_extra);
109 if (nullptr != atoms->pdbinfo)
111 srenew(atoms->pdbinfo, atoms->nr + natom_extra);
113 if (nullptr != atoms->atomtype)
115 srenew(atoms->atomtype, atoms->nr + natom_extra);
117 if (nullptr != atoms->atomtypeB)
119 srenew(atoms->atomtypeB, atoms->nr + natom_extra);
121 for (i = atoms->nr; (i < atoms->nr + natom_extra); i++)
123 atoms->atomname[i] = nullptr;
124 memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
125 if (nullptr != atoms->pdbinfo)
127 std::memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
129 if (nullptr != atoms->atomtype)
131 atoms->atomtype[i] = nullptr;
133 if (nullptr != atoms->atomtypeB)
135 atoms->atomtypeB[i] = nullptr;
138 atoms->nr += natom_extra;
142 srenew(atoms->resinfo, atoms->nres + nres_extra);
143 for (i = atoms->nres; (i < atoms->nres + nres_extra); i++)
145 std::memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
147 atoms->nres += nres_extra;
151 void init_t_atoms(t_atoms* atoms, int natoms, gmx_bool bPdbinfo)
155 snew(atoms->atomname, natoms);
156 atoms->atomtype = nullptr;
157 atoms->atomtypeB = nullptr;
158 snew(atoms->resinfo, natoms);
159 snew(atoms->atom, natoms);
160 atoms->haveMass = FALSE;
161 atoms->haveCharge = FALSE;
162 atoms->haveType = FALSE;
163 atoms->haveBState = FALSE;
164 atoms->havePdbInfo = bPdbinfo;
165 if (atoms->havePdbInfo)
167 snew(atoms->pdbinfo, natoms);
171 atoms->pdbinfo = nullptr;
175 void gmx_pdbinfo_init_default(t_pdbinfo* pdbinfo)
177 pdbinfo->type = epdbATOM;
179 pdbinfo->altloc = ' ';
180 pdbinfo->atomnm[0] = '\0';
181 pdbinfo->occup = 1.0;
183 pdbinfo->bAnisotropic = FALSE;
184 std::fill(pdbinfo->uij, pdbinfo->uij + 6, 0.0);
187 t_atoms* copy_t_atoms(const t_atoms* src)
193 init_t_atoms(dst, src->nr, (nullptr != src->pdbinfo));
195 if (nullptr != src->atomname)
197 snew(dst->atomname, src->nr);
199 if (nullptr != src->atomtype)
201 snew(dst->atomtype, src->nr);
203 if (nullptr != src->atomtypeB)
205 snew(dst->atomtypeB, src->nr);
207 for (i = 0; (i < src->nr); i++)
209 dst->atom[i] = src->atom[i];
210 if (nullptr != src->pdbinfo)
212 dst->pdbinfo[i] = src->pdbinfo[i];
214 if (nullptr != src->atomname)
216 dst->atomname[i] = src->atomname[i];
218 if (nullptr != src->atomtype)
220 dst->atomtype[i] = src->atomtype[i];
222 if (nullptr != src->atomtypeB)
224 dst->atomtypeB[i] = src->atomtypeB[i];
227 dst->haveBState = src->haveBState;
228 dst->haveCharge = src->haveCharge;
229 dst->haveMass = src->haveMass;
230 dst->havePdbInfo = src->havePdbInfo;
231 dst->haveType = src->haveType;
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,
251 ri = &atoms->resinfo[atoms->atom[atom_ind].resind];
252 ri->name = put_symtab(symtab, resname);
256 ri->chainnum = chainnum;
257 ri->chainid = chainid;
260 static void pr_atom(FILE* fp, int indent, const char* title, const t_atom* atom, int n)
264 if (available(fp, atom, indent, title))
266 indent = pr_title_n(fp, indent, title, n);
267 for (i = 0; i < n; i++)
269 pr_indent(fp, indent);
271 "%s[%6d]={type=%3hu, typeB=%3hu, ptype=%8s, m=%12.5e, "
272 "q=%12.5e, mB=%12.5e, qB=%12.5e, resind=%5d, atomnumber=%3d}\n",
273 title, i, atom[i].type, atom[i].typeB, ptype_str[atom[i].ptype], atom[i].m,
274 atom[i].q, atom[i].mB, atom[i].qB, atom[i].resind, atom[i].atomnumber);
279 static void pr_strings2(FILE* fp, int indent, const char* title, char*** nm, char*** nmB, int n, gmx_bool bShowNumbers)
283 if (available(fp, nm, indent, title))
285 indent = pr_title_n(fp, indent, title, n);
286 for (i = 0; i < n; i++)
288 pr_indent(fp, indent);
289 fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n", title, bShowNumbers ? i : -1,
290 *(nm[i]), *(nmB[i]));
295 static void pr_resinfo(FILE* fp, int indent, const char* title, const t_resinfo* resinfo, int n, gmx_bool bShowNumbers)
299 if (available(fp, resinfo, indent, title))
301 indent = pr_title_n(fp, indent, title, n);
302 for (i = 0; i < n; i++)
304 pr_indent(fp, indent);
305 fprintf(fp, "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n", title, bShowNumbers ? i : -1,
306 *(resinfo[i].name), resinfo[i].nr, (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
311 void pr_atoms(FILE* fp, int indent, const char* title, const t_atoms* atoms, 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, gmx_bool bShowNumbers)
327 if (available(fp, atomtypes, indent, title))
329 indent = pr_title(fp, indent, title);
330 for (i = 0; i < atomtypes->nr; i++)
332 pr_indent(fp, indent);
333 fprintf(fp, "atomtype[%3d]={atomnumber=%4d}\n", bShowNumbers ? i : -1,
334 atomtypes->atomnumber[i]);
339 static void compareAtom(FILE* fp, int index, const t_atom* a1, const t_atom* a2, real relativeTolerance, real absoluteTolerance)
343 cmp_us(fp, "atom.type", index, a1->type, a2->type);
344 cmp_us(fp, "atom.ptype", index, a1->ptype, a2->ptype);
345 cmp_int(fp, "atom.resind", index, a1->resind, a2->resind);
346 cmp_int(fp, "atom.atomnumber", index, a1->atomnumber, a2->atomnumber);
347 cmp_real(fp, "atom.m", index, a1->m, a2->m, relativeTolerance, absoluteTolerance);
348 cmp_real(fp, "atom.q", index, a1->q, a2->q, relativeTolerance, absoluteTolerance);
349 cmp_us(fp, "atom.typeB", index, a1->typeB, a2->typeB);
350 cmp_real(fp, "atom.mB", index, a1->mB, a2->mB, relativeTolerance, absoluteTolerance);
351 cmp_real(fp, "atom.qB", index, a1->qB, a2->qB, relativeTolerance, absoluteTolerance);
352 cmp_str(fp, "elem", index, a1->elem, a2->elem);
356 cmp_us(fp, "atom.type", index, a1->type, a1->typeB);
357 cmp_real(fp, "atom.m", index, a1->m, a1->mB, relativeTolerance, absoluteTolerance);
358 cmp_real(fp, "atom.q", index, a1->q, a1->qB, relativeTolerance, absoluteTolerance);
362 static void compareResinfo(FILE* fp, int residue, const t_resinfo& r1, const t_resinfo& r2)
364 fprintf(fp, "comparing t_resinfo\n");
365 cmp_str(fp, "name", residue, *r1.name, *r2.name);
366 cmp_int(fp, "nr", residue, r1.nr, r2.nr);
367 cmp_uc(fp, "ic", residue, r1.ic, r2.ic);
368 cmp_int(fp, "chainnum", residue, r1.chainnum, r2.chainnum);
369 cmp_uc(fp, "chainid", residue, r1.chainid, r2.chainid);
370 if ((r1.rtp || r2.rtp) && (!r1.rtp || !r2.rtp))
372 fprintf(fp, "rtp info is present in topology %d but not in the other\n", r1.rtp ? 1 : 2);
374 if (r1.rtp && r2.rtp)
376 cmp_str(fp, "rtp", residue, *r1.rtp, *r2.rtp);
380 static void comparePdbinfo(FILE* fp,
382 const t_pdbinfo& pdb1,
383 const t_pdbinfo& pdb2,
384 real relativeTolerance,
385 real absoluteTolerance)
387 fprintf(fp, "comparing t_pdbinfo\n");
388 cmp_int(fp, "type", pdb, pdb1.type, pdb2.type);
389 cmp_int(fp, "atomnr", pdb, pdb1.atomnr, pdb2.atomnr);
390 cmp_uc(fp, "altloc", pdb, pdb1.altloc, pdb2.altloc);
391 cmp_str(fp, "atomnm", pdb, pdb1.atomnm, pdb2.atomnm);
392 cmp_real(fp, "occup", pdb, pdb1.occup, pdb2.occup, relativeTolerance, absoluteTolerance);
393 cmp_real(fp, "bfac", pdb, pdb1.bfac, pdb2.bfac, relativeTolerance, absoluteTolerance);
394 cmp_bool(fp, "bAnistropic", pdb, pdb1.bAnisotropic, pdb2.bAnisotropic);
395 for (int i = 0; i < 6; i++)
397 std::string buf = gmx::formatString("uij[%d]", i);
398 cmp_int(fp, buf.c_str(), pdb, pdb1.uij[i], pdb2.uij[i]);
403 void compareAtoms(FILE* fp, const t_atoms* a1, const t_atoms* a2, real relativeTolerance, real absoluteTolerance)
405 fprintf(fp, "comparing atoms\n");
409 cmp_int(fp, "atoms->nr", -1, a1->nr, a2->nr);
410 cmp_int(fp, "atoms->nres", -1, a1->nres, a2->nres);
411 cmp_bool(fp, "atoms->haveMass", -1, a1->haveMass, a2->haveMass);
412 cmp_bool(fp, "atoms->haveCharge", -1, a1->haveCharge, a2->haveCharge);
413 cmp_bool(fp, "atoms->haveType", -1, a1->haveType, a2->haveType);
414 cmp_bool(fp, "atoms->haveBState", -1, a1->haveBState, a2->haveBState);
415 cmp_bool(fp, "atoms->havePdbInfo", -1, a1->havePdbInfo, a2->havePdbInfo);
416 for (int i = 0; i < std::min(a1->nr, a2->nr); i++)
418 compareAtom(fp, i, &(a1->atom[i]), &(a2->atom[i]), relativeTolerance, absoluteTolerance);
419 if (a1->atomname && a2->atomname)
421 cmp_str(fp, "atomname", i, *a1->atomname[i], *a2->atomname[i]);
423 if (a1->havePdbInfo && a2->havePdbInfo)
425 comparePdbinfo(fp, i, a1->pdbinfo[i], a2->pdbinfo[i], relativeTolerance, absoluteTolerance);
427 if (a1->haveType && a2->haveType)
429 cmp_str(fp, "atomtype", i, *a1->atomtype[i], *a2->atomtype[i]);
431 if (a1->haveBState && a2->haveBState)
433 cmp_str(fp, "atomtypeB", i, *a1->atomtypeB[i], *a2->atomtypeB[i]);
436 for (int i = 0; i < std::min(a1->nres, a2->nres); i++)
438 compareResinfo(fp, i, a1->resinfo[i], a2->resinfo[i]);
443 for (int i = 0; (i < a1->nr); i++)
445 compareAtom(fp, i, &(a1->atom[i]), nullptr, relativeTolerance, absoluteTolerance);
450 void atomsSetMassesBasedOnNames(t_atoms* atoms, gmx_bool printMissingMasses)
454 /* We could decide to anyhow assign then or generate a fatal error,
455 * but it's probably most useful to keep the masses we have.
460 int maxWarn = (printMissingMasses ? 10 : 0);
465 bool haveMass = true;
466 for (int i = 0; i < atoms->nr; i++)
468 if (!aps.setAtomProperty(epropMass, *atoms->resinfo[atoms->atom[i].resind].name,
469 *atoms->atomname[i], &atoms->atom[i].m))
473 if (numWarn < maxWarn)
475 fprintf(stderr, "Can not find mass in database for atom %s in residue %d %s\n",
476 *atoms->atomname[i], atoms->resinfo[atoms->atom[i].resind].nr,
477 *atoms->resinfo[atoms->atom[i].resind].name);
486 atoms->haveMass = haveMass;