55f67e928838fd526ceacfceadba2900680c6d95
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda_types.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,2013, 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_CUDA_TYPES_H
40 #define NBNXN_CUDA_TYPES_H
41
42 #include "types/nbnxn_pairlist.h"
43 #include "types/nbnxn_cuda_types_ext.h"
44 #include "../../gmxlib/cuda_tools/cudautils.cuh"
45
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
49
50 /*! Types of electrostatics available in the CUDA nonbonded force kernels. */
51 enum {
52     eelCuEWALD, eelCuEWALD_TWIN, eelCuRF, eelCuCUT, eelCuNR
53 };
54
55 enum {
56     eNbnxnCuKDefault, eNbnxnCuKLegacy, eNbnxnCuKOld, eNbnxnCuKNR
57 };
58
59 #define NBNXN_KVER_OLD(k)      (k == eNbnxnCuKOld)
60 #define NBNXN_KVER_LEGACY(k)   (k == eNbnxnCuKLegacy)
61 #define NBNXN_KVER_DEFAULT(k)  (k == eNbnxnCuKDefault)
62
63 /*! Non-bonded kernel versions. */
64
65 /* All structs prefixed with "cu_" hold data used in GPU calculations and
66  * are passed to the kernels, except cu_timers_t. */
67 typedef struct cu_plist     cu_plist_t;
68 typedef struct cu_atomdata  cu_atomdata_t;
69 typedef struct cu_nbparam   cu_nbparam_t;
70 typedef struct cu_timers    cu_timers_t;
71 typedef struct nb_staging   nb_staging_t;
72
73
74 /*! Staging area for temporary data. The energies get downloaded here first,
75  *  before getting added to the CPU-side aggregate values.
76  */
77 struct nb_staging
78 {
79     float   *e_lj;      /* LJ energy            */
80     float   *e_el;      /* electrostatic energy */
81     float3  *fshift;    /* shift forces         */
82 };
83
84 /*! Nonbonded atom data -- both inputs and outputs. */
85 struct cu_atomdata
86 {
87     int      natoms;            /* number of atoms                              */
88     int      natoms_local;      /* number of local atoms                        */
89     int      nalloc;            /* allocation size for the atom data (xq, f)    */
90
91     float4  *xq;                /* atom coordinates + charges, size natoms      */
92     float3  *f;                 /* force output array, size natoms              */
93     /* TODO: try float2 for the energies */
94     float   *e_lj,              /* LJ energy output, size 1                     */
95             *e_el;              /* Electrostatics energy input, size 1          */
96
97     float3  *fshift;            /* shift forces                                 */
98
99     int      ntypes;            /* number of atom types                         */
100     int     *atom_types;        /* atom type indices, size natoms               */
101
102     float3  *shift_vec;         /* shifts                                       */
103     bool     bShiftVecUploaded; /* true if the shift vector has been uploaded   */
104 };
105
106 /*! Parameters required for the CUDA nonbonded calculations. */
107 struct cu_nbparam
108 {
109     int      eeltype;        /* type of electrostatics                       */
110
111     float    epsfac;         /* charge multiplication factor                 */
112     float    c_rf, two_k_rf; /* Reaction-Field constants                     */
113     float    ewald_beta;     /* Ewald/PME parameter                          */
114     float    sh_ewald;       /* Ewald/PME  correction term                   */
115     float    rvdw_sq;        /* VdW cut-off                                  */
116     float    rcoulomb_sq;    /* Coulomb cut-off                              */
117     float    rlist_sq;       /* pair-list cut-off                            */
118     float    sh_invrc6;      /* LJ potential correction term                 */
119
120     float   *nbfp;           /* nonbonded parameter table with C6/C12 pairs  */
121
122     /* Ewald Coulomb force table */
123     int      coulomb_tab_size;
124     float    coulomb_tab_scale;
125     float   *coulomb_tab;
126 };
127
128 /*! Pair list data */
129 struct cu_plist
130 {
131     int              na_c;        /* number of atoms per cluster                  */
132
133     int              nsci;        /* size of sci, # of i clusters in the list     */
134     int              sci_nalloc;  /* allocation size of sci                       */
135     nbnxn_sci_t     *sci;         /* list of i-cluster ("super-clusters")         */
136
137     int              ncj4;        /* total # of 4*j clusters                      */
138     int              cj4_nalloc;  /* allocation size of cj4                       */
139     nbnxn_cj4_t     *cj4;         /* 4*j cluster list, contains j cluster number
140                                      and index into the i cluster list            */
141     nbnxn_excl_t    *excl;        /* atom interaction bits                        */
142     int              nexcl;       /* count for excl                               */
143     int              excl_nalloc; /* allocation size of excl                      */
144
145     bool             bDoPrune;    /* true if pair-list pruning needs to be
146                                      done during the  current step                */
147 };
148
149 /* CUDA events used for timing GPU kernels and H2D/D2H transfers.
150  * The two-sized arrays hold the local and non-local values and should always
151  * be indexed with eintLocal/eintNonlocal.
152  */
153 struct cu_timers
154 {
155     cudaEvent_t start_atdat, stop_atdat;         /* atom data transfer (every PS step)      */
156     cudaEvent_t start_nb_h2d[2], stop_nb_h2d[2]; /* x/q H2D transfer (every step)           */
157     cudaEvent_t start_nb_d2h[2], stop_nb_d2h[2]; /* f D2H transfer (every step)             */
158     cudaEvent_t start_pl_h2d[2], stop_pl_h2d[2]; /* pair-list H2D transfer (every PS step)  */
159     cudaEvent_t start_nb_k[2], stop_nb_k[2];     /* non-bonded kernels (every step)         */
160 };
161
162 /* Main data structure for CUDA nonbonded force calculations. */
163 struct nbnxn_cuda
164 {
165     cuda_dev_info_t *dev_info;       /* CUDA device information                              */
166     int              kernel_ver;     /* The version of the kernel to be executed on the
167                                         device in use, possible values: eNbnxnCuK*           */
168     bool             bUseTwoStreams; /* true if doing both local/non-local NB work on GPU    */
169     bool             bUseStreamSync; /* true if the standard cudaStreamSynchronize is used
170                                         and not memory polling-based waiting                 */
171     cu_atomdata_t   *atdat;          /* atom data                                            */
172     cu_nbparam_t    *nbparam;        /* parameters required for the non-bonded calc.         */
173     cu_plist_t      *plist[2];       /* pair-list data structures (local and non-local)      */
174     nb_staging_t     nbst;           /* staging area where fshift/energies get downloaded    */
175
176     cudaStream_t     stream[2];      /* local and non-local GPU streams                      */
177
178     /* events used for synchronization */
179     cudaEvent_t    nonlocal_done, misc_ops_done;
180
181     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
182      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
183      * Timer init/uninit is still done even with timing off so only the condition
184      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
185     bool             bDoTime;       /* True if event-based timing is enabled.               */
186     cu_timers_t     *timers;        /* CUDA event-based timers.                             */
187     wallclock_gpu_t *timings;       /* Timing data.                                         */
188 };
189
190 #ifdef __cplusplus
191 }
192 #endif
193
194 #endif  /* NBNXN_CUDA_TYPES_H */