Move initiation of local CPU force H2D transfer to producer
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / topdirs.cpp
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,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "topdirs.h"
41
42 #include <cstdarg>
43 #include <cstdio>
44
45 #include <algorithm>
46
47 #include "gromacs/topology/idef.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/enumerationhelpers.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
52
53 /* Must correspond to the Directive enum in grompp_impl.h */
54 static gmx::EnumerationArray<Directive, const char*> directive_names = {
55     { "defaults",
56       "atomtypes",
57       "bondtypes",
58       "constrainttypes",
59       "pairtypes",
60       "angletypes",
61       "dihedraltypes",
62       "nonbond_params",
63       "implicit_genborn_params",
64       "implicit_surface_params",
65       "cmaptypes",
66       /* All the directives above can not appear after moleculetype */
67       "moleculetype",
68       "atoms",
69       "virtual_sites2",
70       "virtual_sites3",
71       "virtual_sites4",
72       "virtual_sitesn",
73       "bonds",
74       "exclusions",
75       "pairs",
76       "pairs_nb",
77       "angles",
78       "dihedrals",
79       "constraints",
80       "settles",
81       "polarization",
82       "water_polarization",
83       "thole_polarization",
84       "system",
85       "molecules",
86       "position_restraints",
87       "angle_restraints",
88       "angle_restraints_z",
89       "distance_restraints",
90       "orientation_restraints",
91       "dihedral_restraints",
92       "cmap",
93       "intermolecular_interactions",
94       "maxdirs",
95       "invalid",
96       "none" }
97 };
98
99 int ifunc_index(Directive d, int type)
100 {
101     switch (d)
102     {
103         case Directive::d_bondtypes:
104         case Directive::d_bonds:
105             switch (type)
106             {
107                 case 1: return F_BONDS;
108                 case 2: return F_G96BONDS;
109                 case 3: return F_MORSE;
110                 case 4: return F_CUBICBONDS;
111                 case 5: return F_CONNBONDS;
112                 case 6: return F_HARMONIC;
113                 case 7: return F_FENEBONDS;
114                 case 8: return F_TABBONDS;
115                 case 9: return F_TABBONDSNC;
116                 case 10: return F_RESTRBONDS;
117                 default: gmx_fatal(FARGS, "Invalid bond type %d", type);
118             }
119         case Directive::d_angles:
120         case Directive::d_angletypes:
121             switch (type)
122             {
123                 case 1: return F_ANGLES;
124                 case 2: return F_G96ANGLES;
125                 case 3: return F_CROSS_BOND_BONDS;
126                 case 4: return F_CROSS_BOND_ANGLES;
127                 case 5: return F_UREY_BRADLEY;
128                 case 6: return F_QUARTIC_ANGLES;
129                 case 8: return F_TABANGLES;
130                 case 9: return F_LINEAR_ANGLES;
131                 case 10: return F_RESTRANGLES;
132                 default: gmx_fatal(FARGS, "Invalid angle type %d", type);
133             }
134         case Directive::d_pairs:
135         case Directive::d_pairtypes:
136             if (type == 1 || (d == Directive::d_pairtypes && type == 2))
137             {
138                 return F_LJ14;
139             }
140             else if (type == 2)
141             {
142                 return F_LJC14_Q;
143             }
144             else
145             {
146                 gmx_fatal(FARGS, "Invalid pairs type %d", type);
147             }
148         case Directive::d_pairs_nb: return F_LJC_PAIRS_NB;
149         case Directive::d_dihedrals:
150         case Directive::d_dihedraltypes:
151             switch (type)
152             {
153                 case 1: return F_PDIHS;
154                 case 2: return F_IDIHS;
155                 case 3: return F_RBDIHS;
156                 case 4: return F_PIDIHS;
157                 case 5: return F_FOURDIHS;
158                 case 8: return F_TABDIHS;
159                 case 9:
160                     return F_PDIHS; /* proper dihedrals where we allow multiple terms over single bond */
161                 case 10: return F_RESTRDIHS;
162                 case 11: return F_CBTDIHS;
163                 default: gmx_fatal(FARGS, "Invalid dihedral type %d", type);
164             }
165         case Directive::d_cmaptypes:
166         case Directive::d_cmap: return F_CMAP;
167
168         case Directive::d_nonbond_params:
169             if (type == 1)
170             {
171                 return F_LJ;
172             }
173             else
174             {
175                 return F_BHAM;
176             }
177         case Directive::d_vsites2:
178             switch (type)
179             {
180                 case 1: return F_VSITE2;
181                 case 2: return F_VSITE2FD;
182                 default: gmx_fatal(FARGS, "Invalid vsites2 type %d", type);
183             }
184         case Directive::d_vsites3:
185             switch (type)
186             {
187                 case 1: return F_VSITE3;
188                 case 2: return F_VSITE3FD;
189                 case 3: return F_VSITE3FAD;
190                 case 4: return F_VSITE3OUT;
191                 default: gmx_fatal(FARGS, "Invalid vsites3 type %d", type);
192             }
193         case Directive::d_vsites4:
194             switch (type)
195             {
196                 case 1: return F_VSITE4FD;
197                 case 2: return F_VSITE4FDN;
198                 default: gmx_fatal(FARGS, "Invalid vsites4 type %d", type);
199             }
200         case Directive::d_vsitesn: return F_VSITEN;
201         case Directive::d_constraints:
202         case Directive::d_constrainttypes:
203             switch (type)
204             {
205                 case 1: return F_CONSTR;
206                 case 2: return F_CONSTRNC;
207                 default: gmx_fatal(FARGS, "Invalid constraints type %d", type);
208             }
209         case Directive::d_settles: return F_SETTLE;
210         case Directive::d_position_restraints:
211             switch (type)
212             {
213                 case 1: return F_POSRES;
214                 case 2: return F_FBPOSRES;
215                 default: gmx_fatal(FARGS, "Invalid position restraint type %d", type);
216             }
217         case Directive::d_polarization:
218             switch (type)
219             {
220                 case 1: return F_POLARIZATION;
221                 case 2: return F_ANHARM_POL;
222                 default: gmx_fatal(FARGS, "Invalid polarization type %d", type);
223             }
224         case Directive::d_thole_polarization: return F_THOLE_POL;
225         case Directive::d_water_polarization: return F_WATER_POL;
226         case Directive::d_angle_restraints: return F_ANGRES;
227         case Directive::d_angle_restraints_z: return F_ANGRESZ;
228         case Directive::d_distance_restraints: return F_DISRES;
229         case Directive::d_orientation_restraints: return F_ORIRES;
230         case Directive::d_dihedral_restraints: return F_DIHRES;
231         default:
232             gmx_fatal(FARGS, "invalid directive %s in ifunc_index (%s:%d)", dir2str(d), __FILE__, __LINE__);
233     }
234 }
235
236 const char* dir2str(Directive d)
237 {
238     int index = static_cast<int>(d);
239     return directive_names[index];
240 }
241
242 Directive str2dir(char* dstr)
243 {
244     char buf[STRLEN], *ptr;
245
246     /* Hack to be able to read old topologies */
247     if (gmx_strncasecmp_min(dstr, "dummies", 7) == 0)
248     {
249         sprintf(buf, "virtual_sites%s", dstr + 7);
250         ptr = buf;
251     }
252     else
253     {
254         ptr = dstr;
255     }
256
257     for (auto d : gmx::EnumerationWrapper<Directive>())
258     {
259         if (gmx_strcasecmp_min(ptr, dir2str(static_cast<Directive>(d))) == 0)
260         {
261             return static_cast<Directive>(d);
262         }
263     }
264
265     return Directive::d_invalid;
266 }
267
268 static gmx::EnumerationArray<Directive, Directive*> necessary = { { nullptr } };
269
270 static void set_nec(Directive** n, ...)
271 /* Must always have at least one extra argument */
272 {
273     va_list   ap;
274     int       ind = 0;
275     Directive d;
276
277     va_start(ap, n);
278     do
279     {
280         d = static_cast<Directive>(va_arg(ap, int));
281         srenew(*n, ++ind);
282         (*n)[ind - 1] = d;
283     } while (d != Directive::d_none);
284     va_end(ap);
285 }
286
287 void DS_Init(DirStack** DS)
288 {
289     if (necessary[0] == nullptr)
290     {
291         set_nec(&(necessary[Directive::d_defaults]), Directive::d_none);
292         set_nec(&(necessary[Directive::d_atomtypes]), Directive::d_defaults, Directive::d_none);
293         set_nec(&(necessary[Directive::d_bondtypes]), Directive::d_atomtypes, Directive::d_none);
294         set_nec(&(necessary[Directive::d_constrainttypes]), Directive::d_atomtypes, Directive::d_none);
295         set_nec(&(necessary[Directive::d_pairtypes]), Directive::d_atomtypes, Directive::d_none);
296         set_nec(&(necessary[Directive::d_angletypes]), Directive::d_atomtypes, Directive::d_none);
297         set_nec(&(necessary[Directive::d_dihedraltypes]), Directive::d_atomtypes, Directive::d_none);
298         set_nec(&(necessary[Directive::d_nonbond_params]), Directive::d_atomtypes, Directive::d_none);
299         // Note that the content of the next two directives are
300         // ignored, but if grompp reads them in old force field files,
301         // it still needs to understand that they are in a valid place
302         // in the .top structure. It doesn't have to require them to
303         // be in the same place that was valid in old versions (ie. child
304         // directive of [atomtypes]) but any relevant case will
305         // satisfy that.
306         set_nec(&(necessary[Directive::d_implicit_genborn_params]), Directive::d_atomtypes, Directive::d_none);
307         set_nec(&(necessary[Directive::d_implicit_surface_params]), Directive::d_atomtypes, Directive::d_none);
308         set_nec(&(necessary[Directive::d_cmaptypes]), Directive::d_atomtypes, Directive::d_none);
309         set_nec(&(necessary[Directive::d_moleculetype]), Directive::d_atomtypes, Directive::d_none);
310         set_nec(&(necessary[Directive::d_atoms]), Directive::d_moleculetype, Directive::d_none);
311         set_nec(&(necessary[Directive::d_vsites2]), Directive::d_atoms, Directive::d_none);
312         set_nec(&(necessary[Directive::d_vsites3]), Directive::d_atoms, Directive::d_none);
313         set_nec(&(necessary[Directive::d_vsites4]), Directive::d_atoms, Directive::d_none);
314         set_nec(&(necessary[Directive::d_vsitesn]), Directive::d_atoms, Directive::d_none);
315         set_nec(&(necessary[Directive::d_bonds]), Directive::d_atoms, Directive::d_none);
316         set_nec(&(necessary[Directive::d_exclusions]),
317                 Directive::d_bonds,
318                 Directive::d_constraints,
319                 Directive::d_settles,
320                 Directive::d_none);
321         set_nec(&(necessary[Directive::d_pairs]), Directive::d_atoms, Directive::d_none);
322         set_nec(&(necessary[Directive::d_pairs_nb]), Directive::d_atoms, Directive::d_none);
323         set_nec(&(necessary[Directive::d_angles]), Directive::d_atoms, Directive::d_none);
324         set_nec(&(necessary[Directive::d_polarization]), Directive::d_atoms, Directive::d_none);
325         set_nec(&(necessary[Directive::d_water_polarization]), Directive::d_atoms, Directive::d_none);
326         set_nec(&(necessary[Directive::d_thole_polarization]), Directive::d_atoms, Directive::d_none);
327         set_nec(&(necessary[Directive::d_dihedrals]), Directive::d_atoms, Directive::d_none);
328         set_nec(&(necessary[Directive::d_constraints]), Directive::d_atoms, Directive::d_none);
329         set_nec(&(necessary[Directive::d_settles]), Directive::d_atoms, Directive::d_none);
330         set_nec(&(necessary[Directive::d_system]), Directive::d_moleculetype, Directive::d_none);
331         set_nec(&(necessary[Directive::d_molecules]), Directive::d_system, Directive::d_none);
332         set_nec(&(necessary[Directive::d_position_restraints]), Directive::d_atoms, Directive::d_none);
333         set_nec(&(necessary[Directive::d_angle_restraints]), Directive::d_atoms, Directive::d_none);
334         set_nec(&(necessary[Directive::d_angle_restraints_z]), Directive::d_atoms, Directive::d_none);
335         set_nec(&(necessary[Directive::d_distance_restraints]), Directive::d_atoms, Directive::d_none);
336         set_nec(&(necessary[Directive::d_orientation_restraints]), Directive::d_atoms, Directive::d_none);
337         set_nec(&(necessary[Directive::d_dihedral_restraints]), Directive::d_atoms, Directive::d_none);
338         set_nec(&(necessary[Directive::d_cmap]), Directive::d_atoms, Directive::d_none);
339         set_nec(&(necessary[Directive::d_intermolecular_interactions]),
340                 Directive::d_molecules,
341                 Directive::d_none);
342     }
343     *DS = nullptr;
344 }
345
346 void DS_Done(DirStack** DS)
347 {
348     DirStack* D;
349
350     while (*DS != nullptr)
351     {
352         D   = *DS;
353         *DS = (*DS)->prev;
354         sfree(D);
355     }
356 }
357
358 void DS_Push(DirStack** DS, Directive d)
359 {
360     DirStack* D;
361
362     snew(D, 1);
363     D->d    = d;
364     D->prev = *DS;
365     *DS     = D;
366 }
367
368 int DS_Search(DirStack* DS, Directive d)
369 {
370     DirStack* D;
371
372     D = DS;
373     while ((D != nullptr) && (D->d != d))
374     {
375         D = D->prev;
376     }
377
378     return static_cast<int>(D != nullptr);
379 }
380
381 int DS_Check_Order(DirStack* DS, Directive d)
382 {
383     Directive d0;
384     int       i = 0;
385
386     /* Check if parameter definitions appear after a moleculetype directive */
387     if (d < Directive::d_moleculetype && DS_Search(DS, Directive::d_moleculetype))
388     {
389         return FALSE;
390     }
391
392     /* Check if all the necessary directives have appeared before directive d */
393     if (necessary[d][0] == Directive::d_none)
394     {
395         return TRUE;
396     }
397     else
398     {
399         do
400         {
401             d0 = necessary[d][i++];
402             if (DS_Search(DS, d0))
403             {
404                 return TRUE;
405             }
406         } while (d0 != Directive::d_none);
407     }
408     return FALSE;
409 }