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.
49 #include "gromacs/fileio/pdbio.h"
55 #include "gmx_fatal.h"
56 #include "mtop_util.h"
61 #include "gromacs/utility/gmxmpi.h"
63 #ifdef GMX_SIMD_X86_SSE2_OR_HIGHER
65 # include "genborn_sse2_double.h"
66 # include "genborn_allvsall_sse2_double.h"
68 # include "genborn_sse2_single.h"
69 # include "genborn_allvsall_sse2_single.h"
70 # endif /* GMX_DOUBLE */
71 #endif /* SSE or AVX present */
73 #include "genborn_allvsall.h"
75 /*#define DISABLE_SSE*/
84 typedef struct gbtmpnbls {
90 /* This function is exactly the same as the one in bondfree.c. The reason
91 * it is copied here is that the bonded gb-interactions are evaluated
92 * not in calc_bonds, but rather in calc_gb_forces
94 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
98 return pbc_dx_aiuc(pbc, xi, xj, dx);
102 rvec_sub(xi, xj, dx);
107 int init_gb_nblist(int natoms, t_nblist *nl)
109 nl->maxnri = natoms*4;
119 /*nl->nltype = nltype;*/
121 srenew(nl->iinr, nl->maxnri);
122 srenew(nl->gid, nl->maxnri);
123 srenew(nl->shift, nl->maxnri);
124 srenew(nl->jindex, nl->maxnri+1);
131 void gb_pd_send(t_commrec gmx_unused *cr, real gmx_unused *send_data, int gmx_unused nr)
135 int *index, *sendc, *disp;
137 snew(sendc, cr->nnodes);
138 snew(disp, cr->nnodes);
140 index = pd_index(cr);
143 /* Setup count/index arrays */
144 for (i = 0; i < cr->nnodes; i++)
146 sendc[i] = index[i+1]-index[i];
150 /* Do communication */
151 MPI_Gatherv(send_data+index[cur], sendc[cur], GMX_MPI_REAL, send_data, sendc,
152 disp, GMX_MPI_REAL, 0, cr->mpi_comm_mygroup);
153 MPI_Bcast(send_data, nr, GMX_MPI_REAL, 0, cr->mpi_comm_mygroup);
159 int init_gb_still(const t_commrec *cr,
160 const t_atomtypes *atype, t_idef *idef, t_atoms *atoms,
161 gmx_genborn_t *born, int natoms)
164 int i, j, i1, i2, k, m, nbond, nang, ia, ib, ic, id, nb, idx, idx2, at;
168 real r, ri, rj, ri2, ri3, rj2, r2, r3, r4, rk, ratio, term, h, doffset;
169 real p1, p2, p3, factor, cosine, rab, rbc;
176 snew(born->gpol_still_work, natoms+3);
182 pd_at_range(cr, &at0, &at1);
184 for (i = 0; i < natoms; i++)
201 doffset = born->gb_doffset;
203 for (i = 0; i < natoms; i++)
205 born->gpol_globalindex[i] = born->vsolv_globalindex[i] =
206 born->gb_radius_globalindex[i] = 0;
209 /* Compute atomic solvation volumes for Still method */
210 for (i = 0; i < natoms; i++)
212 ri = atype->gb_radius[atoms->atom[i].type];
213 born->gb_radius_globalindex[i] = ri;
215 born->vsolv_globalindex[i] = (4*M_PI/3)*r3;
218 for (j = 0; j < idef->il[F_GB12].nr; j += 3)
220 m = idef->il[F_GB12].iatoms[j];
221 ia = idef->il[F_GB12].iatoms[j+1];
222 ib = idef->il[F_GB12].iatoms[j+2];
224 r = 1.01*idef->iparams[m].gb.st;
226 ri = atype->gb_radius[atoms->atom[ia].type];
227 rj = atype->gb_radius[atoms->atom[ib].type];
233 ratio = (rj2-ri2-r*r)/(2*ri*r);
235 term = (M_PI/3.0)*h*h*(3.0*ri-h);
243 born->vsolv_globalindex[ia] -= term;
246 ratio = (ri2-rj2-r*r)/(2*rj*r);
248 term = (M_PI/3.0)*h*h*(3.0*rj-h);
256 born->vsolv_globalindex[ib] -= term;
262 gmx_sum(natoms, vsol, cr);
264 for (i = 0; i < natoms; i++)
266 born->vsolv_globalindex[i] = born->vsolv_globalindex[i]-vsol[i];
270 /* Get the self-, 1-2 and 1-3 polarization energies for analytical Still
273 for (j = 0; j < natoms; j++)
275 if (born->use_globalindex[j] == 1)
277 born->gpol_globalindex[j] = -0.5*ONE_4PI_EPS0/
278 (atype->gb_radius[atoms->atom[j].type]-doffset+STILL_P1);
283 for (j = 0; j < idef->il[F_GB12].nr; j += 3)
285 m = idef->il[F_GB12].iatoms[j];
286 ia = idef->il[F_GB12].iatoms[j+1];
287 ib = idef->il[F_GB12].iatoms[j+2];
289 r = idef->iparams[m].gb.st;
295 gp[ia] += STILL_P2*born->vsolv_globalindex[ib]/r4;
296 gp[ib] += STILL_P2*born->vsolv_globalindex[ia]/r4;
300 born->gpol_globalindex[ia] = born->gpol_globalindex[ia]+
301 STILL_P2*born->vsolv_globalindex[ib]/r4;
302 born->gpol_globalindex[ib] = born->gpol_globalindex[ib]+
303 STILL_P2*born->vsolv_globalindex[ia]/r4;
308 for (j = 0; j < idef->il[F_GB13].nr; j += 3)
310 m = idef->il[F_GB13].iatoms[j];
311 ia = idef->il[F_GB13].iatoms[j+1];
312 ib = idef->il[F_GB13].iatoms[j+2];
314 r = idef->iparams[m].gb.st;
319 gp[ia] += STILL_P3*born->vsolv[ib]/r4;
320 gp[ib] += STILL_P3*born->vsolv[ia]/r4;
324 born->gpol_globalindex[ia] = born->gpol_globalindex[ia]+
325 STILL_P3*born->vsolv_globalindex[ib]/r4;
326 born->gpol_globalindex[ib] = born->gpol_globalindex[ib]+
327 STILL_P3*born->vsolv_globalindex[ia]/r4;
333 gmx_sum(natoms, gp, cr);
335 for (i = 0; i < natoms; i++)
337 born->gpol_globalindex[i] = born->gpol_globalindex[i]+gp[i];
347 /* Initialize all GB datastructs and compute polarization energies */
348 int init_gb(gmx_genborn_t **p_born,
349 const t_commrec *cr, t_forcerec *fr, const t_inputrec *ir,
350 const gmx_mtop_t *mtop, int gb_algorithm)
352 int i, j, m, ai, aj, jj, natoms, nalloc;
353 real rai, sk, p, doffset;
357 gmx_localtop_t *localtop;
359 natoms = mtop->natoms;
361 atoms = gmx_mtop_global_atoms(mtop);
362 localtop = gmx_mtop_generate_local_top(mtop, ir);
369 snew(born->drobc, natoms);
370 snew(born->bRad, natoms);
372 /* Allocate memory for the global data arrays */
373 snew(born->param_globalindex, natoms+3);
374 snew(born->gpol_globalindex, natoms+3);
375 snew(born->vsolv_globalindex, natoms+3);
376 snew(born->gb_radius_globalindex, natoms+3);
377 snew(born->use_globalindex, natoms+3);
379 snew(fr->invsqrta, natoms);
380 snew(fr->dvda, natoms);
383 fr->dadx_rawptr = NULL;
385 born->gpol_still_work = NULL;
386 born->gpol_hct_work = NULL;
388 /* snew(born->asurf,natoms); */
389 /* snew(born->dasurf,natoms); */
391 /* Initialize the gb neighbourlist */
392 init_gb_nblist(natoms, &(fr->gblist));
394 /* Do the Vsites exclusions (if any) */
395 for (i = 0; i < natoms; i++)
397 jj = atoms.atom[i].type;
398 if (mtop->atomtypes.gb_radius[atoms.atom[i].type] > 0)
400 born->use_globalindex[i] = 1;
404 born->use_globalindex[i] = 0;
407 /* If we have a Vsite, put vs_globalindex[i]=0 */
408 if (C6 (fr->nbfp, fr->ntype, jj, jj) == 0 &&
409 C12(fr->nbfp, fr->ntype, jj, jj) == 0 &&
410 atoms.atom[i].q == 0)
412 born->use_globalindex[i] = 0;
416 /* Copy algorithm parameters from inputrecord to local structure */
417 born->obc_alpha = ir->gb_obc_alpha;
418 born->obc_beta = ir->gb_obc_beta;
419 born->obc_gamma = ir->gb_obc_gamma;
420 born->gb_doffset = ir->gb_dielectric_offset;
421 born->gb_epsilon_solvent = ir->gb_epsilon_solvent;
422 born->epsilon_r = ir->epsilon_r;
424 doffset = born->gb_doffset;
426 /* Set the surface tension */
427 born->sa_surface_tension = ir->sa_surface_tension;
429 /* If Still model, initialise the polarisation energies */
430 if (gb_algorithm == egbSTILL)
432 init_gb_still(cr, &(mtop->atomtypes), &(localtop->idef), &atoms,
437 /* If HCT/OBC, precalculate the sk*atype->S_hct factors */
438 else if (gb_algorithm == egbHCT || gb_algorithm == egbOBC)
441 snew(born->gpol_hct_work, natoms+3);
443 for (i = 0; i < natoms; i++)
445 if (born->use_globalindex[i] == 1)
447 rai = mtop->atomtypes.gb_radius[atoms.atom[i].type]-doffset;
448 sk = rai * mtop->atomtypes.S_hct[atoms.atom[i].type];
449 born->param_globalindex[i] = sk;
450 born->gb_radius_globalindex[i] = rai;
454 born->param_globalindex[i] = 0;
455 born->gb_radius_globalindex[i] = 0;
460 /* Allocate memory for work arrays for temporary use */
461 snew(born->work, natoms+4);
462 snew(born->count, natoms);
463 snew(born->nblist_work, natoms);
465 /* Domain decomposition specific stuff */
474 calc_gb_rad_still(t_commrec *cr, t_forcerec *fr, int natoms, gmx_localtop_t *top,
475 rvec x[], t_nblist *nl,
476 gmx_genborn_t *born, t_mdatoms *md)
478 int i, k, n, nj0, nj1, ai, aj, type;
481 real gpi, dr, dr2, dr4, idr4, rvdw, ratio, ccf, theta, term, rai, raj;
482 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
483 real rinv, idr2, idr6, vaj, dccf, cosq, sinq, prod, gpi2;
485 real vai, prod_ai, icf4, icf6;
487 factor = 0.5*ONE_4PI_EPS0;
490 for (i = 0; i < born->nr; i++)
492 born->gpol_still_work[i] = 0;
495 for (i = 0; i < nl->nri; i++)
500 nj1 = nl->jindex[i+1];
502 /* Load shifts for this list */
503 shift = nl->shift[i];
504 shX = fr->shift_vec[shift][0];
505 shY = fr->shift_vec[shift][1];
506 shZ = fr->shift_vec[shift][2];
510 rai = top->atomtypes.gb_radius[md->typeA[ai]];
511 vai = born->vsolv[ai];
512 prod_ai = STILL_P4*vai;
514 /* Load atom i coordinates, add shift vectors */
515 ix1 = shX + x[ai][0];
516 iy1 = shY + x[ai][1];
517 iz1 = shZ + x[ai][2];
519 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
530 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
531 rinv = gmx_invsqrt(dr2);
536 raj = top->atomtypes.gb_radius[md->typeA[aj]];
540 ratio = dr2 / (rvdw * rvdw);
541 vaj = born->vsolv[aj];
543 if (ratio > STILL_P5INV)
550 theta = ratio*STILL_PIP5;
552 term = 0.5*(1.0-cosq);
554 sinq = 1.0 - cosq*cosq;
555 dccf = 2.0*term*sinq*gmx_invsqrt(sinq)*theta;
560 icf6 = (4*ccf-dccf)*idr6;
561 born->gpol_still_work[aj] += prod_ai*icf4;
564 /* Save ai->aj and aj->ai chain rule terms */
565 fr->dadx[n++] = prod*icf6;
566 fr->dadx[n++] = prod_ai*icf6;
568 born->gpol_still_work[ai] += gpi;
571 /* Parallel summations */
574 gmx_sum(natoms, born->gpol_still_work, cr);
576 else if (DOMAINDECOMP(cr))
578 dd_atom_sum_real(cr->dd, born->gpol_still_work);
581 /* Calculate the radii */
582 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
584 if (born->use[i] != 0)
586 gpi = born->gpol[i]+born->gpol_still_work[i];
588 born->bRad[i] = factor*gmx_invsqrt(gpi2);
589 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
593 /* Extra communication required for DD */
594 if (DOMAINDECOMP(cr))
596 dd_atom_spread_real(cr->dd, born->bRad);
597 dd_atom_spread_real(cr->dd, fr->invsqrta);
606 calc_gb_rad_hct(t_commrec *cr, t_forcerec *fr, int natoms, gmx_localtop_t *top,
607 rvec x[], t_nblist *nl,
608 gmx_genborn_t *born, t_mdatoms *md)
610 int i, k, n, ai, aj, nj0, nj1, at0, at1;
613 real rai, raj, gpi, dr2, dr, sk, sk_ai, sk2, sk2_ai, lij, uij, diff2, tmp, sum_ai;
614 real rad, min_rad, rinv, rai_inv;
615 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
616 real lij2, uij2, lij3, uij3, t1, t2, t3;
617 real lij_inv, dlij, duij, sk2_rinv, prod, log_term;
618 real doffset, raj_inv, dadx_val;
621 doffset = born->gb_doffset;
622 gb_radius = born->gb_radius;
624 for (i = 0; i < born->nr; i++)
626 born->gpol_hct_work[i] = 0;
629 /* Keep the compiler happy */
633 for (i = 0; i < nl->nri; i++)
638 nj1 = nl->jindex[i+1];
640 /* Load shifts for this list */
641 shift = nl->shift[i];
642 shX = fr->shift_vec[shift][0];
643 shY = fr->shift_vec[shift][1];
644 shZ = fr->shift_vec[shift][2];
649 sk_ai = born->param[ai];
650 sk2_ai = sk_ai*sk_ai;
652 /* Load atom i coordinates, add shift vectors */
653 ix1 = shX + x[ai][0];
654 iy1 = shY + x[ai][1];
655 iz1 = shZ + x[ai][2];
659 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
671 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
672 rinv = gmx_invsqrt(dr2);
675 sk = born->param[aj];
678 /* aj -> ai interaction */
699 lij_inv = gmx_invsqrt(lij2);
702 prod = 0.25*sk2_rinv;
704 log_term = log(uij*lij_inv);
706 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term +
711 tmp = tmp + 2.0 * (rai_inv-lij);
714 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
715 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
716 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
718 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
719 /* fr->dadx[n++] = (dlij*t1+duij*t2+t3)*rinv; */
720 /* rb2 is moved to chainrule */
728 fr->dadx[n++] = dadx_val;
731 /* ai -> aj interaction */
732 if (raj < dr + sk_ai)
734 lij = 1.0/(dr-sk_ai);
747 uij = 1.0/(dr+sk_ai);
753 lij_inv = gmx_invsqrt(lij2);
754 sk2 = sk2_ai; /* sk2_ai = sk_ai * sk_ai in i loop above */
756 prod = 0.25 * sk2_rinv;
758 /* log_term = table_log(uij*lij_inv,born->log_table,
759 LOG_TABLE_ACCURACY); */
760 log_term = log(uij*lij_inv);
762 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term +
767 tmp = tmp + 2.0 * (raj_inv-lij);
771 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
772 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
773 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
775 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
776 /* fr->dadx[n++] = (dlij*t1+duij*t2+t3)*rinv; */ /* rb2 is moved to chainrule */
778 born->gpol_hct_work[aj] += 0.5*tmp;
784 fr->dadx[n++] = dadx_val;
787 born->gpol_hct_work[ai] += sum_ai;
790 /* Parallel summations */
793 gmx_sum(natoms, born->gpol_hct_work, cr);
795 else if (DOMAINDECOMP(cr))
797 dd_atom_sum_real(cr->dd, born->gpol_hct_work);
800 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
802 if (born->use[i] != 0)
804 rai = top->atomtypes.gb_radius[md->typeA[i]]-doffset;
805 sum_ai = 1.0/rai - born->gpol_hct_work[i];
806 min_rad = rai + doffset;
809 born->bRad[i] = rad > min_rad ? rad : min_rad;
810 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
814 /* Extra communication required for DD */
815 if (DOMAINDECOMP(cr))
817 dd_atom_spread_real(cr->dd, born->bRad);
818 dd_atom_spread_real(cr->dd, fr->invsqrta);
826 calc_gb_rad_obc(t_commrec *cr, t_forcerec *fr, int natoms, gmx_localtop_t *top,
827 rvec x[], t_nblist *nl, gmx_genborn_t *born, t_mdatoms *md)
829 int i, k, ai, aj, nj0, nj1, n, at0, at1;
832 real rai, raj, gpi, dr2, dr, sk, sk2, lij, uij, diff2, tmp, sum_ai;
833 real rad, min_rad, sum_ai2, sum_ai3, tsum, tchain, rinv, rai_inv, lij_inv, rai_inv2;
834 real log_term, prod, sk2_rinv, sk_ai, sk2_ai;
835 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11;
836 real lij2, uij2, lij3, uij3, dlij, duij, t1, t2, t3;
837 real doffset, raj_inv, dadx_val;
840 /* Keep the compiler happy */
845 doffset = born->gb_doffset;
846 gb_radius = born->gb_radius;
848 for (i = 0; i < born->nr; i++)
850 born->gpol_hct_work[i] = 0;
853 for (i = 0; i < nl->nri; i++)
858 nj1 = nl->jindex[i+1];
860 /* Load shifts for this list */
861 shift = nl->shift[i];
862 shX = fr->shift_vec[shift][0];
863 shY = fr->shift_vec[shift][1];
864 shZ = fr->shift_vec[shift][2];
869 sk_ai = born->param[ai];
870 sk2_ai = sk_ai*sk_ai;
872 /* Load atom i coordinates, add shift vectors */
873 ix1 = shX + x[ai][0];
874 iy1 = shY + x[ai][1];
875 iz1 = shZ + x[ai][2];
879 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
891 dr2 = dx11*dx11+dy11*dy11+dz11*dz11;
892 rinv = gmx_invsqrt(dr2);
895 /* sk is precalculated in init_gb() */
896 sk = born->param[aj];
899 /* aj -> ai interaction */
919 lij_inv = gmx_invsqrt(lij2);
922 prod = 0.25*sk2_rinv;
924 log_term = log(uij*lij_inv);
926 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term + prod*(-diff2);
930 tmp = tmp + 2.0 * (rai_inv-lij);
934 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
935 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
936 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
938 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
946 fr->dadx[n++] = dadx_val;
948 /* ai -> aj interaction */
949 if (raj < dr + sk_ai)
951 lij = 1.0/(dr-sk_ai);
964 uij = 1.0/(dr+sk_ai);
970 lij_inv = gmx_invsqrt(lij2);
971 sk2 = sk2_ai; /* sk2_ai = sk_ai * sk_ai in i loop above */
973 prod = 0.25 * sk2_rinv;
975 /* log_term = table_log(uij*lij_inv,born->log_table,LOG_TABLE_ACCURACY); */
976 log_term = log(uij*lij_inv);
978 tmp = lij-uij + 0.25*dr*diff2 + (0.5*rinv)*log_term + prod*(-diff2);
982 tmp = tmp + 2.0 * (raj_inv-lij);
985 t1 = 0.5*lij2 + prod*lij3 - 0.25*(lij*rinv+lij3*dr);
986 t2 = -0.5*uij2 - 0.25*sk2_rinv*uij3 + 0.25*(uij*rinv+uij3*dr);
987 t3 = 0.125*(1.0+sk2_rinv*rinv)*(-diff2)+0.25*log_term*rinv*rinv;
989 dadx_val = (dlij*t1+t2+t3)*rinv; /* rb2 is moved to chainrule */
991 born->gpol_hct_work[aj] += 0.5*tmp;
998 fr->dadx[n++] = dadx_val;
1001 born->gpol_hct_work[ai] += sum_ai;
1005 /* Parallel summations */
1008 gmx_sum(natoms, born->gpol_hct_work, cr);
1010 else if (DOMAINDECOMP(cr))
1012 dd_atom_sum_real(cr->dd, born->gpol_hct_work);
1015 for (i = 0; i < fr->natoms_force; i++) /* PELA born->nr */
1017 if (born->use[i] != 0)
1019 rai = top->atomtypes.gb_radius[md->typeA[i]];
1023 sum_ai = rai * born->gpol_hct_work[i];
1024 sum_ai2 = sum_ai * sum_ai;
1025 sum_ai3 = sum_ai2 * sum_ai;
1027 tsum = tanh(born->obc_alpha*sum_ai-born->obc_beta*sum_ai2+born->obc_gamma*sum_ai3);
1028 born->bRad[i] = rai_inv - tsum*rai_inv2;
1029 born->bRad[i] = 1.0 / born->bRad[i];
1031 fr->invsqrta[i] = gmx_invsqrt(born->bRad[i]);
1033 tchain = rai * (born->obc_alpha-2*born->obc_beta*sum_ai+3*born->obc_gamma*sum_ai2);
1034 born->drobc[i] = (1.0-tsum*tsum)*tchain*rai_inv2;
1038 /* Extra (local) communication required for DD */
1039 if (DOMAINDECOMP(cr))
1041 dd_atom_spread_real(cr->dd, born->bRad);
1042 dd_atom_spread_real(cr->dd, fr->invsqrta);
1043 dd_atom_spread_real(cr->dd, born->drobc);
1052 int calc_gb_rad(t_commrec *cr, t_forcerec *fr, t_inputrec *ir, gmx_localtop_t *top,
1053 rvec x[], t_nblist *nl, gmx_genborn_t *born, t_mdatoms *md, t_nrnb *nrnb)
1059 if (fr->bAllvsAll && fr->dadx == NULL)
1061 /* We might need up to 8 atoms of padding before and after,
1062 * and another 4 units to guarantee SSE alignment.
1064 fr->nalloc_dadx = 2*(md->homenr+12)*(md->nr/2+1+12);
1065 snew(fr->dadx_rawptr, fr->nalloc_dadx);
1066 fr->dadx = (real *) (((size_t) fr->dadx_rawptr + 16) & (~((size_t) 15)));
1070 /* In the SSE-enabled gb-loops, when writing to dadx, we
1071 * always write 2*4 elements at a time, even in the case with only
1072 * 1-3 j particles, where we only really need to write 2*(1-3)
1073 * elements. This is because we want dadx to be aligned to a 16-
1074 * byte boundary, and being able to use _mm_store/load_ps
1076 ndadx = 2 * (nl->nrj + 3*nl->nri);
1078 /* First, reallocate the dadx array, we need 3 extra for SSE */
1079 if (ndadx + 3 > fr->nalloc_dadx)
1081 fr->nalloc_dadx = over_alloc_large(ndadx) + 3;
1082 srenew(fr->dadx_rawptr, fr->nalloc_dadx);
1083 fr->dadx = (real *) (((size_t) fr->dadx_rawptr + 16) & (~((size_t) 15)));
1089 cnt = md->homenr*(md->nr/2+1);
1091 if (ir->gb_algorithm == egbSTILL)
1093 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1094 if (fr->use_simd_kernels)
1097 genborn_allvsall_calc_still_radii_sse2_double(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
1099 genborn_allvsall_calc_still_radii_sse2_single(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
1104 genborn_allvsall_calc_still_radii(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
1107 genborn_allvsall_calc_still_radii(fr, md, born, top, x[0], cr, &fr->AllvsAll_workgb);
1109 /* 13 flops in outer loop, 47 flops in inner loop */
1110 inc_nrnb(nrnb, eNR_BORN_AVA_RADII_STILL, md->homenr*13+cnt*47);
1112 else if (ir->gb_algorithm == egbHCT || ir->gb_algorithm == egbOBC)
1114 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1115 if (fr->use_simd_kernels)
1118 genborn_allvsall_calc_hct_obc_radii_sse2_double(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1120 genborn_allvsall_calc_hct_obc_radii_sse2_single(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1125 genborn_allvsall_calc_hct_obc_radii(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1128 genborn_allvsall_calc_hct_obc_radii(fr, md, born, ir->gb_algorithm, top, x[0], cr, &fr->AllvsAll_workgb);
1130 /* 24 flops in outer loop, 183 in inner */
1131 inc_nrnb(nrnb, eNR_BORN_AVA_RADII_HCT_OBC, md->homenr*24+cnt*183);
1135 gmx_fatal(FARGS, "Bad gb algorithm for all-vs-all interactions");
1140 /* Switch for determining which algorithm to use for Born radii calculation */
1143 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1144 /* x86 or x86-64 with GCC inline assembly and/or SSE intrinsics */
1145 switch (ir->gb_algorithm)
1148 if (fr->use_simd_kernels)
1150 calc_gb_rad_still_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born);
1154 calc_gb_rad_still(cr, fr, born->nr, top, atype, x, nl, born, md);
1158 if (fr->use_simd_kernels)
1160 calc_gb_rad_hct_obc_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born, md, ir->gb_algorithm);
1164 calc_gb_rad_hct(cr, fr, born->nr, top, atype, x, nl, born, md);
1168 if (fr->use_simd_kernels)
1170 calc_gb_rad_hct_obc_sse2_double(cr, fr, born->nr, top, atype, x[0], nl, born, md, ir->gb_algorithm);
1174 calc_gb_rad_obc(cr, fr, born->nr, top, atype, x, nl, born, md);
1179 gmx_fatal(FARGS, "Unknown double precision sse-enabled algorithm for Born radii calculation: %d", ir->gb_algorithm);
1182 switch (ir->gb_algorithm)
1185 calc_gb_rad_still(cr, fr, born->nr, top, x, nl, born, md);
1188 calc_gb_rad_hct(cr, fr, born->nr, top, x, nl, born, md);
1191 calc_gb_rad_obc(cr, fr, born->nr, top, x, nl, born, md);
1195 gmx_fatal(FARGS, "Unknown double precision algorithm for Born radii calculation: %d", ir->gb_algorithm);
1202 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1203 /* x86 or x86-64 with GCC inline assembly and/or SSE intrinsics */
1204 switch (ir->gb_algorithm)
1207 if (fr->use_simd_kernels)
1209 calc_gb_rad_still_sse2_single(cr, fr, born->nr, top, x[0], nl, born);
1213 calc_gb_rad_still(cr, fr, born->nr, top, x, nl, born, md);
1217 if (fr->use_simd_kernels)
1219 calc_gb_rad_hct_obc_sse2_single(cr, fr, born->nr, top, x[0], nl, born, md, ir->gb_algorithm);
1223 calc_gb_rad_hct(cr, fr, born->nr, top, x, nl, born, md);
1228 if (fr->use_simd_kernels)
1230 calc_gb_rad_hct_obc_sse2_single(cr, fr, born->nr, top, x[0], nl, born, md, ir->gb_algorithm);
1234 calc_gb_rad_obc(cr, fr, born->nr, top, x, nl, born, md);
1239 gmx_fatal(FARGS, "Unknown sse-enabled algorithm for Born radii calculation: %d", ir->gb_algorithm);
1243 switch (ir->gb_algorithm)
1246 calc_gb_rad_still(cr, fr, born->nr, top, x, nl, born, md);
1249 calc_gb_rad_hct(cr, fr, born->nr, top, x, nl, born, md);
1252 calc_gb_rad_obc(cr, fr, born->nr, top, x, nl, born, md);
1256 gmx_fatal(FARGS, "Unknown algorithm for Born radii calculation: %d", ir->gb_algorithm);
1259 #endif /* Single precision sse */
1261 #endif /* Double or single precision */
1263 if (fr->bAllvsAll == FALSE)
1265 switch (ir->gb_algorithm)
1268 /* 17 flops per outer loop iteration, 47 flops per inner loop */
1269 inc_nrnb(nrnb, eNR_BORN_RADII_STILL, nl->nri*17+nl->nrj*47);
1273 /* 61 (assuming 10 for tanh) flops for outer loop iteration, 183 flops per inner loop */
1274 inc_nrnb(nrnb, eNR_BORN_RADII_HCT_OBC, nl->nri*61+nl->nrj*183);
1287 real gb_bonds_tab(rvec x[], rvec f[], rvec fshift[], real *charge, real *p_gbtabscale,
1288 real *invsqrta, real *dvda, real *GBtab, t_idef *idef, real epsilon_r,
1289 real gb_epsilon_solvent, real facel, const t_pbc *pbc, const t_graph *graph)
1291 int i, j, n0, m, nnn, type, ai, aj;
1297 real isaprod, qq, gbscale, gbtabscale, Y, F, Geps, Heps2, Fp, VV, FF, rt, eps, eps2;
1298 real vgb, fgb, vcoul, fijC, dvdatmp, fscal, dvdaj;
1304 t_iatom *forceatoms;
1306 /* Scale the electrostatics by gb_epsilon_solvent */
1307 facel = facel * ((1.0/epsilon_r) - 1.0/gb_epsilon_solvent);
1309 gbtabscale = *p_gbtabscale;
1312 for (j = F_GB12; j <= F_GB14; j++)
1314 forceatoms = idef->il[j].iatoms;
1316 for (i = 0; i < idef->il[j].nr; )
1318 /* To avoid reading in the interaction type, we just increment i to pass over
1319 * the types in the forceatoms array, this saves some memory accesses
1322 ai = forceatoms[i++];
1323 aj = forceatoms[i++];
1325 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx);
1326 rsq11 = iprod(dx, dx);
1328 isai = invsqrta[ai];
1329 iq = (-1)*facel*charge[ai];
1331 rinv11 = gmx_invsqrt(rsq11);
1332 isaj = invsqrta[aj];
1333 isaprod = isai*isaj;
1334 qq = isaprod*iq*charge[aj];
1335 gbscale = isaprod*gbtabscale;
1344 Geps = eps*GBtab[nnn+2];
1345 Heps2 = eps2*GBtab[nnn+3];
1348 FF = Fp+Geps+2.0*Heps2;
1350 fijC = qq*FF*gbscale;
1351 dvdatmp = -(vgb+fijC*r)*0.5;
1352 dvda[aj] = dvda[aj] + dvdatmp*isaj*isaj;
1353 dvda[ai] = dvda[ai] + dvdatmp*isai*isai;
1354 vctot = vctot + vgb;
1355 fgb = -(fijC)*rinv11;
1359 ivec_sub(SHIFT_IVEC(graph, ai), SHIFT_IVEC(graph, aj), dt);
1363 for (m = 0; (m < DIM); m++) /* 15 */
1368 fshift[ki][m] += fscal;
1369 fshift[CENTRAL][m] -= fscal;
1377 real calc_gb_selfcorrections(t_commrec *cr, int natoms,
1378 real *charge, gmx_genborn_t *born, real *dvda, double facel)
1380 int i, ai, at0, at1;
1381 real rai, e, derb, q, q2, fi, rai_inv, vtot;
1385 pd_at_range(cr, &at0, &at1);
1387 else if (DOMAINDECOMP(cr))
1390 at1 = cr->dd->nat_home;
1399 /* Scale the electrostatics by gb_epsilon_solvent */
1400 facel = facel * ((1.0/born->epsilon_r) - 1.0/born->gb_epsilon_solvent);
1404 /* Apply self corrections */
1405 for (i = at0; i < at1; i++)
1409 if (born->use[ai] == 1)
1411 rai = born->bRad[ai];
1417 derb = 0.5*e*rai_inv*rai_inv;
1418 dvda[ai] += derb*rai;
1427 real calc_gb_nonpolar(t_commrec *cr, t_forcerec *fr, int natoms, gmx_genborn_t *born, gmx_localtop_t *top,
1428 real *dvda, t_mdatoms *md)
1430 int ai, i, at0, at1;
1431 real e, es, rai, rbi, term, probe, tmp, factor;
1432 real rbi_inv, rbi_inv2;
1434 /* To keep the compiler happy */
1439 pd_at_range(cr, &at0, &at1);
1441 else if (DOMAINDECOMP(cr))
1444 at1 = cr->dd->nat_home;
1452 /* factor is the surface tension */
1453 factor = born->sa_surface_tension;
1456 // The surface tension factor is 0.0049 for Still model, 0.0054 for HCT/OBC
1457 if(gb_algorithm==egbSTILL)
1459 factor=0.0049*100*CAL2JOULE;
1463 factor=0.0054*100*CAL2JOULE;
1466 /* if(gb_algorithm==egbHCT || gb_algorithm==egbOBC) */
1472 for (i = at0; i < at1; i++)
1476 if (born->use[ai] == 1)
1478 rai = top->atomtypes.gb_radius[md->typeA[ai]];
1479 rbi_inv = fr->invsqrta[ai];
1480 rbi_inv2 = rbi_inv * rbi_inv;
1481 tmp = (rai*rbi_inv2)*(rai*rbi_inv2);
1483 e = factor*term*(rai+probe)*(rai+probe)*tmp;
1484 dvda[ai] = dvda[ai] - 6*e*rbi_inv2;
1494 real calc_gb_chainrule(int natoms, t_nblist *nl, real *dadx, real *dvda, rvec x[], rvec t[], rvec fshift[],
1495 rvec shift_vec[], int gb_algorithm, gmx_genborn_t *born)
1497 int i, k, n, ai, aj, nj0, nj1, n0, n1;
1500 real fgb, fij, rb2, rbi, fix1, fiy1, fiz1;
1501 real ix1, iy1, iz1, jx1, jy1, jz1, dx11, dy11, dz11, rsq11;
1502 real rinv11, tx, ty, tz, rbai, rbaj, fgb_ai;
1512 if (gb_algorithm == egbSTILL)
1514 for (i = n0; i < n1; i++)
1516 rbi = born->bRad[i];
1517 rb[i] = (2 * rbi * rbi * dvda[i])/ONE_4PI_EPS0;
1520 else if (gb_algorithm == egbHCT)
1522 for (i = n0; i < n1; i++)
1524 rbi = born->bRad[i];
1525 rb[i] = rbi * rbi * dvda[i];
1528 else if (gb_algorithm == egbOBC)
1530 for (i = n0; i < n1; i++)
1532 rbi = born->bRad[i];
1533 rb[i] = rbi * rbi * born->drobc[i] * dvda[i];
1537 for (i = 0; i < nl->nri; i++)
1541 nj0 = nl->jindex[i];
1542 nj1 = nl->jindex[i+1];
1544 /* Load shifts for this list */
1545 shift = nl->shift[i];
1546 shX = shift_vec[shift][0];
1547 shY = shift_vec[shift][1];
1548 shZ = shift_vec[shift][2];
1550 /* Load atom i coordinates, add shift vectors */
1551 ix1 = shX + x[ai][0];
1552 iy1 = shY + x[ai][1];
1553 iz1 = shZ + x[ai][2];
1561 for (k = nj0; k < nj1 && nl->jjnr[k] >= 0; k++)
1575 fgb = rbai*dadx[n++];
1576 fgb_ai = rbaj*dadx[n++];
1578 /* Total force between ai and aj is the sum of ai->aj and aj->ai */
1589 /* Update force on atom aj */
1590 t[aj][0] = t[aj][0] - tx;
1591 t[aj][1] = t[aj][1] - ty;
1592 t[aj][2] = t[aj][2] - tz;
1595 /* Update force and shift forces on atom ai */
1596 t[ai][0] = t[ai][0] + fix1;
1597 t[ai][1] = t[ai][1] + fiy1;
1598 t[ai][2] = t[ai][2] + fiz1;
1600 fshift[shift][0] = fshift[shift][0] + fix1;
1601 fshift[shift][1] = fshift[shift][1] + fiy1;
1602 fshift[shift][2] = fshift[shift][2] + fiz1;
1611 calc_gb_forces(t_commrec *cr, t_mdatoms *md, gmx_genborn_t *born, gmx_localtop_t *top,
1612 rvec x[], rvec f[], t_forcerec *fr, t_idef *idef, int gb_algorithm, int sa_algorithm, t_nrnb *nrnb,
1613 const t_pbc *pbc, const t_graph *graph, gmx_enerdata_t *enerd)
1620 const t_pbc *pbc_null;
1631 if (sa_algorithm == esaAPPROX)
1633 /* Do a simple ACE type approximation for the non-polar solvation */
1634 enerd->term[F_NPSOLVATION] += calc_gb_nonpolar(cr, fr, born->nr, born, top, fr->dvda, md);
1637 /* Calculate the bonded GB-interactions using either table or analytical formula */
1638 enerd->term[F_GBPOL] += gb_bonds_tab(x, f, fr->fshift, md->chargeA, &(fr->gbtabscale),
1639 fr->invsqrta, fr->dvda, fr->gbtab.data, idef, born->epsilon_r, born->gb_epsilon_solvent, fr->epsfac, pbc_null, graph);
1641 /* Calculate self corrections to the GB energies - currently only A state used! (FIXME) */
1642 enerd->term[F_GBPOL] += calc_gb_selfcorrections(cr, born->nr, md->chargeA, born, fr->dvda, fr->epsfac);
1644 /* If parallel, sum the derivative of the potential w.r.t the born radii */
1647 gmx_sum(md->nr, fr->dvda, cr);
1649 else if (DOMAINDECOMP(cr))
1651 dd_atom_sum_real(cr->dd, fr->dvda);
1652 dd_atom_spread_real(cr->dd, fr->dvda);
1657 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1658 if (fr->use_simd_kernels)
1661 genborn_allvsall_calc_chainrule_sse2_double(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1663 genborn_allvsall_calc_chainrule_sse2_single(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1668 genborn_allvsall_calc_chainrule(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1671 genborn_allvsall_calc_chainrule(fr, md, born, x[0], f[0], gb_algorithm, fr->AllvsAll_workgb);
1673 cnt = md->homenr*(md->nr/2+1);
1674 /* 9 flops for outer loop, 15 for inner */
1675 inc_nrnb(nrnb, eNR_BORN_AVA_CHAINRULE, md->homenr*9+cnt*15);
1679 #if 0 && defined (GMX_SIMD_X86_SSE2_OR_HIGHER)
1680 if (fr->use_simd_kernels)
1683 calc_gb_chainrule_sse2_double(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda, x[0],
1684 f[0], fr->fshift[0], fr->shift_vec[0], gb_algorithm, born, md);
1686 calc_gb_chainrule_sse2_single(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda, x[0],
1687 f[0], fr->fshift[0], fr->shift_vec[0], gb_algorithm, born, md);
1692 calc_gb_chainrule(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda,
1693 x, f, fr->fshift, fr->shift_vec, gb_algorithm, born, md);
1696 calc_gb_chainrule(fr->natoms_force, &(fr->gblist), fr->dadx, fr->dvda,
1697 x, f, fr->fshift, fr->shift_vec, gb_algorithm, born);
1702 /* 9 flops for outer loop, 15 for inner */
1703 inc_nrnb(nrnb, eNR_BORN_CHAINRULE, fr->gblist.nri*9+fr->gblist.nrj*15);
1707 static void add_j_to_gblist(gbtmpnbl_t *list, int aj)
1709 if (list->naj >= list->aj_nalloc)
1711 list->aj_nalloc = over_alloc_large(list->naj+1);
1712 srenew(list->aj, list->aj_nalloc);
1715 list->aj[list->naj++] = aj;
1718 static gbtmpnbl_t *find_gbtmplist(struct gbtmpnbls *lists, int shift)
1722 /* Search the list with the same shift, if there is one */
1724 while (ind < lists->nlist && shift != lists->list[ind].shift)
1728 if (ind == lists->nlist)
1730 if (lists->nlist == lists->list_nalloc)
1732 lists->list_nalloc++;
1733 srenew(lists->list, lists->list_nalloc);
1734 for (i = lists->nlist; i < lists->list_nalloc; i++)
1736 lists->list[i].aj = NULL;
1737 lists->list[i].aj_nalloc = 0;
1742 lists->list[lists->nlist].shift = shift;
1743 lists->list[lists->nlist].naj = 0;
1747 return &lists->list[ind];
1750 static void add_bondeds_to_gblist(t_ilist *il,
1751 gmx_bool bMolPBC, t_pbc *pbc, t_graph *g, rvec *x,
1752 struct gbtmpnbls *nls)
1754 int ind, j, ai, aj, shift, found;
1760 for (ind = 0; ind < il->nr; ind += 3)
1762 ai = il->iatoms[ind+1];
1763 aj = il->iatoms[ind+2];
1768 rvec_sub(x[ai], x[aj], dx);
1769 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
1770 shift = IVEC2IS(dt);
1774 shift = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
1777 /* Find the list for this shift or create one */
1778 list = find_gbtmplist(&nls[ai], shift);
1782 /* So that we do not add the same bond twice.
1783 * This happens with some constraints between 1-3 atoms
1784 * that are in the bond-list but should not be in the GB nb-list */
1785 for (j = 0; j < list->naj; j++)
1787 if (list->aj[j] == aj)
1797 gmx_incons("ai == aj");
1800 add_j_to_gblist(list, aj);
1806 compare_int (const void * a, const void * b)
1808 return ( *(int*)a - *(int*)b );
1813 int make_gb_nblist(t_commrec *cr, int gb_algorithm,
1814 rvec x[], matrix box,
1815 t_forcerec *fr, t_idef *idef, t_graph *graph, gmx_genborn_t *born)
1817 int i, l, ii, j, k, n, nj0, nj1, ai, aj, at0, at1, found, shift, s;
1822 struct gbtmpnbls *nls;
1823 gbtmpnbl_t *list = NULL;
1825 set_pbc(&pbc, fr->ePBC, box);
1826 nls = born->nblist_work;
1828 for (i = 0; i < born->nr; i++)
1835 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
1838 switch (gb_algorithm)
1842 /* Loop over 1-2, 1-3 and 1-4 interactions */
1843 for (j = F_GB12; j <= F_GB14; j++)
1845 add_bondeds_to_gblist(&idef->il[j], fr->bMolPBC, &pbc, graph, x, nls);
1849 /* Loop over 1-4 interactions */
1850 add_bondeds_to_gblist(&idef->il[F_GB14], fr->bMolPBC, &pbc, graph, x, nls);
1853 gmx_incons("Unknown GB algorithm");
1856 /* Loop over the VDWQQ and VDW nblists to set up the nonbonded part of the GB list */
1857 for (n = 0; (n < fr->nnblists); n++)
1859 for (i = 0; (i < eNL_NR); i++)
1861 nblist = &(fr->nblists[n].nlist_sr[i]);
1863 if (nblist->nri > 0 && (i == eNL_VDWQQ || i == eNL_QQ))
1865 for (j = 0; j < nblist->nri; j++)
1867 ai = nblist->iinr[j];
1868 shift = nblist->shift[j];
1870 /* Find the list for this shift or create one */
1871 list = find_gbtmplist(&nls[ai], shift);
1873 nj0 = nblist->jindex[j];
1874 nj1 = nblist->jindex[j+1];
1876 /* Add all the j-atoms in the non-bonded list to the GB list */
1877 for (k = nj0; k < nj1; k++)
1879 add_j_to_gblist(list, nblist->jjnr[k]);
1886 /* Zero out some counters */
1890 fr->gblist.jindex[0] = fr->gblist.nri;
1892 for (i = 0; i < fr->natoms_force; i++)
1894 for (s = 0; s < nls[i].nlist; s++)
1896 list = &nls[i].list[s];
1898 /* Only add those atoms that actually have neighbours */
1899 if (born->use[i] != 0)
1901 fr->gblist.iinr[fr->gblist.nri] = i;
1902 fr->gblist.shift[fr->gblist.nri] = list->shift;
1905 for (k = 0; k < list->naj; k++)
1907 /* Memory allocation for jjnr */
1908 if (fr->gblist.nrj >= fr->gblist.maxnrj)
1910 fr->gblist.maxnrj += over_alloc_large(fr->gblist.maxnrj);
1914 fprintf(debug, "Increasing GB neighbourlist j size to %d\n", fr->gblist.maxnrj);
1917 srenew(fr->gblist.jjnr, fr->gblist.maxnrj);
1921 if (i == list->aj[k])
1923 gmx_incons("i == list->aj[k]");
1925 fr->gblist.jjnr[fr->gblist.nrj++] = list->aj[k];
1928 fr->gblist.jindex[fr->gblist.nri] = fr->gblist.nrj;
1935 for (i = 0; i < fr->gblist.nri; i++)
1937 nj0 = fr->gblist.jindex[i];
1938 nj1 = fr->gblist.jindex[i+1];
1939 ai = fr->gblist.iinr[i];
1942 for (j = nj0; j < nj1; j++)
1944 if (fr->gblist.jjnr[j] < ai)
1946 fr->gblist.jjnr[j] += fr->natoms_force;
1949 qsort(fr->gblist.jjnr+nj0, nj1-nj0, sizeof(int), compare_int);
1951 for (j = nj0; j < nj1; j++)
1953 if (fr->gblist.jjnr[j] >= fr->natoms_force)
1955 fr->gblist.jjnr[j] -= fr->natoms_force;
1965 void make_local_gb(const t_commrec *cr, gmx_genborn_t *born, int gb_algorithm)
1968 gmx_domdec_t *dd = NULL;
1970 if (DOMAINDECOMP(cr))
1978 /* Single node or particle decomp (global==local), just copy pointers and return */
1979 if (gb_algorithm == egbSTILL)
1981 born->gpol = born->gpol_globalindex;
1982 born->vsolv = born->vsolv_globalindex;
1983 born->gb_radius = born->gb_radius_globalindex;
1987 born->param = born->param_globalindex;
1988 born->gb_radius = born->gb_radius_globalindex;
1991 born->use = born->use_globalindex;
1996 /* Reallocation of local arrays if necessary */
1997 /* fr->natoms_force is equal to dd->nat_tot */
1998 if (DOMAINDECOMP(cr) && dd->nat_tot > born->nalloc)
2002 nalloc = dd->nat_tot;
2004 /* Arrays specific to different gb algorithms */
2005 if (gb_algorithm == egbSTILL)
2007 srenew(born->gpol, nalloc+3);
2008 srenew(born->vsolv, nalloc+3);
2009 srenew(born->gb_radius, nalloc+3);
2010 for (i = born->nalloc; (i < nalloc+3); i++)
2014 born->gb_radius[i] = 0;
2019 srenew(born->param, nalloc+3);
2020 srenew(born->gb_radius, nalloc+3);
2021 for (i = born->nalloc; (i < nalloc+3); i++)
2024 born->gb_radius[i] = 0;
2028 /* All gb-algorithms use the array for vsites exclusions */
2029 srenew(born->use, nalloc+3);
2030 for (i = born->nalloc; (i < nalloc+3); i++)
2035 born->nalloc = nalloc;
2038 /* With dd, copy algorithm specific arrays */
2039 if (gb_algorithm == egbSTILL)
2041 for (i = at0; i < at1; i++)
2043 born->gpol[i] = born->gpol_globalindex[dd->gatindex[i]];
2044 born->vsolv[i] = born->vsolv_globalindex[dd->gatindex[i]];
2045 born->gb_radius[i] = born->gb_radius_globalindex[dd->gatindex[i]];
2046 born->use[i] = born->use_globalindex[dd->gatindex[i]];
2051 for (i = at0; i < at1; i++)
2053 born->param[i] = born->param_globalindex[dd->gatindex[i]];
2054 born->gb_radius[i] = born->gb_radius_globalindex[dd->gatindex[i]];
2055 born->use[i] = born->use_globalindex[dd->gatindex[i]];