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-2008, 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.
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.
40 #include "gromacs/legacyheaders/genborn.h"
45 #include "gromacs/fileio/pdbio.h"
46 #include "gromacs/legacyheaders/domdec.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/network.h"
49 #include "gromacs/legacyheaders/nrnb.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/types/commrec.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/ishift.h"
55 #include "gromacs/pbcutil/mshift.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxmpi.h"
60 #include "gromacs/utility/smalloc.h"
62 #ifdef GMX_SIMD_X86_SSE2_OR_HIGHER
64 # include "gromacs/mdlib/genborn_allvsall_sse2_double.h"
65 # include "gromacs/mdlib/genborn_sse2_double.h"
67 # include "gromacs/mdlib/genborn_allvsall_sse2_single.h"
68 # include "gromacs/mdlib/genborn_sse2_single.h"
69 # endif /* GMX_DOUBLE */
70 #endif /* SSE or AVX present */
72 #include "gromacs/mdlib/genborn_allvsall.h"
74 /*#define DISABLE_SSE*/
83 typedef struct gbtmpnbls {
89 /* This function is exactly the same as the one in listed-forces/bonded.cpp. The reason
90 * it is copied here is that the bonded gb-interactions are evaluated
91 * not in calc_bonds, but rather in calc_gb_forces
93 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
97 return pbc_dx_aiuc(pbc, xi, xj, dx);
101 rvec_sub(xi, xj, dx);
106 static int init_gb_nblist(int natoms, t_nblist *nl)
108 nl->maxnri = natoms*4;
117 /*nl->nltype = nltype;*/
119 srenew(nl->iinr, nl->maxnri);
120 srenew(nl->gid, nl->maxnri);
121 srenew(nl->shift, nl->maxnri);
122 srenew(nl->jindex, nl->maxnri+1);
130 static int init_gb_still(const t_atomtypes *atype, t_idef *idef, t_atoms *atoms,
131 gmx_genborn_t *born, int natoms)
134 int i, j, i1, i2, k, m, nbond, nang, ia, ib, ic, id, nb, idx, idx2, at;
138 real r, ri, rj, ri2, ri3, rj2, r2, r3, r4, rk, ratio, term, h, doffset;
139 real p1, p2, p3, factor, cosine, rab, rbc;
146 snew(born->gpol_still_work, natoms+3);
151 doffset = born->gb_doffset;
153 for (i = 0; i < natoms; i++)
155 born->gpol_globalindex[i] = born->vsolv_globalindex[i] =
156 born->gb_radius_globalindex[i] = 0;
159 /* Compute atomic solvation volumes for Still method */
160 for (i = 0; i < natoms; i++)
162 ri = atype->gb_radius[atoms->atom[i].type];
163 born->gb_radius_globalindex[i] = ri;
165 born->vsolv_globalindex[i] = (4*M_PI/3)*r3;
168 for (j = 0; j < idef->il[F_GB12].nr; j += 3)
170 m = idef->il[F_GB12].iatoms[j];
171 ia = idef->il[F_GB12].iatoms[j+1];
172 ib = idef->il[F_GB12].iatoms[j+2];
174 r = 1.01*idef->iparams[m].gb.st;
176 ri = atype->gb_radius[atoms->atom[ia].type];
177 rj = atype->gb_radius[atoms->atom[ib].type];
183 ratio = (rj2-ri2-r*r)/(2*ri*r);
185 term = (M_PI/3.0)*h*h*(3.0*ri-h);
187 born->vsolv_globalindex[ia] -= term;
189 ratio = (ri2-rj2-r*r)/(2*rj*r);
191 term = (M_PI/3.0)*h*h*(3.0*rj-h);
193 born->vsolv_globalindex[ib] -= term;
196 /* Get the self-, 1-2 and 1-3 polarization energies for analytical Still
199 for (j = 0; j < natoms; j++)
201 if (born->use_globalindex[j] == 1)
203 born->gpol_globalindex[j] = -0.5*ONE_4PI_EPS0/
204 (atype->gb_radius[atoms->atom[j].type]-doffset+STILL_P1);
209 for (j = 0; j < idef->il[F_GB12].nr; j += 3)
211 m = idef->il[F_GB12].iatoms[j];
212 ia = idef->il[F_GB12].iatoms[j+1];
213 ib = idef->il[F_GB12].iatoms[j+2];
215 r = idef->iparams[m].gb.st;
219 born->gpol_globalindex[ia] = born->gpol_globalindex[ia]+
220 STILL_P2*born->vsolv_globalindex[ib]/r4;
221 born->gpol_globalindex[ib] = born->gpol_globalindex[ib]+
222 STILL_P2*born->vsolv_globalindex[ia]/r4;
226 for (j = 0; j < idef->il[F_GB13].nr; j += 3)
228 m = idef->il[F_GB13].iatoms[j];
229 ia = idef->il[F_GB13].iatoms[j+1];
230 ib = idef->il[F_GB13].iatoms[j+2];
232 r = idef->iparams[m].gb.st;
235 born->gpol_globalindex[ia] = born->gpol_globalindex[ia]+
236 STILL_P3*born->vsolv_globalindex[ib]/r4;
237 born->gpol_globalindex[ib] = born->gpol_globalindex[ib]+
238 STILL_P3*born->vsolv_globalindex[ia]/r4;
247 /* Initialize all GB datastructs and compute polarization energies */
248 int init_gb(gmx_genborn_t **p_born,
249 t_forcerec *fr, const t_inputrec *ir,
250 const gmx_mtop_t *mtop, int gb_algorithm)
252 int i, j, m, ai, aj, jj, natoms, nalloc;
253 real rai, sk, p, doffset;
257 gmx_localtop_t *localtop;
259 natoms = mtop->natoms;
261 atoms = gmx_mtop_global_atoms(mtop);
262 localtop = gmx_mtop_generate_local_top(mtop, ir);
269 snew(born->drobc, natoms);
270 snew(born->bRad, natoms);
272 /* Allocate memory for the global data arrays */
273 snew(born->param_globalindex, natoms+3);
274 snew(born->gpol_globalindex, natoms+3);
275 snew(born->vsolv_globalindex, natoms+3);
276 snew(born->gb_radius_globalindex, natoms+3);
277 snew(born->use_globalindex, natoms+3);
279 snew(fr->invsqrta, natoms);
280 snew(fr->dvda, natoms);
283 fr->dadx_rawptr = NULL;
285 born->gpol_still_work = NULL;
286 born->gpol_hct_work = NULL;
288 /* snew(born->asurf,natoms); */
289 /* snew(born->dasurf,natoms); */
291 /* Initialize the gb neighbourlist */
292 init_gb_nblist(natoms, &(fr->gblist));
294 /* Do the Vsites exclusions (if any) */
295 for (i = 0; i < natoms; i++)
297 jj = atoms.atom[i].type;
298 if (mtop->atomtypes.gb_radius[atoms.atom[i].type] > 0)
300 born->use_globalindex[i] = 1;
304 born->use_globalindex[i] = 0;
307 /* If we have a Vsite, put vs_globalindex[i]=0 */
308 if (C6 (fr->nbfp, fr->ntype, jj, jj) == 0 &&
309 C12(fr->nbfp, fr->ntype, jj, jj) == 0 &&
310 atoms.atom[i].q == 0)
312 born->use_globalindex[i] = 0;
316 /* Copy algorithm parameters from inputrecord to local structure */
317 born->obc_alpha = ir->gb_obc_alpha;
318 born->obc_beta = ir->gb_obc_beta;
319 born->obc_gamma = ir->gb_obc_gamma;
320 born->gb_doffset = ir->gb_dielectric_offset;
321 born->gb_epsilon_solvent = ir->gb_epsilon_solvent;
322 born->epsilon_r = ir->epsilon_r;
324 doffset = born->gb_doffset;
326 /* Set the surface tension */
327 born->sa_surface_tension = ir->sa_surface_tension;
329 /* If Still model, initialise the polarisation energies */
330 if (gb_algorithm == egbSTILL)
332 init_gb_still(&(mtop->atomtypes), &(localtop->idef), &atoms,
337 /* If HCT/OBC, precalculate the sk*atype->S_hct factors */
338 else if (gb_algorithm == egbHCT || gb_algorithm == egbOBC)
341 snew(born->gpol_hct_work, natoms+3);
343 for (i = 0; i < natoms; i++)
345 if (born->use_globalindex[i] == 1)
347 rai = mtop->atomtypes.gb_radius[atoms.atom[i].type]-doffset;
348 sk = rai * mtop->atomtypes.S_hct[atoms.atom[i].type];
349 born->param_globalindex[i] = sk;
350 born->gb_radius_globalindex[i] = rai;
354 born->param_globalindex[i] = 0;
355 born->gb_radius_globalindex[i] = 0;
360 /* Allocate memory for work arrays for temporary use */
361 snew(born->work, natoms+4);
362 snew(born->count, natoms);
363 snew(born->nblist_work, natoms);
365 /* Domain decomposition specific stuff */
374 calc_gb_rad_still(t_commrec *cr, t_forcerec *fr, gmx_localtop_t *top,
375 rvec x[], t_nblist *nl,
376 gmx_genborn_t *born, t_mdatoms *md)
378 int i, k, n, nj0, nj1, ai, aj, type;
381 real gpi, dr, dr2, dr4, idr4, rvdw, ratio, ccf, theta, term, rai, raj;
382 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
383 real rinv, idr2, idr6, vaj, dccf, cosq, sinq, prod, gpi2;
385 real vai, prod_ai, icf4, icf6;
387 factor = 0.5*ONE_4PI_EPS0;
390 for (i = 0; i < born->nr; i++)
392 born->gpol_still_work[i] = 0;
395 for (i = 0; i < nl->nri; i++)
400 nj1 = nl->jindex[i+1];
402 /* Load shifts for this list */
403 shift = nl->shift[i];
404 shX = fr->shift_vec[shift][0];
405 shY = fr->shift_vec[shift][1];
406 shZ = fr->shift_vec[shift][2];
410 rai = top->atomtypes.gb_radius[md->typeA[ai]];
411 vai = born->vsolv[ai];
412 prod_ai = STILL_P4*vai;
414 /* Load atom i coordinates, add shift vectors */
415 ix1 = shX + x[ai][0];
416 iy1 = shY + x[ai][1];
417 iz1 = shZ + x[ai][2];
419 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
430 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
431 rinv = gmx_invsqrt(dr2);
436 raj = top->atomtypes.gb_radius[md->typeA[aj]];
440 ratio = dr2 / (rvdw * rvdw);
441 vaj = born->vsolv[aj];
443 if (ratio > STILL_P5INV)
450 theta = ratio*STILL_PIP5;
452 term = 0.5*(1.0-cosq);
454 sinq = 1.0 - cosq*cosq;
455 dccf = 2.0*term*sinq*gmx_invsqrt(sinq)*theta;
460 icf6 = (4*ccf-dccf)*idr6;
461 born->gpol_still_work[aj] += prod_ai*icf4;
464 /* Save ai->aj and aj->ai chain rule terms */
465 fr->dadx[n++] = prod*icf6;
466 fr->dadx[n++] = prod_ai*icf6;
468 born->gpol_still_work[ai] += gpi;
471 /* Parallel summations */
472 if (DOMAINDECOMP(cr))
474 dd_atom_sum_real(cr->dd, born->gpol_still_work);
477 /* Calculate the radii */
478 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
480 if (born->use[i] != 0)
482 gpi = born->gpol[i]+born->gpol_still_work[i];
484 born->bRad[i] = factor*gmx_invsqrt(gpi2);
485 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
489 /* Extra communication required for DD */
490 if (DOMAINDECOMP(cr))
492 dd_atom_spread_real(cr->dd, born->bRad);
493 dd_atom_spread_real(cr->dd, fr->invsqrta);
502 calc_gb_rad_hct(t_commrec *cr, t_forcerec *fr, gmx_localtop_t *top,
503 rvec x[], t_nblist *nl,
504 gmx_genborn_t *born, t_mdatoms *md)
506 int i, k, n, ai, aj, nj0, nj1, at0, at1;
509 real rai, raj, gpi, dr2, dr, sk, sk_ai, sk2, sk2_ai, lij, uij, diff2, tmp, sum_ai;
510 real rad, min_rad, rinv, rai_inv;
511 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
512 real lij2, uij2, lij3, uij3, t1, t2, t3;
513 real lij_inv, dlij, duij, sk2_rinv, prod, log_term;
514 real doffset, raj_inv, dadx_val;
517 doffset = born->gb_doffset;
518 gb_radius = born->gb_radius;
520 for (i = 0; i < born->nr; i++)
522 born->gpol_hct_work[i] = 0;
525 /* Keep the compiler happy */
529 for (i = 0; i < nl->nri; i++)
534 nj1 = nl->jindex[i+1];
536 /* Load shifts for this list */
537 shift = nl->shift[i];
538 shX = fr->shift_vec[shift][0];
539 shY = fr->shift_vec[shift][1];
540 shZ = fr->shift_vec[shift][2];
545 sk_ai = born->param[ai];
546 sk2_ai = sk_ai*sk_ai;
548 /* Load atom i coordinates, add shift vectors */
549 ix1 = shX + x[ai][0];
550 iy1 = shY + x[ai][1];
551 iz1 = shZ + x[ai][2];
555 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
567 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
568 rinv = gmx_invsqrt(dr2);
571 sk = born->param[aj];
574 /* aj -> ai interaction */
595 lij_inv = gmx_invsqrt(lij2);
598 prod = 0.25*sk2_rinv;
600 log_term = log(uij*lij_inv);
602 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term +
607 tmp = tmp + 2.0 * (rai_inv-lij);
610 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
611 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
612 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
614 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
615 /* fr->dadx[n++] = (dlij*t1+duij*t2+t3)*rinv; */
616 /* rb2 is moved to chainrule */
624 fr->dadx[n++] = dadx_val;
627 /* ai -> aj interaction */
628 if (raj < dr + sk_ai)
630 lij = 1.0/(dr-sk_ai);
643 uij = 1.0/(dr+sk_ai);
649 lij_inv = gmx_invsqrt(lij2);
650 sk2 = sk2_ai; /* sk2_ai = sk_ai * sk_ai in i loop above */
652 prod = 0.25 * sk2_rinv;
654 /* log_term = table_log(uij*lij_inv,born->log_table,
655 LOG_TABLE_ACCURACY); */
656 log_term = log(uij*lij_inv);
658 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term +
663 tmp = tmp + 2.0 * (raj_inv-lij);
667 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
668 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
669 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
671 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
672 /* fr->dadx[n++] = (dlij*t1+duij*t2+t3)*rinv; */ /* rb2 is moved to chainrule */
674 born->gpol_hct_work[aj] += 0.5*tmp;
680 fr->dadx[n++] = dadx_val;
683 born->gpol_hct_work[ai] += sum_ai;
686 /* Parallel summations */
687 if (DOMAINDECOMP(cr))
689 dd_atom_sum_real(cr->dd, born->gpol_hct_work);
692 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
694 if (born->use[i] != 0)
696 rai = top->atomtypes.gb_radius[md->typeA[i]]-doffset;
697 sum_ai = 1.0/rai - born->gpol_hct_work[i];
698 min_rad = rai + doffset;
701 born->bRad[i] = rad > min_rad ? rad : min_rad;
702 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
706 /* Extra communication required for DD */
707 if (DOMAINDECOMP(cr))
709 dd_atom_spread_real(cr->dd, born->bRad);
710 dd_atom_spread_real(cr->dd, fr->invsqrta);
718 calc_gb_rad_obc(t_commrec *cr, t_forcerec *fr, gmx_localtop_t *top,
719 rvec x[], t_nblist *nl, gmx_genborn_t *born, t_mdatoms *md)
721 int i, k, ai, aj, nj0, nj1, n, at0, at1;
724 real rai, raj, gpi, dr2, dr, sk, sk2, lij, uij, diff2, tmp, sum_ai;
725 real rad, min_rad, sum_ai2, sum_ai3, tsum, tchain, rinv, rai_inv, lij_inv, rai_inv2;
726 real log_term, prod, sk2_rinv, sk_ai, sk2_ai;
727 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
728 real lij2, uij2, lij3, uij3, dlij, duij, t1, t2, t3;
729 real doffset, raj_inv, dadx_val;
732 /* Keep the compiler happy */
737 doffset = born->gb_doffset;
738 gb_radius = born->gb_radius;
740 for (i = 0; i < born->nr; i++)
742 born->gpol_hct_work[i] = 0;
745 for (i = 0; i < nl->nri; i++)
750 nj1 = nl->jindex[i+1];
752 /* Load shifts for this list */
753 shift = nl->shift[i];
754 shX = fr->shift_vec[shift][0];
755 shY = fr->shift_vec[shift][1];
756 shZ = fr->shift_vec[shift][2];
761 sk_ai = born->param[ai];
762 sk2_ai = sk_ai*sk_ai;
764 /* Load atom i coordinates, add shift vectors */
765 ix1 = shX + x[ai][0];
766 iy1 = shY + x[ai][1];
767 iz1 = shZ + x[ai][2];
771 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
783 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
784 rinv = gmx_invsqrt(dr2);
787 /* sk is precalculated in init_gb() */
788 sk = born->param[aj];
791 /* aj -> ai interaction */
811 lij_inv = gmx_invsqrt(lij2);
814 prod = 0.25*sk2_rinv;
816 log_term = log(uij*lij_inv);
818 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term + prod*(-diff2);
822 tmp = tmp + 2.0 * (rai_inv-lij);
826 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
827 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
828 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
830 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
838 fr->dadx[n++] = dadx_val;
840 /* ai -> aj interaction */
841 if (raj < dr + sk_ai)
843 lij = 1.0/(dr-sk_ai);
856 uij = 1.0/(dr+sk_ai);
862 lij_inv = gmx_invsqrt(lij2);
863 sk2 = sk2_ai; /* sk2_ai = sk_ai * sk_ai in i loop above */
865 prod = 0.25 * sk2_rinv;
867 /* log_term = table_log(uij*lij_inv,born->log_table,LOG_TABLE_ACCURACY); */
868 log_term = log(uij*lij_inv);
870 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term + prod*(-diff2);
874 tmp = tmp + 2.0 * (raj_inv-lij);
877 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
878 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
879 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
881 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
883 born->gpol_hct_work[aj] += 0.5*tmp;
890 fr->dadx[n++] = dadx_val;
893 born->gpol_hct_work[ai] += sum_ai;
897 /* Parallel summations */
898 if (DOMAINDECOMP(cr))
900 dd_atom_sum_real(cr->dd, born->gpol_hct_work);
903 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
905 if (born->use[i] != 0)
907 rai = top->atomtypes.gb_radius[md->typeA[i]];
911 sum_ai = rai * born->gpol_hct_work[i];
912 sum_ai2 = sum_ai * sum_ai;
913 sum_ai3 = sum_ai2 * sum_ai;
915 tsum = tanh(born->obc_alpha*sum_ai-born->obc_beta*sum_ai2+born->obc_gamma*sum_ai3);
916 born->bRad[i] = rai_inv - tsum*rai_inv2;
917 born->bRad[i] = 1.0 / born->bRad[i];
919 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
921 tchain = rai * (born->obc_alpha-2*born->obc_beta*sum_ai+3*born->obc_gamma*sum_ai2);
922 born->drobc[i] = (1.0-tsum*tsum)*tchain*rai_inv2;
926 /* Extra (local) communication required for DD */
927 if (DOMAINDECOMP(cr))
929 dd_atom_spread_real(cr->dd, born->bRad);
930 dd_atom_spread_real(cr->dd, fr->invsqrta);
931 dd_atom_spread_real(cr->dd, born->drobc);
940 int calc_gb_rad(t_commrec *cr, t_forcerec *fr, t_inputrec *ir, gmx_localtop_t *top,
941 rvec x[], t_nblist *nl, gmx_genborn_t *born, t_mdatoms *md, t_nrnb *nrnb)
947 if (fr->bAllvsAll && fr->dadx == NULL)
949 /* We might need up to 8 atoms of padding before and after,
950 * and another 4 units to guarantee SSE alignment.
952 fr->nalloc_dadx = 2*(md->homenr+12)*(md->nr/2+1+12);
953 snew(fr->dadx_rawptr, fr->nalloc_dadx);
954 fr->dadx = (real *) (((size_t) fr->dadx_rawptr + 16) & (~((size_t) 15)));
958 /* In the SSE-enabled gb-loops, when writing to dadx, we
959 * always write 2*4 elements at a time, even in the case with only
960 * 1-3 j particles, where we only really need to write 2*(1-3)
961 * elements. This is because we want dadx to be aligned to a 16-
962 * byte boundary, and being able to use _mm_store/load_ps
964 ndadx = 2 * (nl->nrj + 3*nl->nri);
966 /* First, reallocate the dadx array, we need 3 extra for SSE */
967 if (ndadx + 3 > fr->nalloc_dadx)
969 fr->nalloc_dadx = over_alloc_large(ndadx) + 3;
970 srenew(fr->dadx_rawptr, fr->nalloc_dadx);
971 fr->dadx = (real *) (((size_t) fr->dadx_rawptr + 16) & (~((size_t) 15)));
977 cnt = md->homenr*(md->nr/2+1);
979 if (ir->gb_algorithm == egbSTILL)
981 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
982 if (fr->use_simd_kernels)
985 genborn_allvsall_calc_still_radii_sse2_double(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
987 genborn_allvsall_calc_still_radii_sse2_single(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
992 genborn_allvsall_calc_still_radii(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
995 genborn_allvsall_calc_still_radii(fr, md, born, top, x[0], &fr->AllvsAll_workgb);
997 /* 13 flops in outer loop, 47 flops in inner loop */
998 inc_nrnb(nrnb, eNR_BORN_AVA_RADII_STILL, md->homenr*13+cnt*47);
1000 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1002 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1003 if (fr->use_simd_kernels)
1006 genborn_allvsall_calc_hct_obc_radii_sse2_double(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1008 genborn_allvsall_calc_hct_obc_radii_sse2_single(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1013 genborn_allvsall_calc_hct_obc_radii(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1016 genborn_allvsall_calc_hct_obc_radii(fr, md, born, ir->gb_algorithm, top, x[0], &fr->AllvsAll_workgb);
1018 /* 24 flops in outer loop, 183 in inner */
1019 inc_nrnb(nrnb, eNR_BORN_AVA_RADII_HCT_OBC, md->homenr*24+cnt*183);
1023 gmx_fatal(FARGS, "Bad gb algorithm for all-vs-all interactions");
1028 /* Switch for determining which algorithm to use for Born radii calculation */
1031 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1032 /* x86 or x86-64 with GCC inline assembly and/or SSE intrinsics */
1033 switch (ir->gb_algorithm)
1036 if (fr->use_simd_kernels)
1038 calc_gb_rad_still_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born);
1042 calc_gb_rad_still(cr, fr, top, x, nl, born, md);
1046 if (fr->use_simd_kernels)
1048 calc_gb_rad_hct_obc_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born, md, ir->gb_algorithm);
1052 calc_gb_rad_hct(cr, fr, top, x, nl, born, md);
1056 if (fr->use_simd_kernels)
1058 calc_gb_rad_hct_obc_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born, md, ir->gb_algorithm);
1062 calc_gb_rad_obc(cr, fr, born->nr, top, x, nl, born, md);
1067 gmx_fatal(FARGS, "Unknown double precision sse-enabled algorithm for Born radii calculation: %d", ir->gb_algorithm);
1070 switch (ir->gb_algorithm)
1073 calc_gb_rad_still(cr, fr, top, x, nl, born, md);
1076 calc_gb_rad_hct(cr, fr, top, x, nl, born, md);
1079 calc_gb_rad_obc(cr, fr, top, x, nl, born, md);
1083 gmx_fatal(FARGS, "Unknown double precision algorithm for Born radii calculation: %d", ir->gb_algorithm);
1090 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1091 /* x86 or x86-64 with GCC inline assembly and/or SSE intrinsics */
1092 switch (ir->gb_algorithm)
1095 if (fr->use_simd_kernels)
1097 calc_gb_rad_still_sse2_single(cr, fr, born->nr, top, x[0], nl, born);
1101 calc_gb_rad_still(cr, fr, top, x, nl, born, md);
1105 if (fr->use_simd_kernels)
1107 calc_gb_rad_hct_obc_sse2_single(cr, fr, born->nr, top, x[0], nl, born, md, ir->gb_algorithm);
1111 calc_gb_rad_hct(cr, fr, top, x, nl, born, md);
1116 if (fr->use_simd_kernels)
1118 calc_gb_rad_hct_obc_sse2_single(cr, fr, born->nr, top, x[0], nl, born, md, ir->gb_algorithm);
1122 calc_gb_rad_obc(cr, fr, born->nr, top, x, nl, born, md);
1127 gmx_fatal(FARGS, "Unknown sse-enabled algorithm for Born radii calculation: %d", ir->gb_algorithm);
1131 switch (ir->gb_algorithm)
1134 calc_gb_rad_still(cr, fr, top, x, nl, born, md);
1137 calc_gb_rad_hct(cr, fr, top, x, nl, born, md);
1140 calc_gb_rad_obc(cr, fr, top, x, nl, born, md);
1144 gmx_fatal(FARGS, "Unknown algorithm for Born radii calculation: %d", ir->gb_algorithm);
1147 #endif /* Single precision sse */
1149 #endif /* Double or single precision */
1151 if (fr->bAllvsAll == FALSE)
1153 switch (ir->gb_algorithm)
1156 /* 17 flops per outer loop iteration, 47 flops per inner loop */
1157 inc_nrnb(nrnb, eNR_BORN_RADII_STILL, nl->nri*17+nl->nrj*47);
1161 /* 61 (assuming 10 for tanh) flops for outer loop iteration, 183 flops per inner loop */
1162 inc_nrnb(nrnb, eNR_BORN_RADII_HCT_OBC, nl->nri*61+nl->nrj*183);
1175 real gb_bonds_tab(rvec x[], rvec f[], rvec fshift[], real *charge, real *p_gbtabscale,
1176 real *invsqrta, real *dvda, real *GBtab, t_idef *idef, real epsilon_r,
1177 real gb_epsilon_solvent, real facel, const t_pbc *pbc, const t_graph *graph)
1179 int i, j, n0, m, nnn, type, ai, aj;
1185 real isaprod, qq, gbscale, gbtabscale, Y, F, Geps, Heps2, Fp, VV, FF, rt, eps, eps2;
1186 real vgb, fgb, vcoul, fijC, dvdatmp, fscal, dvdaj;
1192 t_iatom *forceatoms;
1194 /* Scale the electrostatics by gb_epsilon_solvent */
1195 facel = facel * ((1.0/epsilon_r) - 1.0/gb_epsilon_solvent);
1197 gbtabscale = *p_gbtabscale;
1200 for (j = F_GB12; j <= F_GB14; j++)
1202 forceatoms = idef->il[j].iatoms;
1204 for (i = 0; i < idef->il[j].nr; )
1206 /* To avoid reading in the interaction type, we just increment i to pass over
1207 * the types in the forceatoms array, this saves some memory accesses
1210 ai = forceatoms[i++];
1211 aj = forceatoms[i++];
1213 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx);
1214 rsq11 = iprod(dx, dx);
1216 isai = invsqrta[ai];
1217 iq = (-1)*facel*charge[ai];
1219 rinv11 = gmx_invsqrt(rsq11);
1220 isaj = invsqrta[aj];
1221 isaprod = isai*isaj;
1222 qq = isaprod*iq*charge[aj];
1223 gbscale = isaprod*gbtabscale;
1232 Geps = eps*GBtab[nnn+2];
1233 Heps2 = eps2*GBtab[nnn+3];
1236 FF = Fp+Geps+2.0*Heps2;
1238 fijC = qq*FF*gbscale;
1239 dvdatmp = -(vgb+fijC*r)*0.5;
1240 dvda[aj] = dvda[aj] + dvdatmp*isaj*isaj;
1241 dvda[ai] = dvda[ai] + dvdatmp*isai*isai;
1242 vctot = vctot + vgb;
1243 fgb = -(fijC)*rinv11;
1247 ivec_sub(SHIFT_IVEC(graph, ai), SHIFT_IVEC(graph, aj), dt);
1251 for (m = 0; (m < DIM); m++) /* 15 */
1256 fshift[ki][m] += fscal;
1257 fshift[CENTRAL][m] -= fscal;
1265 real calc_gb_selfcorrections(t_commrec *cr, int natoms,
1266 real *charge, gmx_genborn_t *born, real *dvda, double facel)
1268 int i, ai, at0, at1;
1269 real rai, e, derb, q, q2, fi, rai_inv, vtot;
1271 if (DOMAINDECOMP(cr))
1274 at1 = cr->dd->nat_home;
1283 /* Scale the electrostatics by gb_epsilon_solvent */
1284 facel = facel * ((1.0/born->epsilon_r) - 1.0/born->gb_epsilon_solvent);
1288 /* Apply self corrections */
1289 for (i = at0; i < at1; i++)
1293 if (born->use[ai] == 1)
1295 rai = born->bRad[ai];
1301 derb = 0.5*e*rai_inv*rai_inv;
1302 dvda[ai] += derb*rai;
1311 real calc_gb_nonpolar(t_commrec *cr, t_forcerec *fr, int natoms, gmx_genborn_t *born, gmx_localtop_t *top,
1312 real *dvda, t_mdatoms *md)
1314 int ai, i, at0, at1;
1315 real e, es, rai, rbi, term, probe, tmp, factor;
1316 real rbi_inv, rbi_inv2;
1318 /* To keep the compiler happy */
1321 if (DOMAINDECOMP(cr))
1324 at1 = cr->dd->nat_home;
1332 /* factor is the surface tension */
1333 factor = born->sa_surface_tension;
1336 // The surface tension factor is 0.0049 for Still model, 0.0054 for HCT/OBC
1337 if(gb_algorithm==egbSTILL)
1339 factor=0.0049*100*CAL2JOULE;
1343 factor=0.0054*100*CAL2JOULE;
1346 /* if(gb_algorithm==egbHCT || gb_algorithm==egbOBC) */
1352 for (i = at0; i < at1; i++)
1356 if (born->use[ai] == 1)
1358 rai = top->atomtypes.gb_radius[md->typeA[ai]];
1359 rbi_inv = fr->invsqrta[ai];
1360 rbi_inv2 = rbi_inv * rbi_inv;
1361 tmp = (rai*rbi_inv2)*(rai*rbi_inv2);
1363 e = factor*term*(rai+probe)*(rai+probe)*tmp;
1364 dvda[ai] = dvda[ai] - 6*e*rbi_inv2;
1374 real calc_gb_chainrule(int natoms, t_nblist *nl, real *dadx, real *dvda, rvec x[], rvec t[], rvec fshift[],
1375 rvec shift_vec[], int gb_algorithm, gmx_genborn_t *born)
1377 int i, k, n, ai, aj, nj0, nj1, n0, n1;
1380 real fgb, fij, rb2, rbi, fix1, fiy1, fiz1;
1381 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11, rsq11;
1382 real rinv11, tx, ty, tz, rbai, rbaj, fgb_ai;
1392 if (gb_algorithm == egbSTILL)
1394 for (i = n0; i < n1; i++)
1396 rbi = born->bRad[i];
1397 rb[i] = (2 * rbi * rbi * dvda[i])/ONE_4PI_EPS0;
1400 else if (gb_algorithm == egbHCT)
1402 for (i = n0; i < n1; i++)
1404 rbi = born->bRad[i];
1405 rb[i] = rbi * rbi * dvda[i];
1408 else if (gb_algorithm == egbOBC)
1410 for (i = n0; i < n1; i++)
1412 rbi = born->bRad[i];
1413 rb[i] = rbi * rbi * born->drobc[i] * dvda[i];
1417 for (i = 0; i < nl->nri; i++)
1421 nj0 = nl->jindex[i];
1422 nj1 = nl->jindex[i+1];
1424 /* Load shifts for this list */
1425 shift = nl->shift[i];
1426 shX = shift_vec[shift][0];
1427 shY = shift_vec[shift][1];
1428 shZ = shift_vec[shift][2];
1430 /* Load atom i coordinates, add shift vectors */
1431 ix1 = shX + x[ai][0];
1432 iy1 = shY + x[ai][1];
1433 iz1 = shZ + x[ai][2];
1441 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
1455 fgb = rbai*dadx[n++];
1456 fgb_ai = rbaj*dadx[n++];
1458 /* Total force between ai and aj is the sum of ai->aj and aj->ai */
1469 /* Update force on atom aj */
1470 t[aj][0] = t[aj][0] - tx;
1471 t[aj][1] = t[aj][1] - ty;
1472 t[aj][2] = t[aj][2] - tz;
1475 /* Update force and shift forces on atom ai */
1476 t[ai][0] = t[ai][0] + fix1;
1477 t[ai][1] = t[ai][1] + fiy1;
1478 t[ai][2] = t[ai][2] + fiz1;
1480 fshift[shift][0] = fshift[shift][0] + fix1;
1481 fshift[shift][1] = fshift[shift][1] + fiy1;
1482 fshift[shift][2] = fshift[shift][2] + fiz1;
1491 calc_gb_forces(t_commrec *cr, t_mdatoms *md, gmx_genborn_t *born, gmx_localtop_t *top,
1492 rvec x[], rvec f[], t_forcerec *fr, t_idef *idef, int gb_algorithm, int sa_algorithm, t_nrnb *nrnb,
1493 const t_pbc *pbc, const t_graph *graph, gmx_enerdata_t *enerd)
1500 const t_pbc *pbc_null;
1511 if (sa_algorithm == esaAPPROX)
1513 /* Do a simple ACE type approximation for the non-polar solvation */
1514 enerd->term[F_NPSOLVATION] += calc_gb_nonpolar(cr, fr, born->nr, born, top, fr->dvda, md);
1517 /* Calculate the bonded GB-interactions using either table or analytical formula */
1518 enerd->term[F_GBPOL] += gb_bonds_tab(x, f, fr->fshift, md->chargeA, &(fr->gbtabscale),
1519 fr->invsqrta, fr->dvda, fr->gbtab.data, idef, born->epsilon_r, born->gb_epsilon_solvent, fr->epsfac, pbc_null, graph);
1521 /* Calculate self corrections to the GB energies - currently only A state used! (FIXME) */
1522 enerd->term[F_GBPOL] += calc_gb_selfcorrections(cr, born->nr, md->chargeA, born, fr->dvda, fr->epsfac);
1524 /* If parallel, sum the derivative of the potential w.r.t the born radii */
1525 if (DOMAINDECOMP(cr))
1527 dd_atom_sum_real(cr->dd, fr->dvda);
1528 dd_atom_spread_real(cr->dd, fr->dvda);
1533 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1534 if (fr->use_simd_kernels)
1537 genborn_allvsall_calc_chainrule_sse2_double(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1539 genborn_allvsall_calc_chainrule_sse2_single(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1544 genborn_allvsall_calc_chainrule(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1547 genborn_allvsall_calc_chainrule(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1549 cnt = md->homenr*(md->nr/2+1);
1550 /* 9 flops for outer loop, 15 for inner */
1551 inc_nrnb(nrnb, eNR_BORN_AVA_CHAINRULE, md->homenr*9+cnt*15);
1555 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1556 if (fr->use_simd_kernels)
1559 calc_gb_chainrule_sse2_double(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda, x[0],
1560 f[0], fr->fshift[0], fr->shift_vec[0], gb_algorithm, born, md);
1562 calc_gb_chainrule_sse2_single(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda, x[0],
1563 f[0], fr->fshift[0], fr->shift_vec[0], gb_algorithm, born, md);
1568 calc_gb_chainrule(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda,
1569 x, f, fr->fshift, fr->shift_vec, gb_algorithm, born, md);
1572 calc_gb_chainrule(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda,
1573 x, f, fr->fshift, fr->shift_vec, gb_algorithm, born);
1578 /* 9 flops for outer loop, 15 for inner */
1579 inc_nrnb(nrnb, eNR_BORN_CHAINRULE, fr->gblist.nri*9+fr->gblist.nrj*15);
1583 static void add_j_to_gblist(gbtmpnbl_t *list, int aj)
1585 if (list->naj >= list->aj_nalloc)
1587 list->aj_nalloc = over_alloc_large(list->naj+1);
1588 srenew(list->aj, list->aj_nalloc);
1591 list->aj[list->naj++] = aj;
1594 static gbtmpnbl_t *find_gbtmplist(struct gbtmpnbls *lists, int shift)
1598 /* Search the list with the same shift, if there is one */
1600 while (ind < lists->nlist && shift != lists->list[ind].shift)
1604 if (ind == lists->nlist)
1606 if (lists->nlist == lists->list_nalloc)
1608 lists->list_nalloc++;
1609 srenew(lists->list, lists->list_nalloc);
1610 for (i = lists->nlist; i < lists->list_nalloc; i++)
1612 lists->list[i].aj = NULL;
1613 lists->list[i].aj_nalloc = 0;
1618 lists->list[lists->nlist].shift = shift;
1619 lists->list[lists->nlist].naj = 0;
1623 return &lists->list[ind];
1626 static void add_bondeds_to_gblist(t_ilist *il,
1627 gmx_bool bMolPBC, t_pbc *pbc, t_graph *g, rvec *x,
1628 struct gbtmpnbls *nls)
1630 int ind, j, ai, aj, shift, found;
1636 for (ind = 0; ind < il->nr; ind += 3)
1638 ai = il->iatoms[ind+1];
1639 aj = il->iatoms[ind+2];
1644 rvec_sub(x[ai], x[aj], dx);
1645 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
1646 shift = IVEC2IS(dt);
1650 shift = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1653 /* Find the list for this shift or create one */
1654 list = find_gbtmplist(&nls[ai], shift);
1658 /* So that we do not add the same bond twice.
1659 * This happens with some constraints between 1-3 atoms
1660 * that are in the bond-list but should not be in the GB nb-list */
1661 for (j = 0; j < list->naj; j++)
1663 if (list->aj[j] == aj)
1673 gmx_incons("ai == aj");
1676 add_j_to_gblist(list, aj);
1682 compare_int (const void * a, const void * b)
1684 return ( *(int*)a - *(int*)b );
1689 int make_gb_nblist(t_commrec *cr, int gb_algorithm,
1690 rvec x[], matrix box,
1691 t_forcerec *fr, t_idef *idef, t_graph *graph, gmx_genborn_t *born)
1693 int i, l, ii, j, k, n, nj0, nj1, ai, aj, at0, at1, found, shift, s;
1698 struct gbtmpnbls *nls;
1699 gbtmpnbl_t *list = NULL;
1701 set_pbc(&pbc, fr->ePBC, box);
1702 nls = born->nblist_work;
1704 for (i = 0; i < born->nr; i++)
1711 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
1714 switch (gb_algorithm)
1718 /* Loop over 1-2, 1-3 and 1-4 interactions */
1719 for (j = F_GB12; j <= F_GB14; j++)
1721 add_bondeds_to_gblist(&idef->il[j], fr->bMolPBC, &pbc, graph, x, nls);
1725 /* Loop over 1-4 interactions */
1726 add_bondeds_to_gblist(&idef->il[F_GB14], fr->bMolPBC, &pbc, graph, x, nls);
1729 gmx_incons("Unknown GB algorithm");
1732 /* Loop over the VDWQQ and VDW nblists to set up the nonbonded part of the GB list */
1733 for (n = 0; (n < fr->nnblists); n++)
1735 for (i = 0; (i < eNL_NR); i++)
1737 nblist = &(fr->nblists[n].nlist_sr[i]);
1739 if (nblist->nri > 0 && (i == eNL_VDWQQ || i == eNL_QQ))
1741 for (j = 0; j < nblist->nri; j++)
1743 ai = nblist->iinr[j];
1744 shift = nblist->shift[j];
1746 /* Find the list for this shift or create one */
1747 list = find_gbtmplist(&nls[ai], shift);
1749 nj0 = nblist->jindex[j];
1750 nj1 = nblist->jindex[j+1];
1752 /* Add all the j-atoms in the non-bonded list to the GB list */
1753 for (k = nj0; k < nj1; k++)
1755 add_j_to_gblist(list, nblist->jjnr[k]);
1762 /* Zero out some counters */
1766 fr->gblist.jindex[0] = fr->gblist.nri;
1768 for (i = 0; i < fr->natoms_force; i++)
1770 for (s = 0; s < nls[i].nlist; s++)
1772 list = &nls[i].list[s];
1774 /* Only add those atoms that actually have neighbours */
1775 if (born->use[i] != 0)
1777 fr->gblist.iinr[fr->gblist.nri] = i;
1778 fr->gblist.shift[fr->gblist.nri] = list->shift;
1781 for (k = 0; k < list->naj; k++)
1783 /* Memory allocation for jjnr */
1784 if (fr->gblist.nrj >= fr->gblist.maxnrj)
1786 fr->gblist.maxnrj += over_alloc_large(fr->gblist.maxnrj);
1790 fprintf(debug, "Increasing GB neighbourlist j size to %d\n", fr->gblist.maxnrj);
1793 srenew(fr->gblist.jjnr, fr->gblist.maxnrj);
1797 if (i == list->aj[k])
1799 gmx_incons("i == list->aj[k]");
1801 fr->gblist.jjnr[fr->gblist.nrj++] = list->aj[k];
1804 fr->gblist.jindex[fr->gblist.nri] = fr->gblist.nrj;
1811 for (i = 0; i < fr->gblist.nri; i++)
1813 nj0 = fr->gblist.jindex[i];
1814 nj1 = fr->gblist.jindex[i+1];
1815 ai = fr->gblist.iinr[i];
1818 for (j = nj0; j < nj1; j++)
1820 if (fr->gblist.jjnr[j] < ai)
1822 fr->gblist.jjnr[j] += fr->natoms_force;
1825 qsort(fr->gblist.jjnr+nj0, nj1-nj0, sizeof(int), compare_int);
1827 for (j = nj0; j < nj1; j++)
1829 if (fr->gblist.jjnr[j] >= fr->natoms_force)
1831 fr->gblist.jjnr[j] -= fr->natoms_force;
1841 void make_local_gb(const t_commrec *cr, gmx_genborn_t *born, int gb_algorithm)
1844 gmx_domdec_t *dd = NULL;
1846 if (DOMAINDECOMP(cr))
1854 /* Single node, just copy pointers and return */
1855 if (gb_algorithm == egbSTILL)
1857 born->gpol = born->gpol_globalindex;
1858 born->vsolv = born->vsolv_globalindex;
1859 born->gb_radius = born->gb_radius_globalindex;
1863 born->param = born->param_globalindex;
1864 born->gb_radius = born->gb_radius_globalindex;
1867 born->use = born->use_globalindex;
1872 /* Reallocation of local arrays if necessary */
1873 /* fr->natoms_force is equal to dd->nat_tot */
1874 if (DOMAINDECOMP(cr) && dd->nat_tot > born->nalloc)
1878 nalloc = dd->nat_tot;
1880 /* Arrays specific to different gb algorithms */
1881 if (gb_algorithm == egbSTILL)
1883 srenew(born->gpol, nalloc+3);
1884 srenew(born->vsolv, nalloc+3);
1885 srenew(born->gb_radius, nalloc+3);
1886 for (i = born->nalloc; (i < nalloc+3); i++)
1890 born->gb_radius[i] = 0;
1895 srenew(born->param, nalloc+3);
1896 srenew(born->gb_radius, nalloc+3);
1897 for (i = born->nalloc; (i < nalloc+3); i++)
1900 born->gb_radius[i] = 0;
1904 /* All gb-algorithms use the array for vsites exclusions */
1905 srenew(born->use, nalloc+3);
1906 for (i = born->nalloc; (i < nalloc+3); i++)
1911 born->nalloc = nalloc;
1914 /* With dd, copy algorithm specific arrays */
1915 if (gb_algorithm == egbSTILL)
1917 for (i = at0; i < at1; i++)
1919 born->gpol[i] = born->gpol_globalindex[dd->gatindex[i]];
1920 born->vsolv[i] = born->vsolv_globalindex[dd->gatindex[i]];
1921 born->gb_radius[i] = born->gb_radius_globalindex[dd->gatindex[i]];
1922 born->use[i] = born->use_globalindex[dd->gatindex[i]];
1927 for (i = at0; i < at1; i++)
1929 born->param[i] = born->param_globalindex[dd->gatindex[i]];
1930 born->gb_radius[i] = born->gb_radius_globalindex[dd->gatindex[i]];
1931 born->use[i] = born->use_globalindex[dd->gatindex[i]];