Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / topdirs.cpp
blobb955ba2d91ad3dfd5dedff082b2af3d4e3c5ed9f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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,2017, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "topdirs.h"
41 #include <stdarg.h>
42 #include <stdio.h>
44 #include "gromacs/utility/cstringutil.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
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 "intermolecular_interactions",
89 "invalid"
92 int ifunc_index(directive d, int type)
94 switch (d)
96 case d_bondtypes:
97 case d_bonds:
98 switch (type)
100 case 1:
101 return F_BONDS;
102 case 2:
103 return F_G96BONDS;
104 case 3:
105 return F_MORSE;
106 case 4:
107 return F_CUBICBONDS;
108 case 5:
109 return F_CONNBONDS;
110 case 6:
111 return F_HARMONIC;
112 case 7:
113 return F_FENEBONDS;
114 case 8:
115 return F_TABBONDS;
116 case 9:
117 return F_TABBONDSNC;
118 case 10:
119 return F_RESTRBONDS;
120 default:
121 gmx_fatal(FARGS, "Invalid bond type %d", type);
122 break;
124 break;
125 case d_angles:
126 case d_angletypes:
127 switch (type)
129 case 1:
130 return F_ANGLES;
131 case 2:
132 return F_G96ANGLES;
133 case 3:
134 return F_CROSS_BOND_BONDS;
135 case 4:
136 return F_CROSS_BOND_ANGLES;
137 case 5:
138 return F_UREY_BRADLEY;
139 case 6:
140 return F_QUARTIC_ANGLES;
141 case 8:
142 return F_TABANGLES;
143 case 9:
144 return F_LINEAR_ANGLES;
145 case 10:
146 return F_RESTRANGLES;
147 default:
148 gmx_fatal(FARGS, "Invalid angle type %d", type);
149 break;
151 break;
152 case d_pairs:
153 case d_pairtypes:
154 if (type == 1 || (d == d_pairtypes && type == 2))
156 return F_LJ14;
158 else if (type == 2)
160 return F_LJC14_Q;
162 else
164 gmx_fatal(FARGS, "Invalid pairs type %d", type);
166 break;
167 case d_pairs_nb:
168 return F_LJC_PAIRS_NB;
169 case d_dihedrals:
170 case d_dihedraltypes:
171 switch (type)
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 case 10:
188 return F_RESTRDIHS;
189 case 11:
190 return F_CBTDIHS;
191 default:
192 gmx_fatal(FARGS, "Invalid dihedral type %d", type);
194 break;
195 case d_cmaptypes:
196 case d_cmap:
197 return F_CMAP;
199 case d_nonbond_params:
200 if (type == 1)
202 return F_LJ;
204 else
206 return F_BHAM;
208 case d_vsites2:
209 return F_VSITE2;
210 case d_vsites3:
211 switch (type)
213 case 1:
214 return F_VSITE3;
215 case 2:
216 return F_VSITE3FD;
217 case 3:
218 return F_VSITE3FAD;
219 case 4:
220 return F_VSITE3OUT;
221 default:
222 gmx_fatal(FARGS, "Invalid vsites3 type %d", type);
224 break;
225 case d_vsites4:
226 switch (type)
228 case 1:
229 return F_VSITE4FD;
230 case 2:
231 return F_VSITE4FDN;
232 default:
233 gmx_fatal(FARGS, "Invalid vsites4 type %d", type);
235 break;
236 case d_vsitesn:
237 return F_VSITEN;
238 case d_constraints:
239 case d_constrainttypes:
240 switch (type)
242 case 1:
243 return F_CONSTR;
244 case 2:
245 return F_CONSTRNC;
246 default:
247 gmx_fatal(FARGS, "Invalid constraints type %d", type);
249 break;
250 case d_settles:
251 return F_SETTLE;
252 case d_position_restraints:
253 switch (type)
255 case 1:
256 return F_POSRES;
257 case 2:
258 return F_FBPOSRES;
259 default:
260 gmx_fatal(FARGS, "Invalid position restraint type %d", type);
262 break;
263 case d_polarization:
264 switch (type)
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);
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__);
292 return -1;
295 const char *dir2str (directive d)
297 if (d < d_maxdir)
299 return directive_names[d];
301 else
303 return directive_names[d_maxdir];
307 directive str2dir (char *dstr)
309 int d;
310 char buf[STRLEN], *ptr;
312 /* Hack to be able to read old topologies */
313 if (gmx_strncasecmp_min(dstr, "dummies", 7) == 0)
315 sprintf(buf, "virtual_sites%s", dstr+7);
316 ptr = buf;
318 else
320 ptr = dstr;
323 for (d = 0; (d < d_maxdir); d++)
325 if (gmx_strcasecmp_min(ptr, dir2str((directive)d)) == 0)
327 return (directive)d;
331 return d_invalid;
334 static directive **necessary = nullptr;
336 static void set_nec(directive **n, ...)
337 /* Must always have at least one extra argument */
339 va_list ap;
340 int ind = 0;
341 directive d;
343 va_start(ap, n);
346 d = (directive)va_arg(ap, int);
347 srenew(*n, ++ind);
348 (*n)[ind-1] = d;
350 while (d != d_none);
351 va_end(ap);
354 void DS_Init(DirStack **DS)
356 if (necessary == nullptr)
358 int i;
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 set_nec(&(necessary[d_intermolecular_interactions]), d_molecules, d_none);
400 for (i = 0; (i < d_maxdir); i++)
402 if (debug)
404 fprintf(debug, "%20s: ", dir2str((directive)i));
406 if (necessary[i])
408 directive d;
409 int j = 0;
412 d = necessary[i][j++];
413 if (debug)
415 fprintf(debug, "%20s ", dir2str(d));
418 while (d != d_none);
420 if (debug)
422 fprintf(debug, "\n");
426 *DS = nullptr;
430 void DS_Done (DirStack **DS)
432 DirStack *D;
434 while (*DS != nullptr)
436 D = *DS;
437 *DS = (*DS)->prev;
438 sfree (D);
442 void DS_Push (DirStack **DS, directive d)
444 DirStack *D;
446 snew(D, 1);
447 D->d = d;
448 D->prev = *DS;
449 *DS = D;
452 int DS_Search(DirStack *DS, directive d)
454 DirStack *D;
456 D = DS;
457 while ((D != nullptr) && (D->d != d))
459 D = D->prev;
462 return (D != nullptr);
465 int DS_Check_Order(DirStack *DS, directive d)
467 directive d0;
468 int i = 0;
470 /* Check if parameter definitions appear after a moleculetype directive */
471 if (d < d_moleculetype && DS_Search(DS, d_moleculetype))
473 return FALSE;
476 /* Check if all the necessary directives have appeared before directive d */
477 if (necessary[d][0] == d_none)
479 return TRUE;
481 else
485 d0 = necessary[d][i++];
486 if (DS_Search(DS, d0))
488 return TRUE;
491 while (d0 != d_none);
493 return FALSE;