9f951b344a40212fdd145247c5796b5ffee66c45
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_data_mgmt.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \file
36  *  \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
37  *
38  *  \author Szilard Pall <pall.szilard@gmail.com>
39  */
40 #include "gmxpre.h"
41
42 #include <assert.h>
43 #include <stdarg.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46
47 #include "gromacs/gpu_utils/cudautils.cuh"
48 #include "gromacs/gpu_utils/gpu_utils.h"
49 #include "gromacs/gpu_utils/pmalloc_cuda.h"
50 #include "gromacs/hardware/gpu_hw_info.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/mdlib/force_flags.h"
53 #include "gromacs/mdlib/nb_verlet.h"
54 #include "gromacs/mdlib/nbnxn_consts.h"
55 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/real.h"
64 #include "gromacs/utility/smalloc.h"
65
66 #include "nbnxn_cuda_types.h"
67
68 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
69
70 /* This is a heuristically determined parameter for the Fermi, Kepler
71  * and Maxwell architectures for the minimum size of ci lists by multiplying
72  * this constant with the # of multiprocessors on the current device.
73  * Since the maximum number of blocks per multiprocessor is 16, the ideal
74  * count for small systems is 32 or 48 blocks per multiprocessor. Because
75  * there is a bit of fluctuations in the generated block counts, we use
76  * a target of 44 instead of the ideal value of 48.
77  */
78 static unsigned int gpu_min_ci_balanced_factor = 44;
79
80 /* Functions from nbnxn_cuda.cu */
81 extern void nbnxn_cuda_set_cacheconfig(const gmx_device_info_t *devinfo);
82 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref();
83 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref();
84 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref();
85
86
87 /* Fw. decl. */
88 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
89
90 /* Fw. decl, */
91 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
92                                           const gmx_device_info_t *dev_info);
93
94
95 /*! \brief Return whether texture objects are used on this device.
96  *
97  * \param[in]   pointer to the GPU device info structure to inspect for texture objects support
98  * \return      true if texture objects are used on this device
99  */
100 static bool use_texobj(const gmx_device_info_t *dev_info)
101 {
102     assert(!c_disableCudaTextures);
103     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
104     return (dev_info->prop.major >= 3);
105 }
106
107 /*! \brief Return whether combination rules are used.
108  *
109  * \param[in]   pointer to nonbonded paramter struct
110  * \return      true if combination rules are used in this run, false otherwise
111  */
112 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
113 {
114     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
115             nbparam->vdwtype == evdwCuCUTCOMBLB);
116 }
117
118 /*! \brief Set up float texture object.
119  *
120  * Set up texture object for float data and bind it to the device memory
121  * \p devPtr points to.
122  *
123  * \param[out] texObj   texture object to initialize
124  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
125  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
126  */
127 static void setup1DFloatTexture(cudaTextureObject_t &texObj,
128                                 void                *devPtr,
129                                 size_t               sizeInBytes)
130 {
131     assert(!c_disableCudaTextures);
132
133     cudaError_t      stat;
134     cudaResourceDesc rd;
135     cudaTextureDesc  td;
136
137     memset(&rd, 0, sizeof(rd));
138     rd.resType                  = cudaResourceTypeLinear;
139     rd.res.linear.devPtr        = devPtr;
140     rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
141     rd.res.linear.desc.x        = 32;
142     rd.res.linear.sizeInBytes   = sizeInBytes;
143
144     memset(&td, 0, sizeof(td));
145     td.readMode                 = cudaReadModeElementType;
146     stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
147     CU_RET_ERR(stat, "cudaCreateTextureObject failed");
148 }
149
150 /*! \brief Set up float texture reference.
151  *
152  * Set up texture object for float data and bind it to the device memory
153  * \p devPtr points to.
154  *
155  * \param[out] texObj   texture reference to initialize
156  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
157  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
158  */
159 static void setup1DFloatTexture(const struct texture<float, 1, cudaReadModeElementType> *texRef,
160                                 const void                                              *devPtr,
161                                 size_t                                                   sizeInBytes)
162 {
163     assert(!c_disableCudaTextures);
164
165     cudaError_t           stat;
166     cudaChannelFormatDesc cd;
167
168     cd   = cudaCreateChannelDesc<float>();
169     stat = cudaBindTexture(NULL, texRef, devPtr, &cd, sizeInBytes);
170     CU_RET_ERR(stat, "cudaBindTexture failed");
171 }
172
173
174 /*! \brief Initialized the Ewald Coulomb correction GPU table.
175
176     Tabulates the Ewald Coulomb force and initializes the size/scale
177     and the table GPU array. If called with an already allocated table,
178     it just re-uploads the table.
179  */
180 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
181                                            cu_nbparam_t              *nbp,
182                                            const gmx_device_info_t   *dev_info)
183 {
184     float       *coul_tab;
185     cudaError_t  stat;
186
187     if (nbp->coulomb_tab != NULL)
188     {
189         nbnxn_cuda_free_nbparam_table(nbp, dev_info);
190     }
191
192     /* initialize table data in nbp and crete/copy into in global mem */
193     stat = cudaMalloc((void **)&coul_tab, ic->tabq_size*sizeof(*coul_tab));
194     CU_RET_ERR(stat, "cudaMalloc failed on coulumb_tab");
195     cu_copy_H2D(coul_tab, ic->tabq_coul_F, ic->tabq_size*sizeof(*coul_tab));
196
197     nbp->coulomb_tab       = coul_tab;
198     nbp->coulomb_tab_size  = ic->tabq_size;
199     nbp->coulomb_tab_scale = ic->tabq_scale;
200
201     if (!c_disableCudaTextures)
202     {
203         if (use_texobj(dev_info))
204         {
205             setup1DFloatTexture(nbp->coulomb_tab_texobj, nbp->coulomb_tab,
206                                 nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
207         }
208         else
209         {
210             setup1DFloatTexture(&nbnxn_cuda_get_coulomb_tab_texref(), nbp->coulomb_tab,
211                                 nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
212         }
213     }
214 }
215
216
217 /*! Initializes the atomdata structure first time, it only gets filled at
218     pair-search. */
219 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
220 {
221     cudaError_t stat;
222
223     ad->ntypes  = ntypes;
224     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
225     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
226     ad->bShiftVecUploaded = false;
227
228     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
229     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
230
231     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
232     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
233     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
234     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
235
236     /* initialize to NULL poiters to data that is not allocated here and will
237        need reallocation in nbnxn_cuda_init_atomdata */
238     ad->xq = NULL;
239     ad->f  = NULL;
240
241     /* size -1 indicates that the respective array hasn't been initialized yet */
242     ad->natoms = -1;
243     ad->nalloc = -1;
244 }
245
246 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
247     earlier GPUs, single or twin cut-off. */
248 static int pick_ewald_kernel_type(bool                     bTwinCut,
249                                   const gmx_device_info_t *dev_info)
250 {
251     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
252     int  kernel_type;
253
254     /* Benchmarking/development environment variables to force the use of
255        analytical or tabulated Ewald kernel. */
256     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
257     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
258
259     if (bForceAnalyticalEwald && bForceTabulatedEwald)
260     {
261         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
262                    "requested through environment variables.");
263     }
264
265     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
266     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
267     {
268         bUseAnalyticalEwald = true;
269
270         if (debug)
271         {
272             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
273         }
274     }
275     else
276     {
277         bUseAnalyticalEwald = false;
278
279         if (debug)
280         {
281             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
282         }
283     }
284
285     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
286        forces it (use it for debugging/benchmarking only). */
287     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
288     {
289         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
290     }
291     else
292     {
293         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
294     }
295
296     return kernel_type;
297 }
298
299 /*! Copies all parameters related to the cut-off from ic to nbp */
300 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
301                                   const interaction_const_t *ic,
302                                   const NbnxnListParameters *listParams)
303 {
304     nbp->ewald_beta        = ic->ewaldcoeff_q;
305     nbp->sh_ewald          = ic->sh_ewald;
306     nbp->epsfac            = ic->epsfac;
307     nbp->two_k_rf          = 2.0 * ic->k_rf;
308     nbp->c_rf              = ic->c_rf;
309     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
310     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
311     nbp->rlistOuter_sq     = listParams->rlistOuter * listParams->rlistOuter;
312     nbp->rlistInner_sq     = listParams->rlistInner * listParams->rlistInner;
313     nbp->useDynamicPruning = listParams->useDynamicPruning;
314
315     nbp->sh_lj_ewald       = ic->sh_lj_ewald;
316     nbp->ewaldcoeff_lj     = ic->ewaldcoeff_lj;
317
318     nbp->rvdw_switch       = ic->rvdw_switch;
319     nbp->dispersion_shift  = ic->dispersion_shift;
320     nbp->repulsion_shift   = ic->repulsion_shift;
321     nbp->vdw_switch        = ic->vdw_switch;
322 }
323
324 /*! \brief Initialize LJ parameter lookup table.
325  *
326  * Initializes device memory and copies data from host an binds
327  * a texture to allocated device memory to be used for LJ parameter
328  * lookup.
329  *
330  * \param[out] devPtr    device pointer to the memory to be allocated
331  * \param[out] texObj    texture object to be initialized
332  * \param[out] texRef    texture reference to be initialized
333  * \param[in]  hostPtr   pointer to the host memory to be uploaded to the device
334  * \param[in]  numElem   number of elements in the hostPtr
335  * \param[in]  devInfo   pointer to the info struct of the device in use
336  */
337 static void initParamLookupTable(float                    * &devPtr,
338                                  cudaTextureObject_t       &texObj,
339                                  const struct texture<float, 1, cudaReadModeElementType> *texRef,
340                                  const float               *hostPtr,
341                                  int                        numElem,
342                                  const gmx_device_info_t   *devInfo)
343 {
344     cudaError_t stat;
345
346     size_t      sizeInBytes = numElem*sizeof(*devPtr);
347
348     stat  = cudaMalloc((void **)&devPtr, sizeInBytes);
349     CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
350     cu_copy_H2D(devPtr, (void *)hostPtr, sizeInBytes);
351
352     if (!c_disableCudaTextures)
353     {
354         if (use_texobj(devInfo))
355         {
356             setup1DFloatTexture(texObj, devPtr, sizeInBytes);
357         }
358         else
359         {
360             setup1DFloatTexture(texRef, devPtr, sizeInBytes);
361         }
362     }
363 }
364
365 /*! Initializes the nonbonded parameter data structure. */
366 static void init_nbparam(cu_nbparam_t              *nbp,
367                          const interaction_const_t *ic,
368                          const NbnxnListParameters *listParams,
369                          const nbnxn_atomdata_t    *nbat,
370                          const gmx_device_info_t   *dev_info)
371 {
372     int         ntypes;
373
374     ntypes  = nbat->ntype;
375
376     set_cutoff_parameters(nbp, ic, listParams);
377
378     /* The kernel code supports LJ combination rules (geometric and LB) for
379      * all kernel types, but we only generate useful combination rule kernels.
380      * We currently only use LJ combination rule (geometric and LB) kernels
381      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
382      * with PME and 20% with RF, the other kernels speed up about half as much.
383      * For LJ force-switch the geometric rule would give 7% speed-up, but this
384      * combination is rarely used. LJ force-switch with LB rule is more common,
385      * but gives only 1% speed-up.
386      */
387     if (ic->vdwtype == evdwCUT)
388     {
389         switch (ic->vdw_modifier)
390         {
391             case eintmodNONE:
392             case eintmodPOTSHIFT:
393                 switch (nbat->comb_rule)
394                 {
395                     case ljcrNONE:
396                         nbp->vdwtype = evdwCuCUT;
397                         break;
398                     case ljcrGEOM:
399                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
400                         break;
401                     case ljcrLB:
402                         nbp->vdwtype = evdwCuCUTCOMBLB;
403                         break;
404                     default:
405                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
406                         break;
407                 }
408                 break;
409             case eintmodFORCESWITCH:
410                 nbp->vdwtype = evdwCuFSWITCH;
411                 break;
412             case eintmodPOTSWITCH:
413                 nbp->vdwtype = evdwCuPSWITCH;
414                 break;
415             default:
416                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
417                 break;
418         }
419     }
420     else if (ic->vdwtype == evdwPME)
421     {
422         if (ic->ljpme_comb_rule == ljcrGEOM)
423         {
424             assert(nbat->comb_rule == ljcrGEOM);
425             nbp->vdwtype = evdwCuEWALDGEOM;
426         }
427         else
428         {
429             assert(nbat->comb_rule == ljcrLB);
430             nbp->vdwtype = evdwCuEWALDLB;
431         }
432     }
433     else
434     {
435         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
436     }
437
438     if (ic->eeltype == eelCUT)
439     {
440         nbp->eeltype = eelCuCUT;
441     }
442     else if (EEL_RF(ic->eeltype))
443     {
444         nbp->eeltype = eelCuRF;
445     }
446     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
447     {
448         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
449         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
450     }
451     else
452     {
453         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
454         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
455     }
456
457     /* generate table for PME */
458     nbp->coulomb_tab = NULL;
459     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
460     {
461         init_ewald_coulomb_force_table(ic, nbp, dev_info);
462     }
463
464     /* set up LJ parameter lookup table */
465     if (!useLjCombRule(nbp))
466     {
467         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
468                              &nbnxn_cuda_get_nbfp_texref(),
469                              nbat->nbfp, 2*ntypes*ntypes, dev_info);
470     }
471
472     /* set up LJ-PME parameter lookup table */
473     if (ic->vdwtype == evdwPME)
474     {
475         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
476                              &nbnxn_cuda_get_nbfp_comb_texref(),
477                              nbat->nbfp_comb, 2*ntypes, dev_info);
478     }
479 }
480
481 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
482  *  electrostatic type switching to twin cut-off (or back) if needed. */
483 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
484                                         const interaction_const_t   *ic,
485                                         const NbnxnListParameters   *listParams)
486 {
487     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
488     {
489         return;
490     }
491     gmx_nbnxn_cuda_t *nb    = nbv->gpu_nbv;
492     cu_nbparam_t     *nbp   = nb->nbparam;
493
494     set_cutoff_parameters(nbp, ic, listParams);
495
496     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
497                                                  nb->dev_info);
498
499     init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
500 }
501
502 /*! Initializes the pair list data structure. */
503 static void init_plist(cu_plist_t *pl)
504 {
505     /* initialize to NULL pointers to data that is not allocated here and will
506        need reallocation in nbnxn_gpu_init_pairlist */
507     pl->sci      = NULL;
508     pl->cj4      = NULL;
509     pl->imask    = NULL;
510     pl->excl     = NULL;
511
512     /* size -1 indicates that the respective array hasn't been initialized yet */
513     pl->na_c           = -1;
514     pl->nsci           = -1;
515     pl->sci_nalloc     = -1;
516     pl->ncj4           = -1;
517     pl->cj4_nalloc     = -1;
518     pl->nimask         = -1;
519     pl->imask_nalloc   = -1;
520     pl->nexcl          = -1;
521     pl->excl_nalloc    = -1;
522     pl->haveFreshList  = false;
523 }
524
525 /*! Initializes the timer data structure. */
526 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
527 {
528     cudaError_t stat;
529     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
530
531     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
532     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
533     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
534     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
535
536     /* The non-local counters/stream (second in the array) are needed only with DD. */
537     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
538     {
539         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
540         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
541         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
542         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
543
544         stat = cudaEventCreateWithFlags(&(t->start_prune_k[i]), eventflags);
545         CU_RET_ERR(stat, "cudaEventCreate on start_prune_k failed");
546         stat = cudaEventCreateWithFlags(&(t->stop_prune_k[i]), eventflags);
547         CU_RET_ERR(stat, "cudaEventCreate on stop_prune_k failed");
548
549         stat = cudaEventCreateWithFlags(&(t->start_rollingPrune_k[i]), eventflags);
550         CU_RET_ERR(stat, "cudaEventCreate on start_rollingPrune_k failed");
551         stat = cudaEventCreateWithFlags(&(t->stop_rollingPrune_k[i]), eventflags);
552         CU_RET_ERR(stat, "cudaEventCreate on stop_rollingPrune_k failed");
553
554         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
555         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
556         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
557         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
558
559         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
560         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
561         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
562         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
563
564         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
565         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
566         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
567         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
568
569         t->didPairlistH2D[i]  = false;
570         t->didPrune[i]        = false;
571         t->didRollingPrune[i] = false;
572     }
573 }
574
575 /*! Initializes the timings data structure. */
576 static void init_timings(gmx_wallclock_gpu_t *t)
577 {
578     int i, j;
579
580     t->nb_h2d_t = 0.0;
581     t->nb_d2h_t = 0.0;
582     t->nb_c     = 0;
583     t->pl_h2d_t = 0.0;
584     t->pl_h2d_c = 0;
585     for (i = 0; i < 2; i++)
586     {
587         for (j = 0; j < 2; j++)
588         {
589             t->ktime[i][j].t = 0.0;
590             t->ktime[i][j].c = 0;
591         }
592     }
593     t->pruneTime.c        = 0;
594     t->pruneTime.t        = 0.0;
595     t->dynamicPruneTime.c = 0;
596     t->dynamicPruneTime.t = 0.0;
597 }
598
599 /*! Initializes simulation constant data. */
600 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t               *nb,
601                                   const interaction_const_t      *ic,
602                                   const NbnxnListParameters      *listParams,
603                                   const nonbonded_verlet_group_t *nbv_group)
604 {
605     init_atomdata_first(nb->atdat, nbv_group[0].nbat->ntype);
606     init_nbparam(nb->nbparam, ic, listParams, nbv_group[0].nbat, nb->dev_info);
607
608     /* clear energy and shift force outputs */
609     nbnxn_cuda_clear_e_fshift(nb);
610 }
611
612 void nbnxn_gpu_init(gmx_nbnxn_cuda_t         **p_nb,
613                     const gmx_device_info_t   *deviceInfo,
614                     const interaction_const_t *ic,
615                     const NbnxnListParameters *listParams,
616                     nonbonded_verlet_group_t  *nbv_grp,
617                     int                        /*rank*/,
618                     gmx_bool                   bLocalAndNonlocal)
619 {
620     cudaError_t       stat;
621     gmx_nbnxn_cuda_t *nb;
622
623     if (p_nb == NULL)
624     {
625         return;
626     }
627
628     snew(nb, 1);
629     snew(nb->atdat, 1);
630     snew(nb->nbparam, 1);
631     snew(nb->plist[eintLocal], 1);
632     if (bLocalAndNonlocal)
633     {
634         snew(nb->plist[eintNonlocal], 1);
635     }
636
637     nb->bUseTwoStreams = bLocalAndNonlocal;
638
639     snew(nb->timers, 1);
640     snew(nb->timings, 1);
641
642     /* init nbst */
643     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
644     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
645     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
646
647     init_plist(nb->plist[eintLocal]);
648
649     /* set device info, just point it to the right GPU among the detected ones */
650     nb->dev_info = deviceInfo;
651
652     /* local/non-local GPU streams */
653     stat = cudaStreamCreate(&nb->stream[eintLocal]);
654     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
655     if (nb->bUseTwoStreams)
656     {
657         init_plist(nb->plist[eintNonlocal]);
658
659         /* Note that the device we're running on does not have to support
660          * priorities, because we are querying the priority range which in this
661          * case will be a single value.
662          */
663         int highest_priority;
664         stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
665         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
666
667         stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
668                                             cudaStreamDefault,
669                                             highest_priority);
670         CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
671     }
672
673     /* init events for sychronization (timing disabled for performance reasons!) */
674     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
675     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
676     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
677     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
678
679     /* CUDA timing disabled as event timers don't work:
680        - with multiple streams = domain-decomposition;
681        - when turned off by GMX_DISABLE_CUDA_TIMING/GMX_DISABLE_GPU_TIMING.
682      */
683     nb->bDoTime = (!nb->bUseTwoStreams &&
684                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL) &&
685                    (getenv("GMX_DISABLE_GPU_TIMING") == NULL));
686
687     if (nb->bDoTime)
688     {
689         init_timers(nb->timers, nb->bUseTwoStreams);
690         init_timings(nb->timings);
691     }
692
693     /* set the kernel type for the current GPU */
694     /* pick L1 cache configuration */
695     nbnxn_cuda_set_cacheconfig(nb->dev_info);
696
697     nbnxn_cuda_init_const(nb, ic, listParams, nbv_grp);
698
699     *p_nb = nb;
700
701     if (debug)
702     {
703         fprintf(debug, "Initialized CUDA data structures.\n");
704     }
705 }
706
707 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t       *nb,
708                              const nbnxn_pairlist_t *h_plist,
709                              int                     iloc)
710 {
711     char          sbuf[STRLEN];
712     cudaError_t   stat;
713     bool          bDoTime    = nb->bDoTime;
714     cudaStream_t  stream     = nb->stream[iloc];
715     cu_plist_t   *d_plist    = nb->plist[iloc];
716
717     if (d_plist->na_c < 0)
718     {
719         d_plist->na_c = h_plist->na_ci;
720     }
721     else
722     {
723         if (d_plist->na_c != h_plist->na_ci)
724         {
725             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
726                     d_plist->na_c, h_plist->na_ci);
727             gmx_incons(sbuf);
728         }
729     }
730
731     if (bDoTime)
732     {
733         stat = cudaEventRecord(nb->timers->start_pl_h2d[iloc], stream);
734         CU_RET_ERR(stat, "cudaEventRecord failed");
735         nb->timers->didPairlistH2D[iloc] = true;
736     }
737
738     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
739                         &d_plist->nsci, &d_plist->sci_nalloc,
740                         h_plist->nsci,
741                         stream, true);
742
743     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
744                         &d_plist->ncj4, &d_plist->cj4_nalloc,
745                         h_plist->ncj4,
746                         stream, true);
747
748     /* this call only allocates space on the device (no data is transferred) */
749     cu_realloc_buffered((void **)&d_plist->imask, NULL, sizeof(*d_plist->imask),
750                         &d_plist->nimask, &d_plist->imask_nalloc,
751                         h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
752                         stream, true);
753
754     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
755                         &d_plist->nexcl, &d_plist->excl_nalloc,
756                         h_plist->nexcl,
757                         stream, true);
758
759     if (bDoTime)
760     {
761         stat = cudaEventRecord(nb->timers->stop_pl_h2d[iloc], stream);
762         CU_RET_ERR(stat, "cudaEventRecord failed");
763     }
764
765     /* the next use of thist list we be the first one, so we need to prune */
766     d_plist->haveFreshList = true;
767 }
768
769 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
770                                const nbnxn_atomdata_t *nbatom)
771 {
772     cu_atomdata_t *adat  = nb->atdat;
773     cudaStream_t   ls    = nb->stream[eintLocal];
774
775     /* only if we have a dynamic box */
776     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
777     {
778         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
779                           SHIFTS * sizeof(*adat->shift_vec), ls);
780         adat->bShiftVecUploaded = true;
781     }
782 }
783
784 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
785 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
786 {
787     cudaError_t    stat;
788     cu_atomdata_t *adat  = nb->atdat;
789     cudaStream_t   ls    = nb->stream[eintLocal];
790
791     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
792     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
793 }
794
795 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
796 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
797 {
798     cudaError_t    stat;
799     cu_atomdata_t *adat  = nb->atdat;
800     cudaStream_t   ls    = nb->stream[eintLocal];
801
802     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
803     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
804     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
805     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
806     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
807     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
808 }
809
810 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
811 {
812     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
813     /* clear shift force array and energies if the outputs were
814        used in the current step */
815     if (flags & GMX_FORCE_VIRIAL)
816     {
817         nbnxn_cuda_clear_e_fshift(nb);
818     }
819 }
820
821 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t              *nb,
822                              const struct nbnxn_atomdata_t *nbat)
823 {
824     cudaError_t    stat;
825     int            nalloc, natoms;
826     bool           realloced;
827     bool           bDoTime   = nb->bDoTime;
828     cu_timers_t   *timers    = nb->timers;
829     cu_atomdata_t *d_atdat   = nb->atdat;
830     cudaStream_t   ls        = nb->stream[eintLocal];
831
832     natoms    = nbat->natoms;
833     realloced = false;
834
835     if (bDoTime)
836     {
837         /* time async copy */
838         stat = cudaEventRecord(timers->start_atdat, ls);
839         CU_RET_ERR(stat, "cudaEventRecord failed");
840     }
841
842     /* need to reallocate if we have to copy more atoms than the amount of space
843        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
844     if (natoms > d_atdat->nalloc)
845     {
846         nalloc = over_alloc_small(natoms);
847
848         /* free up first if the arrays have already been initialized */
849         if (d_atdat->nalloc != -1)
850         {
851             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
852             cu_free_buffered(d_atdat->xq);
853             cu_free_buffered(d_atdat->atom_types);
854             cu_free_buffered(d_atdat->lj_comb);
855         }
856
857         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
858         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
859         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
860         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
861         if (useLjCombRule(nb->nbparam))
862         {
863             stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
864             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
865         }
866         else
867         {
868             stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
869             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
870         }
871
872         d_atdat->nalloc = nalloc;
873         realloced       = true;
874     }
875
876     d_atdat->natoms       = natoms;
877     d_atdat->natoms_local = nbat->natoms_local;
878
879     /* need to clear GPU f output if realloc happened */
880     if (realloced)
881     {
882         nbnxn_cuda_clear_f(nb, nalloc);
883     }
884
885     if (useLjCombRule(nb->nbparam))
886     {
887         cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
888                           natoms*sizeof(*d_atdat->lj_comb), ls);
889     }
890     else
891     {
892         cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
893                           natoms*sizeof(*d_atdat->atom_types), ls);
894     }
895
896     if (bDoTime)
897     {
898         stat = cudaEventRecord(timers->stop_atdat, ls);
899         CU_RET_ERR(stat, "cudaEventRecord failed");
900     }
901 }
902
903 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
904                                           const gmx_device_info_t *dev_info)
905 {
906     cudaError_t stat;
907
908     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
909     {
910         if (!c_disableCudaTextures)
911         {
912             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
913             if (use_texobj(dev_info))
914             {
915                 stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
916                 CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
917             }
918             else
919             {
920                 GMX_UNUSED_VALUE(dev_info);
921                 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
922                 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
923             }
924         }
925         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
926     }
927 }
928
929 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
930 {
931     cudaError_t      stat;
932     cu_atomdata_t   *atdat;
933     cu_nbparam_t    *nbparam;
934     cu_plist_t      *plist, *plist_nl;
935     cu_timers_t     *timers;
936
937     if (nb == NULL)
938     {
939         return;
940     }
941
942     atdat       = nb->atdat;
943     nbparam     = nb->nbparam;
944     plist       = nb->plist[eintLocal];
945     plist_nl    = nb->plist[eintNonlocal];
946     timers      = nb->timers;
947
948     nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
949
950     stat = cudaEventDestroy(nb->nonlocal_done);
951     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
952     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
953     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
954
955     if (nb->bDoTime)
956     {
957         stat = cudaEventDestroy(timers->start_atdat);
958         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
959         stat = cudaEventDestroy(timers->stop_atdat);
960         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
961
962         /* The non-local counters/stream (second in the array) are needed only with DD. */
963         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
964         {
965             stat = cudaEventDestroy(timers->start_nb_k[i]);
966             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
967             stat = cudaEventDestroy(timers->stop_nb_k[i]);
968             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
969
970             stat = cudaEventDestroy(timers->start_prune_k[i]);
971             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_prune_k");
972             stat = cudaEventDestroy(timers->stop_prune_k[i]);
973             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_prune_k");
974
975             stat = cudaEventDestroy(timers->start_rollingPrune_k[i]);
976             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_rollingPrune_k");
977             stat = cudaEventDestroy(timers->stop_rollingPrune_k[i]);
978             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_rollingPrune_k");
979
980             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
981             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
982             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
983             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
984
985             stat = cudaStreamDestroy(nb->stream[i]);
986             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
987
988             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
989             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
990             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
991             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
992
993             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
994             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
995             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
996             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
997         }
998     }
999
1000     if (!useLjCombRule(nb->nbparam))
1001     {
1002         if (!c_disableCudaTextures)
1003         {
1004             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
1005             if (use_texobj(nb->dev_info))
1006             {
1007                 stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
1008                 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
1009             }
1010             else
1011             {
1012                 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
1013                 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
1014             }
1015         }
1016         cu_free_buffered(nbparam->nbfp);
1017     }
1018
1019     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
1020     {
1021         if (!c_disableCudaTextures)
1022         {
1023             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
1024             if (use_texobj(nb->dev_info))
1025             {
1026                 stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
1027                 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
1028             }
1029             else
1030             {
1031                 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
1032                 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
1033             }
1034         }
1035         cu_free_buffered(nbparam->nbfp_comb);
1036     }
1037
1038     stat = cudaFree(atdat->shift_vec);
1039     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
1040     stat = cudaFree(atdat->fshift);
1041     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
1042
1043     stat = cudaFree(atdat->e_lj);
1044     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
1045     stat = cudaFree(atdat->e_el);
1046     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
1047
1048     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1049     cu_free_buffered(atdat->xq);
1050     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1051     cu_free_buffered(atdat->lj_comb);
1052
1053     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1054     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1055     cu_free_buffered(plist->imask, &plist->nimask, &plist->imask_nalloc);
1056     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1057     if (nb->bUseTwoStreams)
1058     {
1059         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1060         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1061         cu_free_buffered(plist_nl->imask, &plist_nl->nimask, &plist_nl->imask_nalloc);
1062         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1063     }
1064
1065     sfree(atdat);
1066     sfree(nbparam);
1067     sfree(plist);
1068     if (nb->bUseTwoStreams)
1069     {
1070         sfree(plist_nl);
1071     }
1072     sfree(timers);
1073     sfree(nb->timings);
1074     sfree(nb);
1075
1076     if (debug)
1077     {
1078         fprintf(debug, "Cleaned up CUDA data structures.\n");
1079     }
1080 }
1081
1082 gmx_wallclock_gpu_t * nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
1083 {
1084     return (nb != NULL && nb->bDoTime) ? nb->timings : NULL;
1085 }
1086
1087 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
1088 {
1089     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1090     {
1091         init_timings(nbv->gpu_nbv->timings);
1092     }
1093 }
1094
1095 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
1096 {
1097     return nb != NULL ?
1098            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
1099
1100 }
1101
1102 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
1103 {
1104     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1105             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1106 }