Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / chargegroup.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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "vec.h"
45 #include "pbc.h"
46 #include "smalloc.h"
47 #include "gmx_fatal.h"
48 #include "chargegroup.h"
49
50
51 void calc_chargegroup_radii(const gmx_mtop_t *mtop, rvec *x,
52                             real *rvdw1, real *rvdw2,
53                             real *rcoul1, real *rcoul2)
54 {
55     real            r2v1, r2v2, r2c1, r2c2, r2;
56     int             ntype, i, j, mb, m, cg, a_mol, a0, a1, a;
57     gmx_bool       *bLJ;
58     gmx_molblock_t *molb;
59     gmx_moltype_t  *molt;
60     t_block        *cgs;
61     t_atom         *atom;
62     rvec            cen;
63
64     r2v1 = 0;
65     r2v2 = 0;
66     r2c1 = 0;
67     r2c2 = 0;
68
69     ntype = mtop->ffparams.atnr;
70     snew(bLJ, ntype);
71     for (i = 0; i < ntype; i++)
72     {
73         bLJ[i] = FALSE;
74         for (j = 0; j < ntype; j++)
75         {
76             if (mtop->ffparams.iparams[i*ntype+j].lj.c6  != 0 ||
77                 mtop->ffparams.iparams[i*ntype+j].lj.c12 != 0)
78             {
79                 bLJ[i] = TRUE;
80             }
81         }
82     }
83
84     a_mol = 0;
85     for (mb = 0; mb < mtop->nmolblock; mb++)
86     {
87         molb = &mtop->molblock[mb];
88         molt = &mtop->moltype[molb->type];
89         cgs  = &molt->cgs;
90         atom = molt->atoms.atom;
91         for (m = 0; m < molb->nmol; m++)
92         {
93             for (cg = 0; cg < cgs->nr; cg++)
94             {
95                 a0 = cgs->index[cg];
96                 a1 = cgs->index[cg+1];
97                 if (a1 - a0 > 1)
98                 {
99                     clear_rvec(cen);
100                     for (a = a0; a < a1; a++)
101                     {
102                         rvec_inc(cen, x[a_mol+a]);
103                     }
104                     svmul(1.0/(a1-a0), cen, cen);
105                     for (a = a0; a < a1; a++)
106                     {
107                         r2 = distance2(cen, x[a_mol+a]);
108                         if (r2 > r2v2 && (bLJ[atom[a].type ] ||
109                                           bLJ[atom[a].typeB]))
110                         {
111                             if (r2 > r2v1)
112                             {
113                                 r2v2 = r2v1;
114                                 r2v1 = r2;
115                             }
116                             else
117                             {
118                                 r2v2 = r2;
119                             }
120                         }
121                         if (r2 > r2c2 &&
122                             (atom[a].q != 0 || atom[a].qB != 0))
123                         {
124                             if (r2 > r2c1)
125                             {
126                                 r2c2 = r2c1;
127                                 r2c1 = r2;
128                             }
129                             else
130                             {
131                                 r2c2 = r2;
132                             }
133                         }
134                     }
135                 }
136             }
137             a_mol += molb->natoms_mol;
138         }
139     }
140
141     sfree(bLJ);
142
143     *rvdw1  = sqrt(r2v1);
144     *rvdw2  = sqrt(r2v2);
145     *rcoul1 = sqrt(r2c1);
146     *rcoul2 = sqrt(r2c2);
147 }
148
149 void calc_cgcm(FILE gmx_unused *fplog, int cg0, int cg1, t_block *cgs,
150                rvec pos[], rvec cg_cm[])
151 {
152     int      icg, k, k0, k1, d;
153     rvec     cg;
154     real     nrcg, inv_ncg;
155     atom_id *cgindex;
156
157 #ifdef DEBUG
158     fprintf(fplog, "Calculating centre of geometry for charge groups %d to %d\n",
159             cg0, cg1);
160 #endif
161     cgindex = cgs->index;
162
163     /* Compute the center of geometry for all charge groups */
164     for (icg = cg0; (icg < cg1); icg++)
165     {
166         k0      = cgindex[icg];
167         k1      = cgindex[icg+1];
168         nrcg    = k1-k0;
169         if (nrcg == 1)
170         {
171             copy_rvec(pos[k0], cg_cm[icg]);
172         }
173         else
174         {
175             inv_ncg = 1.0/nrcg;
176
177             clear_rvec(cg);
178             for (k = k0; (k < k1); k++)
179             {
180                 for (d = 0; (d < DIM); d++)
181                 {
182                     cg[d] += pos[k][d];
183                 }
184             }
185             for (d = 0; (d < DIM); d++)
186             {
187                 cg_cm[icg][d] = inv_ncg*cg[d];
188             }
189         }
190     }
191 }
192
193 void put_charge_groups_in_box(FILE gmx_unused *fplog, int cg0, int cg1,
194                               int ePBC, matrix box, t_block *cgs,
195                               rvec pos[], rvec cg_cm[])
196
197 {
198     int      npbcdim, icg, k, k0, k1, d, e;
199     rvec     cg;
200     real     nrcg, inv_ncg;
201     atom_id *cgindex;
202     gmx_bool bTric;
203
204     if (ePBC == epbcNONE)
205     {
206         gmx_incons("Calling put_charge_groups_in_box for a system without PBC");
207     }
208
209 #ifdef DEBUG
210     fprintf(fplog, "Putting cgs %d to %d in box\n", cg0, cg1);
211 #endif
212     cgindex = cgs->index;
213
214     if (ePBC == epbcXY)
215     {
216         npbcdim = 2;
217     }
218     else
219     {
220         npbcdim = 3;
221     }
222
223     bTric = TRICLINIC(box);
224
225     for (icg = cg0; (icg < cg1); icg++)
226     {
227         /* First compute the center of geometry for this charge group */
228         k0      = cgindex[icg];
229         k1      = cgindex[icg+1];
230         nrcg    = k1-k0;
231
232         if (nrcg == 1)
233         {
234             copy_rvec(pos[k0], cg_cm[icg]);
235         }
236         else
237         {
238             inv_ncg = 1.0/nrcg;
239
240             clear_rvec(cg);
241             for (k = k0; (k < k1); k++)
242             {
243                 for (d = 0; d < DIM; d++)
244                 {
245                     cg[d] += pos[k][d];
246                 }
247             }
248             for (d = 0; d < DIM; d++)
249             {
250                 cg_cm[icg][d] = inv_ncg*cg[d];
251             }
252         }
253         /* Now check pbc for this cg */
254         if (bTric)
255         {
256             for (d = npbcdim-1; d >= 0; d--)
257             {
258                 while (cg_cm[icg][d] < 0)
259                 {
260                     for (e = d; e >= 0; e--)
261                     {
262                         cg_cm[icg][e] += box[d][e];
263                         for (k = k0; (k < k1); k++)
264                         {
265                             pos[k][e] += box[d][e];
266                         }
267                     }
268                 }
269                 while (cg_cm[icg][d] >= box[d][d])
270                 {
271                     for (e = d; e >= 0; e--)
272                     {
273                         cg_cm[icg][e] -= box[d][e];
274                         for (k = k0; (k < k1); k++)
275                         {
276                             pos[k][e] -= box[d][e];
277                         }
278                     }
279                 }
280             }
281         }
282         else
283         {
284             for (d = 0; d < npbcdim; d++)
285             {
286                 while (cg_cm[icg][d] < 0)
287                 {
288                     cg_cm[icg][d] += box[d][d];
289                     for (k = k0; (k < k1); k++)
290                     {
291                         pos[k][d] += box[d][d];
292                     }
293                 }
294                 while (cg_cm[icg][d] >= box[d][d])
295                 {
296                     cg_cm[icg][d] -= box[d][d];
297                     for (k = k0; (k < k1); k++)
298                     {
299                         pos[k][d] -= box[d][d];
300                     }
301                 }
302             }
303         }
304 #ifdef DEBUG_PBC
305         for (d = 0; (d < npbcdim); d++)
306         {
307             if ((cg_cm[icg][d] < 0) || (cg_cm[icg][d] >= box[d][d]))
308             {
309                 gmx_fatal(FARGS, "cg_cm[%d] = %15f  %15f  %15f\n"
310                           "box = %15f  %15f  %15f\n",
311                           icg, cg_cm[icg][XX], cg_cm[icg][YY], cg_cm[icg][ZZ],
312                           box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
313             }
314         }
315 #endif
316     }
317 }