3615ed3eabeaf79957adc6f841f52f43de324931
[alexxy/gromacs.git] / src / gromacs / topology / atoms.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "atoms.h"
41
42 #include <cstdio>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/topology/atomprop.h"
48 #include "gromacs/topology/symtab.h"
49 #include "gromacs/utility/compare.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/utility/txtdump.h"
53
54 const char* ptype_str[eptNR + 1] = { "Atom", "Nucleus", "Shell", "Bond", "VSite", nullptr };
55
56 void init_atom(t_atoms* at)
57 {
58     at->nr          = 0;
59     at->nres        = 0;
60     at->atom        = nullptr;
61     at->resinfo     = nullptr;
62     at->atomname    = nullptr;
63     at->atomtype    = nullptr;
64     at->atomtypeB   = nullptr;
65     at->pdbinfo     = nullptr;
66     at->haveMass    = FALSE;
67     at->haveCharge  = FALSE;
68     at->haveType    = FALSE;
69     at->haveBState  = FALSE;
70     at->havePdbInfo = FALSE;
71 }
72
73 void init_atomtypes(t_atomtypes* at)
74 {
75     at->nr         = 0;
76     at->atomnumber = nullptr;
77 }
78
79 void done_atom(t_atoms* at)
80 {
81     sfree(at->atom);
82     sfree(at->resinfo);
83     sfree(at->atomname);
84     sfree(at->atomtype);
85     sfree(at->atomtypeB);
86     sfree(at->pdbinfo);
87     init_atom(at);
88 }
89
90 void done_and_delete_atoms(t_atoms* atoms)
91 {
92     done_atom(atoms);
93     delete atoms;
94 }
95
96 void done_atomtypes(t_atomtypes* atype)
97 {
98     atype->nr = 0;
99     sfree(atype->atomnumber);
100 }
101
102 void add_t_atoms(t_atoms* atoms, int natom_extra, int nres_extra)
103 {
104     if (natom_extra > 0)
105     {
106         srenew(atoms->atomname, atoms->nr + natom_extra);
107         srenew(atoms->atom, atoms->nr + natom_extra);
108         if (nullptr != atoms->pdbinfo)
109         {
110             srenew(atoms->pdbinfo, atoms->nr + natom_extra);
111         }
112         if (nullptr != atoms->atomtype)
113         {
114             srenew(atoms->atomtype, atoms->nr + natom_extra);
115         }
116         if (nullptr != atoms->atomtypeB)
117         {
118             srenew(atoms->atomtypeB, atoms->nr + natom_extra);
119         }
120         for (int i = atoms->nr; (i < atoms->nr + natom_extra); i++)
121         {
122             atoms->atomname[i] = nullptr;
123             memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
124             if (nullptr != atoms->pdbinfo)
125             {
126                 std::memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
127             }
128             if (nullptr != atoms->atomtype)
129             {
130                 atoms->atomtype[i] = nullptr;
131             }
132             if (nullptr != atoms->atomtypeB)
133             {
134                 atoms->atomtypeB[i] = nullptr;
135             }
136         }
137         atoms->nr += natom_extra;
138     }
139     if (nres_extra > 0)
140     {
141         srenew(atoms->resinfo, atoms->nres + nres_extra);
142         for (int i = atoms->nres; (i < atoms->nres + nres_extra); i++)
143         {
144             std::memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
145         }
146         atoms->nres += nres_extra;
147     }
148 }
149
150 void init_t_atoms(t_atoms* atoms, int natoms, gmx_bool bPdbinfo)
151 {
152     atoms->nr   = natoms;
153     atoms->nres = 0;
154     snew(atoms->atomname, natoms);
155     atoms->atomtype  = nullptr;
156     atoms->atomtypeB = nullptr;
157     snew(atoms->resinfo, natoms);
158     snew(atoms->atom, natoms);
159     atoms->haveMass    = FALSE;
160     atoms->haveCharge  = FALSE;
161     atoms->haveType    = FALSE;
162     atoms->haveBState  = FALSE;
163     atoms->havePdbInfo = bPdbinfo;
164     if (atoms->havePdbInfo)
165     {
166         snew(atoms->pdbinfo, natoms);
167     }
168     else
169     {
170         atoms->pdbinfo = nullptr;
171     }
172 }
173
174 void gmx_pdbinfo_init_default(t_pdbinfo* pdbinfo)
175 {
176     pdbinfo->type         = epdbATOM;
177     pdbinfo->atomnr       = 0;
178     pdbinfo->altloc       = ' ';
179     pdbinfo->atomnm[0]    = '\0';
180     pdbinfo->occup        = 1.0;
181     pdbinfo->bfac         = 0.0;
182     pdbinfo->bAnisotropic = FALSE;
183     std::fill(pdbinfo->uij, pdbinfo->uij + 6, 0.0);
184 }
185
186 t_atoms* copy_t_atoms(const t_atoms* src)
187 {
188     t_atoms* dst = nullptr;
189
190     snew(dst, 1);
191     init_t_atoms(dst, src->nr, (nullptr != src->pdbinfo));
192     dst->nr = src->nr;
193     if (nullptr != src->atomname)
194     {
195         snew(dst->atomname, src->nr);
196     }
197     if (nullptr != src->atomtype)
198     {
199         snew(dst->atomtype, src->nr);
200     }
201     if (nullptr != src->atomtypeB)
202     {
203         snew(dst->atomtypeB, src->nr);
204     }
205     for (int i = 0; (i < src->nr); i++)
206     {
207         dst->atom[i] = src->atom[i];
208         if (nullptr != src->pdbinfo)
209         {
210             dst->pdbinfo[i] = src->pdbinfo[i];
211         }
212         if (nullptr != src->atomname)
213         {
214             dst->atomname[i] = src->atomname[i];
215         }
216         if (nullptr != src->atomtype)
217         {
218             dst->atomtype[i] = src->atomtype[i];
219         }
220         if (nullptr != src->atomtypeB)
221         {
222             dst->atomtypeB[i] = src->atomtypeB[i];
223         }
224     }
225     dst->haveBState  = src->haveBState;
226     dst->haveCharge  = src->haveCharge;
227     dst->haveMass    = src->haveMass;
228     dst->havePdbInfo = src->havePdbInfo;
229     dst->haveType    = src->haveType;
230     dst->nres        = src->nres;
231     for (int i = 0; (i < src->nres); i++)
232     {
233         dst->resinfo[i] = src->resinfo[i];
234     }
235     return dst;
236 }
237
238 void t_atoms_set_resinfo(t_atoms*      atoms,
239                          int           atom_ind,
240                          t_symtab*     symtab,
241                          const char*   resname,
242                          int           resnr,
243                          unsigned char ic,
244                          int           chainnum,
245                          char          chainid)
246 {
247     t_resinfo* ri = &atoms->resinfo[atoms->atom[atom_ind].resind];
248     ri->name      = put_symtab(symtab, resname);
249     ri->rtp       = nullptr;
250     ri->nr        = resnr;
251     ri->ic        = ic;
252     ri->chainnum  = chainnum;
253     ri->chainid   = chainid;
254 }
255
256 static void pr_atom(FILE* fp, int indent, const char* title, const t_atom* atom, int n)
257 {
258     if (available(fp, atom, indent, title))
259     {
260         indent = pr_title_n(fp, indent, title, n);
261         for (int i = 0; i < n; i++)
262         {
263             pr_indent(fp, indent);
264             fprintf(fp,
265                     "%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,
268                     i,
269                     atom[i].type,
270                     atom[i].typeB,
271                     ptype_str[atom[i].ptype],
272                     atom[i].m,
273                     atom[i].q,
274                     atom[i].mB,
275                     atom[i].qB,
276                     atom[i].resind,
277                     atom[i].atomnumber);
278         }
279     }
280 }
281
282 static void pr_strings2(FILE* fp, int indent, const char* title, char*** nm, char*** nmB, int n, gmx_bool bShowNumbers)
283 {
284     if (available(fp, nm, indent, title))
285     {
286         indent = pr_title_n(fp, indent, title, n);
287         for (int i = 0; i < n; i++)
288         {
289             pr_indent(fp, indent);
290             fprintf(fp, "%s[%d]={name=\"%s\",nameB=\"%s\"}\n", title, bShowNumbers ? i : -1, *(nm[i]), *(nmB[i]));
291         }
292     }
293 }
294
295 static void pr_resinfo(FILE* fp, int indent, const char* title, const t_resinfo* resinfo, int n, gmx_bool bShowNumbers)
296 {
297     if (available(fp, resinfo, indent, title))
298     {
299         indent = pr_title_n(fp, indent, title, n);
300         for (int i = 0; i < n; i++)
301         {
302             pr_indent(fp, indent);
303             fprintf(fp,
304                     "%s[%d]={name=\"%s\", nr=%d, ic='%c'}\n",
305                     title,
306                     bShowNumbers ? i : -1,
307                     *(resinfo[i].name),
308                     resinfo[i].nr,
309                     (resinfo[i].ic == '\0') ? ' ' : resinfo[i].ic);
310         }
311     }
312 }
313
314 void pr_atoms(FILE* fp, int indent, const char* title, const t_atoms* atoms, gmx_bool bShownumbers)
315 {
316     if (available(fp, atoms, indent, title))
317     {
318         indent = pr_title(fp, indent, title);
319         pr_atom(fp, indent, "atom", atoms->atom, atoms->nr);
320         pr_strings(fp, indent, "atom", atoms->atomname, atoms->nr, bShownumbers);
321         pr_strings2(fp, indent, "type", atoms->atomtype, atoms->atomtypeB, atoms->nr, bShownumbers);
322         pr_resinfo(fp, indent, "residue", atoms->resinfo, atoms->nres, bShownumbers);
323     }
324 }
325
326
327 void pr_atomtypes(FILE* fp, int indent, const char* title, const t_atomtypes* atomtypes, gmx_bool bShowNumbers)
328 {
329     if (available(fp, atomtypes, indent, title))
330     {
331         indent = pr_title(fp, indent, title);
332         for (int i = 0; i < atomtypes->nr; i++)
333         {
334             pr_indent(fp, indent);
335             fprintf(fp, "atomtype[%3d]={atomnumber=%4d}\n", bShowNumbers ? i : -1, atomtypes->atomnumber[i]);
336         }
337     }
338 }
339
340 static void compareAtom(FILE* fp, int index, const t_atom* a1, const t_atom* a2, real relativeTolerance, real absoluteTolerance)
341 {
342     if (a2)
343     {
344         cmp_us(fp, "atom.type", index, a1->type, a2->type);
345         cmp_us(fp, "atom.ptype", index, a1->ptype, a2->ptype);
346         cmp_int(fp, "atom.resind", index, a1->resind, a2->resind);
347         cmp_int(fp, "atom.atomnumber", index, a1->atomnumber, a2->atomnumber);
348         cmp_real(fp, "atom.m", index, a1->m, a2->m, relativeTolerance, absoluteTolerance);
349         cmp_real(fp, "atom.q", index, a1->q, a2->q, relativeTolerance, absoluteTolerance);
350         cmp_us(fp, "atom.typeB", index, a1->typeB, a2->typeB);
351         cmp_real(fp, "atom.mB", index, a1->mB, a2->mB, relativeTolerance, absoluteTolerance);
352         cmp_real(fp, "atom.qB", index, a1->qB, a2->qB, relativeTolerance, absoluteTolerance);
353         cmp_str(fp, "elem", index, a1->elem, a2->elem);
354     }
355     else
356     {
357         cmp_us(fp, "atom.type", index, a1->type, a1->typeB);
358         cmp_real(fp, "atom.m", index, a1->m, a1->mB, relativeTolerance, absoluteTolerance);
359         cmp_real(fp, "atom.q", index, a1->q, a1->qB, relativeTolerance, absoluteTolerance);
360     }
361 }
362
363 static void compareResinfo(FILE* fp, int residue, const t_resinfo& r1, const t_resinfo& r2)
364 {
365     fprintf(fp, "comparing t_resinfo\n");
366     cmp_str(fp, "name", residue, *r1.name, *r2.name);
367     cmp_int(fp, "nr", residue, r1.nr, r2.nr);
368     cmp_uc(fp, "ic", residue, r1.ic, r2.ic);
369     cmp_int(fp, "chainnum", residue, r1.chainnum, r2.chainnum);
370     cmp_uc(fp, "chainid", residue, r1.chainid, r2.chainid);
371     if ((r1.rtp || r2.rtp) && (!r1.rtp || !r2.rtp))
372     {
373         fprintf(fp, "rtp info is present in topology %d but not in the other\n", r1.rtp ? 1 : 2);
374     }
375     if (r1.rtp && r2.rtp)
376     {
377         cmp_str(fp, "rtp", residue, *r1.rtp, *r2.rtp);
378     }
379 }
380
381 static void comparePdbinfo(FILE*            fp,
382                            int              pdb,
383                            const t_pdbinfo& pdb1,
384                            const t_pdbinfo& pdb2,
385                            real             relativeTolerance,
386                            real             absoluteTolerance)
387 {
388     fprintf(fp, "comparing t_pdbinfo\n");
389     cmp_int(fp, "type", pdb, pdb1.type, pdb2.type);
390     cmp_int(fp, "atomnr", pdb, pdb1.atomnr, pdb2.atomnr);
391     cmp_uc(fp, "altloc", pdb, pdb1.altloc, pdb2.altloc);
392     cmp_str(fp, "atomnm", pdb, pdb1.atomnm, pdb2.atomnm);
393     cmp_real(fp, "occup", pdb, pdb1.occup, pdb2.occup, relativeTolerance, absoluteTolerance);
394     cmp_real(fp, "bfac", pdb, pdb1.bfac, pdb2.bfac, relativeTolerance, absoluteTolerance);
395     cmp_bool(fp, "bAnistropic", pdb, pdb1.bAnisotropic, pdb2.bAnisotropic);
396     for (int i = 0; i < 6; i++)
397     {
398         std::string buf = gmx::formatString("uij[%d]", i);
399         cmp_int(fp, buf.c_str(), pdb, pdb1.uij[i], pdb2.uij[i]);
400     }
401 }
402
403
404 void compareAtoms(FILE* fp, const t_atoms* a1, const t_atoms* a2, real relativeTolerance, real absoluteTolerance)
405 {
406     fprintf(fp, "comparing atoms\n");
407
408     if (a2)
409     {
410         cmp_int(fp, "atoms->nr", -1, a1->nr, a2->nr);
411         cmp_int(fp, "atoms->nres", -1, a1->nres, a2->nres);
412         cmp_bool(fp, "atoms->haveMass", -1, a1->haveMass, a2->haveMass);
413         cmp_bool(fp, "atoms->haveCharge", -1, a1->haveCharge, a2->haveCharge);
414         cmp_bool(fp, "atoms->haveType", -1, a1->haveType, a2->haveType);
415         cmp_bool(fp, "atoms->haveBState", -1, a1->haveBState, a2->haveBState);
416         cmp_bool(fp, "atoms->havePdbInfo", -1, a1->havePdbInfo, a2->havePdbInfo);
417         for (int i = 0; i < std::min(a1->nr, a2->nr); i++)
418         {
419             compareAtom(fp, i, &(a1->atom[i]), &(a2->atom[i]), relativeTolerance, absoluteTolerance);
420             if (a1->atomname && a2->atomname)
421             {
422                 cmp_str(fp, "atomname", i, *a1->atomname[i], *a2->atomname[i]);
423             }
424             if (a1->havePdbInfo && a2->havePdbInfo)
425             {
426                 comparePdbinfo(fp, i, a1->pdbinfo[i], a2->pdbinfo[i], relativeTolerance, absoluteTolerance);
427             }
428             if (a1->haveType && a2->haveType)
429             {
430                 cmp_str(fp, "atomtype", i, *a1->atomtype[i], *a2->atomtype[i]);
431             }
432             if (a1->haveBState && a2->haveBState)
433             {
434                 cmp_str(fp, "atomtypeB", i, *a1->atomtypeB[i], *a2->atomtypeB[i]);
435             }
436         }
437         for (int i = 0; i < std::min(a1->nres, a2->nres); i++)
438         {
439             compareResinfo(fp, i, a1->resinfo[i], a2->resinfo[i]);
440         }
441     }
442     else
443     {
444         for (int i = 0; (i < a1->nr); i++)
445         {
446             compareAtom(fp, i, &(a1->atom[i]), nullptr, relativeTolerance, absoluteTolerance);
447         }
448     }
449 }
450
451 void atomsSetMassesBasedOnNames(t_atoms* atoms, gmx_bool printMissingMasses)
452 {
453     if (atoms->haveMass)
454     {
455         /* We could decide to anyhow assign then or generate a fatal error,
456          * but it's probably most useful to keep the masses we have.
457          */
458         return;
459     }
460
461     int maxWarn = (printMissingMasses ? 10 : 0);
462     int numWarn = 0;
463
464     AtomProperties aps;
465
466     bool haveMass = true;
467     for (int i = 0; i < atoms->nr; i++)
468     {
469         if (!aps.setAtomProperty(epropMass,
470                                  *atoms->resinfo[atoms->atom[i].resind].name,
471                                  *atoms->atomname[i],
472                                  &atoms->atom[i].m))
473         {
474             haveMass = false;
475
476             if (numWarn < maxWarn)
477             {
478                 fprintf(stderr,
479                         "Can not find mass in database for atom %s in residue %d %s\n",
480                         *atoms->atomname[i],
481                         atoms->resinfo[atoms->atom[i].resind].nr,
482                         *atoms->resinfo[atoms->atom[i].resind].name);
483                 numWarn++;
484             }
485             else
486             {
487                 break;
488             }
489         }
490     }
491     atoms->haveMass = haveMass;
492 }