static void constr_dum3OUT(rvec xi,rvec xj,rvec xk,rvec x,real a,real b,real c) { rvec xij,xik,temp; rvec_sub(xj,xi,xij); rvec_sub(xk,xi,xik); oprod(xij,xik,temp); /* 15 Flops */ x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX]; x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY]; x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ]; /* 18 Flops */ /* TOTAL: 33 flops */ }
void calc_h_pos(int nht, rvec xa[], rvec xh[]) { #define alfaH (DEG2RAD*109.5) #define distH 0.1 #define alfaCOM (DEG2RAD*117) #define alfaCO (DEG2RAD*121) #define alfaCOA (DEG2RAD*115) #define distO 0.123 #define distOA 0.125 #define distOM 0.136 rvec sa,sb,sij; real s6,rij,ra,rb,xd; int d; s6=0.5*sqrt(3.e0); /* common work for constructing one, two or three dihedral hydrogens */ switch (nht) { case 2: case 3: case 4: case 8: case 9: rij = 0.e0; for(d=0; (d<DIM); d++) { xd = xAJ[d]; sij[d] = xAI[d]-xd; sb[d] = xd-xAK[d]; rij += sqr(sij[d]); } rij = sqrt(rij); sa[XX] = sij[YY]*sb[ZZ]-sij[ZZ]*sb[YY]; sa[YY] = sij[ZZ]*sb[XX]-sij[XX]*sb[ZZ]; sa[ZZ] = sij[XX]*sb[YY]-sij[YY]*sb[XX]; ra = 0.e0; for(d=0; (d<DIM); d++) { sij[d] = sij[d]/rij; ra += sqr(sa[d]); } ra = sqrt(ra); for(d=0; (d<DIM); d++) sa[d] = sa[d]/ra; sb[XX] = sa[YY]*sij[ZZ]-sa[ZZ]*sij[YY]; sb[YY] = sa[ZZ]*sij[XX]-sa[XX]*sij[ZZ]; sb[ZZ] = sa[XX]*sij[YY]-sa[YY]*sij[XX]; break; }/* end switch */ switch (nht) { case 1: /* construct one planar hydrogen (peptide,rings) */ rij = 0.e0; rb = 0.e0; for(d=0; (d<DIM); d++) { sij[d] = xAI[d]-xAJ[d]; sb[d] = xAI[d]-xAK[d]; rij += sqr(sij[d]); rb += sqr(sb[d]); } rij = sqrt(rij); rb = sqrt(rb); ra = 0.e0; for(d=0; (d<DIM); d++) { sa[d] = sij[d]/rij+sb[d]/rb; ra += sqr(sa[d]); } ra = sqrt(ra); for(d=0; (d<DIM); d++) xH1[d] = xAI[d]+distH*sa[d]/ra; break; case 2: /* one single hydrogen, e.g. hydroxyl */ for(d=0; (d<DIM); d++) { xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d]; } break; case 3: /* two planar hydrogens, e.g. -NH2 */ for(d=0; (d<DIM); d++) { xH1[d] = xAI[d]-distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d]; xH2[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d]; } break; case 4: /* two or three tetrahedral hydrogens, e.g. -CH3 */ for(d=0; (d<DIM); d++) { xH1[d] = xAI[d]+distH*sin(alfaH)*sb[d]-distH*cos(alfaH)*sij[d]; xH2[d] = ( xAI[d] - distH*sin(alfaH)*0.5*sb[d] + distH*sin(alfaH)*s6*sa[d] - distH*cos(alfaH)*sij[d] ); if ( fabs(xH3[XX]-NOTSET)>GMX_REAL_MIN && fabs(xH3[YY]-NOTSET)>GMX_REAL_MIN && fabs(xH3[ZZ]-NOTSET)>GMX_REAL_MIN ) xH3[d] = ( xAI[d] - distH*sin(alfaH)*0.5*sb[d] - distH*sin(alfaH)*s6*sa[d] - distH*cos(alfaH)*sij[d] ); } break; case 5: { /* one tetrahedral hydrogen, e.g. C3CH */ real center; rvec dxc; for(d=0; (d<DIM); d++) { center=(xAJ[d]+xAK[d]+xAL[d])/3.0; dxc[d]=xAI[d]-center; } center=norm(dxc); for(d=0; (d<DIM); d++) xH1[d]=xAI[d]+dxc[d]*distH/center; break; } case 6: { /* two tetrahedral hydrogens, e.g. C-CH2-C */ rvec rBB,rCC1,rCC2,rNN; real bb,nn; for(d=0; (d<DIM); d++) rBB[d]=xAI[d]-0.5*(xAJ[d]+xAK[d]); bb=norm(rBB); rvec_sub(xAI,xAJ,rCC1); rvec_sub(xAI,xAK,rCC2); oprod(rCC1,rCC2,rNN); nn=norm(rNN); for(d=0; (d<DIM); d++) { xH1[d]=xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb+ sin(alfaH/2.0)*rNN[d]/nn); xH2[d]=xAI[d]+distH*(cos(alfaH/2.0)*rBB[d]/bb- sin(alfaH/2.0)*rNN[d]/nn); } break; } case 7: /* two water hydrogens */ gen_waterhydrogen(xa, xh); break; case 8: /* two carboxyl oxygens, -COO- */ for(d=0; (d<DIM); d++) { xH1[d] = xAI[d]-distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d]; xH2[d] = xAI[d]+distOM*sin(alfaCOM)*sb[d]-distOM*cos(alfaCOM)*sij[d]; } break; case 9: { /* carboxyl oxygens and hydrogen, -COOH */ rvec xa2[4]; /* i,j,k,l */ /* first add two oxygens */ for(d=0; (d<DIM); d++) { xH1[d] = xAI[d]-distO *sin(alfaCO )*sb[d]-distO *cos(alfaCO )*sij[d]; xH2[d] = xAI[d]+distOA*sin(alfaCOA)*sb[d]-distOA*cos(alfaCOA)*sij[d]; } /* now use rule 2 to add hydrogen to 2nd oxygen */ copy_rvec(xH2, xa2[0]); /* new i = n' */ copy_rvec(xAI, xa2[1]); /* new j = i */ copy_rvec(xAJ, xa2[2]); /* new k = j */ copy_rvec(xAK, xa2[3]); /* new l = k, not used */ calc_h_pos(2, xa2, (xh+2)); break; } default: fatal_error(0,"Invalid argument (%d) for nht in routine genh\n",nht); } /* end switch */ }
void PREFIX angle_restraint(int i_c, struct mtd_data_s *mtd_data) { int i, j, iat, k, ix; real myp[3][3], totmasse[3], d[3][3]; rvec rij, sij, tij; real mod_rij, mod_sij; real eps, sign0; real caa, cbb, cab, ccc, sab; real ac, Vac, bc, dVac, cos_psi, sin_psi; totmasse[0] = totmasse[1] = totmasse[2] = 0.; myp[0][0] = myp[0][1] = myp[0][2] = 0.; myp[1][0] = myp[1][1] = myp[1][2] = 0.; myp[2][0] = myp[2][1] = myp[2][2] = 0.; k = 0; for(j=0;j<3;j++) { for(i=0;i<colvar.list[i_c][j];i++){ iat = colvar.cvatoms[i_c][k]; myp[j][0] += mtd_data->mass[iat]*mtd_data->pos[iat][0]; myp[j][1] += mtd_data->mass[iat]*mtd_data->pos[iat][1]; myp[j][2] += mtd_data->mass[iat]*mtd_data->pos[iat][2]; totmasse[j] += mtd_data->mass[iat]; k++; } myp[j][0] /= totmasse[j]; myp[j][1] /= totmasse[j]; myp[j][2] /= totmasse[j]; } minimal_image(myp[1], myp[0], &mod_rij, rij); minimal_image(myp[1], myp[2], &mod_sij, sij); caa = norm2(rij); cab = iprod(rij,sij); cbb = norm2(sij); oprod(rij,sij,tij); sab = norm(tij); ccc = 1.0/sqrt(caa*cbb); ac = cab*ccc; // cosine of teta Vac = acos(ac); dVac = -ccc/sqrt(1.-ac*ac); for(ix=0;ix<3;ix++) { d[0][ix] = -dVac*(-cab/caa*rij[ix]+sij[ix]); d[1][ix] = dVac*(-cab/caa*rij[ix]-cab/cbb*sij[ix]+rij[ix]+sij[ix]); d[2][ix] = dVac*(cab/cbb*sij[ix]-rij[ix]); } // Now we do appropriate multiplication of the derivatives to get what we are interested in if( colvar.doTrig[i_c]==0 ){ colvar.ss0[i_c] = Vac; } else if( colvar.doTrig[i_c]==1 ){ cos_psi=cos( Vac ); for(i=0;i<3;++i){d[0][i]*= cos_psi; d[1][i]*= cos_psi; d[2][i]*= cos_psi;} colvar.ss0[i_c] = sin( Vac ); } else if( colvar.doTrig[i_c]==2 ){ sin_psi=-sin( Vac ); for(i=0;i<3;++i){d[0][i]*= sin_psi; d[1][i]*= sin_psi; d[2][i]*= sin_psi;} colvar.ss0[i_c] = cos( Vac ); } else{ plumed_error("No trigonometric mode defined in torsion restraint"); } k=0; for(j=0;j<3;j++) { for(i=0;i<colvar.list[i_c][j];i++){ iat = colvar.cvatoms[i_c][k]; for(ix=0;ix<3;ix++) colvar.myder[i_c][k][ix] = d[j][ix]*mtd_data->mass[iat]/totmasse[j]; k++; } } //colvar.ss0[i_c] = Vac; }