Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / gmxlib / chargegroup.cpp
blobdcb2b656466e8156f673b48d8305e4130a259642
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2018, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "chargegroup.h"
41 #include <cmath>
43 #include "gromacs/math/vec.h"
44 #include "gromacs/pbcutil/pbc.h"
45 #include "gromacs/topology/topology.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
50 void calc_chargegroup_radii(const gmx_mtop_t *mtop, rvec *x,
51 real *rvdw1, real *rvdw2,
52 real *rcoul1, real *rcoul2)
54 real r2v1, r2v2, r2c1, r2c2, r2;
55 int ntype, i, j, m, cg, a_mol, a0, a1, a;
56 gmx_bool *bLJ;
57 rvec cen;
59 r2v1 = 0;
60 r2v2 = 0;
61 r2c1 = 0;
62 r2c2 = 0;
64 ntype = mtop->ffparams.atnr;
65 snew(bLJ, ntype);
66 for (i = 0; i < ntype; i++)
68 bLJ[i] = FALSE;
69 for (j = 0; j < ntype; j++)
71 if (mtop->ffparams.iparams[i*ntype+j].lj.c6 != 0 ||
72 mtop->ffparams.iparams[i*ntype+j].lj.c12 != 0)
74 bLJ[i] = TRUE;
79 a_mol = 0;
80 for (const gmx_molblock_t &molb : mtop->molblock)
82 const gmx_moltype_t *molt = &mtop->moltype[molb.type];
83 const t_block *cgs = &molt->cgs;
84 const t_atom *atom = molt->atoms.atom;
85 for (m = 0; m < molb.nmol; m++)
87 for (cg = 0; cg < cgs->nr; cg++)
89 a0 = cgs->index[cg];
90 a1 = cgs->index[cg+1];
91 if (a1 - a0 > 1)
93 clear_rvec(cen);
94 for (a = a0; a < a1; a++)
96 rvec_inc(cen, x[a_mol+a]);
98 svmul(1.0/(a1-a0), cen, cen);
99 for (a = a0; a < a1; a++)
101 r2 = distance2(cen, x[a_mol+a]);
102 if (r2 > r2v2 && (bLJ[atom[a].type ] ||
103 bLJ[atom[a].typeB]))
105 if (r2 > r2v1)
107 r2v2 = r2v1;
108 r2v1 = r2;
110 else
112 r2v2 = r2;
115 if (r2 > r2c2 &&
116 (atom[a].q != 0 || atom[a].qB != 0))
118 if (r2 > r2c1)
120 r2c2 = r2c1;
121 r2c1 = r2;
123 else
125 r2c2 = r2;
131 a_mol += molt->atoms.nr;
135 sfree(bLJ);
137 *rvdw1 = std::sqrt(r2v1);
138 *rvdw2 = std::sqrt(r2v2);
139 *rcoul1 = std::sqrt(r2c1);
140 *rcoul2 = std::sqrt(r2c2);
143 void calc_cgcm(FILE gmx_unused *fplog, int cg0, int cg1, const t_block *cgs,
144 rvec pos[], rvec cg_cm[])
146 int icg, k, k0, k1, d;
147 rvec cg;
148 real nrcg, inv_ncg;
149 int *cgindex;
151 #ifdef DEBUG
152 fprintf(fplog, "Calculating centre of geometry for charge groups %d to %d\n",
153 cg0, cg1);
154 #endif
155 cgindex = cgs->index;
157 /* Compute the center of geometry for all charge groups */
158 for (icg = cg0; (icg < cg1); icg++)
160 k0 = cgindex[icg];
161 k1 = cgindex[icg+1];
162 nrcg = k1-k0;
163 if (nrcg == 1)
165 copy_rvec(pos[k0], cg_cm[icg]);
167 else
169 inv_ncg = 1.0/nrcg;
171 clear_rvec(cg);
172 for (k = k0; (k < k1); k++)
174 for (d = 0; (d < DIM); d++)
176 cg[d] += pos[k][d];
179 for (d = 0; (d < DIM); d++)
181 cg_cm[icg][d] = inv_ncg*cg[d];
187 void put_charge_groups_in_box(FILE gmx_unused *fplog, int cg0, int cg1,
188 int ePBC, matrix box, const t_block *cgs,
189 rvec pos[], rvec cg_cm[])
192 int npbcdim, icg, k, k0, k1, d, e;
193 rvec cg;
194 real nrcg, inv_ncg;
195 int *cgindex;
196 gmx_bool bTric;
198 if (ePBC == epbcNONE)
200 gmx_incons("Calling put_charge_groups_in_box for a system without PBC");
203 #ifdef DEBUG
204 fprintf(fplog, "Putting cgs %d to %d in box\n", cg0, cg1);
205 #endif
206 cgindex = cgs->index;
208 if (ePBC == epbcXY)
210 npbcdim = 2;
212 else
214 npbcdim = 3;
217 bTric = TRICLINIC(box);
219 for (icg = cg0; (icg < cg1); icg++)
221 /* First compute the center of geometry for this charge group */
222 k0 = cgindex[icg];
223 k1 = cgindex[icg+1];
224 nrcg = k1-k0;
226 if (nrcg == 1)
228 copy_rvec(pos[k0], cg_cm[icg]);
230 else
232 inv_ncg = 1.0/nrcg;
234 clear_rvec(cg);
235 for (k = k0; (k < k1); k++)
237 for (d = 0; d < DIM; d++)
239 cg[d] += pos[k][d];
242 for (d = 0; d < DIM; d++)
244 cg_cm[icg][d] = inv_ncg*cg[d];
247 /* Now check pbc for this cg */
248 if (bTric)
250 for (d = npbcdim-1; d >= 0; d--)
252 while (cg_cm[icg][d] < 0)
254 for (e = d; e >= 0; e--)
256 cg_cm[icg][e] += box[d][e];
257 for (k = k0; (k < k1); k++)
259 pos[k][e] += box[d][e];
263 while (cg_cm[icg][d] >= box[d][d])
265 for (e = d; e >= 0; e--)
267 cg_cm[icg][e] -= box[d][e];
268 for (k = k0; (k < k1); k++)
270 pos[k][e] -= box[d][e];
276 else
278 for (d = 0; d < npbcdim; d++)
280 while (cg_cm[icg][d] < 0)
282 cg_cm[icg][d] += box[d][d];
283 for (k = k0; (k < k1); k++)
285 pos[k][d] += box[d][d];
288 while (cg_cm[icg][d] >= box[d][d])
290 cg_cm[icg][d] -= box[d][d];
291 for (k = k0; (k < k1); k++)
293 pos[k][d] -= box[d][d];
298 #ifdef DEBUG_PBC
299 for (d = 0; (d < npbcdim); d++)
301 if ((cg_cm[icg][d] < 0) || (cg_cm[icg][d] >= box[d][d]))
303 gmx_fatal(FARGS, "cg_cm[%d] = %15f %15f %15f\n"
304 "box = %15f %15f %15f\n",
305 icg, cg_cm[icg][XX], cg_cm[icg][YY], cg_cm[icg][ZZ],
306 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
309 #endif