Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / include / types / nbnxn_pairlist.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-2012, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * 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
39 #ifndef _nbnxn_pairlist_h
40 #define _nbnxn_pairlist_h
41
42 #ifdef __cplusplus
43 extern "C" {
44 #endif
45
46 /* A buffer data structure of 64 bytes
47  * to be placed at the beginning and end of structs
48  * to avoid cache invalidation of the real contents
49  * of the struct by writes to neighboring memory.
50  */
51 typedef struct {
52     int dummy[16];
53 } gmx_cache_protect_t;
54
55 /* Abstract type for pair searching data */
56 typedef struct nbnxn_search * nbnxn_search_t;
57
58 /* Function that should return a pointer *ptr to memory
59  * of size nbytes.
60  * Error handling should be done within this function.
61  */
62 typedef void nbnxn_alloc_t(void **ptr,size_t nbytes);
63
64 /* Function that should free the memory pointed to by *ptr.
65  * NULL should not be passed to this function.
66  */
67 typedef void nbnxn_free_t(void *ptr);
68
69 typedef struct {
70     int      cj;    /* The j-cluster                    */
71     unsigned excl;  /* The exclusion (interaction) bits */
72 } nbnxn_cj_t;
73
74 #define NBNXN_CI_SHIFT          127
75 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
76 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
77 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
78
79 /* Simple pair-list i-unit */
80 typedef struct {
81     int ci;             /* i-cluster             */
82     int shift;          /* Shift vector index plus possible flags */
83     int cj_ind_start;   /* Start index into cj   */
84     int cj_ind_end;     /* End index into cj     */
85 } nbnxn_ci_t;
86
87 /* Grouped pair-list i-unit */
88 typedef struct {
89     int sci;            /* i-super-cluster       */
90     int shift;          /* Shift vector index plus possible flags */
91     int cj4_ind_start;  /* Start index into cj4  */
92     int cj4_ind_end;    /* End index into cj4    */
93 } nbnxn_sci_t;
94
95 typedef struct {
96     unsigned imask;        /* The i-cluster interactions mask for 1 warp  */
97     int excl_ind;          /* Index into the exclusion array for 1 warp   */
98 } nbnxn_im_ei_t;
99
100 typedef struct {
101     int cj[4];             /* The 4 j-clusters                            */
102     nbnxn_im_ei_t imei[2]; /* The i-cluster mask data       for 2 warps   */
103 } nbnxn_cj4_t;
104
105 typedef struct {
106     unsigned pair[32];     /* Exclusion bits for one warp,                *
107                             * each unsigned has bit for 4*8 i clusters    */
108 } nbnxn_excl_t;
109
110 typedef struct {
111     gmx_cache_protect_t cp0;
112
113     nbnxn_alloc_t *alloc;
114     nbnxn_free_t  *free;
115
116     gmx_bool bSimple;      /* Simple list has na_sc=na_s and uses cj   *
117                             * Complex list uses cj4                    */
118
119     int      na_ci;        /* The number of atoms per i-cluster        */
120     int      na_cj;        /* The number of atoms per j-cluster        */
121     int      na_sc;        /* The number of atoms per super cluster    */
122     real     rlist;        /* The radius for constructing the list     */
123     int      nci;          /* The number of i-clusters in the list     */
124     nbnxn_ci_t *ci;        /* The i-cluster list, size nci             */
125     int      ci_nalloc;    /* The allocation size of ci                */
126     int      nsci;         /* The number of i-super-clusters in the list */
127     nbnxn_sci_t *sci;      /* The i-super-cluster list                 */
128     int      sci_nalloc;   /* The allocation size of sci               */
129
130     int      ncj;          /* The number of j-clusters in the list     */
131     nbnxn_cj_t *cj;        /* The j-cluster list, size ncj             */
132     int      cj_nalloc;    /* The allocation size of cj                */
133
134     int      ncj4;         /* The total number of 4*j clusters         */
135     nbnxn_cj4_t *cj4;      /* The 4*j cluster list, size ncj4          */
136     int      cj4_nalloc;   /* The allocation size of cj4               */
137     int      nexcl;        /* The count for excl                       */
138     nbnxn_excl_t *excl;    /* Atom interaction bits (non-exclusions)   */
139     int      excl_nalloc;  /* The allocation size for excl             */
140     int      nci_tot;      /* The total number of i clusters           */
141
142     struct nbnxn_list_work *work;
143
144     gmx_cache_protect_t cp1;
145 } nbnxn_pairlist_t;
146
147 typedef struct {
148     int          nnbl;      /* number of lists */
149     nbnxn_pairlist_t **nbl; /* lists */
150     gmx_bool     bCombined; /* TRUE if lists get combined into one (the 1st) */
151     gmx_bool     bSimple;   /* TRUE if the list of of type "simple"
152                                (na_sc=na_s, no super-clusters used) */
153     int          natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
154     int          natpair_lj;  /* Total number of atom pairs for LJ kernel   */
155     int          natpair_q;   /* Total number of atom pairs for Q kernel    */
156 } nbnxn_pairlist_set_t;
157
158 enum { nbatXYZ, nbatXYZQ, nbatX4, nbatX8 };
159
160 typedef struct {
161     real *f;      /* f, size natoms*fstride                             */
162     real *fshift; /* Shift force array, size SHIFTS*DIM                 */
163     int  nV;      /* The size of *Vvdw and *Vc                          */
164     real *Vvdw;   /* Temporary Van der Waals group energy storage       */
165     real *Vc;     /* Temporary Coulomb group energy storage             */
166     int  nVS;     /* The size of *VSvdw and *VSc                        */
167     real *VSvdw;  /* Temporary SIMD Van der Waals group energy storage  */
168     real *VSc;    /* Temporary SIMD Coulomb group energy storage        */
169 } nbnxn_atomdata_output_t;
170
171 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
172 enum { ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR };
173
174 typedef struct {
175     nbnxn_alloc_t *alloc;
176     nbnxn_free_t  *free;
177     int  ntype;      /* The number of different atom types                 */
178     real *nbfp;      /* Lennard-Jones 6*C6 and 12*C12 params, size ntype^2*2 */
179     int  comb_rule;  /* Combination rule, see enum above                   */
180     real *nbfp_comb; /* LJ parameter per atom type, size ntype*2           */
181     real *nbfp_s4;   /* As nbfp, but with stride 4, size ntype^2*4         */
182     int  natoms;     /* Number of atoms                                    */
183     int  natoms_local;  /* Number of local atoms                           */
184     int  *type;      /* Atom types                                         */
185     real *lj_comb;   /* LJ parameters per atom for combining for pairs     */
186     int  XFormat;    /* The format of x (and q), enum                      */
187     int  FFormat;    /* The format of f, enum                              */
188     real *q;         /* Charges, can be NULL if incorporated in x          */
189     int  na_c;       /* The number of atoms per cluster                    */
190     int  nenergrp;   /* The number of energy groups                        */
191     int  neg_2log;   /* Log2 of nenergrp                                   */
192     int  *energrp;   /* The energy groups per cluster, can be NULL         */
193     gmx_bool bDynamicBox; /* Do we need to update shift_vec every step?    */
194     rvec *shift_vec; /* Shift vectors, copied from t_forcerec              */
195     int  xstride;    /* stride for a coordinate in x (usually 3 or 4)      */
196     int  fstride;    /* stride for a coordinate in f (usually 3 or 4)      */
197     real *x;         /* x and possibly q, size natoms*xstride              */
198     int  nout;       /* The number of force arrays                         */
199     nbnxn_atomdata_output_t *out;  /* Output data structures               */
200     int  nalloc;     /* Allocation size of all arrays (for x/f *x/fstride) */
201 } nbnxn_atomdata_t;
202
203 #ifdef __cplusplus
204 }
205 #endif
206
207 #endif