added Verlet scheme and NxN non-bonded functionality
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda_types.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2012, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14  *
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35
36 #ifndef NBNXN_CUDA_TYPES_H
37 #define NBNXN_CUDA_TYPES_H
38
39 #include "types/nbnxn_pairlist.h"
40 #include "types/nbnxn_cuda_types_ext.h"
41 #include "../../gmxlib/cuda_tools/cudautils.cuh"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47 /*! Types of electrostatics available in the CUDA nonbonded force kernels. */
48 enum { eelCuEWALD, eelCuEWALD_TWIN, eelCuRF, eelCuCUT, eelCuNR };
49
50 enum { eNbnxnCuKDefault, eNbnxnCuKLegacy, eNbnxnCuKOld, eNbnxnCuKNR };
51
52 #define NBNXN_KVER_OLD(k)      (k == eNbnxnCuKOld)
53 #define NBNXN_KVER_LEGACY(k)   (k == eNbnxnCuKLegacy)
54 #define NBNXN_KVER_DEFAULT(k)  (k == eNbnxnCuKDefault)
55
56 /*! Non-bonded kernel versions. */
57
58 /* All structs prefixed with "cu_" hold data used in GPU calculations and
59  * are passed to the kernels, except cu_timers_t. */
60 typedef struct cu_plist     cu_plist_t;
61 typedef struct cu_atomdata  cu_atomdata_t;
62 typedef struct cu_nbparam   cu_nbparam_t;
63 typedef struct cu_timers    cu_timers_t;
64 typedef struct nb_staging   nb_staging_t;
65
66
67 /*! Staging area for temporary data. The energies get downloaded here first,
68  *  before getting added to the CPU-side aggregate values.
69  */
70 struct nb_staging
71 {
72     float   *e_lj;      /* LJ energy            */
73     float   *e_el;      /* electrostatic energy */
74     float3  *fshift;    /* shift forces         */
75 };
76
77 /*! Nonbonded atom data -- both inputs and outputs. */
78 struct cu_atomdata
79 {
80     int     natoms;             /* number of atoms                              */
81     int     natoms_local;       /* number of local atoms                        */
82     int     nalloc;             /* allocation size for the atom data (xq, f)    */
83
84     float4  *xq;                /* atom coordinates + charges, size natoms      */
85     float3  *f;                 /* force output array, size natoms              */
86     /* TODO: try float2 for the energies */
87     float   *e_lj,              /* LJ energy output, size 1                     */
88             *e_el;              /* Electrostatics energy input, size 1          */
89
90     float3  *fshift;            /* shift forces                                 */
91
92     int     ntypes;             /* number of atom types                         */
93     int     *atom_types;        /* atom type indices, size natoms               */
94
95     float3  *shift_vec;         /* shifts                                       */
96     bool    bShiftVecUploaded;  /* true if the shift vector has been uploaded   */
97 };
98
99 /*! Parameters required for the CUDA nonbonded calculations. */
100 struct cu_nbparam
101 {
102     int     eeltype;        /* type of electrostatics                       */
103
104     float   epsfac;         /* charge multiplication factor                 */
105     float   c_rf, two_k_rf; /* Reaction-Field constants                     */
106     float   ewald_beta;     /* Ewald/PME parameter                          */
107     float   sh_ewald;       /* Ewald/PME  correction term                   */
108     float   rvdw_sq;        /* VdW cut-off                                  */
109     float   rcoulomb_sq;    /* Coulomb cut-off                              */
110     float   rlist_sq;       /* pair-list cut-off                            */
111     float   sh_invrc6;      /* LJ potential correction term                 */
112
113     float   *nbfp;          /* nonbonded parameter table with C6/C12 pairs  */
114
115     /* Ewald Coulomb force table */
116     int     coulomb_tab_size;
117     float   coulomb_tab_scale;
118     float   *coulomb_tab;
119 };
120
121 /*! Pair list data */
122 struct cu_plist
123 {
124     int             na_c;       /* number of atoms per cluster                  */
125
126     int             nsci;       /* size of sci, # of i clusters in the list     */
127     int             sci_nalloc; /* allocation size of sci                       */
128     nbnxn_sci_t     *sci;       /* list of i-cluster ("super-clusters")         */
129
130     int             ncj4;       /* total # of 4*j clusters                      */
131     int             cj4_nalloc; /* allocation size of cj4                       */
132     nbnxn_cj4_t     *cj4;       /* 4*j cluster list, contains j cluster number
133                                    and index into the i cluster list            */
134     nbnxn_excl_t    *excl;      /* atom interaction bits                        */
135     int             nexcl;      /* count for excl                               */
136     int             excl_nalloc;/* allocation size of excl                      */
137
138     bool            bDoPrune;   /* true if pair-list pruning needs to be
139                                    done during the  current step                */
140 };
141
142 /* CUDA events used for timing GPU kernels and H2D/D2H transfers.
143  * The two-sized arrays hold the local and non-local values and should always
144  * be indexed with eintLocal/eintNonlocal.
145  */
146 struct cu_timers
147 {
148     cudaEvent_t start_atdat, stop_atdat;         /* atom data transfer (every PS step)      */
149     cudaEvent_t start_nb_h2d[2], stop_nb_h2d[2]; /* x/q H2D transfer (every step)           */
150     cudaEvent_t start_nb_d2h[2], stop_nb_d2h[2]; /* f D2H transfer (every step)             */
151     cudaEvent_t start_pl_h2d[2], stop_pl_h2d[2]; /* pair-list H2D transfer (every PS step)  */
152     cudaEvent_t start_nb_k[2], stop_nb_k[2];     /* non-bonded kernels (every step)         */
153 };
154
155 /* Main data structure for CUDA nonbonded force calculations. */
156 struct nbnxn_cuda
157 {
158     cuda_dev_info_t *dev_info;      /* CUDA device information                              */
159     int             kernel_ver;     /* The version of the kernel to be executed on the
160                                        device in use, possible values: eNbnxnCuK*           */
161     bool            bUseTwoStreams; /* true if doing both local/non-local NB work on GPU    */
162     bool            bUseStreamSync; /* true if the standard cudaStreamSynchronize is used
163                                        and not memory polling-based waiting                 */
164     cu_atomdata_t   *atdat;         /* atom data                                            */
165     cu_nbparam_t    *nbparam;       /* parameters required for the non-bonded calc.         */
166     cu_plist_t      *plist[2];      /* pair-list data structures (local and non-local)      */
167     nb_staging_t    nbst;           /* staging area where fshift/energies get downloaded    */
168
169     cudaStream_t    stream[2];      /* local and non-local GPU streams                      */
170
171     /* events used for synchronization */
172     cudaEvent_t    nonlocal_done, misc_ops_done;
173
174     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
175      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
176      * Timer init/uninit is still done even with timing off so only the condition
177      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
178     bool            bDoTime;        /* True if event-based timing is enabled.               */
179     cu_timers_t     *timers;        /* CUDA event-based timers.                             */
180     wallclock_gpu_t *timings;       /* Timing data.                                         */
181 };
182
183 #ifdef __cplusplus
184 }
185 #endif
186
187 #endif  /* NBNXN_CUDA_TYPES_H */