Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / topdirs.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) 2013,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 <stdio.h>
42 #include <stdarg.h>
43
44 #include "sysstuff.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "string2.h"
48 #include "gmx_fatal.h"
49 #include "topdirs.h"
50
51 /* Must correspond to the directive enum in grompp-impl.h */
52 static const char *directive_names[d_maxdir+1] = {
53     "defaults",
54     "atomtypes",
55     "bondtypes",
56     "constrainttypes",
57     "pairtypes",
58     "angletypes",
59     "dihedraltypes",
60     "nonbond_params",
61     "implicit_genborn_params",
62     "implicit_surface_params",
63     "cmaptypes",
64     /* All the directives above can not appear after moleculetype */
65     "moleculetype",
66     "atoms",
67     "virtual_sites2",
68     "virtual_sites3",
69     "virtual_sites4",
70     "virtual_sitesn",
71     "bonds",
72     "exclusions",
73     "pairs",
74     "pairs_nb",
75     "angles",
76     "dihedrals",
77     "constraints",
78     "settles",
79     "polarization",
80     "water_polarization",
81     "thole_polarization",
82     "system",
83     "molecules",
84     "position_restraints",
85     "angle_restraints",
86     "angle_restraints_z",
87     "distance_restraints",
88     "orientation_restraints",
89     "dihedral_restraints",
90     "cmap",
91     "invalid"
92 };
93
94 int ifunc_index(directive d, int type)
95 {
96     switch (d)
97     {
98         case d_bondtypes:
99         case d_bonds:
100             switch (type)
101             {
102                 case 1:
103                     return F_BONDS;
104                 case 2:
105                     return F_G96BONDS;
106                 case 3:
107                     return F_MORSE;
108                 case 4:
109                     return F_CUBICBONDS;
110                 case 5:
111                     return F_CONNBONDS;
112                 case 6:
113                     return F_HARMONIC;
114                 case 7:
115                     return F_FENEBONDS;
116                 case 8:
117                     return F_TABBONDS;
118                 case 9:
119                     return F_TABBONDSNC;
120                 case 10:
121                     return F_RESTRBONDS;
122                 default:
123                     gmx_fatal(FARGS, "Invalid bond type %d", type);
124                     break;
125             }
126             break;
127         case d_angles:
128         case d_angletypes:
129             switch (type)
130             {
131                 case 1:
132                     return F_ANGLES;
133                 case 2:
134                     return F_G96ANGLES;
135                 case 3:
136                     return F_CROSS_BOND_BONDS;
137                 case 4:
138                     return F_CROSS_BOND_ANGLES;
139                 case 5:
140                     return F_UREY_BRADLEY;
141                 case 6:
142                     return F_QUARTIC_ANGLES;
143                 case 8:
144                     return F_TABANGLES;
145                 case 9:
146                     return F_LINEAR_ANGLES;
147                 default:
148                     gmx_fatal(FARGS, "Invalid angle type %d", type);
149                     break;
150             }
151             break;
152         case d_pairs:
153         case d_pairtypes:
154             if (type == 1 || (d == d_pairtypes && type == 2))
155             {
156                 return F_LJ14;
157             }
158             else if (type == 2)
159             {
160                 return F_LJC14_Q;
161             }
162             else
163             {
164                 gmx_fatal(FARGS, "Invalid pairs type %d", type);
165             }
166             break;
167         case d_pairs_nb:
168             return F_LJC_PAIRS_NB;
169         case d_dihedrals:
170         case d_dihedraltypes:
171             switch (type)
172             {
173                 case 1:
174                     return F_PDIHS;
175                 case 2:
176                     return F_IDIHS;
177                 case 3:
178                     return F_RBDIHS;
179                 case 4:
180                     return F_PIDIHS;
181                 case 5:
182                     return F_FOURDIHS;
183                 case 8:
184                     return F_TABDIHS;
185                 case 9:
186                     return F_PDIHS; /* proper dihedrals where we allow multiple terms over single bond */
187                 default:
188                     gmx_fatal(FARGS, "Invalid dihedral type %d", type);
189             }
190             break;
191         case d_cmaptypes:
192         case d_cmap:
193             return F_CMAP;
194
195         case d_nonbond_params:
196             if (type == 1)
197             {
198                 return F_LJ;
199             }
200             else
201             {
202                 return F_BHAM;
203             }
204         case d_vsites2:
205             return F_VSITE2;
206         case d_vsites3:
207             switch (type)
208             {
209                 case 1:
210                     return F_VSITE3;
211                 case 2:
212                     return F_VSITE3FD;
213                 case 3:
214                     return F_VSITE3FAD;
215                 case 4:
216                     return F_VSITE3OUT;
217                 default:
218                     gmx_fatal(FARGS, "Invalid vsites3 type %d", type);
219             }
220             break;
221         case d_vsites4:
222             switch (type)
223             {
224                 case 1:
225                     return F_VSITE4FD;
226                 case 2:
227                     return F_VSITE4FDN;
228                 default:
229                     gmx_fatal(FARGS, "Invalid vsites4 type %d", type);
230             }
231             break;
232         case d_vsitesn:
233             return F_VSITEN;
234         case d_constraints:
235         case d_constrainttypes:
236             switch (type)
237             {
238                 case 1:
239                     return F_CONSTR;
240                 case 2:
241                     return F_CONSTRNC;
242                 default:
243                     gmx_fatal(FARGS, "Invalid constraints type %d", type);
244             }
245             break;
246         case d_settles:
247             return F_SETTLE;
248         case d_position_restraints:
249             switch (type)
250             {
251                 case 1:
252                     return F_POSRES;
253                 case 2:
254                     return F_FBPOSRES;
255                     break;
256                 default:
257                     gmx_fatal(FARGS, "Invalid position restraint type %d", type);
258             }
259             break;
260         case d_polarization:
261             switch (type)
262             {
263                 case 1:
264                     return F_POLARIZATION;
265                 case 2:
266                     return F_ANHARM_POL;
267                 default:
268                     gmx_fatal(FARGS, "Invalid polarization type %d", type);
269             }
270             break;
271         case d_thole_polarization:
272             return F_THOLE_POL;
273         case d_water_polarization:
274             return F_WATER_POL;
275         case d_angle_restraints:
276             return F_ANGRES;
277         case d_angle_restraints_z:
278             return F_ANGRESZ;
279         case d_distance_restraints:
280             return F_DISRES;
281         case d_orientation_restraints:
282             return F_ORIRES;
283         case d_dihedral_restraints:
284             return F_DIHRES;
285         default:
286             gmx_fatal(FARGS, "invalid directive %s in ifunc_index (%s:%s)",
287                       dir2str(d), __FILE__, __LINE__);
288     }
289     return -1;
290 }
291
292 const char *dir2str (directive d)
293 {
294     if (d < d_maxdir)
295     {
296         return directive_names[d];
297     }
298     else
299     {
300         return directive_names[d_maxdir];
301     }
302 }
303
304 directive str2dir (char *dstr)
305 {
306     int  d;
307     char buf[STRLEN], *ptr;
308
309     /* Hack to be able to read old topologies */
310     if (gmx_strncasecmp_min(dstr, "dummies", 7) == 0)
311     {
312         sprintf(buf, "virtual_sites%s", dstr+7);
313         ptr = buf;
314     }
315     else
316     {
317         ptr = dstr;
318     }
319
320     for (d = 0; (d < d_maxdir); d++)
321     {
322         if (gmx_strcasecmp_min(ptr, dir2str((directive)d)) == 0)
323         {
324             return (directive)d;
325         }
326     }
327
328     return d_invalid;
329 }
330
331 static directive **necessary = NULL;
332
333 static void set_nec(directive **n, ...)
334 /* Must always have at least one extra argument */
335 {
336     va_list   ap;
337     int       ind = 0;
338     directive d;
339
340     va_start(ap, n);
341     do
342     {
343         d = (directive)va_arg(ap, int);
344         srenew(*n, ++ind);
345         (*n)[ind-1] = d;
346     }
347     while (d != d_none);
348     va_end(ap);
349 }
350
351 void DS_Init(DirStack **DS)
352 {
353     if (necessary == NULL)
354     {
355         int i;
356
357         snew(necessary, d_maxdir);
358         set_nec(&(necessary[d_defaults]), d_none);
359         set_nec(&(necessary[d_atomtypes]), d_defaults, d_none);
360         set_nec(&(necessary[d_bondtypes]), d_atomtypes, d_none);
361         set_nec(&(necessary[d_constrainttypes]), d_atomtypes, d_none);
362         set_nec(&(necessary[d_pairtypes]), d_atomtypes, d_none);
363         set_nec(&(necessary[d_angletypes]), d_atomtypes, d_none);
364         set_nec(&(necessary[d_dihedraltypes]), d_atomtypes, d_none);
365         set_nec(&(necessary[d_nonbond_params]), d_atomtypes, d_none);
366         set_nec(&(necessary[d_implicit_genborn_params]), d_atomtypes, d_none);
367         set_nec(&(necessary[d_implicit_surface_params]), d_atomtypes, d_none);
368         set_nec(&(necessary[d_cmaptypes]), d_atomtypes, d_none);
369         set_nec(&(necessary[d_moleculetype]), d_atomtypes, d_none);
370         set_nec(&(necessary[d_atoms]), d_moleculetype, d_none);
371         set_nec(&(necessary[d_vsites2]), d_atoms, d_none);
372         set_nec(&(necessary[d_vsites3]), d_atoms, d_none);
373         set_nec(&(necessary[d_vsites4]), d_atoms, d_none);
374         set_nec(&(necessary[d_vsitesn]), d_atoms, d_none);
375         set_nec(&(necessary[d_bonds]), d_atoms, d_none);
376         set_nec(&(necessary[d_exclusions]), d_bonds, d_constraints, d_settles, d_none);
377         set_nec(&(necessary[d_pairs]), d_atoms, d_none);
378         set_nec(&(necessary[d_pairs_nb]), d_atoms, d_none);
379         set_nec(&(necessary[d_angles]), d_atoms, d_none);
380         set_nec(&(necessary[d_polarization]), d_atoms, d_none);
381         set_nec(&(necessary[d_water_polarization]), d_atoms, d_none);
382         set_nec(&(necessary[d_thole_polarization]), d_atoms, d_none);
383         set_nec(&(necessary[d_dihedrals]), d_atoms, d_none);
384         set_nec(&(necessary[d_constraints]), d_atoms, d_none);
385         set_nec(&(necessary[d_settles]), d_atoms, d_none);
386         set_nec(&(necessary[d_system]), d_moleculetype, d_none);
387         set_nec(&(necessary[d_molecules]), d_system, d_none);
388         set_nec(&(necessary[d_position_restraints]), d_atoms, d_none);
389         set_nec(&(necessary[d_angle_restraints]), d_atoms, d_none);
390         set_nec(&(necessary[d_angle_restraints_z]), d_atoms, d_none);
391         set_nec(&(necessary[d_distance_restraints]), d_atoms, d_none);
392         set_nec(&(necessary[d_orientation_restraints]), d_atoms, d_none);
393         set_nec(&(necessary[d_dihedral_restraints]), d_atoms, d_none);
394         set_nec(&(necessary[d_cmap]), d_atoms, d_none);
395
396         for (i = 0; (i < d_maxdir); i++)
397         {
398             if (debug)
399             {
400                 fprintf(debug, "%20s:  ", dir2str((directive)i));
401             }
402             if (necessary[i])
403             {
404                 directive d;
405                 int       j = 0;
406                 do
407                 {
408                     d = necessary[i][j++];
409                     if (debug)
410                     {
411                         fprintf(debug, "%20s  ", dir2str(d));
412                     }
413                 }
414                 while (d != d_none);
415             }
416             if (debug)
417             {
418                 fprintf(debug, "\n");
419             }
420         }
421     }
422     *DS = NULL;
423
424 }
425
426 void DS_Done (DirStack **DS)
427 {
428     DirStack *D;
429
430     while (*DS != NULL)
431     {
432         D   = *DS;
433         *DS = (*DS)->prev;
434         sfree (D);
435     }
436 }
437
438 void DS_Push (DirStack **DS, directive d)
439 {
440     DirStack *D;
441
442     snew(D, 1);
443     D->d    = d;
444     D->prev = *DS;
445     *DS     = D;
446 }
447
448 int DS_Search(DirStack *DS, directive d)
449 {
450     DirStack *D;
451
452     D = DS;
453     while ((D != NULL) && (D->d != d))
454     {
455         D = D->prev;
456     }
457
458     return (D != NULL);
459 }
460
461 int DS_Check_Order(DirStack *DS, directive d)
462 {
463     directive d0;
464     int       i = 0;
465
466     /* Check if parameter definitions appear after a moleculetype directive */
467     if (d < d_moleculetype && DS_Search(DS, d_moleculetype))
468     {
469         return FALSE;
470     }
471
472     /* Check if all the necessary directives have appeared before directive d */
473     if (necessary[d][0] == d_none)
474     {
475         return TRUE;
476     }
477     else
478     {
479         do
480         {
481             d0 = necessary[d][i++];
482             if (DS_Search(DS, d0))
483             {
484                 return TRUE;
485             }
486         }
487         while (d0 != d_none);
488     }
489     return FALSE;
490 }