Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / toputil.c
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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <string.h>
43
44 #include "macros.h"
45 #include "topdirs.h"
46 #include "toputil.h"
47 #include "gpp_atomtype.h"
48
49 #include "gromacs/topology/block.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
53
54 /* UTILITIES */
55
56 void set_p_string(t_param *p, const char *s)
57 {
58     if (s)
59     {
60         if (strlen(s) < sizeof(p->s)-1)
61         {
62             strncpy(p->s, s, sizeof(p->s));
63         }
64         else
65         {
66             gmx_fatal(FARGS, "Increase MAXSLEN in include/grompp-impl.h to at least %d,"
67                       " or shorten your definition of bonds like %s to at most %d",
68                       strlen(s)+1, s, MAXSLEN-1);
69         }
70     }
71     else
72     {
73         strcpy(p->s, "");
74     }
75 }
76
77 void pr_alloc (int extra, t_params *pr)
78 {
79     int i, j;
80
81     /* get new space for arrays */
82     if (extra < 0)
83     {
84         gmx_fatal(FARGS, "Trying to make array smaller.\n");
85     }
86     if (extra == 0)
87     {
88         return;
89     }
90     if ((pr->nr == 0) && (pr->param != NULL))
91     {
92         fprintf(stderr, "Warning: dangling pointer at %lx\n",
93                 (unsigned long)pr->param);
94         pr->param = NULL;
95     }
96     if (pr->nr+extra > pr->maxnr)
97     {
98         pr->maxnr = max(1.2*pr->maxnr, pr->maxnr + extra);
99         srenew(pr->param, pr->maxnr);
100         for (i = pr->nr; (i < pr->maxnr); i++)
101         {
102             for (j = 0; (j < MAXATOMLIST); j++)
103             {
104                 pr->param[i].a[j] = 0;
105             }
106             for (j = 0; (j < MAXFORCEPARAM); j++)
107             {
108                 pr->param[i].c[j] = 0;
109             }
110             set_p_string(&(pr->param[i]), "");
111         }
112     }
113 }
114
115 void init_plist(t_params plist[])
116 {
117     int i;
118
119     for (i = 0; (i < F_NRE); i++)
120     {
121         plist[i].nr    = 0;
122         plist[i].maxnr = 0;
123         plist[i].param = NULL;
124
125         /* CMAP */
126         plist[i].ncmap        = 0;
127         plist[i].cmap         = NULL;
128         plist[i].grid_spacing = 0;
129         plist[i].nc           = 0;
130         plist[i].nct          = 0;
131         plist[i].cmap_types   = NULL;
132     }
133 }
134
135 void cp_param(t_param *dest, t_param *src)
136 {
137     int j;
138
139     for (j = 0; (j < MAXATOMLIST); j++)
140     {
141         dest->a[j] = src->a[j];
142     }
143     for (j = 0; (j < MAXFORCEPARAM); j++)
144     {
145         dest->c[j] = src->c[j];
146     }
147     strncpy(dest->s, src->s, sizeof(dest->s));
148 }
149
150 void add_param_to_list(t_params *list, t_param *b)
151 {
152     int j;
153
154     /* allocate one position extra */
155     pr_alloc (1, list);
156
157     /* fill the arrays */
158     for (j = 0; (j < MAXFORCEPARAM); j++)
159     {
160         list->param[list->nr].c[j]   = b->c[j];
161     }
162     for (j = 0; (j < MAXATOMLIST); j++)
163     {
164         list->param[list->nr].a[j]   = b->a[j];
165     }
166     memset(list->param[list->nr].s, 0, sizeof(list->param[list->nr].s));
167
168     list->nr++;
169 }
170
171
172 void init_molinfo(t_molinfo *mol)
173 {
174     mol->nrexcl     = 0;
175     mol->excl_set   = FALSE;
176     mol->bProcessed = FALSE;
177     init_plist(mol->plist);
178     init_block(&mol->cgs);
179     init_block(&mol->mols);
180     init_blocka(&mol->excls);
181     init_atom(&mol->atoms);
182 }
183
184 /* FREEING MEMORY */
185
186 void done_bt (t_params *pl)
187 {
188     sfree(pl->param);
189 }
190
191 void done_mi(t_molinfo *mi)
192 {
193     int i;
194
195     done_atom (&(mi->atoms));
196     done_block(&(mi->cgs));
197     done_block(&(mi->mols));
198     for (i = 0; (i < F_NRE); i++)
199     {
200         done_bt(&(mi->plist[i]));
201     }
202 }
203
204 /* PRINTING STRUCTURES */
205
206 void print_bt(FILE *out, directive d, gpp_atomtype_t at,
207               int ftype, int fsubtype, t_params plist[],
208               gmx_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     gmx_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 < ((int)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     atom_id     i;
388     gmx_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(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_t atype, t_atoms *at, int *cgnr,
435                  gmx_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(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 (debug)
451     {
452         fprintf(debug, "This molecule has %d atoms and %d residues\n",
453                 at->nr, at->nres);
454     }
455
456     if (at->nres)
457     {
458         /* if the information is present... */
459         for (i = 0; (i < at->nr); i++)
460         {
461             ri = at->atom[i].resind;
462             if ((i == 0 || ri != at->atom[i-1].resind) &&
463                 at->resinfo[ri].rtp != NULL)
464             {
465                 qres = get_residue_charge(at, i);
466                 fprintf(out, "; residue %3d %-3s rtp %-4s q ",
467                         at->resinfo[ri].nr,
468                         *at->resinfo[ri].name,
469                         *at->resinfo[ri].rtp);
470                 if (fabs(qres) < 0.001)
471                 {
472                     fprintf(out, " %s", "0.0");
473                 }
474                 else
475                 {
476                     fprintf(out, "%+3.1f", qres);
477                 }
478                 fprintf(out, "\n");
479             }
480             tpA = at->atom[i].type;
481             if ((tpnmA = get_atomtype_name(tpA, atype)) == NULL)
482             {
483                 gmx_fatal(FARGS, "tpA = %d, i= %d in print_atoms", tpA, i);
484             }
485
486             fprintf(out, "%6d %10s %6d%c %5s %6s %6d %10g %10g",
487                     i+1, tpnmA,
488                     at->resinfo[ri].nr,
489                     at->resinfo[ri].ic,
490                     bRTPresname ?
491                     *(at->resinfo[at->atom[i].resind].rtp) :
492                     *(at->resinfo[at->atom[i].resind].name),
493                     *(at->atomname[i]), cgnr[i],
494                     at->atom[i].q, at->atom[i].m);
495             if (PERTURBED(at->atom[i]))
496             {
497                 tpB = at->atom[i].typeB;
498                 if ((tpnmB = get_atomtype_name(tpB, atype)) == NULL)
499                 {
500                     gmx_fatal(FARGS, "tpB = %d, i= %d in print_atoms", tpB, i);
501                 }
502                 fprintf(out, " %6s %10g %10g",
503                         tpnmB, at->atom[i].qB, at->atom[i].mB);
504             }
505             qtot += (double)at->atom[i].q;
506             if (fabs(qtot) < 4*GMX_REAL_EPS)
507             {
508                 qtot = 0;
509             }
510             fprintf(out, "   ; qtot %.4g\n", qtot);
511         }
512     }
513     fprintf(out, "\n");
514     fflush(out);
515 }
516
517 void print_bondeds(FILE *out, int natoms, directive d,
518                    int ftype, int fsubtype, t_params plist[])
519 {
520     t_symtab       stab;
521     gpp_atomtype_t atype;
522     t_param       *param;
523     t_atom        *a;
524     int            i;
525
526     atype = init_atomtype();
527     snew(a, 1);
528     snew(param, 1);
529     open_symtab(&stab);
530     for (i = 0; (i < natoms); i++)
531     {
532         char buf[12];
533         sprintf(buf, "%4d", (i+1));
534         add_atomtype(atype, &stab, a, buf, param, 0, 0, 0, 0, 0, 0, 0);
535     }
536     print_bt(out, d, atype, ftype, fsubtype, plist, TRUE);
537
538     done_symtab(&stab);
539     sfree(a);
540     sfree(param);
541     done_atomtype(atype);
542 }