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