Remove unnecessary config.h includes
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic_cg.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include <math.h>
40
41 #include "gromacs/legacyheaders/types/simple.h"
42 #include "gromacs/math/vec.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "nb_generic_cg.h"
45 #include "gromacs/legacyheaders/nonbonded.h"
46 #include "nb_kernel.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/utility/fatalerror.h"
50
51 void
52 gmx_nb_generic_cg_kernel(t_nblist *                nlist,
53                          rvec *                    xx,
54                          rvec *                    ff,
55                          t_forcerec *              fr,
56                          t_mdatoms *               mdatoms,
57                          nb_kernel_data_t *        kernel_data,
58                          t_nrnb *                  nrnb)
59 {
60     int           nri, ntype, table_nelements, ielec, ivdw;
61     real          facel, gbtabscale;
62     int           n, is3, i3, k, nj0, nj1, j3, ggid, nnn, n0;
63     int           ai0, ai1, ai, aj0, aj1, aj;
64     real          shX, shY, shZ;
65     real          fscal, tx, ty, tz;
66     real          rinvsq;
67     real          iq;
68     real          qq, vcoul, krsq, vctot;
69     int           nti, nvdwparam;
70     int           tj;
71     real          rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
72     real          rinvsix;
73     real          Vvdwtot;
74     real          Vvdw_rep, Vvdw_disp;
75     real          ix, iy, iz, fix, fiy, fiz;
76     real          jx, jy, jz;
77     real          dx, dy, dz, rsq, rinv;
78     real          c6, c12, cexp1, cexp2, br;
79     real *        charge;
80     real *        shiftvec;
81     real *        vdwparam;
82     int *         shift;
83     int *         type;
84     t_excl *      excl;
85     real *        fshift;
86     real *        Vc;
87     real *        Vvdw;
88     real          tabscale;
89     real *        VFtab;
90     real *        x;
91     real *        f;
92
93     x                   = xx[0];
94     f                   = ff[0];
95     ielec               = nlist->ielec;
96     ivdw                = nlist->ivdw;
97
98     fshift              = fr->fshift[0];
99     Vc                  = kernel_data->energygrp_elec;
100     Vvdw                = kernel_data->energygrp_vdw;
101     tabscale            = kernel_data->table_elec_vdw->scale;
102     VFtab               = kernel_data->table_elec_vdw->data;
103
104     /* avoid compiler warnings for cases that cannot happen */
105     nnn                 = 0;
106     vcoul               = 0.0;
107     eps                 = 0.0;
108     eps2                = 0.0;
109
110     /* 3 VdW parameters for buckingham, otherwise 2 */
111     nvdwparam           = (nlist->ivdw == 2) ? 3 : 2;
112     table_nelements     = (ielec == 3) ? 4 : 0;
113     table_nelements    += (ivdw == 3) ? 8 : 0;
114
115     charge              = mdatoms->chargeA;
116     type                = mdatoms->typeA;
117     facel               = fr->epsfac;
118     shiftvec            = fr->shift_vec[0];
119     vdwparam            = fr->nbfp;
120     ntype               = fr->ntype;
121
122     for (n = 0; (n < nlist->nri); n++)
123     {
124         is3              = 3*nlist->shift[n];
125         shX              = shiftvec[is3];
126         shY              = shiftvec[is3+1];
127         shZ              = shiftvec[is3+2];
128         nj0              = nlist->jindex[n];
129         nj1              = nlist->jindex[n+1];
130         ai0              = nlist->iinr[n];
131         ai1              = nlist->iinr_end[n];
132         vctot            = 0;
133         Vvdwtot          = 0;
134         fix              = 0;
135         fiy              = 0;
136         fiz              = 0;
137
138         for (k = nj0; (k < nj1); k++)
139         {
140             aj0              = nlist->jjnr[k];
141             aj1              = nlist->jjnr_end[k];
142             excl             = &nlist->excl[k*MAX_CGCGSIZE];
143
144             for (ai = ai0; (ai < ai1); ai++)
145             {
146                 i3               = ai*3;
147                 ix               = shX + x[i3+0];
148                 iy               = shY + x[i3+1];
149                 iz               = shZ + x[i3+2];
150                 iq               = facel*charge[ai];
151                 nti              = nvdwparam*ntype*type[ai];
152
153                 /* Note that this code currently calculates
154                  * all LJ and Coulomb interactions,
155                  * even if the LJ parameters or charges are zero.
156                  * If required, this can be optimized.
157                  */
158
159                 for (aj = aj0; (aj < aj1); aj++)
160                 {
161                     /* Check if this interaction is excluded */
162                     if (excl[aj-aj0] & (1<<(ai-ai0)))
163                     {
164                         continue;
165                     }
166
167                     j3               = aj*3;
168                     jx               = x[j3+0];
169                     jy               = x[j3+1];
170                     jz               = x[j3+2];
171                     dx               = ix - jx;
172                     dy               = iy - jy;
173                     dz               = iz - jz;
174                     rsq              = dx*dx+dy*dy+dz*dz;
175                     rinv             = gmx_invsqrt(rsq);
176                     rinvsq           = rinv*rinv;
177                     fscal            = 0;
178
179                     if (ielec == 3 || ivdw == 3)
180                     {
181                         r                = rsq*rinv;
182                         rt               = r*tabscale;
183                         n0               = rt;
184                         eps              = rt-n0;
185                         eps2             = eps*eps;
186                         nnn              = table_nelements*n0;
187                     }
188
189                     /* Coulomb interaction. ielec==0 means no interaction */
190                     if (ielec > 0)
191                     {
192                         qq               = iq*charge[aj];
193
194                         switch (ielec)
195                         {
196                             case 1:
197                                 /* Vanilla cutoff coulomb */
198                                 vcoul            = qq*rinv;
199                                 fscal            = vcoul*rinvsq;
200                                 break;
201
202                             case 2:
203                                 /* Reaction-field */
204                                 krsq             = fr->k_rf*rsq;
205                                 vcoul            = qq*(rinv+krsq-fr->c_rf);
206                                 fscal            = qq*(rinv-2.0*krsq)*rinvsq;
207                                 break;
208
209                             case 3:
210                                 /* Tabulated coulomb */
211                                 Y                = VFtab[nnn];
212                                 F                = VFtab[nnn+1];
213                                 Geps             = eps*VFtab[nnn+2];
214                                 Heps2            = eps2*VFtab[nnn+3];
215                                 nnn             += 4;
216                                 Fp               = F+Geps+Heps2;
217                                 VV               = Y+eps*Fp;
218                                 FF               = Fp+Geps+2.0*Heps2;
219                                 vcoul            = qq*VV;
220                                 fscal            = -qq*FF*tabscale*rinv;
221                                 break;
222
223                             case 4:
224                                 /* GB */
225                                 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
226                                 break;
227
228                             default:
229                                 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
230                                 break;
231                         }
232                         vctot            = vctot+vcoul;
233                     }  /* End of coulomb interactions */
234
235
236                     /* VdW interaction. ivdw==0 means no interaction */
237                     if (ivdw > 0)
238                     {
239                         tj               = nti+nvdwparam*type[aj];
240
241                         switch (ivdw)
242                         {
243                             case 1:
244                                 /* Vanilla Lennard-Jones cutoff */
245                                 c6               = vdwparam[tj];
246                                 c12              = vdwparam[tj+1];
247
248                                 rinvsix          = rinvsq*rinvsq*rinvsq;
249                                 Vvdw_disp        = c6*rinvsix;
250                                 Vvdw_rep         = c12*rinvsix*rinvsix;
251                                 fscal           += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
252                                 Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
253                                 break;
254
255                             case 2:
256                                 /* Buckingham */
257                                 c6               = vdwparam[tj];
258                                 cexp1            = vdwparam[tj+1];
259                                 cexp2            = vdwparam[tj+2];
260
261                                 rinvsix          = rinvsq*rinvsq*rinvsq;
262                                 Vvdw_disp        = c6*rinvsix;
263                                 br               = cexp2*rsq*rinv;
264                                 Vvdw_rep         = cexp1*exp(-br);
265                                 fscal           += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
266                                 Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
267                                 break;
268
269                             case 3:
270                                 /* Tabulated VdW */
271                                 c6               = vdwparam[tj];
272                                 c12              = vdwparam[tj+1];
273
274                                 Y                = VFtab[nnn];
275                                 F                = VFtab[nnn+1];
276                                 Geps             = eps*VFtab[nnn+2];
277                                 Heps2            = eps2*VFtab[nnn+3];
278                                 Fp               = F+Geps+Heps2;
279                                 VV               = Y+eps*Fp;
280                                 FF               = Fp+Geps+2.0*Heps2;
281                                 Vvdw_disp        = c6*VV;
282                                 fijD             = c6*FF;
283                                 nnn             += 4;
284                                 Y                = VFtab[nnn];
285                                 F                = VFtab[nnn+1];
286                                 Geps             = eps*VFtab[nnn+2];
287                                 Heps2            = eps2*VFtab[nnn+3];
288                                 Fp               = F+Geps+Heps2;
289                                 VV               = Y+eps*Fp;
290                                 FF               = Fp+Geps+2.0*Heps2;
291                                 Vvdw_rep         = c12*VV;
292                                 fijR             = c12*FF;
293                                 fscal           += -(fijD+fijR)*tabscale*rinv;
294                                 Vvdwtot          = Vvdwtot + Vvdw_disp + Vvdw_rep;
295                                 break;
296
297                             default:
298                                 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
299                                 break;
300                         }
301                     }  /* end VdW interactions */
302
303
304                     tx               = fscal*dx;
305                     ty               = fscal*dy;
306                     tz               = fscal*dz;
307                     f[i3+0]         += tx;
308                     f[i3+1]         += ty;
309                     f[i3+2]         += tz;
310                     f[j3+0]         -= tx;
311                     f[j3+1]         -= ty;
312                     f[j3+2]         -= tz;
313                     fix             += tx;
314                     fiy             += ty;
315                     fiz             += tz;
316                 }
317             }
318         }
319
320         fshift[is3]     += fix;
321         fshift[is3+1]   += fiy;
322         fshift[is3+2]   += fiz;
323         ggid             = nlist->gid[n];
324         Vc[ggid]        += vctot;
325         Vvdw[ggid]      += Vvdwtot;
326     }
327     /* Estimate flops, average for generic cg kernel:
328      * 12  flops per outer iteration
329      * 100 flops per inner iteration
330      */
331     inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);
332 }