0c2808c832699ca1ae545dccd7fe999602507a14
[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, 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 #include "config.h"
36
37 #include <assert.h>
38 #include <stdarg.h>
39 #include <stdlib.h>
40 #include <stdio.h>
41
42 #include <cuda.h>
43
44 #include "gromacs/legacyheaders/tables.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/legacyheaders/types/enums.h"
47 #include "gromacs/mdlib/nb_verlet.h"
48 #include "gromacs/legacyheaders/types/interaction_const.h"
49 #include "gromacs/legacyheaders/types/force_flags.h"
50 #include "../nbnxn_consts.h"
51 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
52
53 #include "nbnxn_cuda_types.h"
54 #include "../../gmxlib/cuda_tools/cudautils.cuh"
55 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
56 #include "gromacs/legacyheaders/pmalloc_cuda.h"
57 #include "gromacs/legacyheaders/gpu_utils.h"
58
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/utility/common.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/smalloc.h"
64
65 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
66
67 /* This is a heuristically determined parameter for the Fermi architecture for
68  * the minimum size of ci lists by multiplying this constant with the # of
69  * multiprocessors on the current device.
70  */
71 static unsigned int gpu_min_ci_balanced_factor = 40;
72
73 /* Functions from nbnxn_cuda.cu */
74 extern void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo);
75 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref();
76 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref();
77 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref();
78
79 /* We should actually be using md_print_warn in md_logging.c,
80  * but we can't include mpi.h in CUDA code.
81  */
82 static void md_print_warn(FILE       *fplog,
83                           const char *fmt, ...)
84 {
85     va_list ap;
86
87     if (fplog != NULL)
88     {
89         /* We should only print to stderr on the master node,
90          * in most cases fplog is only set on the master node, so this works.
91          */
92         va_start(ap, fmt);
93         fprintf(stderr, "\n");
94         vfprintf(stderr, fmt, ap);
95         fprintf(stderr, "\n");
96         va_end(ap);
97
98         va_start(ap, fmt);
99         fprintf(fplog, "\n");
100         vfprintf(fplog, fmt, ap);
101         fprintf(fplog, "\n");
102         va_end(ap);
103     }
104 }
105
106
107 /* Fw. decl. */
108 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
109
110
111 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
112     and the table GPU array. If called with an already allocated table,
113     it just re-uploads the table.
114  */
115 static void init_ewald_coulomb_force_table(cu_nbparam_t          *nbp,
116                                            const cuda_dev_info_t *dev_info)
117 {
118     float       *ftmp, *coul_tab;
119     int          tabsize;
120     double       tabscale;
121     cudaError_t  stat;
122
123     tabsize     = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
124     /* Subtract 2 iso 1 to avoid access out of range due to rounding */
125     tabscale    = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
126
127     pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
128
129     table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
130                                 1/tabscale, nbp->ewald_beta, v_q_ewald_lr);
131
132     /* If the table pointer == NULL the table is generated the first time =>
133        the array pointer will be saved to nbparam and the texture is bound.
134      */
135     coul_tab = nbp->coulomb_tab;
136     if (coul_tab == NULL)
137     {
138         stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
139         CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
140
141         nbp->coulomb_tab = coul_tab;
142
143 #ifdef TEXOBJ_SUPPORTED
144         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
145         if (dev_info->prop.major >= 3)
146         {
147             cudaResourceDesc rd;
148             memset(&rd, 0, sizeof(rd));
149             rd.resType                  = cudaResourceTypeLinear;
150             rd.res.linear.devPtr        = nbp->coulomb_tab;
151             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
152             rd.res.linear.desc.x        = 32;
153             rd.res.linear.sizeInBytes   = tabsize*sizeof(*coul_tab);
154
155             cudaTextureDesc td;
156             memset(&td, 0, sizeof(td));
157             td.readMode                 = cudaReadModeElementType;
158             stat = cudaCreateTextureObject(&nbp->coulomb_tab_texobj, &rd, &td, NULL);
159             CU_RET_ERR(stat, "cudaCreateTextureObject on coulomb_tab_texobj failed");
160         }
161         else
162 #endif
163         {
164             GMX_UNUSED_VALUE(dev_info);
165             cudaChannelFormatDesc cd   = cudaCreateChannelDesc<float>();
166             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
167                                    coul_tab, &cd, tabsize*sizeof(*coul_tab));
168             CU_RET_ERR(stat, "cudaBindTexture on coulomb_tab_texref failed");
169         }
170     }
171
172     cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
173
174     nbp->coulomb_tab_size     = tabsize;
175     nbp->coulomb_tab_scale    = tabscale;
176
177     pfree(ftmp);
178 }
179
180
181 /*! Initializes the atomdata structure first time, it only gets filled at
182     pair-search. */
183 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
184 {
185     cudaError_t stat;
186
187     ad->ntypes  = ntypes;
188     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
189     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
190     ad->bShiftVecUploaded = false;
191
192     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
193     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
194
195     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
196     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
197     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
198     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
199
200     /* initialize to NULL poiters to data that is not allocated here and will
201        need reallocation in nbnxn_cuda_init_atomdata */
202     ad->xq = NULL;
203     ad->f  = NULL;
204
205     /* size -1 indicates that the respective array hasn't been initialized yet */
206     ad->natoms = -1;
207     ad->nalloc = -1;
208 }
209
210 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
211     earlier GPUs, single or twin cut-off. */
212 static int pick_ewald_kernel_type(bool                   bTwinCut,
213                                   const cuda_dev_info_t *dev_info)
214 {
215     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
216     int  kernel_type;
217
218     /* Benchmarking/development environment variables to force the use of
219        analytical or tabulated Ewald kernel. */
220     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
221     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
222
223     if (bForceAnalyticalEwald && bForceTabulatedEwald)
224     {
225         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
226                    "requested through environment variables.");
227     }
228
229     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
230     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
231     {
232         bUseAnalyticalEwald = true;
233
234         if (debug)
235         {
236             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
237         }
238     }
239     else
240     {
241         bUseAnalyticalEwald = false;
242
243         if (debug)
244         {
245             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
246         }
247     }
248
249     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
250        forces it (use it for debugging/benchmarking only). */
251     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
252     {
253         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
254     }
255     else
256     {
257         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
258     }
259
260     return kernel_type;
261 }
262
263 /*! Copies all parameters related to the cut-off from ic to nbp */
264 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
265                                   const interaction_const_t *ic)
266 {
267     nbp->ewald_beta       = ic->ewaldcoeff_q;
268     nbp->sh_ewald         = ic->sh_ewald;
269     nbp->epsfac           = ic->epsfac;
270     nbp->two_k_rf         = 2.0 * ic->k_rf;
271     nbp->c_rf             = ic->c_rf;
272     nbp->rvdw_sq          = ic->rvdw * ic->rvdw;
273     nbp->rcoulomb_sq      = ic->rcoulomb * ic->rcoulomb;
274     nbp->rlist_sq         = ic->rlist * ic->rlist;
275
276     nbp->sh_lj_ewald      = ic->sh_lj_ewald;
277     nbp->ewaldcoeff_lj    = ic->ewaldcoeff_lj;
278
279     nbp->rvdw_switch      = ic->rvdw_switch;
280     nbp->dispersion_shift = ic->dispersion_shift;
281     nbp->repulsion_shift  = ic->repulsion_shift;
282     nbp->vdw_switch       = ic->vdw_switch;
283 }
284
285 /*! Initializes the nonbonded parameter data structure. */
286 static void init_nbparam(cu_nbparam_t              *nbp,
287                          const interaction_const_t *ic,
288                          const nbnxn_atomdata_t    *nbat,
289                          const cuda_dev_info_t     *dev_info)
290 {
291     cudaError_t stat;
292     int         ntypes, nnbfp, nnbfp_comb;
293
294     ntypes  = nbat->ntype;
295
296     set_cutoff_parameters(nbp, ic);
297
298     if (ic->vdwtype == evdwCUT)
299     {
300         switch (ic->vdw_modifier)
301         {
302             case eintmodNONE:
303             case eintmodPOTSHIFT:
304                 nbp->vdwtype = evdwCuCUT;
305                 break;
306             case eintmodFORCESWITCH:
307                 nbp->vdwtype = evdwCuFSWITCH;
308                 break;
309             case eintmodPOTSWITCH:
310                 nbp->vdwtype = evdwCuPSWITCH;
311                 break;
312             default:
313                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
314                 break;
315         }
316     }
317     else if (ic->vdwtype == evdwPME)
318     {
319         if (ic->ljpme_comb_rule == ljcrGEOM)
320         {
321             assert(nbat->comb_rule == ljcrGEOM);
322             nbp->vdwtype = evdwCuEWALDGEOM;
323         }
324         else
325         {
326             assert(nbat->comb_rule == ljcrLB);
327             nbp->vdwtype = evdwCuEWALDLB;
328         }
329     }
330     else
331     {
332         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
333     }
334
335     if (ic->eeltype == eelCUT)
336     {
337         nbp->eeltype = eelCuCUT;
338     }
339     else if (EEL_RF(ic->eeltype))
340     {
341         nbp->eeltype = eelCuRF;
342     }
343     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
344     {
345         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
346         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
347     }
348     else
349     {
350         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
351         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
352     }
353
354     /* generate table for PME */
355     nbp->coulomb_tab = NULL;
356     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
357     {
358         init_ewald_coulomb_force_table(nbp, dev_info);
359     }
360
361     nnbfp      = 2*ntypes*ntypes;
362     nnbfp_comb = 2*ntypes;
363
364     stat  = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
365     CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
366     cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
367
368
369     if (ic->vdwtype == evdwPME)
370     {
371         stat  = cudaMalloc((void **)&nbp->nbfp_comb, nnbfp_comb*sizeof(*nbp->nbfp_comb));
372         CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp_comb");
373         cu_copy_H2D(nbp->nbfp_comb, nbat->nbfp_comb, nnbfp_comb*sizeof(*nbp->nbfp_comb));
374     }
375
376 #ifdef TEXOBJ_SUPPORTED
377     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
378     if (dev_info->prop.major >= 3)
379     {
380         cudaResourceDesc rd;
381         cudaTextureDesc  td;
382
383         memset(&rd, 0, sizeof(rd));
384         rd.resType                  = cudaResourceTypeLinear;
385         rd.res.linear.devPtr        = nbp->nbfp;
386         rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
387         rd.res.linear.desc.x        = 32;
388         rd.res.linear.sizeInBytes   = nnbfp*sizeof(*nbp->nbfp);
389
390         memset(&td, 0, sizeof(td));
391         td.readMode                 = cudaReadModeElementType;
392         stat = cudaCreateTextureObject(&nbp->nbfp_texobj, &rd, &td, NULL);
393         CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_texobj failed");
394
395         if (ic->vdwtype == evdwPME)
396         {
397             memset(&rd, 0, sizeof(rd));
398             rd.resType                  = cudaResourceTypeLinear;
399             rd.res.linear.devPtr        = nbp->nbfp_comb;
400             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
401             rd.res.linear.desc.x        = 32;
402             rd.res.linear.sizeInBytes   = nnbfp_comb*sizeof(*nbp->nbfp_comb);
403
404             memset(&td, 0, sizeof(td));
405             td.readMode = cudaReadModeElementType;
406             stat        = cudaCreateTextureObject(&nbp->nbfp_comb_texobj, &rd, &td, NULL);
407             CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_comb_texobj failed");
408         }
409     }
410     else
411 #endif
412     {
413         cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
414         stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
415                                nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
416         CU_RET_ERR(stat, "cudaBindTexture on nbfp_texref failed");
417
418         if (ic->vdwtype == evdwPME)
419         {
420             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_comb_texref(),
421                                    nbp->nbfp_comb, &cd, nnbfp_comb*sizeof(*nbp->nbfp_comb));
422             CU_RET_ERR(stat, "cudaBindTexture on nbfp_comb_texref failed");
423         }
424     }
425 }
426
427 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
428  *  electrostatic type switching to twin cut-off (or back) if needed. */
429 void nbnxn_cuda_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
430                                          const interaction_const_t   *ic)
431 {
432     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_CUDA)
433     {
434         return;
435     }
436     nbnxn_cuda_ptr_t cu_nb = nbv->cu_nbv;
437     cu_nbparam_t    *nbp   = cu_nb->nbparam;
438
439     set_cutoff_parameters(nbp, ic);
440
441     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
442                                                  cu_nb->dev_info);
443
444     init_ewald_coulomb_force_table(cu_nb->nbparam, cu_nb->dev_info);
445 }
446
447 /*! Initializes the pair list data structure. */
448 static void init_plist(cu_plist_t *pl)
449 {
450     /* initialize to NULL pointers to data that is not allocated here and will
451        need reallocation in nbnxn_cuda_init_pairlist */
452     pl->sci     = NULL;
453     pl->cj4     = NULL;
454     pl->excl    = NULL;
455
456     /* size -1 indicates that the respective array hasn't been initialized yet */
457     pl->na_c        = -1;
458     pl->nsci        = -1;
459     pl->sci_nalloc  = -1;
460     pl->ncj4        = -1;
461     pl->cj4_nalloc  = -1;
462     pl->nexcl       = -1;
463     pl->excl_nalloc = -1;
464     pl->bDoPrune    = false;
465 }
466
467 /*! Initializes the timer data structure. */
468 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
469 {
470     cudaError_t stat;
471     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
472
473     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
474     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
475     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
476     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
477
478     /* The non-local counters/stream (second in the array) are needed only with DD. */
479     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
480     {
481         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
482         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
483         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
484         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
485
486
487         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
488         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
489         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
490         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
491
492         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
493         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
494         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
495         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
496
497         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
498         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
499         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
500         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
501     }
502 }
503
504 /*! Initializes the timings data structure. */
505 static void init_timings(wallclock_gpu_t *t)
506 {
507     int i, j;
508
509     t->nb_h2d_t = 0.0;
510     t->nb_d2h_t = 0.0;
511     t->nb_c     = 0;
512     t->pl_h2d_t = 0.0;
513     t->pl_h2d_c = 0;
514     for (i = 0; i < 2; i++)
515     {
516         for (j = 0; j < 2; j++)
517         {
518             t->ktime[i][j].t = 0.0;
519             t->ktime[i][j].c = 0;
520         }
521     }
522 }
523
524 void nbnxn_cuda_init(FILE                 *fplog,
525                      nbnxn_cuda_ptr_t     *p_cu_nb,
526                      const gmx_gpu_info_t *gpu_info,
527                      const gmx_gpu_opt_t  *gpu_opt,
528                      int                   my_gpu_index,
529                      gmx_bool              bLocalAndNonlocal)
530 {
531     cudaError_t       stat;
532     nbnxn_cuda_ptr_t  nb;
533     char              sbuf[STRLEN];
534     bool              bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
535     int               cuda_drv_ver;
536
537     assert(gpu_info);
538
539     if (p_cu_nb == NULL)
540     {
541         return;
542     }
543
544     snew(nb, 1);
545     snew(nb->atdat, 1);
546     snew(nb->nbparam, 1);
547     snew(nb->plist[eintLocal], 1);
548     if (bLocalAndNonlocal)
549     {
550         snew(nb->plist[eintNonlocal], 1);
551     }
552
553     nb->bUseTwoStreams = bLocalAndNonlocal;
554
555     snew(nb->timers, 1);
556     snew(nb->timings, 1);
557
558     /* init nbst */
559     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
560     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
561     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
562
563     init_plist(nb->plist[eintLocal]);
564
565     /* set device info, just point it to the right GPU among the detected ones */
566     nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
567
568     /* local/non-local GPU streams */
569     stat = cudaStreamCreate(&nb->stream[eintLocal]);
570     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
571     if (nb->bUseTwoStreams)
572     {
573         init_plist(nb->plist[eintNonlocal]);
574
575         /* CUDA stream priority available in the CUDA RT 5.5 API.
576          * Note that the device we're running on does not have to support
577          * priorities, because we are querying the priority range which in this
578          * case will be a single value.
579          */
580 #if CUDA_VERSION >= 5500
581         {
582             int highest_priority;
583             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
584             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
585
586             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
587                                                 cudaStreamDefault,
588                                                 highest_priority);
589             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
590         }
591 #else
592         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
593         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
594 #endif
595     }
596
597     /* init events for sychronization (timing disabled for performance reasons!) */
598     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
599     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
600     stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
601     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
602
603     /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
604      * (which increases with shorter time/step) caused by a known CUDA driver bug.
605      * To work around the issue we'll use an (admittedly fragile) memory polling
606      * waiting to preserve performance. This requires support for atomic
607      * operations and only works on x86/x86_64.
608      * With polling wait event-timing also needs to be disabled.
609      *
610      * The overhead is greatly reduced in API v5.0 drivers and the improvement
611      * is independent of runtime version. Hence, with API v5.0 drivers and later
612      * we won't switch to polling.
613      *
614      * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
615      * ranks so we will also disable it in that case.
616      */
617
618     bStreamSync    = getenv("GMX_CUDA_STREAMSYNC") != NULL;
619     bNoStreamSync  = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
620
621 #ifdef TMPI_ATOMICS
622     bTMPIAtomics = true;
623 #else
624     bTMPIAtomics = false;
625 #endif
626
627 #ifdef GMX_TARGET_X86
628     bX86 = true;
629 #else
630     bX86 = false;
631 #endif
632
633     if (bStreamSync && bNoStreamSync)
634     {
635         gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
636     }
637
638     stat = cudaDriverGetVersion(&cuda_drv_ver);
639     CU_RET_ERR(stat, "cudaDriverGetVersion failed");
640
641     bOldDriver = (cuda_drv_ver < 5000);
642
643     if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
644     {
645         /* Polling wait should be used instead of cudaStreamSynchronize only if:
646          *   - ECC is ON & driver is old (checked above),
647          *   - we're on x86/x86_64,
648          *   - atomics are available, and
649          *   - GPUs are not being shared.
650          */
651         bool bShouldUsePollSync = (bX86 && bTMPIAtomics &&
652                                    (gmx_count_gpu_dev_shared(gpu_opt) < 1));
653
654         if (bStreamSync)
655         {
656             nb->bUseStreamSync = true;
657
658             /* only warn if polling should be used */
659             if (bShouldUsePollSync)
660             {
661                 md_print_warn(fplog,
662                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
663                               "      cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
664             }
665         }
666         else
667         {
668             nb->bUseStreamSync = !bShouldUsePollSync;
669
670             if (bShouldUsePollSync)
671             {
672                 md_print_warn(fplog,
673                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
674                               "      cause performance loss. Switching to the alternative polling GPU wait.\n"
675                               "      If you encounter issues, switch back to standard GPU waiting by setting\n"
676                               "      the GMX_CUDA_STREAMSYNC environment variable.\n");
677             }
678             else
679             {
680                 /* Tell the user that the ECC+old driver combination can be bad */
681                 sprintf(sbuf,
682                         "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
683                         "      A known bug in this driver version can cause performance loss.\n"
684                         "      However, the polling wait workaround can not be used because\n%s\n"
685                         "      Consider updating the driver or turning ECC off.",
686                         (bX86 && bTMPIAtomics) ?
687                         "      GPU(s) are being oversubscribed." :
688                         "      atomic operations are not supported by the platform/CPU+compiler.");
689                 md_print_warn(fplog, sbuf);
690             }
691         }
692     }
693     else
694     {
695         if (bNoStreamSync)
696         {
697             nb->bUseStreamSync = false;
698
699             md_print_warn(fplog,
700                           "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
701         }
702         else
703         {
704             /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
705             nb->bUseStreamSync = true;
706         }
707     }
708
709     /* CUDA timing disabled as event timers don't work:
710        - with multiple streams = domain-decomposition;
711        - with the polling waiting hack (without cudaStreamSynchronize);
712        - when turned off by GMX_DISABLE_CUDA_TIMING.
713      */
714     nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
715                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
716
717     if (nb->bDoTime)
718     {
719         init_timers(nb->timers, nb->bUseTwoStreams);
720         init_timings(nb->timings);
721     }
722
723     /* set the kernel type for the current GPU */
724     /* pick L1 cache configuration */
725     nbnxn_cuda_set_cacheconfig(nb->dev_info);
726
727     *p_cu_nb = nb;
728
729     if (debug)
730     {
731         fprintf(debug, "Initialized CUDA data structures.\n");
732     }
733 }
734
735 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
736                            const interaction_const_t      *ic,
737                            const nonbonded_verlet_group_t *nbv_group)
738 {
739     init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
740     init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
741
742     /* clear energy and shift force outputs */
743     nbnxn_cuda_clear_e_fshift(cu_nb);
744 }
745
746 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t        cu_nb,
747                               const nbnxn_pairlist_t *h_plist,
748                               int                     iloc)
749 {
750     char          sbuf[STRLEN];
751     cudaError_t   stat;
752     bool          bDoTime    = cu_nb->bDoTime;
753     cudaStream_t  stream     = cu_nb->stream[iloc];
754     cu_plist_t   *d_plist    = cu_nb->plist[iloc];
755
756     if (d_plist->na_c < 0)
757     {
758         d_plist->na_c = h_plist->na_ci;
759     }
760     else
761     {
762         if (d_plist->na_c != h_plist->na_ci)
763         {
764             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
765                     d_plist->na_c, h_plist->na_ci);
766             gmx_incons(sbuf);
767         }
768     }
769
770     if (bDoTime)
771     {
772         stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
773         CU_RET_ERR(stat, "cudaEventRecord failed");
774     }
775
776     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
777                         &d_plist->nsci, &d_plist->sci_nalloc,
778                         h_plist->nsci,
779                         stream, true);
780
781     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
782                         &d_plist->ncj4, &d_plist->cj4_nalloc,
783                         h_plist->ncj4,
784                         stream, true);
785
786     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
787                         &d_plist->nexcl, &d_plist->excl_nalloc,
788                         h_plist->nexcl,
789                         stream, true);
790
791     if (bDoTime)
792     {
793         stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
794         CU_RET_ERR(stat, "cudaEventRecord failed");
795     }
796
797     /* need to prune the pair list during the next step */
798     d_plist->bDoPrune = true;
799 }
800
801 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t        cu_nb,
802                                 const nbnxn_atomdata_t *nbatom)
803 {
804     cu_atomdata_t *adat  = cu_nb->atdat;
805     cudaStream_t   ls    = cu_nb->stream[eintLocal];
806
807     /* only if we have a dynamic box */
808     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
809     {
810         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
811                           SHIFTS * sizeof(*adat->shift_vec), ls);
812         adat->bShiftVecUploaded = true;
813     }
814 }
815
816 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
817 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
818 {
819     cudaError_t    stat;
820     cu_atomdata_t *adat  = cu_nb->atdat;
821     cudaStream_t   ls    = cu_nb->stream[eintLocal];
822
823     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
824     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
825 }
826
827 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
828 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
829 {
830     cudaError_t    stat;
831     cu_atomdata_t *adat  = cu_nb->atdat;
832     cudaStream_t   ls    = cu_nb->stream[eintLocal];
833
834     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
835     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
836     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
837     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
838     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
839     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
840 }
841
842 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
843 {
844     nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
845     /* clear shift force array and energies if the outputs were
846        used in the current step */
847     if (flags & GMX_FORCE_VIRIAL)
848     {
849         nbnxn_cuda_clear_e_fshift(cu_nb);
850     }
851 }
852
853 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t        cu_nb,
854                               const nbnxn_atomdata_t *nbat)
855 {
856     cudaError_t    stat;
857     int            nalloc, natoms;
858     bool           realloced;
859     bool           bDoTime   = cu_nb->bDoTime;
860     cu_timers_t   *timers    = cu_nb->timers;
861     cu_atomdata_t *d_atdat   = cu_nb->atdat;
862     cudaStream_t   ls        = cu_nb->stream[eintLocal];
863
864     natoms    = nbat->natoms;
865     realloced = false;
866
867     if (bDoTime)
868     {
869         /* time async copy */
870         stat = cudaEventRecord(timers->start_atdat, ls);
871         CU_RET_ERR(stat, "cudaEventRecord failed");
872     }
873
874     /* need to reallocate if we have to copy more atoms than the amount of space
875        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
876     if (natoms > d_atdat->nalloc)
877     {
878         nalloc = over_alloc_small(natoms);
879
880         /* free up first if the arrays have already been initialized */
881         if (d_atdat->nalloc != -1)
882         {
883             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
884             cu_free_buffered(d_atdat->xq);
885             cu_free_buffered(d_atdat->atom_types);
886         }
887
888         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
889         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
890         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
891         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
892
893         stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
894         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
895
896         d_atdat->nalloc = nalloc;
897         realloced       = true;
898     }
899
900     d_atdat->natoms       = natoms;
901     d_atdat->natoms_local = nbat->natoms_local;
902
903     /* need to clear GPU f output if realloc happened */
904     if (realloced)
905     {
906         nbnxn_cuda_clear_f(cu_nb, nalloc);
907     }
908
909     cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
910                       natoms*sizeof(*d_atdat->atom_types), ls);
911
912     if (bDoTime)
913     {
914         stat = cudaEventRecord(timers->stop_atdat, ls);
915         CU_RET_ERR(stat, "cudaEventRecord failed");
916     }
917 }
918
919 void nbnxn_cuda_free(nbnxn_cuda_ptr_t cu_nb)
920 {
921     cudaError_t      stat;
922     cu_atomdata_t   *atdat;
923     cu_nbparam_t    *nbparam;
924     cu_plist_t      *plist, *plist_nl;
925     cu_timers_t     *timers;
926
927     if (cu_nb == NULL)
928     {
929         return;
930     }
931
932     atdat       = cu_nb->atdat;
933     nbparam     = cu_nb->nbparam;
934     plist       = cu_nb->plist[eintLocal];
935     plist_nl    = cu_nb->plist[eintNonlocal];
936     timers      = cu_nb->timers;
937
938     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
939     {
940
941 #ifdef TEXOBJ_SUPPORTED
942         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
943         if (cu_nb->dev_info->prop.major >= 3)
944         {
945             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
946             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
947         }
948         else
949 #endif
950         {
951             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
952             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
953         }
954         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
955     }
956
957     stat = cudaEventDestroy(cu_nb->nonlocal_done);
958     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
959     stat = cudaEventDestroy(cu_nb->misc_ops_done);
960     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
961
962     if (cu_nb->bDoTime)
963     {
964         stat = cudaEventDestroy(timers->start_atdat);
965         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
966         stat = cudaEventDestroy(timers->stop_atdat);
967         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
968
969         /* The non-local counters/stream (second in the array) are needed only with DD. */
970         for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
971         {
972             stat = cudaEventDestroy(timers->start_nb_k[i]);
973             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
974             stat = cudaEventDestroy(timers->stop_nb_k[i]);
975             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
976
977             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
978             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
979             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
980             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
981
982             stat = cudaStreamDestroy(cu_nb->stream[i]);
983             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
984
985             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
986             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
987             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
988             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
989
990             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
991             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
992             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
993             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
994         }
995     }
996
997 #ifdef TEXOBJ_SUPPORTED
998     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
999     if (cu_nb->dev_info->prop.major >= 3)
1000     {
1001         stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
1002         CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
1003     }
1004     else
1005 #endif
1006     {
1007         stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
1008         CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
1009     }
1010     cu_free_buffered(nbparam->nbfp);
1011
1012     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
1013     {
1014 #ifdef TEXOBJ_SUPPORTED
1015         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
1016         if (cu_nb->dev_info->prop.major >= 3)
1017         {
1018             stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
1019             CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
1020         }
1021         else
1022 #endif
1023         {
1024             stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
1025             CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
1026         }
1027         cu_free_buffered(nbparam->nbfp_comb);
1028     }
1029
1030     stat = cudaFree(atdat->shift_vec);
1031     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
1032     stat = cudaFree(atdat->fshift);
1033     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
1034
1035     stat = cudaFree(atdat->e_lj);
1036     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
1037     stat = cudaFree(atdat->e_el);
1038     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
1039
1040     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1041     cu_free_buffered(atdat->xq);
1042     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1043
1044     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1045     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1046     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1047     if (cu_nb->bUseTwoStreams)
1048     {
1049         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1050         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1051         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1052     }
1053
1054     sfree(atdat);
1055     sfree(nbparam);
1056     sfree(plist);
1057     if (cu_nb->bUseTwoStreams)
1058     {
1059         sfree(plist_nl);
1060     }
1061     sfree(timers);
1062     sfree(cu_nb->timings);
1063     sfree(cu_nb);
1064
1065     if (debug)
1066     {
1067         fprintf(debug, "Cleaned up CUDA data structures.\n");
1068     }
1069 }
1070
1071 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
1072 {
1073     cudaError_t  stat;
1074     cudaStream_t stream = cu_nb->stream[iloc];
1075
1076     stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
1077     CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
1078 }
1079
1080 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
1081 {
1082     return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
1083 }
1084
1085 void nbnxn_cuda_reset_timings(nonbonded_verlet_t* nbv)
1086 {
1087     if (nbv->cu_nbv && nbv->cu_nbv->bDoTime)
1088     {
1089         init_timings(nbv->cu_nbv->timings);
1090     }
1091 }
1092
1093 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
1094 {
1095     return cu_nb != NULL ?
1096            gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
1097
1098 }
1099
1100 gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
1101 {
1102     return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1103             (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1104 }