real do_listed_vdw_q(int ftype,int nbonds, const t_iatom iatoms[],const t_iparams iparams[], const rvec x[],rvec f[],rvec fshift[], const t_pbc *pbc,const t_graph *g, real lambda,real *dvdlambda, const t_mdatoms *md, const t_forcerec *fr,gmx_grppairener_t *grppener, int *global_atom_index,gmx_mc_move *mc_move) { static bool bWarn=FALSE; real eps,r2,*tab,rtab2=0; rvec dx,x14[2],f14[2]; int i,ai,aj,itype; int typeA[2]={0,0},typeB[2]={0,1}; real chargeA[2]={0,0},chargeB[2]; int gid,shift_vir,shift_f; int j_index[] = { 0, 1 }; int i0=0,i1=1,i2=2; ivec dt; int outeriter,inneriter; int nthreads = 1; int count; real krf,crf,tabscale; int ntype=0; real *nbfp=NULL; real *egnb=NULL,*egcoul=NULL; t_nblist tmplist; int icoul,ivdw; bool bMolPBC,bFreeEnergy; #if GMX_THREADS pthread_mutex_t mtx; #else void * mtx = NULL; #endif #if GMX_THREADS pthread_mutex_initialize(&mtx); #endif bMolPBC = fr->bMolPBC; switch (ftype) { case F_LJ14: case F_LJC14_Q: eps = fr->epsfac*fr->fudgeQQ; ntype = 1; egnb = grppener->ener[egLJ14]; egcoul = grppener->ener[egCOUL14]; break; case F_LJC_PAIRS_NB: eps = fr->epsfac; ntype = 1; egnb = grppener->ener[egLJSR]; egcoul = grppener->ener[egCOULSR]; break; default: gmx_fatal(FARGS,"Unknown function type %d in do_nonbonded14", ftype); } tab = fr->tab14.tab; rtab2 = sqr(fr->tab14.r); tabscale = fr->tab14.scale; krf = fr->k_rf; crf = fr->c_rf; /* Determine the values for icoul/ivdw. */ if (fr->bEwald) { icoul = 1; } else if(fr->bcoultab) { icoul = 3; } else if(fr->eeltype == eelRF_NEC) { icoul = 2; } else { icoul = 1; } if(fr->bvdwtab) { ivdw = 3; } else if(fr->bBHAM) { ivdw = 2; } else { ivdw = 1; } /* We don't do SSE or altivec here, due to large overhead for 4-fold * unrolling on short lists */ bFreeEnergy = FALSE; for(i=0; (i<nbonds); ) { itype = iatoms[i++]; ai = iatoms[i++]; aj = iatoms[i++]; gid = GID(md->cENER[ai],md->cENER[aj],md->nenergrp); switch (ftype) { case F_LJ14: bFreeEnergy = (fr->efep != efepNO && ((md->nPerturbed && (md->bPerturbed[ai] || md->bPerturbed[aj])) || iparams[itype].lj14.c6A != iparams[itype].lj14.c6B || iparams[itype].lj14.c12A != iparams[itype].lj14.c12B)); chargeA[0] = md->chargeA[ai]; chargeA[1] = md->chargeA[aj]; nbfp = (real *)&(iparams[itype].lj14.c6A); break; case F_LJC14_Q: eps = fr->epsfac*iparams[itype].ljc14.fqq; chargeA[0] = iparams[itype].ljc14.qi; chargeA[1] = iparams[itype].ljc14.qj; nbfp = (real *)&(iparams[itype].ljc14.c6); break; case F_LJC_PAIRS_NB: chargeA[0] = iparams[itype].ljcnb.qi; chargeA[1] = iparams[itype].ljcnb.qj; nbfp = (real *)&(iparams[itype].ljcnb.c6); break; } if (!bMolPBC) { /* This is a bonded interaction, atoms are in the same box */ shift_f = CENTRAL; r2 = distance2(x[ai],x[aj]); } else { /* Apply full periodic boundary conditions */ shift_f = pbc_dx_aiuc(pbc,x[ai],x[aj],dx); r2 = norm2(dx); } if (r2 >= rtab2) { if (!bWarn) { fprintf(stderr,"Warning: 1-4 interaction between %d and %d " "at distance %.3f which is larger than the 1-4 table size %.3f nm\n", glatnr(global_atom_index,ai), glatnr(global_atom_index,aj), sqrt(r2), sqrt(rtab2)); fprintf(stderr,"These are ignored for the rest of the simulation\n"); fprintf(stderr,"This usually means your system is exploding,\n" "if not, you should increase table-extension in your mdp file\n" "or with user tables increase the table size\n"); bWarn = TRUE; } if (debug) fprintf(debug,"%8f %8f %8f\n%8f %8f %8f\n1-4 (%d,%d) interaction not within cut-off! r=%g. Ignored\n", x[ai][XX],x[ai][YY],x[ai][ZZ], x[aj][XX],x[aj][YY],x[aj][ZZ], glatnr(global_atom_index,ai), glatnr(global_atom_index,aj), sqrt(r2)); } else { copy_rvec(x[ai],x14[0]); copy_rvec(x[aj],x14[1]); clear_rvec(f14[0]); clear_rvec(f14[1]); #ifdef DEBUG fprintf(debug,"LJ14: grp-i=%2d, grp-j=%2d, ngrp=%2d, GID=%d\n", md->cENER[ai],md->cENER[aj],md->nenergrp,gid); #endif outeriter = inneriter = count = 0; if (bFreeEnergy) { chargeB[0] = md->chargeB[ai]; chargeB[1] = md->chargeB[aj]; /* We pass &(iparams[itype].lj14.c6A) as LJ parameter matrix * to the innerloops. * Here we use that the LJ-14 parameters are stored in iparams * as c6A,c12A,c6B,c12B, which are referenced correctly * in the innerloops if we assign type combinations 0-0 and 0-1 * to atom pair ai-aj in topologies A and B respectively. */ if(ivdw==2) { gmx_fatal(FARGS,"Cannot do free energy Buckingham interactions."); } count = 0; gmx_nb_free_energy_kernel(icoul, ivdw, i1, &i0, j_index, &i1, &shift_f, fr->shift_vec[0], fshift[0], &gid, x14[0], f14[0], chargeA, chargeB, eps, krf, crf, fr->ewaldcoeff, egcoul, typeA, typeB, ntype, nbfp, egnb, tabscale, tab, lambda, dvdlambda, fr->sc_alpha, fr->sc_power, fr->sc_sigma6, TRUE, &outeriter, &inneriter); } else { /* Not perturbed - call kernel 330 */ nb_kernel330 ( &i1, &i0, j_index, &i1, &shift_f, fr->shift_vec[0], fshift[0], &gid, x14[0], f14[0], chargeA, &eps, &krf, &crf, egcoul, typeA, &ntype, nbfp, egnb, &tabscale, tab, mc_move ? mc_move->enerd[F_COUL14] : NULL, mc_move ? mc_move->enerd[F_LJ14] : NULL, (mc_move) ? mc_move->enerd_prev[F_COUL14] : NULL, (mc_move) ? mc_move->enerd_prev[F_LJ14] : NULL, (mc_move) ? &mc_move->start : NULL, (mc_move) ? &mc_move->end : NULL, (mc_move) ? &mc_move->homenr : NULL, (mc_move) ? mc_move->sum_index : NULL, NULL, NULL, NULL, NULL, &nthreads, &count, (void *)&mtx, &outeriter, &inneriter, NULL); } /* Add the forces */ rvec_inc(f[ai],f14[0]); rvec_dec(f[aj],f14[0]); if (g) { /* Correct the shift forces using the graph */ ivec_sub(SHIFT_IVEC(g,ai),SHIFT_IVEC(g,aj),dt); shift_vir = IVEC2IS(dt); rvec_inc(fshift[shift_vir],f14[0]); rvec_dec(fshift[CENTRAL],f14[0]); } /* flops: eNR_KERNEL_OUTER + eNR_KERNEL330 + 12 */ } } return 0.0; }
void dd_move_f_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac, rvec *f, rvec *fshift) { gmx_specatsend_t *spas; rvec *vbuf; int n, n0, n1, d, dim, dir, i; ivec vis; int is; gmx_bool bPBC, bScrew; n = spac->at_end; for (d = dd->ndim-1; d >= 0; d--) { dim = dd->dim[d]; if (dd->nc[dim] > 2) { /* Pulse the grid forward and backward */ spas = spac->spas[d]; n0 = spas[0].nrecv; n1 = spas[1].nrecv; n -= n1 + n0; vbuf = spac->vbuf; /* Send and receive the coordinates */ dd_sendrecv2_rvec(dd, d, f+n+n1, n0, vbuf, spas[0].nsend, f+n, n1, vbuf+spas[0].nsend, spas[1].nsend); for (dir = 0; dir < 2; dir++) { bPBC = ((dir == 0 && dd->ci[dim] == 0) || (dir == 1 && dd->ci[dim] == dd->nc[dim]-1)); bScrew = (bPBC && dd->bScrewPBC && dim == XX); spas = &spac->spas[d][dir]; /* Sum the buffer into the required forces */ if (!bPBC || (!bScrew && fshift == NULL)) { for (i = 0; i < spas->nsend; i++) { rvec_inc(f[spas->a[i]], *vbuf); vbuf++; } } else { clear_ivec(vis); vis[dim] = (dir == 0 ? 1 : -1); is = IVEC2IS(vis); if (!bScrew) { /* Sum and add to shift forces */ for (i = 0; i < spas->nsend; i++) { rvec_inc(f[spas->a[i]], *vbuf); rvec_inc(fshift[is], *vbuf); vbuf++; } } else { /* Rotate the forces */ for (i = 0; i < spas->nsend; i++) { f[spas->a[i]][XX] += (*vbuf)[XX]; f[spas->a[i]][YY] -= (*vbuf)[YY]; f[spas->a[i]][ZZ] -= (*vbuf)[ZZ]; if (fshift) { rvec_inc(fshift[is], *vbuf); } vbuf++; } } } } } else { /* Two cells, so we only need to communicate one way */ spas = &spac->spas[d][0]; n -= spas->nrecv; /* Send and receive the coordinates */ dd_sendrecv_rvec(dd, d, dddirForward, f+n, spas->nrecv, spac->vbuf, spas->nsend); /* Sum the buffer into the required forces */ if (dd->bScrewPBC && dim == XX && (dd->ci[dim] == 0 || dd->ci[dim] == dd->nc[dim]-1)) { for (i = 0; i < spas->nsend; i++) { /* Rotate the force */ f[spas->a[i]][XX] += spac->vbuf[i][XX]; f[spas->a[i]][YY] -= spac->vbuf[i][YY]; f[spas->a[i]][ZZ] -= spac->vbuf[i][ZZ]; } } else { for (i = 0; i < spas->nsend; i++) { rvec_inc(f[spas->a[i]], spac->vbuf[i]); } } } } }
real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[], const rvec x[], rvec f[], rvec fshift[], const t_pbc *pbc, const t_graph *g, real gmx_unused lambda, real gmx_unused *dvdlambda, const t_mdatoms gmx_unused *md, t_fcdata *fcd, int gmx_unused *global_atom_index) { atom_id ai, aj; int fa, d, i, type, ex, power, ki = CENTRAL; ivec dt; real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac; rvec r, Sr, fij; real vtot; const t_oriresdata *od; gmx_bool bTAV; vtot = 0; od = &(fcd->orires); if (od->fc != 0) { bTAV = (od->edt != 0); smooth_fc = od->fc; if (bTAV) { /* Smoothly switch on the restraining when time averaging is used */ smooth_fc *= (1.0 - od->exp_min_t_tau); } d = 0; for (fa = 0; fa < nfa; fa += 3) { type = forceatoms[fa]; ai = forceatoms[fa+1]; aj = forceatoms[fa+2]; if (pbc) { ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r); } else { rvec_sub(x[ai], x[aj], r); } r2 = norm2(r); invr = gmx_invsqrt(r2); invr2 = invr*invr; ex = ip[type].orires.ex; power = ip[type].orires.power; fc = smooth_fc*ip[type].orires.kfac; dev = od->otav[d] - ip[type].orires.obs; /* NOTE: * there is no real potential when time averaging is applied */ vtot += 0.5*fc*sqr(dev); if (bTAV) { /* Calculate the force as the sqrt of tav times instantaneous */ devins = od->oins[d] - ip[type].orires.obs; if (dev*devins <= 0) { dev = 0; } else { dev = std::sqrt(dev*devins); if (devins < 0) { dev = -dev; } } } pfac = fc*ip[type].orires.c*invr2; for (i = 0; i < power; i++) { pfac *= invr; } mvmul(od->S[ex], r, Sr); for (i = 0; i < DIM; i++) { fij[i] = -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]); } if (g) { ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt); ki = IVEC2IS(dt); } for (i = 0; i < DIM; i++) { f[ai][i] += fij[i]; f[aj][i] -= fij[i]; fshift[ki][i] += fij[i]; fshift[CENTRAL][i] -= fij[i]; } d++; } } return vtot; /* Approx. 80*nfa/3 flops */ }
int pbc_dx_aiuc(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx) { int i, j, is; rvec dx_start, trial; real d2min, d2trial; ivec ishift, ishift_start; rvec_sub(x1, x2, dx); clear_ivec(ishift); switch (pbc->ePBCDX) { case epbcdxRECTANGULAR: for (i = 0; i < DIM; i++) { if (dx[i] > pbc->hbox_diag[i]) { dx[i] -= pbc->fbox_diag[i]; ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { dx[i] += pbc->fbox_diag[i]; ishift[i]++; } } break; case epbcdxTRICLINIC: /* For triclinic boxes the performance difference between * if/else and two while loops is negligible. * However, the while version can cause extreme delays * before a simulation crashes due to large forces which * can cause unlimited displacements. * Also allowing multiple shifts would index fshift beyond bounds. */ for (i = DIM-1; i >= 1; i--) { if (dx[i] > pbc->hbox_diag[i]) { for (j = i; j >= 0; j--) { dx[j] -= pbc->box[i][j]; } ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { for (j = i; j >= 0; j--) { dx[j] += pbc->box[i][j]; } ishift[i]++; } } /* Allow 2 shifts in x */ if (dx[XX] > pbc->hbox_diag[XX]) { dx[XX] -= pbc->fbox_diag[XX]; ishift[XX]--; if (dx[XX] > pbc->hbox_diag[XX]) { dx[XX] -= pbc->fbox_diag[XX]; ishift[XX]--; } } else if (dx[XX] <= pbc->mhbox_diag[XX]) { dx[XX] += pbc->fbox_diag[XX]; ishift[XX]++; if (dx[XX] <= pbc->mhbox_diag[XX]) { dx[XX] += pbc->fbox_diag[XX]; ishift[XX]++; } } /* dx is the distance in a rectangular box */ d2min = norm2(dx); if (d2min > pbc->max_cutoff2) { copy_rvec(dx, dx_start); copy_ivec(ishift, ishift_start); d2min = norm2(dx); /* Now try all possible shifts, when the distance is within max_cutoff * it must be the shortest possible distance. */ i = 0; while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec)) { rvec_add(dx_start, pbc->tric_vec[i], trial); d2trial = norm2(trial); if (d2trial < d2min) { copy_rvec(trial, dx); ivec_add(ishift_start, pbc->tric_shift[i], ishift); d2min = d2trial; } i++; } } break; case epbcdx2D_RECT: for (i = 0; i < DIM; i++) { if (i != pbc->dim) { if (dx[i] > pbc->hbox_diag[i]) { dx[i] -= pbc->fbox_diag[i]; ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { dx[i] += pbc->fbox_diag[i]; ishift[i]++; } } } break; case epbcdx2D_TRIC: d2min = 0; for (i = DIM-1; i >= 1; i--) { if (i != pbc->dim) { if (dx[i] > pbc->hbox_diag[i]) { for (j = i; j >= 0; j--) { dx[j] -= pbc->box[i][j]; } ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { for (j = i; j >= 0; j--) { dx[j] += pbc->box[i][j]; } ishift[i]++; } d2min += dx[i]*dx[i]; } } if (pbc->dim != XX) { /* Allow 2 shifts in x */ if (dx[XX] > pbc->hbox_diag[XX]) { dx[XX] -= pbc->fbox_diag[XX]; ishift[XX]--; if (dx[XX] > pbc->hbox_diag[XX]) { dx[XX] -= pbc->fbox_diag[XX]; ishift[XX]--; } } else if (dx[XX] <= pbc->mhbox_diag[XX]) { dx[XX] += pbc->fbox_diag[XX]; ishift[XX]++; if (dx[XX] <= pbc->mhbox_diag[XX]) { dx[XX] += pbc->fbox_diag[XX]; ishift[XX]++; } } d2min += dx[XX]*dx[XX]; } if (d2min > pbc->max_cutoff2) { copy_rvec(dx, dx_start); copy_ivec(ishift, ishift_start); /* Now try all possible shifts, when the distance is within max_cutoff * it must be the shortest possible distance. */ i = 0; while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec)) { rvec_add(dx_start, pbc->tric_vec[i], trial); d2trial = 0; for (j = 0; j < DIM; j++) { if (j != pbc->dim) { d2trial += trial[j]*trial[j]; } } if (d2trial < d2min) { copy_rvec(trial, dx); ivec_add(ishift_start, pbc->tric_shift[i], ishift); d2min = d2trial; } i++; } } break; case epbcdx1D_RECT: i = pbc->dim; if (dx[i] > pbc->hbox_diag[i]) { dx[i] -= pbc->fbox_diag[i]; ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { dx[i] += pbc->fbox_diag[i]; ishift[i]++; } break; case epbcdx1D_TRIC: i = pbc->dim; if (dx[i] > pbc->hbox_diag[i]) { rvec_dec(dx, pbc->box[i]); ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { rvec_inc(dx, pbc->box[i]); ishift[i]++; } break; case epbcdxSCREW_RECT: /* The shift definition requires x first */ if (dx[XX] > pbc->hbox_diag[XX]) { dx[XX] -= pbc->fbox_diag[XX]; ishift[XX]--; } else if (dx[XX] <= pbc->mhbox_diag[XX]) { dx[XX] += pbc->fbox_diag[XX]; ishift[XX]++; } if (ishift[XX] == 1 || ishift[XX] == -1) { /* Rotate around the x-axis in the middle of the box */ dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY]; dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ]; } /* Normal pbc for y and z */ for (i = YY; i <= ZZ; i++) { if (dx[i] > pbc->hbox_diag[i]) { dx[i] -= pbc->fbox_diag[i]; ishift[i]--; } else if (dx[i] <= pbc->mhbox_diag[i]) { dx[i] += pbc->fbox_diag[i]; ishift[i]++; } } break; case epbcdxNOPBC: case epbcdxUNSUPPORTED: break; default: gmx_fatal(FARGS, "Internal error in pbc_dx_aiuc, set_pbc_dd or set_pbc has not been called"); break; } is = IVEC2IS(ishift); if (debug) { range_check_mesg(is, 0, SHIFTS, "PBC shift vector index range check."); } return is; }
real ta_disres(int nfa, const t_iatom forceatoms[], const t_iparams ip[], const rvec x[], rvec f[], rvec fshift[], const t_pbc *pbc, const t_graph *g, real gmx_unused lambda, real gmx_unused *dvdlambda, const t_mdatoms gmx_unused *md, t_fcdata *fcd, int gmx_unused *global_atom_index) { const real sixth = 1.0/6.0; const real seven_three = 7.0/3.0; atom_id ai, aj; int fa, res, npair, p, pair, ki = CENTRAL, m; int type; rvec dx; real weight_rt_1; real smooth_fc, Rt, Rtav, rt2, *Rtl_6, *Rt_6, *Rtav_6; real k0, f_scal = 0, fmax_scal, fk_scal, fij; real tav_viol, instant_viol, mixed_viol, violtot, vtot; real tav_viol_Rtav7, instant_viol_Rtav7; real up1, up2, low; gmx_bool bConservative, bMixed, bViolation; ivec it, jt, dt; t_disresdata *dd; int dr_weighting; gmx_bool dr_bMixed; dd = &(fcd->disres); dr_weighting = dd->dr_weighting; dr_bMixed = dd->dr_bMixed; Rtl_6 = dd->Rtl_6; Rt_6 = dd->Rt_6; Rtav_6 = dd->Rtav_6; tav_viol = instant_viol = mixed_viol = tav_viol_Rtav7 = instant_viol_Rtav7 = 0; smooth_fc = dd->dr_fc; if (dd->dr_tau != 0) { /* scaling factor to smoothly turn on the restraint forces * * when using time averaging */ smooth_fc *= (1.0 - dd->exp_min_t_tau); } violtot = 0; vtot = 0; /* 'loop' over all atom pairs (pair_nr=fa/3) involved in restraints, * * the total number of atoms pairs is nfa/3 */ res = 0; fa = 0; while (fa < nfa) { type = forceatoms[fa]; /* Take action depending on restraint, calculate scalar force */ npair = ip[type].disres.npair; up1 = ip[type].disres.up1; up2 = ip[type].disres.up2; low = ip[type].disres.low; k0 = smooth_fc*ip[type].disres.kfac; /* save some flops when there is only one pair */ if (ip[type].disres.type != 2) { bConservative = (dr_weighting == edrwConservative) && (npair > 1); bMixed = dr_bMixed; Rt = pow(Rt_6[res], -sixth); Rtav = pow(Rtav_6[res], -sixth); } else { /* When rtype=2 use instantaneous not ensemble avereged distance */ bConservative = (npair > 1); bMixed = FALSE; Rt = pow(Rtl_6[res], -sixth); Rtav = Rt; } if (Rtav > up1) { bViolation = TRUE; tav_viol = Rtav - up1; } else if (Rtav < low) { bViolation = TRUE; tav_viol = Rtav - low; } else { bViolation = FALSE; } if (bViolation) { /* NOTE: * there is no real potential when time averaging is applied */ vtot += 0.5*k0*sqr(tav_viol); if (1/vtot == 0) { printf("vtot is inf: %f\n", vtot); } if (!bMixed) { f_scal = -k0*tav_viol; violtot += fabs(tav_viol); } else { if (Rt > up1) { if (tav_viol > 0) { instant_viol = Rt - up1; } else { bViolation = FALSE; } } else if (Rt < low) { if (tav_viol < 0) { instant_viol = Rt - low; } else { bViolation = FALSE; } } else { bViolation = FALSE; } if (bViolation) { mixed_viol = sqrt(tav_viol*instant_viol); f_scal = -k0*mixed_viol; violtot += mixed_viol; } } } if (bViolation) { fmax_scal = -k0*(up2-up1); /* Correct the force for the number of restraints */ if (bConservative) { f_scal = max(f_scal, fmax_scal); if (!bMixed) { f_scal *= Rtav/Rtav_6[res]; } else { f_scal /= 2*mixed_viol; tav_viol_Rtav7 = tav_viol*Rtav/Rtav_6[res]; instant_viol_Rtav7 = instant_viol*Rt/Rt_6[res]; } } else { f_scal /= (real)npair; f_scal = max(f_scal, fmax_scal); } /* Exert the force ... */ /* Loop over the atom pairs of 'this' restraint */ for (p = 0; p < npair; p++) { pair = fa/3; ai = forceatoms[fa+1]; aj = forceatoms[fa+2]; if (pbc) { ki = pbc_dx_aiuc(pbc, x[ai], x[aj], dx); } else { rvec_sub(x[ai], x[aj], dx); } rt2 = iprod(dx, dx); weight_rt_1 = gmx_invsqrt(rt2); if (bConservative) { if (!dr_bMixed) { weight_rt_1 *= pow(dd->rm3tav[pair], seven_three); } else { weight_rt_1 *= tav_viol_Rtav7*pow(dd->rm3tav[pair], seven_three)+ instant_viol_Rtav7*pow(dd->rt[pair], -7); } } fk_scal = f_scal*weight_rt_1; if (g) { ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt); ki = IVEC2IS(dt); } for (m = 0; m < DIM; m++) { fij = fk_scal*dx[m]; f[ai][m] += fij; f[aj][m] -= fij; fshift[ki][m] += fij; fshift[CENTRAL][m] -= fij; } fa += 3; } } else { /* No violation so force and potential contributions */ fa += 3*npair; } res++; } dd->sumviol = violtot; /* Return energy */ return vtot; }
real orires(int nfa,t_iatom forceatoms[],t_iparams ip[], rvec x[],rvec f[],t_forcerec *fr,t_graph *g, matrix box,real lambda,real *dvdlambda, t_mdatoms *md,int ngrp,real egnb[],real egcoul[], t_fcdata *fcd) { atom_id ai,aj; int fa,d,i,type,ex,power,ki; ivec dt; real r2,invr,invr2,fc,smooth_fc,dev,devins,pfac; rvec r,Sr,fij; real vtot; t_oriresdata *od; bool bTAV; vtot = 0; od = &(fcd->orires); if (fabs(od->fc) > GMX_REAL_MIN) { bTAV = (fabs(od->edt) > GMX_REAL_MIN); /* Smoothly switch on the restraining when time averaging is used */ smooth_fc = od->fc*(1.0 - od->exp_min_t_tau); d = 0; for(fa=0; fa<nfa; fa+=3) { type = forceatoms[fa]; ai = forceatoms[fa+1]; aj = forceatoms[fa+2]; rvec_sub(x[ai],x[aj],r); r2 = norm2(r); invr = invsqrt(r2); invr2 = invr*invr; ex = ip[type].orires.ex; power = ip[type].orires.pow; fc = smooth_fc*ip[type].orires.kfac; dev = od->otav[d] - ip[type].orires.obs; /* NOTE: there is no real potential when time averaging is applied */ vtot += 0.5*fc*sqr(dev); if (bTAV) { /* Calculate the force as the sqrt of tav times instantaneous */ devins = od->oins[d] - ip[type].orires.obs; if (dev*devins <= 0) dev = 0; else { dev = sqrt(dev*devins); if (devins < 0) dev = -dev; } } pfac = fc*ip[type].orires.c*invr2; for(i=0; i<power; i++) pfac *= invr; mvmul(od->S[ex],r,Sr); for(i=0; i<DIM; i++) fij[i] = -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr,r)*r[i]); ivec_sub(SHIFT_IVEC(g,ai),SHIFT_IVEC(g,aj),dt); ki=IVEC2IS(dt); for(i=0; i<DIM; i++) { f[ai][i] += fij[i]; f[aj][i] -= fij[i]; fr->fshift[ki][i] += fij[i]; fr->fshift[CENTRAL][i] -= fij[i]; } d++; } } return vtot; /* Approx. 80*nfa/3 flops */ }