6427c9d26c64bc83b45bf707d34aae55706c9721
[alexxy/gromacs.git] / src / gromacs / mdtypes / nblist.h
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) 2012,2014,2015,2019,2020,2021, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_MDTYPES_NBLIST_H
38 #define GMX_MDTYPES_NBLIST_H
39
40 #include <vector>
41
42 #include "gromacs/utility/alignedallocator.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
45
46 typedef unsigned long t_excl;
47
48 /* The interactions contained in a (possibly merged) table
49  * for computing electrostatic, VDW repulsion and/or VDW dispersion
50  * contributions.
51  */
52 enum gmx_table_interaction
53 {
54     GMX_TABLE_INTERACTION_ELEC,
55     GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
56     GMX_TABLE_INTERACTION_VDWEXPREP_VDWDISP,
57     GMX_TABLE_INTERACTION_VDWDISP,
58     GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP,
59     GMX_TABLE_INTERACTION_ELEC_VDWEXPREP_VDWDISP,
60     GMX_TABLE_INTERACTION_ELEC_VDWDISP,
61     GMX_TABLE_INTERACTION_NR
62 };
63
64 /* Different formats for table data. Cubic spline tables are typically stored
65  * with the four Y,F,G,H intermediate values (check tables.c for format), which
66  * makes it easy to load with a single 4-way SIMD instruction too.
67  * Linear tables only need one value per table point, or two if both V and F
68  * are calculated. However, with SIMD instructions this makes the loads unaligned,
69  * and in that case we store the data as F, D=F(i+1)-F(i), V, and then a blank value,
70  * which again makes it possible to load as a single instruction.
71  */
72 enum gmx_table_format
73 {
74     GMX_TABLE_FORMAT_CUBICSPLINE_YFGH,
75     GMX_TABLE_FORMAT_LINEAR_VF,
76     GMX_TABLE_FORMAT_LINEAR_V,
77     GMX_TABLE_FORMAT_LINEAR_F,
78     GMX_TABLE_FORMAT_LINEAR_FDV0,
79     GMX_TABLE_FORMAT_NR
80 };
81
82 enum
83 {
84     eNL_VDWQQ,
85     eNL_VDW,
86     eNL_QQ,
87     eNL_VDWQQ_FREE,
88     eNL_VDW_FREE,
89     eNL_QQ_FREE,
90     eNL_VDWQQ_WATER,
91     eNL_QQ_WATER,
92     eNL_VDWQQ_WATERWATER,
93     eNL_QQ_WATERWATER,
94     eNL_NR
95 };
96
97 #define MAX_CG 1024
98
99 typedef struct
100 {
101     int ncg;
102     int nj;
103     int jcg[MAX_CG];
104 } t_ns_buf;
105
106
107 /* The maximum charge group size because of minimum size of t_excl
108  * could be 32 bits.
109  */
110 #define MAX_CHARGEGROUP_SIZE 32
111
112 /* The maximum charge group size for CG-CG nblists.
113  * The excl entry in t_nblist uses blocks of this size.
114  */
115 #define MAX_CGCGSIZE 32
116
117 typedef struct t_nblist
118 {
119     int igeometry; /* The type of list (atom, water, etc.)  */
120     int ielec;     /* Coulomb loop type index for kernels   */
121     int ielecmod;  /* Coulomb modifier (e.g. switch/shift)  */
122     int ivdw;      /* VdW loop type index for kernels       */
123     int ivdwmod;   /* VdW modifier (e.g. switch/shift)      */
124     int type;      /* Type of interaction, listed in
125                       gmx_nblist_interaction_type           */
126
127     int     nri, maxnri; /* Current/max number of i particles      */
128     int     nrj, maxnrj; /* Current/max number of j particles      */
129     int*    iinr;        /* The i-elements                        */
130     int*    iinr_end;    /* The end atom, only with enlistCG      */
131     int*    gid;         /* Index in energy arrays                */
132     int*    shift;       /* Shift vector index                    */
133     int*    jindex;      /* Index in jjnr                         */
134     int*    jjnr;        /* The j-atom list                       */
135     int*    jjnr_end;    /* The end atom, only with enltypeCG     */
136     char*   excl_fep;    /* Exclusions for FEP with Verlet scheme */
137     t_excl* excl;        /* Exclusions, only with enltypeCG       */
138
139     /* We use separate pointers for kernels that compute both potential
140      * and force (vf suffix), only potential (v) or only force (f)
141      */
142     void* kernelptr_vf;
143     void* kernelptr_v;
144     void* kernelptr_f;
145
146     /* Pad the list of neighbors for each i atom with "-1" entries up to the
147      * simd_padding_width, if it is larger than 0. This is necessary for many
148      * accelerated kernels using single-instruction multiple-data operations
149      * internally.
150      */
151     int simd_padding_width;
152
153 } t_nblist;
154
155 /* For atom I =  nblist->iinr[N] (0 <= N < nblist->nri) there can be
156  * several neighborlists (N's), for different energy groups (gid) and
157  * different shifts (shift).
158  * For corresponding J atoms for each list start at:
159  * nblist->jjnr[JI]
160  * with nblist->jindex[N] <= JI < nblist->jindex[N+1]
161  *
162  * enlist is of the form enlistUNIT1_UNIT2:
163  * UNIT ATOM:  there is one atom: iinr[N] or jjnr[JI]
164  * UNIT SPC:   there are 3 atoms: iinr[N],iinr[N]+1,iinr[N]+2, jjnr analog.
165  * UNIT TIP4P: there are 4 atoms: iinr[N],...,iinr[N]+3, jjnr analog.
166  * UNIT CG:    there are N atoms: iinr[N],...,iinr_end[N]-1, jjnr analog.
167  *
168  * Clear?
169  */
170
171 /* Structure describing the data in a single table */
172 struct t_forcetable
173 {
174     t_forcetable(enum gmx_table_interaction interaction, enum gmx_table_format format);
175
176     ~t_forcetable();
177
178     enum gmx_table_interaction interaction; /* Types of interactions stored in this table */
179     enum gmx_table_format      format;      /* Interpolation type and data format */
180
181     real r;     /* range of the table */
182     int  n;     /* n+1 is the number of table points */
183     real scale; /* distance (nm) between two table points */
184     std::vector<real, gmx::AlignedAllocator<real>> data; /* the actual table data */
185
186     /* Some information about the table layout. This can also be derived from the interpolation
187      * type and the table interactions, but it is convenient to have here for sanity checks, and it
188      * makes it much easier to access the tables in the nonbonded kernels when we can set the data
189      * from variables. It is always true that stride = formatsize*ninteractions
190      */
191     int formatsize; /* Number of fp variables for each table point (1 for F, 2 for VF, 4 for YFGH, etc.) */
192     int ninteractions; /* Number of interactions in table, 1 for coul-only, 3 for coul+rep+disp. */
193     int stride; /* Distance to next table point (number of fp variables per table point in total) */
194 };
195
196 struct gmx_ns_t
197 {
198     gmx_bool   bCGlist;
199     int*       simple_aaj;
200     t_excl*    bexcl;
201     gmx_bool*  bHaveVdW;
202     t_ns_buf** ns_buf;
203     gmx_bool*  bExcludeAlleg;
204     int        nra_alloc;
205     int        cg_alloc;
206     int**      nl_sr;
207     int*       nsr;
208     int**      nl_lr_ljc;
209     int**      nl_lr_one;
210     int*       nlr_ljc;
211     int*       nlr_one;
212     /* the nblists should probably go in here */
213     gmx_bool nblist_initialized; /* has the nblist been initialized?  */
214     int      dump_nl;            /* neighbour list dump level (from env. var. GMX_DUMP_NL)*/
215 };
216
217 #endif /* GMX_MDTYPES_NBLIST_H */