Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_atomdata.c
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
33 #ifdef HAVE_CONFIG_H
34 #include <config.h>
35 #endif
36
37 #include <math.h>
38 #include <string.h>
39 #include "smalloc.h"
40 #include "macros.h"
41 #include "vec.h"
42 #include "nbnxn_consts.h"
43 #include "nbnxn_internal.h"
44 #include "nbnxn_search.h"
45 #include "nbnxn_atomdata.h"
46 #include "gmx_omp_nthreads.h"
47
48 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
49 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
50 {
51     *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
52 }
53
54 /* Free function for memory allocated with nbnxn_alloc_aligned */
55 void nbnxn_free_aligned(void *ptr)
56 {
57     sfree_aligned(ptr);
58 }
59
60 /* Reallocation wrapper function for nbnxn data structures */
61 void nbnxn_realloc_void(void **ptr,
62                         int nbytes_copy, int nbytes_new,
63                         nbnxn_alloc_t *ma,
64                         nbnxn_free_t  *mf)
65 {
66     void *ptr_new;
67
68     ma(&ptr_new, nbytes_new);
69
70     if (nbytes_new > 0 && ptr_new == NULL)
71     {
72         gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
73     }
74
75     if (nbytes_copy > 0)
76     {
77         if (nbytes_new < nbytes_copy)
78         {
79             gmx_incons("In nbnxn_realloc_void: new size less than copy size");
80         }
81         memcpy(ptr_new, *ptr, nbytes_copy);
82     }
83     if (*ptr != NULL)
84     {
85         mf(*ptr);
86     }
87     *ptr = ptr_new;
88 }
89
90 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
91 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
92 {
93     int t;
94
95     nbnxn_realloc_void((void **)&nbat->type,
96                        nbat->natoms*sizeof(*nbat->type),
97                        n*sizeof(*nbat->type),
98                        nbat->alloc, nbat->free);
99     nbnxn_realloc_void((void **)&nbat->lj_comb,
100                        nbat->natoms*2*sizeof(*nbat->lj_comb),
101                        n*2*sizeof(*nbat->lj_comb),
102                        nbat->alloc, nbat->free);
103     if (nbat->XFormat != nbatXYZQ)
104     {
105         nbnxn_realloc_void((void **)&nbat->q,
106                            nbat->natoms*sizeof(*nbat->q),
107                            n*sizeof(*nbat->q),
108                            nbat->alloc, nbat->free);
109     }
110     if (nbat->nenergrp > 1)
111     {
112         nbnxn_realloc_void((void **)&nbat->energrp,
113                            nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
114                            n/nbat->na_c*sizeof(*nbat->energrp),
115                            nbat->alloc, nbat->free);
116     }
117     nbnxn_realloc_void((void **)&nbat->x,
118                        nbat->natoms*nbat->xstride*sizeof(*nbat->x),
119                        n*nbat->xstride*sizeof(*nbat->x),
120                        nbat->alloc, nbat->free);
121     for (t = 0; t < nbat->nout; t++)
122     {
123         /* Allocate one element extra for possible signaling with CUDA */
124         nbnxn_realloc_void((void **)&nbat->out[t].f,
125                            nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
126                            n*nbat->fstride*sizeof(*nbat->out[t].f),
127                            nbat->alloc, nbat->free);
128     }
129     nbat->nalloc = n;
130 }
131
132 /* Initializes an nbnxn_atomdata_output_t data structure */
133 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
134                                        int nb_kernel_type,
135                                        int nenergrp, int stride,
136                                        nbnxn_alloc_t *ma)
137 {
138     int cj_size;
139
140     out->f = NULL;
141     ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
142     out->nV = nenergrp*nenergrp;
143     ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
144     ma((void **)&out->Vc, out->nV*sizeof(*out->Vc  ));
145
146     if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
147         nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
148     {
149         cj_size  = nbnxn_kernel_to_cj_size(nb_kernel_type);
150         out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
151         ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
152         ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc  ));
153     }
154     else
155     {
156         out->nVS = 0;
157     }
158 }
159
160 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
161                                  const int *in, int fill, int *innb)
162 {
163     int i, j;
164
165     j = 0;
166     for (i = 0; i < na; i++)
167     {
168         innb[j++] = in[a[i]];
169     }
170     /* Complete the partially filled last cell with fill */
171     for (; i < na_round; i++)
172     {
173         innb[j++] = fill;
174     }
175 }
176
177 static void clear_nbat_real(int na, int nbatFormat, real *xnb, int a0)
178 {
179     int a, d, j, c;
180
181     switch (nbatFormat)
182     {
183         case nbatXYZ:
184             for (a = 0; a < na; a++)
185             {
186                 for (d = 0; d < DIM; d++)
187                 {
188                     xnb[(a0+a)*STRIDE_XYZ+d] = 0;
189                 }
190             }
191             break;
192         case nbatXYZQ:
193             for (a = 0; a < na; a++)
194             {
195                 for (d = 0; d < DIM; d++)
196                 {
197                     xnb[(a0+a)*STRIDE_XYZQ+d] = 0;
198                 }
199             }
200             break;
201         case nbatX4:
202             j = X4_IND_A(a0);
203             c = a0 & (PACK_X4-1);
204             for (a = 0; a < na; a++)
205             {
206                 xnb[j+XX*PACK_X4] = 0;
207                 xnb[j+YY*PACK_X4] = 0;
208                 xnb[j+ZZ*PACK_X4] = 0;
209                 j++;
210                 c++;
211                 if (c == PACK_X4)
212                 {
213                     j += (DIM-1)*PACK_X4;
214                     c  = 0;
215                 }
216             }
217             break;
218         case nbatX8:
219             j = X8_IND_A(a0);
220             c = a0 & (PACK_X8-1);
221             for (a = 0; a < na; a++)
222             {
223                 xnb[j+XX*PACK_X8] = 0;
224                 xnb[j+YY*PACK_X8] = 0;
225                 xnb[j+ZZ*PACK_X8] = 0;
226                 j++;
227                 c++;
228                 if (c == PACK_X8)
229                 {
230                     j += (DIM-1)*PACK_X8;
231                     c  = 0;
232                 }
233             }
234             break;
235     }
236 }
237
238 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
239                             rvec *x, int nbatFormat, real *xnb, int a0,
240                             int cx, int cy, int cz)
241 {
242     int i, j, c;
243
244 /* We might need to place filler particles to fill up the cell to na_round.
245  * The coefficients (LJ and q) for such particles are zero.
246  * But we might still get NaN as 0*NaN when distances are too small.
247  * We hope that -107 nm is far away enough from to zero
248  * to avoid accidental short distances to particles shifted down for pbc.
249  */
250 #define NBAT_FAR_AWAY 107
251
252     switch (nbatFormat)
253     {
254         case nbatXYZ:
255             j = a0*STRIDE_XYZ;
256             for (i = 0; i < na; i++)
257             {
258                 xnb[j++] = x[a[i]][XX];
259                 xnb[j++] = x[a[i]][YY];
260                 xnb[j++] = x[a[i]][ZZ];
261             }
262             /* Complete the partially filled last cell with copies of the last element.
263              * This simplifies the bounding box calculation and avoid
264              * numerical issues with atoms that are coincidentally close.
265              */
266             for (; i < na_round; i++)
267             {
268                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
269                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
270                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
271             }
272             break;
273         case nbatXYZQ:
274             j = a0*STRIDE_XYZQ;
275             for (i = 0; i < na; i++)
276             {
277                 xnb[j++] = x[a[i]][XX];
278                 xnb[j++] = x[a[i]][YY];
279                 xnb[j++] = x[a[i]][ZZ];
280                 j++;
281             }
282             /* Complete the partially filled last cell with particles far apart */
283             for (; i < na_round; i++)
284             {
285                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cx);
286                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cy);
287                 xnb[j++] = -NBAT_FAR_AWAY*(1 + cz + i);
288                 j++;
289             }
290             break;
291         case nbatX4:
292             j = X4_IND_A(a0);
293             c = a0 & (PACK_X4-1);
294             for (i = 0; i < na; i++)
295             {
296                 xnb[j+XX*PACK_X4] = x[a[i]][XX];
297                 xnb[j+YY*PACK_X4] = x[a[i]][YY];
298                 xnb[j+ZZ*PACK_X4] = x[a[i]][ZZ];
299                 j++;
300                 c++;
301                 if (c == PACK_X4)
302                 {
303                     j += (DIM-1)*PACK_X4;
304                     c  = 0;
305                 }
306             }
307             /* Complete the partially filled last cell with particles far apart */
308             for (; i < na_round; i++)
309             {
310                 xnb[j+XX*PACK_X4] = -NBAT_FAR_AWAY*(1 + cx);
311                 xnb[j+YY*PACK_X4] = -NBAT_FAR_AWAY*(1 + cy);
312                 xnb[j+ZZ*PACK_X4] = -NBAT_FAR_AWAY*(1 + cz + i);
313                 j++;
314                 c++;
315                 if (c == PACK_X4)
316                 {
317                     j += (DIM-1)*PACK_X4;
318                     c  = 0;
319                 }
320             }
321             break;
322         case nbatX8:
323             j = X8_IND_A(a0);
324             c = a0 & (PACK_X8 - 1);
325             for (i = 0; i < na; i++)
326             {
327                 xnb[j+XX*PACK_X8] = x[a[i]][XX];
328                 xnb[j+YY*PACK_X8] = x[a[i]][YY];
329                 xnb[j+ZZ*PACK_X8] = x[a[i]][ZZ];
330                 j++;
331                 c++;
332                 if (c == PACK_X8)
333                 {
334                     j += (DIM-1)*PACK_X8;
335                     c  = 0;
336                 }
337             }
338             /* Complete the partially filled last cell with particles far apart */
339             for (; i < na_round; i++)
340             {
341                 xnb[j+XX*PACK_X8] = -NBAT_FAR_AWAY*(1 + cx);
342                 xnb[j+YY*PACK_X8] = -NBAT_FAR_AWAY*(1 + cy);
343                 xnb[j+ZZ*PACK_X8] = -NBAT_FAR_AWAY*(1 + cz + i);
344                 j++;
345                 c++;
346                 if (c == PACK_X8)
347                 {
348                     j += (DIM-1)*PACK_X8;
349                     c  = 0;
350                 }
351             }
352             break;
353         default:
354             gmx_incons("Unsupported nbnxn_atomdata_t format");
355     }
356 }
357
358 /* Determines the combination rule (or none) to be used, stores it,
359  * and sets the LJ parameters required with the rule.
360  */
361 static void set_combination_rule_data(nbnxn_atomdata_t *nbat)
362 {
363     int  nt, i, j;
364     real c6, c12;
365
366     nt = nbat->ntype;
367
368     switch (nbat->comb_rule)
369     {
370         case  ljcrGEOM:
371             nbat->comb_rule = ljcrGEOM;
372
373             for (i = 0; i < nt; i++)
374             {
375                 /* Copy the diagonal from the nbfp matrix */
376                 nbat->nbfp_comb[i*2  ] = sqrt(nbat->nbfp[(i*nt+i)*2  ]);
377                 nbat->nbfp_comb[i*2+1] = sqrt(nbat->nbfp[(i*nt+i)*2+1]);
378             }
379             break;
380         case ljcrLB:
381             for (i = 0; i < nt; i++)
382             {
383                 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
384                 c6  = nbat->nbfp[(i*nt+i)*2  ];
385                 c12 = nbat->nbfp[(i*nt+i)*2+1];
386                 if (c6 > 0 && c12 > 0)
387                 {
388                     /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
389                      * so we get 6*C6 and 12*C12 after combining.
390                      */
391                     nbat->nbfp_comb[i*2  ] = 0.5*pow(c12/c6, 1.0/6.0);
392                     nbat->nbfp_comb[i*2+1] = sqrt(c6*c6/c12);
393                 }
394                 else
395                 {
396                     nbat->nbfp_comb[i*2  ] = 0;
397                     nbat->nbfp_comb[i*2+1] = 0;
398                 }
399             }
400             break;
401         case ljcrNONE:
402             /* nbfp_s4 stores two parameters using a stride of 4,
403              * because this would suit x86 SIMD single-precision
404              * quad-load intrinsics. There's a slight inefficiency in
405              * allocating and initializing nbfp_s4 when it might not
406              * be used, but introducing the conditional code is not
407              * really worth it. */
408             nbat->alloc((void **)&nbat->nbfp_s4, nt*nt*4*sizeof(*nbat->nbfp_s4));
409             for (i = 0; i < nt; i++)
410             {
411                 for (j = 0; j < nt; j++)
412                 {
413                     nbat->nbfp_s4[(i*nt+j)*4+0] = nbat->nbfp[(i*nt+j)*2+0];
414                     nbat->nbfp_s4[(i*nt+j)*4+1] = nbat->nbfp[(i*nt+j)*2+1];
415                     nbat->nbfp_s4[(i*nt+j)*4+2] = 0;
416                     nbat->nbfp_s4[(i*nt+j)*4+3] = 0;
417                 }
418             }
419             break;
420         default:
421             gmx_incons("Unknown combination rule");
422             break;
423     }
424 }
425
426 /* Initializes an nbnxn_atomdata_t data structure */
427 void nbnxn_atomdata_init(FILE *fp,
428                          nbnxn_atomdata_t *nbat,
429                          int nb_kernel_type,
430                          int ntype, const real *nbfp,
431                          int n_energygroups,
432                          int nout,
433                          nbnxn_alloc_t *alloc,
434                          nbnxn_free_t  *free)
435 {
436     int      i, j;
437     real     c6, c12, tol;
438     char    *ptr;
439     gmx_bool simple, bCombGeom, bCombLB;
440
441     if (alloc == NULL)
442     {
443         nbat->alloc = nbnxn_alloc_aligned;
444     }
445     else
446     {
447         nbat->alloc = alloc;
448     }
449     if (free == NULL)
450     {
451         nbat->free = nbnxn_free_aligned;
452     }
453     else
454     {
455         nbat->free = free;
456     }
457
458     if (debug)
459     {
460         fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
461     }
462     nbat->ntype = ntype + 1;
463     nbat->alloc((void **)&nbat->nbfp,
464                 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
465     nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
466
467     /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
468      * force-field floating point parameters.
469      */
470     tol = 1e-5;
471     ptr = getenv("GMX_LJCOMB_TOL");
472     if (ptr != NULL)
473     {
474         double dbl;
475
476         sscanf(ptr, "%lf", &dbl);
477         tol = dbl;
478     }
479     bCombGeom = TRUE;
480     bCombLB   = TRUE;
481
482     /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
483      * to check for the LB rule.
484      */
485     for (i = 0; i < ntype; i++)
486     {
487         c6  = nbfp[(i*ntype+i)*2  ]/6.0;
488         c12 = nbfp[(i*ntype+i)*2+1]/12.0;
489         if (c6 > 0 && c12 > 0)
490         {
491             nbat->nbfp_comb[i*2  ] = pow(c12/c6, 1.0/6.0);
492             nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
493         }
494         else if (c6 == 0 && c12 == 0)
495         {
496             nbat->nbfp_comb[i*2  ] = 0;
497             nbat->nbfp_comb[i*2+1] = 0;
498         }
499         else
500         {
501             /* Can not use LB rule with only dispersion or repulsion */
502             bCombLB = FALSE;
503         }
504     }
505
506     for (i = 0; i < nbat->ntype; i++)
507     {
508         for (j = 0; j < nbat->ntype; j++)
509         {
510             if (i < ntype && j < ntype)
511             {
512                 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
513                  * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
514                  */
515                 c6  = nbfp[(i*ntype+j)*2  ];
516                 c12 = nbfp[(i*ntype+j)*2+1];
517                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = c6;
518                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
519
520                 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
521                 bCombGeom = bCombGeom &&
522                     gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2  ]*nbfp[(j*ntype+j)*2  ], tol) &&
523                     gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
524
525                 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
526                 c6     /= 6.0;
527                 c12    /= 12.0;
528                 bCombLB = bCombLB &&
529                     ((c6 == 0 && c12 == 0 &&
530                       (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
531                      (c6 > 0 && c12 > 0 &&
532                       gmx_within_tol(pow(c12/c6, 1.0/6.0), 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
533                       gmx_within_tol(0.25*c6*c6/c12, sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
534             }
535             else
536             {
537                 /* Add zero parameters for the additional dummy atom type */
538                 nbat->nbfp[(i*nbat->ntype+j)*2  ] = 0;
539                 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
540             }
541         }
542     }
543     if (debug)
544     {
545         fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
546                 bCombGeom, bCombLB);
547     }
548
549     simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
550
551     if (simple)
552     {
553         /* We prefer the geometic combination rule,
554          * as that gives a slightly faster kernel than the LB rule.
555          */
556         if (bCombGeom)
557         {
558             nbat->comb_rule = ljcrGEOM;
559         }
560         else if (bCombLB)
561         {
562             nbat->comb_rule = ljcrLB;
563         }
564         else
565         {
566             nbat->comb_rule = ljcrNONE;
567
568             nbat->free(nbat->nbfp_comb);
569         }
570
571         if (fp)
572         {
573             if (nbat->comb_rule == ljcrNONE)
574             {
575                 fprintf(fp, "Using full Lennard-Jones parameter combination matrix\n\n");
576             }
577             else
578             {
579                 fprintf(fp, "Using %s Lennard-Jones combination rule\n\n",
580                         nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
581             }
582         }
583
584         set_combination_rule_data(nbat);
585     }
586     else
587     {
588         nbat->comb_rule = ljcrNONE;
589
590         nbat->free(nbat->nbfp_comb);
591     }
592
593     nbat->natoms  = 0;
594     nbat->type    = NULL;
595     nbat->lj_comb = NULL;
596     if (simple)
597     {
598         int pack_x;
599
600         switch (nb_kernel_type)
601         {
602             case nbnxnk4xN_SIMD_4xN:
603             case nbnxnk4xN_SIMD_2xNN:
604                 pack_x = max(NBNXN_CPU_CLUSTER_I_SIZE,
605                              nbnxn_kernel_to_cj_size(nb_kernel_type));
606                 switch (pack_x)
607                 {
608                     case 4:
609                         nbat->XFormat = nbatX4;
610                         break;
611                     case 8:
612                         nbat->XFormat = nbatX8;
613                         break;
614                     default:
615                         gmx_incons("Unsupported packing width");
616                 }
617                 break;
618             default:
619                 nbat->XFormat = nbatXYZ;
620                 break;
621         }
622
623         nbat->FFormat = nbat->XFormat;
624     }
625     else
626     {
627         nbat->XFormat = nbatXYZQ;
628         nbat->FFormat = nbatXYZ;
629     }
630     nbat->q        = NULL;
631     nbat->nenergrp = n_energygroups;
632     if (!simple)
633     {
634         /* Energy groups not supported yet for super-sub lists */
635         if (n_energygroups > 1 && fp != NULL)
636         {
637             fprintf(fp, "\nNOTE: With GPUs, reporting energy group contributions is not supported\n\n");
638         }
639         nbat->nenergrp = 1;
640     }
641     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
642     if (nbat->nenergrp > 64)
643     {
644         gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
645     }
646     nbat->neg_2log = 1;
647     while (nbat->nenergrp > (1<<nbat->neg_2log))
648     {
649         nbat->neg_2log++;
650     }
651     nbat->energrp = NULL;
652     nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
653     nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
654     nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
655     nbat->x       = NULL;
656
657 #ifdef GMX_NBNXN_SIMD
658     if (simple)
659     {
660         /* Set the diagonal cluster pair interaction mask setup data.
661          * In the kernel we check 0 < j - i to generate the masks.
662          * Here we store j - i for generating the mask for the first i (i=0);
663          * we substract 0.5 to avoid rounding issues.
664          * In the kernel we can subtract 1 to generate the mask for the next i.
665          */
666         const int simd_width = GMX_SIMD_WIDTH_HERE;
667         int       simd_4xn_diag_ind_size, simd_interaction_size, j;
668
669         simd_4xn_diag_ind_size = max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
670         snew_aligned(nbat->simd_4xn_diagonal_j_minus_i,
671                      simd_4xn_diag_ind_size, NBNXN_MEM_ALIGN);
672         for (j = 0; j < simd_4xn_diag_ind_size; j++)
673         {
674             nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
675         }
676
677         snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i,
678                      simd_width, NBNXN_MEM_ALIGN);
679         for (j = 0; j < simd_width/2; j++)
680         {
681             /* The j-cluster size is half the SIMD width */
682             nbat->simd_2xnn_diagonal_j_minus_i[j]              = j - 0.5;
683             /* The next half of the SIMD width is for i + 1 */
684             nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
685         }
686
687         /* We use up to 32 bits for exclusion masking.
688          * The same masks are used for the 4xN and 2x(N+N) kernels.
689          * The masks are read either into epi32 SIMD registers or into
690          * real SIMD registers (together with a cast).
691          * In single precision this means the real and epi32 SIMD registers
692          * are of equal size.
693          * In double precision the epi32 registers can be smaller than
694          * the real registers, so depending on the architecture, we might
695          * need to use two, identical, 32-bit masks per real.
696          */
697         simd_interaction_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
698         snew_aligned(nbat->simd_exclusion_filter1, simd_interaction_size,   NBNXN_MEM_ALIGN);
699         snew_aligned(nbat->simd_exclusion_filter2, simd_interaction_size*2, NBNXN_MEM_ALIGN);
700         
701         for (j = 0; j < simd_interaction_size; j++)
702         {
703             /* Set the consecutive bits for filters pair exclusions masks */
704             nbat->simd_exclusion_filter1[j]       = (1U << j);
705             nbat->simd_exclusion_filter2[j*2 + 0] = (1U << j);
706             nbat->simd_exclusion_filter2[j*2 + 1] = (1U << j);
707         }
708     }
709 #endif
710
711     /* Initialize the output data structures */
712     nbat->nout    = nout;
713     snew(nbat->out, nbat->nout);
714     nbat->nalloc  = 0;
715     for (i = 0; i < nbat->nout; i++)
716     {
717         nbnxn_atomdata_output_init(&nbat->out[i],
718                                    nb_kernel_type,
719                                    nbat->nenergrp, 1<<nbat->neg_2log,
720                                    nbat->alloc);
721     }
722     nbat->buffer_flags.flag        = NULL;
723     nbat->buffer_flags.flag_nalloc = 0;
724 }
725
726 static void copy_lj_to_nbat_lj_comb_x4(const real *ljparam_type,
727                                        const int *type, int na,
728                                        real *ljparam_at)
729 {
730     int is, k, i;
731
732     /* The LJ params follow the combination rule:
733      * copy the params for the type array to the atom array.
734      */
735     for (is = 0; is < na; is += PACK_X4)
736     {
737         for (k = 0; k < PACK_X4; k++)
738         {
739             i = is + k;
740             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
741             ljparam_at[is*2+PACK_X4+k] = ljparam_type[type[i]*2+1];
742         }
743     }
744 }
745
746 static void copy_lj_to_nbat_lj_comb_x8(const real *ljparam_type,
747                                        const int *type, int na,
748                                        real *ljparam_at)
749 {
750     int is, k, i;
751
752     /* The LJ params follow the combination rule:
753      * copy the params for the type array to the atom array.
754      */
755     for (is = 0; is < na; is += PACK_X8)
756     {
757         for (k = 0; k < PACK_X8; k++)
758         {
759             i = is + k;
760             ljparam_at[is*2        +k] = ljparam_type[type[i]*2  ];
761             ljparam_at[is*2+PACK_X8+k] = ljparam_type[type[i]*2+1];
762         }
763     }
764 }
765
766 /* Sets the atom type and LJ data in nbnxn_atomdata_t */
767 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t    *nbat,
768                                          int                  ngrid,
769                                          const nbnxn_search_t nbs,
770                                          const int           *type)
771 {
772     int                 g, i, ncz, ash;
773     const nbnxn_grid_t *grid;
774
775     for (g = 0; g < ngrid; g++)
776     {
777         grid = &nbs->grid[g];
778
779         /* Loop over all columns and copy and fill */
780         for (i = 0; i < grid->ncx*grid->ncy; i++)
781         {
782             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
783             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
784
785             copy_int_to_nbat_int(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
786                                  type, nbat->ntype-1, nbat->type+ash);
787
788             if (nbat->comb_rule != ljcrNONE)
789             {
790                 if (nbat->XFormat == nbatX4)
791                 {
792                     copy_lj_to_nbat_lj_comb_x4(nbat->nbfp_comb,
793                                                nbat->type+ash, ncz*grid->na_sc,
794                                                nbat->lj_comb+ash*2);
795                 }
796                 else if (nbat->XFormat == nbatX8)
797                 {
798                     copy_lj_to_nbat_lj_comb_x8(nbat->nbfp_comb,
799                                                nbat->type+ash, ncz*grid->na_sc,
800                                                nbat->lj_comb+ash*2);
801                 }
802             }
803         }
804     }
805 }
806
807 /* Sets the charges in nbnxn_atomdata_t *nbat */
808 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t    *nbat,
809                                        int                  ngrid,
810                                        const nbnxn_search_t nbs,
811                                        const real          *charge)
812 {
813     int                 g, cxy, ncz, ash, na, na_round, i, j;
814     real               *q;
815     const nbnxn_grid_t *grid;
816
817     for (g = 0; g < ngrid; g++)
818     {
819         grid = &nbs->grid[g];
820
821         /* Loop over all columns and copy and fill */
822         for (cxy = 0; cxy < grid->ncx*grid->ncy; cxy++)
823         {
824             ash      = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
825             na       = grid->cxy_na[cxy];
826             na_round = (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
827
828             if (nbat->XFormat == nbatXYZQ)
829             {
830                 q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
831                 for (i = 0; i < na; i++)
832                 {
833                     *q = charge[nbs->a[ash+i]];
834                     q += STRIDE_XYZQ;
835                 }
836                 /* Complete the partially filled last cell with zeros */
837                 for (; i < na_round; i++)
838                 {
839                     *q = 0;
840                     q += STRIDE_XYZQ;
841                 }
842             }
843             else
844             {
845                 q = nbat->q + ash;
846                 for (i = 0; i < na; i++)
847                 {
848                     *q = charge[nbs->a[ash+i]];
849                     q++;
850                 }
851                 /* Complete the partially filled last cell with zeros */
852                 for (; i < na_round; i++)
853                 {
854                     *q = 0;
855                     q++;
856                 }
857             }
858         }
859     }
860 }
861
862 /* Copies the energy group indices to a reordered and packed array */
863 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
864                                   int na_c, int bit_shift,
865                                   const int *in, int *innb)
866 {
867     int i, j, sa, at;
868     int comb;
869
870     j = 0;
871     for (i = 0; i < na; i += na_c)
872     {
873         /* Store na_c energy group numbers into one int */
874         comb = 0;
875         for (sa = 0; sa < na_c; sa++)
876         {
877             at = a[i+sa];
878             if (at >= 0)
879             {
880                 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
881             }
882         }
883         innb[j++] = comb;
884     }
885     /* Complete the partially filled last cell with fill */
886     for (; i < na_round; i += na_c)
887     {
888         innb[j++] = 0;
889     }
890 }
891
892 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
893 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t    *nbat,
894                                             int                  ngrid,
895                                             const nbnxn_search_t nbs,
896                                             const int           *atinfo)
897 {
898     int                 g, i, ncz, ash;
899     const nbnxn_grid_t *grid;
900
901     for (g = 0; g < ngrid; g++)
902     {
903         grid = &nbs->grid[g];
904
905         /* Loop over all columns and copy and fill */
906         for (i = 0; i < grid->ncx*grid->ncy; i++)
907         {
908             ncz = grid->cxy_ind[i+1] - grid->cxy_ind[i];
909             ash = (grid->cell0 + grid->cxy_ind[i])*grid->na_sc;
910
911             copy_egp_to_nbat_egps(nbs->a+ash, grid->cxy_na[i], ncz*grid->na_sc,
912                                   nbat->na_c, nbat->neg_2log,
913                                   atinfo, nbat->energrp+(ash>>grid->na_c_2log));
914         }
915     }
916 }
917
918 /* Sets all required atom parameter data in nbnxn_atomdata_t */
919 void nbnxn_atomdata_set(nbnxn_atomdata_t    *nbat,
920                         int                  locality,
921                         const nbnxn_search_t nbs,
922                         const t_mdatoms     *mdatoms,
923                         const int           *atinfo)
924 {
925     int ngrid;
926
927     if (locality == eatLocal)
928     {
929         ngrid = 1;
930     }
931     else
932     {
933         ngrid = nbs->ngrid;
934     }
935
936     nbnxn_atomdata_set_atomtypes(nbat, ngrid, nbs, mdatoms->typeA);
937
938     nbnxn_atomdata_set_charges(nbat, ngrid, nbs, mdatoms->chargeA);
939
940     if (nbat->nenergrp > 1)
941     {
942         nbnxn_atomdata_set_energygroups(nbat, ngrid, nbs, atinfo);
943     }
944 }
945
946 /* Copies the shift vector array to nbnxn_atomdata_t */
947 void nbnxn_atomdata_copy_shiftvec(gmx_bool          bDynamicBox,
948                                   rvec             *shift_vec,
949                                   nbnxn_atomdata_t *nbat)
950 {
951     int i;
952
953     nbat->bDynamicBox = bDynamicBox;
954     for (i = 0; i < SHIFTS; i++)
955     {
956         copy_rvec(shift_vec[i], nbat->shift_vec[i]);
957     }
958 }
959
960 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
961 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
962                                      int                  locality,
963                                      gmx_bool             FillLocal,
964                                      rvec                *x,
965                                      nbnxn_atomdata_t    *nbat)
966 {
967     int g0 = 0, g1 = 0;
968     int nth, th;
969
970     switch (locality)
971     {
972         case eatAll:
973             g0 = 0;
974             g1 = nbs->ngrid;
975             break;
976         case eatLocal:
977             g0 = 0;
978             g1 = 1;
979             break;
980         case eatNonlocal:
981             g0 = 1;
982             g1 = nbs->ngrid;
983             break;
984     }
985
986     if (FillLocal)
987     {
988         nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
989     }
990
991     nth = gmx_omp_nthreads_get(emntPairsearch);
992
993 #pragma omp parallel for num_threads(nth) schedule(static)
994     for (th = 0; th < nth; th++)
995     {
996         int g;
997
998         for (g = g0; g < g1; g++)
999         {
1000             const nbnxn_grid_t *grid;
1001             int                 cxy0, cxy1, cxy;
1002
1003             grid = &nbs->grid[g];
1004
1005             cxy0 = (grid->ncx*grid->ncy* th   +nth-1)/nth;
1006             cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1007
1008             for (cxy = cxy0; cxy < cxy1; cxy++)
1009             {
1010                 int na, ash, na_fill;
1011
1012                 na  = grid->cxy_na[cxy];
1013                 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1014
1015                 if (g == 0 && FillLocal)
1016                 {
1017                     na_fill =
1018                         (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1019                 }
1020                 else
1021                 {
1022                     /* We fill only the real particle locations.
1023                      * We assume the filling entries at the end have been
1024                      * properly set before during ns.
1025                      */
1026                     na_fill = na;
1027                 }
1028                 copy_rvec_to_nbat_real(nbs->a+ash, na, na_fill, x,
1029                                        nbat->XFormat, nbat->x, ash,
1030                                        0, 0, 0);
1031             }
1032         }
1033     }
1034 }
1035
1036 static void
1037 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1038                            int i0, int i1)
1039 {
1040     int i;
1041
1042     for (i = i0; i < i1; i++)
1043     {
1044         dest[i] = 0;
1045     }
1046 }
1047
1048 static void
1049 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1050                             gmx_bool bDestSet,
1051                             real ** gmx_restrict src,
1052                             int nsrc,
1053                             int i0, int i1)
1054 {
1055     int i, s;
1056
1057     if (bDestSet)
1058     {
1059         /* The destination buffer contains data, add to it */
1060         for (i = i0; i < i1; i++)
1061         {
1062             for (s = 0; s < nsrc; s++)
1063             {
1064                 dest[i] += src[s][i];
1065             }
1066         }
1067     }
1068     else
1069     {
1070         /* The destination buffer is unitialized, set it first */
1071         for (i = i0; i < i1; i++)
1072         {
1073             dest[i] = src[0][i];
1074             for (s = 1; s < nsrc; s++)
1075             {
1076                 dest[i] += src[s][i];
1077             }
1078         }
1079     }
1080 }
1081
1082 static void
1083 nbnxn_atomdata_reduce_reals_simd(real * gmx_restrict dest,
1084                                  gmx_bool bDestSet,
1085                                  real ** gmx_restrict src,
1086                                  int nsrc,
1087                                  int i0, int i1)
1088 {
1089 #ifdef GMX_NBNXN_SIMD
1090 /* The SIMD width here is actually independent of that in the kernels,
1091  * but we use the same width for simplicity (usually optimal anyhow).
1092  */
1093     int       i, s;
1094     gmx_mm_pr dest_SSE, src_SSE;
1095
1096     if (bDestSet)
1097     {
1098         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1099         {
1100             dest_SSE = gmx_load_pr(dest+i);
1101             for (s = 0; s < nsrc; s++)
1102             {
1103                 src_SSE  = gmx_load_pr(src[s]+i);
1104                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1105             }
1106             gmx_store_pr(dest+i, dest_SSE);
1107         }
1108     }
1109     else
1110     {
1111         for (i = i0; i < i1; i += GMX_SIMD_WIDTH_HERE)
1112         {
1113             dest_SSE = gmx_load_pr(src[0]+i);
1114             for (s = 1; s < nsrc; s++)
1115             {
1116                 src_SSE  = gmx_load_pr(src[s]+i);
1117                 dest_SSE = gmx_add_pr(dest_SSE, src_SSE);
1118             }
1119             gmx_store_pr(dest+i, dest_SSE);
1120         }
1121     }
1122 #endif
1123 }
1124
1125 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1126 static void
1127 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1128                                     const nbnxn_atomdata_t *nbat,
1129                                     nbnxn_atomdata_output_t *out,
1130                                     int nfa,
1131                                     int a0, int a1,
1132                                     rvec *f)
1133 {
1134     int         a, i, fa;
1135     const int  *cell;
1136     const real *fnb;
1137
1138     cell = nbs->cell;
1139
1140     /* Loop over all columns and copy and fill */
1141     switch (nbat->FFormat)
1142     {
1143         case nbatXYZ:
1144         case nbatXYZQ:
1145             if (nfa == 1)
1146             {
1147                 fnb = out[0].f;
1148
1149                 for (a = a0; a < a1; a++)
1150                 {
1151                     i = cell[a]*nbat->fstride;
1152
1153                     f[a][XX] += fnb[i];
1154                     f[a][YY] += fnb[i+1];
1155                     f[a][ZZ] += fnb[i+2];
1156                 }
1157             }
1158             else
1159             {
1160                 for (a = a0; a < a1; a++)
1161                 {
1162                     i = cell[a]*nbat->fstride;
1163
1164                     for (fa = 0; fa < nfa; fa++)
1165                     {
1166                         f[a][XX] += out[fa].f[i];
1167                         f[a][YY] += out[fa].f[i+1];
1168                         f[a][ZZ] += out[fa].f[i+2];
1169                     }
1170                 }
1171             }
1172             break;
1173         case nbatX4:
1174             if (nfa == 1)
1175             {
1176                 fnb = out[0].f;
1177
1178                 for (a = a0; a < a1; a++)
1179                 {
1180                     i = X4_IND_A(cell[a]);
1181
1182                     f[a][XX] += fnb[i+XX*PACK_X4];
1183                     f[a][YY] += fnb[i+YY*PACK_X4];
1184                     f[a][ZZ] += fnb[i+ZZ*PACK_X4];
1185                 }
1186             }
1187             else
1188             {
1189                 for (a = a0; a < a1; a++)
1190                 {
1191                     i = X4_IND_A(cell[a]);
1192
1193                     for (fa = 0; fa < nfa; fa++)
1194                     {
1195                         f[a][XX] += out[fa].f[i+XX*PACK_X4];
1196                         f[a][YY] += out[fa].f[i+YY*PACK_X4];
1197                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X4];
1198                     }
1199                 }
1200             }
1201             break;
1202         case nbatX8:
1203             if (nfa == 1)
1204             {
1205                 fnb = out[0].f;
1206
1207                 for (a = a0; a < a1; a++)
1208                 {
1209                     i = X8_IND_A(cell[a]);
1210
1211                     f[a][XX] += fnb[i+XX*PACK_X8];
1212                     f[a][YY] += fnb[i+YY*PACK_X8];
1213                     f[a][ZZ] += fnb[i+ZZ*PACK_X8];
1214                 }
1215             }
1216             else
1217             {
1218                 for (a = a0; a < a1; a++)
1219                 {
1220                     i = X8_IND_A(cell[a]);
1221
1222                     for (fa = 0; fa < nfa; fa++)
1223                     {
1224                         f[a][XX] += out[fa].f[i+XX*PACK_X8];
1225                         f[a][YY] += out[fa].f[i+YY*PACK_X8];
1226                         f[a][ZZ] += out[fa].f[i+ZZ*PACK_X8];
1227                     }
1228                 }
1229             }
1230             break;
1231         default:
1232             gmx_incons("Unsupported nbnxn_atomdata_t format");
1233     }
1234 }
1235
1236 /* Add the force array(s) from nbnxn_atomdata_t to f */
1237 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t    nbs,
1238                                     int                     locality,
1239                                     const nbnxn_atomdata_t *nbat,
1240                                     rvec                   *f)
1241 {
1242     int a0 = 0, na = 0;
1243     int nth, th;
1244
1245     nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1246
1247     switch (locality)
1248     {
1249         case eatAll:
1250             a0 = 0;
1251             na = nbs->natoms_nonlocal;
1252             break;
1253         case eatLocal:
1254             a0 = 0;
1255             na = nbs->natoms_local;
1256             break;
1257         case eatNonlocal:
1258             a0 = nbs->natoms_local;
1259             na = nbs->natoms_nonlocal - nbs->natoms_local;
1260             break;
1261     }
1262
1263     nth = gmx_omp_nthreads_get(emntNonbonded);
1264
1265     if (nbat->nout > 1)
1266     {
1267         if (locality != eatAll)
1268         {
1269             gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1270         }
1271
1272         /* Reduce the force thread output buffers into buffer 0, before adding
1273          * them to the, differently ordered, "real" force buffer.
1274          */
1275 #pragma omp parallel for num_threads(nth) schedule(static)
1276         for (th = 0; th < nth; th++)
1277         {
1278             const nbnxn_buffer_flags_t *flags;
1279             int   b0, b1, b;
1280             int   i0, i1;
1281             int   nfptr;
1282             real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1283             int   out;
1284
1285             flags = &nbat->buffer_flags;
1286
1287             /* Calculate the cell-block range for our thread */
1288             b0 = (flags->nflag* th   )/nth;
1289             b1 = (flags->nflag*(th+1))/nth;
1290
1291             for (b = b0; b < b1; b++)
1292             {
1293                 i0 =  b   *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1294                 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1295
1296                 nfptr = 0;
1297                 for (out = 1; out < nbat->nout; out++)
1298                 {
1299                     if (flags->flag[b] & (1U<<out))
1300                     {
1301                         fptr[nfptr++] = nbat->out[out].f;
1302                     }
1303                 }
1304                 if (nfptr > 0)
1305                 {
1306 #ifdef GMX_NBNXN_SIMD
1307                     nbnxn_atomdata_reduce_reals_simd
1308 #else
1309                     nbnxn_atomdata_reduce_reals
1310 #endif
1311                         (nbat->out[0].f,
1312                         flags->flag[b] & (1U<<0),
1313                         fptr, nfptr,
1314                         i0, i1);
1315                 }
1316                 else if (!(flags->flag[b] & (1U<<0)))
1317                 {
1318                     nbnxn_atomdata_clear_reals(nbat->out[0].f,
1319                                                i0, i1);
1320                 }
1321             }
1322         }
1323     }
1324
1325 #pragma omp parallel for num_threads(nth) schedule(static)
1326     for (th = 0; th < nth; th++)
1327     {
1328         nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1329                                             nbat->out,
1330                                             1,
1331                                             a0+((th+0)*na)/nth,
1332                                             a0+((th+1)*na)/nth,
1333                                             f);
1334     }
1335
1336     nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1337 }
1338
1339 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1340 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1341                                               rvec                   *fshift)
1342 {
1343     const nbnxn_atomdata_output_t *out;
1344     int  th;
1345     int  s;
1346     rvec sum;
1347
1348     out = nbat->out;
1349
1350     for (s = 0; s < SHIFTS; s++)
1351     {
1352         clear_rvec(sum);
1353         for (th = 0; th < nbat->nout; th++)
1354         {
1355             sum[XX] += out[th].fshift[s*DIM+XX];
1356             sum[YY] += out[th].fshift[s*DIM+YY];
1357             sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1358         }
1359         rvec_inc(fshift[s], sum);
1360     }
1361 }