Rename all source files from - to _ for consistency.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / toputil.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) 2012,2014,2015,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "toputil.h"
40
41 #include <climits>
42 #include <cmath>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/gmxpreprocess/notset.h"
50 #include "gromacs/gmxpreprocess/topdirs.h"
51 #include "gromacs/topology/block.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/topology/symtab.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/smalloc.h"
57
58 /* UTILITIES */
59
60 void set_p_string(t_param *p, const char *s)
61 {
62     if (s)
63     {
64         if (strlen(s) < sizeof(p->s)-1)
65         {
66             strncpy(p->s, s, sizeof(p->s));
67         }
68         else
69         {
70             gmx_fatal(FARGS, "Increase MAXSLEN in the grompp code to at least %zu,"
71                       " or shorten your definition of bonds like %s to at most %d",
72                       strlen(s)+1, s, MAXSLEN-1);
73         }
74     }
75     else
76     {
77         strcpy(p->s, "");
78     }
79 }
80
81 void pr_alloc (int extra, t_params *pr)
82 {
83     int i, j;
84
85     /* get new space for arrays */
86     if (extra < 0)
87     {
88         gmx_fatal(FARGS, "Trying to make array smaller.\n");
89     }
90     if (extra == 0)
91     {
92         return;
93     }
94     GMX_ASSERT(pr->nr != 0 || pr->param == nullptr, "Invalid t_params object");
95     if (pr->nr+extra > pr->maxnr)
96     {
97         pr->maxnr = std::max(static_cast<int>(1.2*pr->maxnr), pr->maxnr + extra);
98         srenew(pr->param, pr->maxnr);
99         for (i = pr->nr; (i < pr->maxnr); i++)
100         {
101             for (j = 0; (j < MAXATOMLIST); j++)
102             {
103                 pr->param[i].a[j] = 0;
104             }
105             for (j = 0; (j < MAXFORCEPARAM); j++)
106             {
107                 pr->param[i].c[j] = 0;
108             }
109             set_p_string(&(pr->param[i]), "");
110         }
111     }
112 }
113
114 void init_plist(t_params plist[])
115 {
116     int i;
117
118     for (i = 0; (i < F_NRE); i++)
119     {
120         plist[i].nr    = 0;
121         plist[i].maxnr = 0;
122         plist[i].param = nullptr;
123
124         /* CMAP */
125         plist[i].ncmap        = 0;
126         plist[i].cmap         = nullptr;
127         plist[i].grid_spacing = 0;
128         plist[i].nc           = 0;
129         plist[i].nct          = 0;
130         plist[i].cmap_types   = nullptr;
131     }
132 }
133
134 void done_plist(t_params *plist)
135 {
136     for (int i = 0; i < F_NRE; i++)
137     {
138         t_params *pl = &plist[i];
139         sfree(pl->param);
140         sfree(pl->cmap);
141         sfree(pl->cmap_types);
142     }
143 }
144
145 void cp_param(t_param *dest, t_param *src)
146 {
147     int j;
148
149     for (j = 0; (j < MAXATOMLIST); j++)
150     {
151         dest->a[j] = src->a[j];
152     }
153     for (j = 0; (j < MAXFORCEPARAM); j++)
154     {
155         dest->c[j] = src->c[j];
156     }
157     strncpy(dest->s, src->s, sizeof(dest->s));
158 }
159
160 void add_param_to_list(t_params *list, t_param *b)
161 {
162     int j;
163
164     /* allocate one position extra */
165     pr_alloc (1, list);
166
167     /* fill the arrays */
168     for (j = 0; (j < MAXFORCEPARAM); j++)
169     {
170         list->param[list->nr].c[j]   = b->c[j];
171     }
172     for (j = 0; (j < MAXATOMLIST); j++)
173     {
174         list->param[list->nr].a[j]   = b->a[j];
175     }
176     memset(list->param[list->nr].s, 0, sizeof(list->param[list->nr].s));
177
178     list->nr++;
179 }
180
181
182 void init_molinfo(t_molinfo *mol)
183 {
184     mol->nrexcl     = 0;
185     mol->excl_set   = FALSE;
186     mol->bProcessed = FALSE;
187     init_plist(mol->plist);
188     init_block(&mol->cgs);
189     init_block(&mol->mols);
190     init_blocka(&mol->excls);
191     init_t_atoms(&mol->atoms, 0, FALSE);
192 }
193
194 /* FREEING MEMORY */
195
196 void done_mi(t_molinfo *mi)
197 {
198     done_atom (&(mi->atoms));
199     done_block(&(mi->cgs));
200     done_block(&(mi->mols));
201     done_plist(mi->plist);
202 }
203
204 /* PRINTING STRUCTURES */
205
206 static void print_bt(FILE *out, Directive d, gpp_atomtype *at,
207                      int ftype, int fsubtype, t_params plist[],
208                      bool bFullDih)
209 {
210     /* This dihp is a DIRTY patch because the dih-types do not use
211      * all four atoms to determine the type.
212      */
213     const int    dihp[2][2] = { { 1, 2 }, { 0, 3 } };
214     t_params    *bt;
215     int          i, j, f, nral, nrfp;
216     bool         bDih = FALSE, bSwapParity;
217
218     bt = &(plist[ftype]);
219
220     if (!bt->nr)
221     {
222         return;
223     }
224
225     f = 0;
226     switch (ftype)
227     {
228         case F_G96ANGLES:
229             f = 1;
230             break;
231         case F_G96BONDS:
232             f = 1;
233             break;
234         case F_MORSE:
235             f = 2;
236             break;
237         case F_CUBICBONDS:
238             f = 3;
239             break;
240         case F_CONNBONDS:
241             f = 4;
242             break;
243         case F_HARMONIC:
244             f = 5;
245             break;
246         case F_CROSS_BOND_ANGLES:
247             f = 2;
248             break;
249         case F_CROSS_BOND_BONDS:
250             f = 3;
251             break;
252         case F_UREY_BRADLEY:
253             f = 4;
254             break;
255         case F_PDIHS:
256         case F_RBDIHS:
257         case F_FOURDIHS:
258             bDih = TRUE;
259             break;
260         case F_IDIHS:
261             f    = 1;
262             bDih = TRUE;
263             break;
264         case F_CONSTRNC:
265             f = 1;
266             break;
267         case F_VSITE3FD:
268             f = 1;
269             break;
270         case F_VSITE3FAD:
271             f = 2;
272             break;
273         case F_VSITE3OUT:
274             f = 3;
275             break;
276         case F_VSITE4FDN:
277             f = 1;
278             break;
279         case F_CMAP:
280             f = 1;
281             break;
282
283         default:
284             bDih = FALSE;
285     }
286     if (bFullDih)
287     {
288         bDih = FALSE;
289     }
290     if (fsubtype)
291     {
292         f = fsubtype-1;
293     }
294
295     nral = NRAL(ftype);
296     nrfp = NRFP(ftype);
297
298     /* header */
299     fprintf(out, "[ %s ]\n", dir2str(d));
300     fprintf(out, "; ");
301     if (!bDih)
302     {
303         fprintf (out, "%3s  %4s", "ai", "aj");
304         for (j = 2; (j < nral); j++)
305         {
306             fprintf (out, "  %3c%c", 'a', 'i'+j);
307         }
308     }
309     else
310     {
311         for (j = 0; (j < 2); j++)
312         {
313             fprintf (out, "%3c%c", 'a', 'i'+dihp[f][j]);
314         }
315     }
316
317     fprintf (out, " funct");
318     for (j = 0; (j < nrfp); j++)
319     {
320         fprintf (out, " %12c%1d", 'c', j);
321     }
322     fprintf (out, "\n");
323
324     /* print bondtypes */
325     for (i = 0; (i < bt->nr); i++)
326     {
327         bSwapParity = (bt->param[i].c0() == NOTSET) && (bt->param[i].c1() == -1);
328         if (!bDih)
329         {
330             for (j = 0; (j < nral); j++)
331             {
332                 fprintf (out, "%5s ", get_atomtype_name(bt->param[i].a[j], at));
333             }
334         }
335         else
336         {
337             for (j = 0; (j < 2); j++)
338             {
339                 fprintf (out, "%5s ", get_atomtype_name(bt->param[i].a[dihp[f][j]], at));
340             }
341         }
342         fprintf (out, "%5d ", bSwapParity ? -f-1 : f+1);
343
344         if (bt->param[i].s[0])
345         {
346             fprintf(out, "   %s", bt->param[i].s);
347         }
348         else
349         {
350             for (j = 0; (j < nrfp && (bt->param[i].c[j] != NOTSET)); j++)
351             {
352                 fprintf (out, "%13.6e ", bt->param[i].c[j]);
353             }
354         }
355
356         fprintf (out, "\n");
357     }
358     fprintf (out, "\n");
359     fflush (out);
360 }
361
362 void print_blocka(FILE *out, const char *szName,
363                   const char *szIndex, const char *szA,
364                   t_blocka *block)
365 {
366     int i, j;
367
368     fprintf (out, "; %s\n", szName);
369     fprintf (out, "; %4s    %s\n", szIndex, szA);
370     for (i = 0; (i < block->nr); i++)
371     {
372         for (i = 0; (i < block->nr); i++)
373         {
374             fprintf (out, "%6d", i+1);
375             for (j = block->index[i]; (j < (block->index[i+1])); j++)
376             {
377                 fprintf (out, "%5d", block->a[j]+1);
378             }
379             fprintf (out, "\n");
380         }
381         fprintf (out, "\n");
382     }
383 }
384
385 void print_excl(FILE *out, int natoms, t_excls excls[])
386 {
387     int         i;
388     bool        have_excl;
389     int         j;
390
391     have_excl = FALSE;
392     for (i = 0; i < natoms && !have_excl; i++)
393     {
394         have_excl = (excls[i].nr > 0);
395     }
396
397     if (have_excl)
398     {
399         fprintf (out, "[ %s ]\n", dir2str(Directive::d_exclusions));
400         fprintf (out, "; %4s    %s\n", "i", "excluded from i");
401         for (i = 0; i < natoms; i++)
402         {
403             if (excls[i].nr > 0)
404             {
405                 fprintf (out, "%6d ", i+1);
406                 for (j = 0; j < excls[i].nr; j++)
407                 {
408                     fprintf (out, " %5d", excls[i].e[j]+1);
409                 }
410                 fprintf (out, "\n");
411             }
412         }
413         fprintf (out, "\n");
414         fflush(out);
415     }
416 }
417
418 static double get_residue_charge(const t_atoms *atoms, int at)
419 {
420     int    ri;
421     double q;
422
423     ri = atoms->atom[at].resind;
424     q  = 0;
425     while (at < atoms->nr && atoms->atom[at].resind == ri)
426     {
427         q += atoms->atom[at].q;
428         at++;
429     }
430
431     return q;
432 }
433
434 void print_atoms(FILE *out, gpp_atomtype *atype, t_atoms *at, int *cgnr,
435                  bool bRTPresname)
436 {
437     int         i, ri;
438     int         tpA, tpB;
439     const char *as;
440     char       *tpnmA, *tpnmB;
441     double      qres, qtot;
442
443     as = dir2str(Directive::d_atoms);
444     fprintf(out, "[ %s ]\n", as);
445     fprintf(out, "; %4s %10s %6s %7s%6s %6s %10s %10s %6s %10s %10s\n",
446             "nr", "type", "resnr", "residue", "atom", "cgnr", "charge", "mass", "typeB", "chargeB", "massB");
447
448     qtot  = 0;
449
450     if (at->nres)
451     {
452         /* if the information is present... */
453         for (i = 0; (i < at->nr); i++)
454         {
455             ri = at->atom[i].resind;
456             if ((i == 0 || ri != at->atom[i-1].resind) &&
457                 at->resinfo[ri].rtp != nullptr)
458             {
459                 qres = get_residue_charge(at, i);
460                 fprintf(out, "; residue %3d %-3s rtp %-4s q ",
461                         at->resinfo[ri].nr,
462                         *at->resinfo[ri].name,
463                         *at->resinfo[ri].rtp);
464                 if (fabs(qres) < 0.001)
465                 {
466                     fprintf(out, " %s", "0.0");
467                 }
468                 else
469                 {
470                     fprintf(out, "%+3.1f", qres);
471                 }
472                 fprintf(out, "\n");
473             }
474             tpA = at->atom[i].type;
475             if ((tpnmA = get_atomtype_name(tpA, atype)) == nullptr)
476             {
477                 gmx_fatal(FARGS, "tpA = %d, i= %d in print_atoms", tpA, i);
478             }
479
480             /* This is true by construction, but static analysers don't know */
481             GMX_ASSERT(!bRTPresname || at->resinfo[at->atom[i].resind].rtp, "-rtpres did not have residue name available");
482             fprintf(out, "%6d %10s %6d%c %5s %6s %6d %10g %10g",
483                     i+1, tpnmA,
484                     at->resinfo[ri].nr,
485                     at->resinfo[ri].ic,
486                     bRTPresname ?
487                     *(at->resinfo[at->atom[i].resind].rtp) :
488                     *(at->resinfo[at->atom[i].resind].name),
489                     *(at->atomname[i]), cgnr[i],
490                     at->atom[i].q, at->atom[i].m);
491             if (PERTURBED(at->atom[i]))
492             {
493                 tpB = at->atom[i].typeB;
494                 if ((tpnmB = get_atomtype_name(tpB, atype)) == nullptr)
495                 {
496                     gmx_fatal(FARGS, "tpB = %d, i= %d in print_atoms", tpB, i);
497                 }
498                 fprintf(out, " %6s %10g %10g",
499                         tpnmB, at->atom[i].qB, at->atom[i].mB);
500             }
501             // Accumulate the total charge to help troubleshoot issues.
502             qtot += static_cast<double>(at->atom[i].q);
503             // Round it to zero if it is close to zero, because
504             // printing -9.34e-5 confuses users.
505             if (fabs(qtot) < 0.0001)
506             {
507                 qtot = 0;
508             }
509             // Write the total charge for the last atom of the system
510             // and/or residue, because generally that's where it is
511             // expected to be an integer.
512             if (i == at->nr-1 || ri != at->atom[i+1].resind)
513             {
514                 fprintf(out, "   ; qtot %.4g\n", qtot);
515             }
516             else
517             {
518                 fputs("\n", out);
519             }
520         }
521     }
522     fprintf(out, "\n");
523     fflush(out);
524 }
525
526 void print_bondeds(FILE *out, int natoms, Directive d,
527                    int ftype, int fsubtype, t_params plist[])
528 {
529     t_symtab       stab;
530     gpp_atomtype  *atype;
531     t_param       *param;
532     t_atom        *a;
533     int            i;
534
535     atype = init_atomtype();
536     snew(a, 1);
537     snew(param, 1);
538     open_symtab(&stab);
539     for (i = 0; (i < natoms); i++)
540     {
541         char buf[12];
542         sprintf(buf, "%4d", (i+1));
543         add_atomtype(atype, &stab, a, buf, param, 0, 0);
544     }
545     print_bt(out, d, atype, ftype, fsubtype, plist, TRUE);
546
547     done_symtab(&stab);
548     sfree(a);
549     sfree(param);
550     done_atomtype(atype);
551 }