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