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