3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
46 #include "gmx_fatal.h"
49 /* Must correspond to the directive enum in grompp.h */
50 static const char *directive_names[d_maxdir+1] = {
59 "implicit_genborn_params",
60 "implicit_surface_params",
62 /* All the directives above can not appear after moleculetype */
82 "position_restraints",
85 "distance_restraints",
86 "orientation_restraints",
87 "dihedral_restraints",
92 int ifunc_index(directive d, int type)
121 gmx_fatal(FARGS, "Invalid bond type %d", type);
134 return F_CROSS_BOND_BONDS;
136 return F_CROSS_BOND_ANGLES;
138 return F_UREY_BRADLEY;
140 return F_QUARTIC_ANGLES;
144 return F_LINEAR_ANGLES;
146 gmx_fatal(FARGS, "Invalid angle type %d", type);
152 if (type == 1 || (d == d_pairtypes && type == 2))
162 gmx_fatal(FARGS, "Invalid pairs type %d", type);
166 return F_LJC_PAIRS_NB;
168 case d_dihedraltypes:
184 return F_PDIHS; /* proper dihedrals where we allow multiple terms over single bond */
186 gmx_fatal(FARGS, "Invalid dihedral type %d", type);
193 case d_nonbond_params:
216 gmx_fatal(FARGS, "Invalid vsites3 type %d", type);
227 gmx_fatal(FARGS, "Invalid vsites4 type %d", type);
233 case d_constrainttypes:
241 gmx_fatal(FARGS, "Invalid constraints type %d", type);
246 case d_position_restraints:
255 gmx_fatal(FARGS, "Invalid position restraint type %d", type);
262 return F_POLARIZATION;
266 gmx_fatal(FARGS, "Invalid polarization type %d", type);
269 case d_thole_polarization:
271 case d_water_polarization:
273 case d_angle_restraints:
275 case d_angle_restraints_z:
277 case d_distance_restraints:
279 case d_orientation_restraints:
281 case d_dihedral_restraints:
284 gmx_fatal(FARGS, "invalid directive %s in ifunc_index (%s:%s)",
285 dir2str(d), __FILE__, __LINE__);
290 const char *dir2str (directive d)
294 return directive_names[d];
298 return directive_names[d_maxdir];
302 directive str2dir (char *dstr)
305 char buf[STRLEN], *ptr;
307 /* Hack to be able to read old topologies */
308 if (gmx_strncasecmp_min(dstr, "dummies", 7) == 0)
310 sprintf(buf, "virtual_sites%s", dstr+7);
318 for (d = 0; (d < d_maxdir); d++)
320 if (gmx_strcasecmp_min(ptr, dir2str((directive)d)) == 0)
329 static directive **necessary = NULL;
331 static void set_nec(directive **n, ...)
332 /* Must always have at least one extra argument */
341 d = (directive)va_arg(ap, int);
349 void DS_Init(DirStack **DS)
351 if (necessary == NULL)
355 snew(necessary, d_maxdir);
356 set_nec(&(necessary[d_defaults]), d_none);
357 set_nec(&(necessary[d_atomtypes]), d_defaults, d_none);
358 set_nec(&(necessary[d_bondtypes]), d_atomtypes, d_none);
359 set_nec(&(necessary[d_constrainttypes]), d_atomtypes, d_none);
360 set_nec(&(necessary[d_pairtypes]), d_atomtypes, d_none);
361 set_nec(&(necessary[d_angletypes]), d_atomtypes, d_none);
362 set_nec(&(necessary[d_dihedraltypes]), d_atomtypes, d_none);
363 set_nec(&(necessary[d_nonbond_params]), d_atomtypes, d_none);
364 set_nec(&(necessary[d_implicit_genborn_params]), d_atomtypes, d_none);
365 set_nec(&(necessary[d_implicit_surface_params]), d_atomtypes, d_none);
366 set_nec(&(necessary[d_cmaptypes]), d_atomtypes, d_none);
367 set_nec(&(necessary[d_moleculetype]), d_atomtypes, d_none);
368 set_nec(&(necessary[d_atoms]), d_moleculetype, d_none);
369 set_nec(&(necessary[d_vsites2]), d_atoms, d_none);
370 set_nec(&(necessary[d_vsites3]), d_atoms, d_none);
371 set_nec(&(necessary[d_vsites4]), d_atoms, d_none);
372 set_nec(&(necessary[d_vsitesn]), d_atoms, d_none);
373 set_nec(&(necessary[d_bonds]), d_atoms, d_none);
374 set_nec(&(necessary[d_exclusions]), d_bonds, d_constraints, d_settles, d_none);
375 set_nec(&(necessary[d_pairs]), d_atoms, d_none);
376 set_nec(&(necessary[d_pairs_nb]), d_atoms, d_none);
377 set_nec(&(necessary[d_angles]), d_atoms, d_none);
378 set_nec(&(necessary[d_polarization]), d_atoms, d_none);
379 set_nec(&(necessary[d_water_polarization]), d_atoms, d_none);
380 set_nec(&(necessary[d_thole_polarization]), d_atoms, d_none);
381 set_nec(&(necessary[d_dihedrals]), d_atoms, d_none);
382 set_nec(&(necessary[d_constraints]), d_atoms, d_none);
383 set_nec(&(necessary[d_settles]), d_atoms, d_none);
384 set_nec(&(necessary[d_system]), d_moleculetype, d_none);
385 set_nec(&(necessary[d_molecules]), d_system, d_none);
386 set_nec(&(necessary[d_position_restraints]), d_atoms, d_none);
387 set_nec(&(necessary[d_angle_restraints]), d_atoms, d_none);
388 set_nec(&(necessary[d_angle_restraints_z]), d_atoms, d_none);
389 set_nec(&(necessary[d_distance_restraints]), d_atoms, d_none);
390 set_nec(&(necessary[d_orientation_restraints]), d_atoms, d_none);
391 set_nec(&(necessary[d_dihedral_restraints]), d_atoms, d_none);
392 set_nec(&(necessary[d_cmap]), d_atoms, d_none);
394 for (i = 0; (i < d_maxdir); i++)
398 fprintf(debug, "%20s: ", dir2str((directive)i));
406 d = necessary[i][j++];
409 fprintf(debug, "%20s ", dir2str(d));
416 fprintf(debug, "\n");
424 void DS_Done (DirStack **DS)
436 void DS_Push (DirStack **DS, directive d)
446 int DS_Search(DirStack *DS, directive d)
451 while ((D != NULL) && (D->d != d))
459 int DS_Check_Order(DirStack *DS, directive d)
464 /* Check if parameter definitions appear after a moleculetype directive */
465 if (d < d_moleculetype && DS_Search(DS, d_moleculetype))
470 /* Check if all the necessary directives have appeared before directive d */
471 if (necessary[d][0] == d_none)
479 d0 = necessary[d][i++];
480 if (DS_Search(DS, d0))
485 while (d0 != d_none);