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, 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.
45 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gmx_fatal.h"
51 /* Must correspond to the directive enum in grompp-impl.h */
52 static const char *directive_names[d_maxdir+1] = {
61 "implicit_genborn_params",
62 "implicit_surface_params",
64 /* All the directives above can not appear after moleculetype */
84 "position_restraints",
87 "distance_restraints",
88 "orientation_restraints",
89 "dihedral_restraints",
94 int ifunc_index(directive d, int type)
123 gmx_fatal(FARGS, "Invalid bond type %d", type);
136 return F_CROSS_BOND_BONDS;
138 return F_CROSS_BOND_ANGLES;
140 return F_UREY_BRADLEY;
142 return F_QUARTIC_ANGLES;
146 return F_LINEAR_ANGLES;
148 return F_RESTRANGLES;
150 gmx_fatal(FARGS, "Invalid angle type %d", type);
156 if (type == 1 || (d == d_pairtypes && type == 2))
166 gmx_fatal(FARGS, "Invalid pairs type %d", type);
170 return F_LJC_PAIRS_NB;
172 case d_dihedraltypes:
188 return F_PDIHS; /* proper dihedrals where we allow multiple terms over single bond */
194 gmx_fatal(FARGS, "Invalid dihedral type %d", type);
201 case d_nonbond_params:
224 gmx_fatal(FARGS, "Invalid vsites3 type %d", type);
235 gmx_fatal(FARGS, "Invalid vsites4 type %d", type);
241 case d_constrainttypes:
249 gmx_fatal(FARGS, "Invalid constraints type %d", type);
254 case d_position_restraints:
263 gmx_fatal(FARGS, "Invalid position restraint type %d", type);
270 return F_POLARIZATION;
274 gmx_fatal(FARGS, "Invalid polarization type %d", type);
277 case d_thole_polarization:
279 case d_water_polarization:
281 case d_angle_restraints:
283 case d_angle_restraints_z:
285 case d_distance_restraints:
287 case d_orientation_restraints:
289 case d_dihedral_restraints:
292 gmx_fatal(FARGS, "invalid directive %s in ifunc_index (%s:%s)",
293 dir2str(d), __FILE__, __LINE__);
298 const char *dir2str (directive d)
302 return directive_names[d];
306 return directive_names[d_maxdir];
310 directive str2dir (char *dstr)
313 char buf[STRLEN], *ptr;
315 /* Hack to be able to read old topologies */
316 if (gmx_strncasecmp_min(dstr, "dummies", 7) == 0)
318 sprintf(buf, "virtual_sites%s", dstr+7);
326 for (d = 0; (d < d_maxdir); d++)
328 if (gmx_strcasecmp_min(ptr, dir2str((directive)d)) == 0)
337 static directive **necessary = NULL;
339 static void set_nec(directive **n, ...)
340 /* Must always have at least one extra argument */
349 d = (directive)va_arg(ap, int);
357 void DS_Init(DirStack **DS)
359 if (necessary == NULL)
363 snew(necessary, d_maxdir);
364 set_nec(&(necessary[d_defaults]), d_none);
365 set_nec(&(necessary[d_atomtypes]), d_defaults, d_none);
366 set_nec(&(necessary[d_bondtypes]), d_atomtypes, d_none);
367 set_nec(&(necessary[d_constrainttypes]), d_atomtypes, d_none);
368 set_nec(&(necessary[d_pairtypes]), d_atomtypes, d_none);
369 set_nec(&(necessary[d_angletypes]), d_atomtypes, d_none);
370 set_nec(&(necessary[d_dihedraltypes]), d_atomtypes, d_none);
371 set_nec(&(necessary[d_nonbond_params]), d_atomtypes, d_none);
372 set_nec(&(necessary[d_implicit_genborn_params]), d_atomtypes, d_none);
373 set_nec(&(necessary[d_implicit_surface_params]), d_atomtypes, d_none);
374 set_nec(&(necessary[d_cmaptypes]), d_atomtypes, d_none);
375 set_nec(&(necessary[d_moleculetype]), d_atomtypes, d_none);
376 set_nec(&(necessary[d_atoms]), d_moleculetype, d_none);
377 set_nec(&(necessary[d_vsites2]), d_atoms, d_none);
378 set_nec(&(necessary[d_vsites3]), d_atoms, d_none);
379 set_nec(&(necessary[d_vsites4]), d_atoms, d_none);
380 set_nec(&(necessary[d_vsitesn]), d_atoms, d_none);
381 set_nec(&(necessary[d_bonds]), d_atoms, d_none);
382 set_nec(&(necessary[d_exclusions]), d_bonds, d_constraints, d_settles, d_none);
383 set_nec(&(necessary[d_pairs]), d_atoms, d_none);
384 set_nec(&(necessary[d_pairs_nb]), d_atoms, d_none);
385 set_nec(&(necessary[d_angles]), d_atoms, d_none);
386 set_nec(&(necessary[d_polarization]), d_atoms, d_none);
387 set_nec(&(necessary[d_water_polarization]), d_atoms, d_none);
388 set_nec(&(necessary[d_thole_polarization]), d_atoms, d_none);
389 set_nec(&(necessary[d_dihedrals]), d_atoms, d_none);
390 set_nec(&(necessary[d_constraints]), d_atoms, d_none);
391 set_nec(&(necessary[d_settles]), d_atoms, d_none);
392 set_nec(&(necessary[d_system]), d_moleculetype, d_none);
393 set_nec(&(necessary[d_molecules]), d_system, d_none);
394 set_nec(&(necessary[d_position_restraints]), d_atoms, d_none);
395 set_nec(&(necessary[d_angle_restraints]), d_atoms, d_none);
396 set_nec(&(necessary[d_angle_restraints_z]), d_atoms, d_none);
397 set_nec(&(necessary[d_distance_restraints]), d_atoms, d_none);
398 set_nec(&(necessary[d_orientation_restraints]), d_atoms, d_none);
399 set_nec(&(necessary[d_dihedral_restraints]), d_atoms, d_none);
400 set_nec(&(necessary[d_cmap]), d_atoms, d_none);
402 for (i = 0; (i < d_maxdir); i++)
406 fprintf(debug, "%20s: ", dir2str((directive)i));
414 d = necessary[i][j++];
417 fprintf(debug, "%20s ", dir2str(d));
424 fprintf(debug, "\n");
432 void DS_Done (DirStack **DS)
444 void DS_Push (DirStack **DS, directive d)
454 int DS_Search(DirStack *DS, directive d)
459 while ((D != NULL) && (D->d != d))
467 int DS_Check_Order(DirStack *DS, directive d)
472 /* Check if parameter definitions appear after a moleculetype directive */
473 if (d < d_moleculetype && DS_Search(DS, d_moleculetype))
478 /* Check if all the necessary directives have appeared before directive d */
479 if (necessary[d][0] == d_none)
487 d0 = necessary[d][i++];
488 if (DS_Search(DS, d0))
493 while (d0 != d_none);