Replaced NOTSET from typedefs.h by local solutions.
[gromacs.git] / src / gromacs / gmxpreprocess / readadress.cpp
blob459dddb0728bdfcec978a281563efb049502fccf
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009 Christoph Junghans, Brad Lambeth.
5 * Copyright (c) 2012,2014,2015, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #include "gmxpre.h"
38 #include <stdlib.h>
39 #include <string.h>
41 #include "gromacs/gmxpreprocess/readir.h"
42 #include "gromacs/legacyheaders/names.h"
43 #include "gromacs/legacyheaders/types/inputrec.h"
44 #include "gromacs/utility/cstringutil.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
48 #define MAXPTR 254
50 static char adress_refs[STRLEN], adress_tf_grp_names[STRLEN], adress_cg_grp_names[STRLEN];
52 void read_adressparams(int *ninp_p, t_inpfile **inp_p, t_adress *adress, warninp_t wi)
54 int nadress_refs, i;
55 const char *tmp;
56 char *ptr1[MAXPTR];
59 int ninp;
60 t_inpfile *inp;
62 ninp = *ninp_p;
63 inp = *inp_p;
65 EETYPE("adress_type", adress->type, eAdresstype_names);
66 RTYPE ("adress_const_wf", adress->const_wf, 1);
67 RTYPE ("adress_ex_width", adress->ex_width, 0);
68 RTYPE ("adress_hy_width", adress->hy_width, 0);
69 RTYPE ("adress_ex_forcecap", adress->ex_forcecap, 0);
70 EETYPE("adress_interface_correction", adress->icor, eAdressICtype_names);
71 EETYPE("adress_site", adress->site, eAdressSITEtype_names);
72 STYPE ("adress_reference_coords", adress_refs, NULL);
73 STYPE ("adress_tf_grp_names", adress_tf_grp_names, NULL);
74 STYPE ("adress_cg_grp_names", adress_cg_grp_names, NULL);
75 EETYPE("adress_do_hybridpairs", adress->do_hybridpairs, yesno_names);
77 nadress_refs = str_nelem(adress_refs, MAXPTR, ptr1);
79 for (i = 0; (i < nadress_refs); i++) /*read vector components*/
81 adress->refs[i] = strtod(ptr1[i], NULL);
83 for (; (i < DIM); i++) /*remaining undefined components of the vector set to zero*/
85 adress->refs[i] = 0;
88 *ninp_p = ninp;
89 *inp_p = inp;
92 void do_adress_index(t_adress *adress, gmx_groups_t *groups, char **gnames, t_grpopts *opts, warninp_t wi)
94 int nr, i, j, k;
95 char *ptr1[MAXPTR];
96 int nadress_cg_grp_names, nadress_tf_grp_names;
98 /* AdResS coarse grained groups input */
100 nr = groups->grps[egcENER].nr;
101 snew(adress->group_explicit, nr);
102 for (i = 0; i < nr; i++)
104 adress->group_explicit[i] = TRUE;
106 adress->n_energy_grps = nr;
108 nadress_cg_grp_names = str_nelem(adress_cg_grp_names, MAXPTR, ptr1);
110 if (nadress_cg_grp_names > 0)
112 for (i = 0; i < nadress_cg_grp_names; i++)
114 /* search for the group name mathching the tf group name */
115 k = 0;
116 while ((k < nr) &&
117 gmx_strcasecmp(ptr1[i], (char*)(gnames[groups->grps[egcENER].nm_ind[k]])))
119 k++;
121 if (k == nr)
123 gmx_fatal(FARGS, "Adress cg energy group %s not found\n", ptr1[i]);
125 adress->group_explicit[k] = FALSE;
126 printf ("AdResS: Energy group %s is treated as coarse-grained \n",
127 (char*)(gnames[groups->grps[egcENER].nm_ind[k]]));
129 /* set energy group exclusions between all coarse-grained and explicit groups */
130 for (j = 0; j < nr; j++)
132 for (k = 0; k < nr; k++)
134 if (!(adress->group_explicit[k] == adress->group_explicit[j]))
136 opts->egp_flags[nr * j + k] |= EGP_EXCL;
137 if (debug)
139 fprintf(debug, "AdResS excl %s %s \n",
140 (char*)(gnames[groups->grps[egcENER].nm_ind[j]]),
141 (char*)(gnames[groups->grps[egcENER].nm_ind[k]]));
147 else
149 warning(wi, "For an AdResS simulation at least one coarse-grained energy group has to be specified in adress_cg_grp_names");
153 /* AdResS multiple tf tables input */
154 nadress_tf_grp_names = str_nelem(adress_tf_grp_names, MAXPTR, ptr1);
155 adress->n_tf_grps = nadress_tf_grp_names;
156 snew(adress->tf_table_index, nadress_tf_grp_names);
158 nr = groups->grps[egcENER].nr;
160 if (nadress_tf_grp_names > 0)
162 for (i = 0; i < nadress_tf_grp_names; i++)
164 /* search for the group name mathching the tf group name */
165 k = 0;
166 while ((k < nr) &&
167 gmx_strcasecmp(ptr1[i], (char*)(gnames[groups->grps[egcENER].nm_ind[k]])))
169 k++;
171 if (k == nr)
173 gmx_fatal(FARGS, "Adress tf energy group %s not found\n", ptr1[i]);
176 adress->tf_table_index[i] = k;
177 if (debug)
179 fprintf(debug, "found tf group %s id %d \n", ptr1[i], k);
181 if (adress->group_explicit[k])
183 gmx_fatal(FARGS, "Thermodynamic force group %s is not a coarse-grained group in adress_cg_grp_names. The thermodynamic force has to act on the coarse-grained vsite of a molecule.\n", ptr1[i]);
188 /* end AdResS multiple tf tables input */