2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, 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.
51 #include "mtop_util.h"
53 void init_orires(FILE *fplog, const gmx_mtop_t *mtop,
56 const t_commrec *cr, t_oriresdata *od,
59 int i, j, d, ex, nmol, *nr_ex;
62 gmx_mtop_ilistloop_t iloop;
64 gmx_mtop_atomloop_all_t aloop;
66 const gmx_multisim_t *ms;
68 od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
71 /* Not doing orientation restraints */
77 gmx_fatal(FARGS, "Orientation restraints do not work with more than one domain (ie. MPI rank).");
79 /* Orientation restraints */
87 od->fc = ir->orires_fc;
96 iloop = gmx_mtop_ilistloop_init(mtop);
97 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
99 for (i = 0; i < il[F_ORIRES].nr; i += 3)
101 ex = mtop->ffparams.iparams[il[F_ORIRES].iatoms[i]].orires.ex;
105 for (j = od->nex; j < ex+1; j++)
114 snew(od->S, od->nex);
115 /* When not doing time averaging, the instaneous and time averaged data
116 * are indentical and the pointers can point to the same memory.
118 snew(od->Dinsl, od->nr);
121 snew(od->Dins, od->nr);
125 od->Dins = od->Dinsl;
128 if (ir->orires_tau == 0)
136 snew(od->Dtav, od->nr);
137 od->edt = exp(-ir->delta_t/ir->orires_tau);
138 od->edt_1 = 1.0 - od->edt;
140 /* Extend the state with the orires history */
141 state->flags |= (1<<estORIRE_INITF);
142 state->hist.orire_initf = 1;
143 state->flags |= (1<<estORIRE_DTAV);
144 state->hist.norire_Dtav = od->nr*5;
145 snew(state->hist.orire_Dtav, state->hist.norire_Dtav);
148 snew(od->oinsl, od->nr);
151 snew(od->oins, od->nr);
155 od->oins = od->oinsl;
157 if (ir->orires_tau == 0)
163 snew(od->otav, od->nr);
165 snew(od->tmp, od->nex);
166 snew(od->TMP, od->nex);
167 for (ex = 0; ex < od->nex; ex++)
169 snew(od->TMP[ex], 5);
170 for (i = 0; i < 5; i++)
172 snew(od->TMP[ex][i], 5);
177 for (i = 0; i < mtop->natoms; i++)
179 if (ggrpnr(&mtop->groups, egcORFIT, i) == 0)
184 snew(od->mref, od->nref);
185 snew(od->xref, od->nref);
186 snew(od->xtmp, od->nref);
188 snew(od->eig, od->nex*12);
190 /* Determine the reference structure on the master node.
191 * Copy it to the other nodes after checking multi compatibility,
192 * so we are sure the subsystems match before copying.
197 aloop = gmx_mtop_atomloop_all_init(mtop);
198 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
200 if (mtop->groups.grpnr[egcORFIT] == NULL ||
201 mtop->groups.grpnr[egcORFIT][i] == 0)
203 /* Not correct for free-energy with changing masses */
204 od->mref[j] = atom->m;
205 if (ms == NULL || MASTERSIM(ms))
207 copy_rvec(xref[i], od->xref[j]);
208 for (d = 0; d < DIM; d++)
210 com[d] += od->mref[j]*xref[i][d];
217 svmul(1.0/mtot, com, com);
218 if (ms == NULL || MASTERSIM(ms))
220 for (j = 0; j < od->nref; j++)
222 rvec_dec(od->xref[j], com);
226 fprintf(fplog, "Found %d orientation experiments\n", od->nex);
227 for (i = 0; i < od->nex; i++)
229 fprintf(fplog, " experiment %d has %d restraints\n", i+1, nr_ex[i]);
234 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n",
239 fprintf(fplog, " the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
241 check_multi_int(fplog, ms, od->nr,
242 "the number of orientation restraints",
244 check_multi_int(fplog, ms, od->nref,
245 "the number of fit atoms for orientation restraining",
247 check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
248 /* Copy the reference coordinates from the master to the other nodes */
249 gmx_sum_sim(DIM*od->nref, od->xref[0], ms);
252 please_cite(fplog, "Hess2003");
255 void diagonalize_orires_tensors(t_oriresdata *od)
257 int ex, i, j, nrot, ord[DIM], t;
263 for (i = 0; i < DIM; i++)
267 snew(od->eig_diag, DIM);
269 for (i = 0; i < DIM; i++)
275 for (ex = 0; ex < od->nex; ex++)
277 /* Rotate the S tensor back to the reference frame */
278 mmul(od->R, od->S[ex], TMP);
279 mtmul(TMP, od->R, S);
280 for (i = 0; i < DIM; i++)
282 for (j = 0; j < DIM; j++)
284 od->M[i][j] = S[i][j];
288 jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
290 for (i = 0; i < DIM; i++)
294 for (i = 0; i < DIM; i++)
296 for (j = i+1; j < DIM; j++)
298 if (sqr(od->eig_diag[ord[j]]) > sqr(od->eig_diag[ord[i]]))
307 for (i = 0; i < DIM; i++)
309 od->eig[ex*12 + i] = od->eig_diag[ord[i]];
311 for (i = 0; i < DIM; i++)
313 for (j = 0; j < DIM; j++)
315 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
321 void print_orires_log(FILE *log, t_oriresdata *od)
326 diagonalize_orires_tensors(od);
328 for (ex = 0; ex < od->nex; ex++)
330 eig = od->eig + ex*12;
331 fprintf(log, " Orientation experiment %d:\n", ex+1);
332 fprintf(log, " order parameter: %g\n", eig[0]);
333 for (i = 0; i < DIM; i++)
335 fprintf(log, " eig: %6.3f %6.3f %6.3f %6.3f\n",
336 (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
345 real calc_orires_dev(const gmx_multisim_t *ms,
346 int nfa, const t_iatom forceatoms[], const t_iparams ip[],
347 const t_mdatoms *md, const rvec x[], const t_pbc *pbc,
348 t_fcdata *fcd, history_t *hist)
350 int fa, d, i, j, type, ex, nref;
351 real edt, edt_1, invn, pfac, r2, invr, corrfac, weight, wsv2, sw, dev;
353 rvec5 *Dinsl, *Dins, *Dtav, *rhs;
356 rvec *xref, *xtmp, com, r_unrot, r;
359 const real two_thr = 2.0/3.0;
365 /* This means that this is not the master node */
366 gmx_fatal(FARGS, "Orientation restraints are only supported on the master node, use less processors");
369 bTAV = (od->edt != 0);
385 od->exp_min_t_tau = hist->orire_initf*edt;
387 /* Correction factor to correct for the lack of history
390 corrfac = 1.0/(1.0 - od->exp_min_t_tau);
409 for (i = 0; i < md->nr; i++)
411 if (md->cORF[i] == 0)
413 copy_rvec(x[i], xtmp[j]);
414 mref[j] = md->massT[i];
415 for (d = 0; d < DIM; d++)
417 com[d] += mref[j]*xref[j][d];
423 svmul(1.0/mtot, com, com);
424 for (j = 0; j < nref; j++)
426 rvec_dec(xtmp[j], com);
428 /* Calculate the rotation matrix to rotate x to the reference orientation */
429 calc_fit_R(DIM, nref, mref, xref, xtmp, R);
433 for (fa = 0; fa < nfa; fa += 3)
435 type = forceatoms[fa];
438 pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
442 rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
444 mvmul(R, r_unrot, r);
446 invr = gmx_invsqrt(r2);
447 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
448 pfac = ip[type].orires.c*invr*invr*3;
449 for (i = 0; i < ip[type].orires.power; i++)
453 Dinsl[d][0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
454 Dinsl[d][1] = pfac*(2*r[0]*r[1]);
455 Dinsl[d][2] = pfac*(2*r[0]*r[2]);
456 Dinsl[d][3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
457 Dinsl[d][4] = pfac*(2*r[1]*r[2]);
461 for (i = 0; i < 5; i++)
463 Dins[d][i] = Dinsl[d][i]*invn;
472 gmx_sum_sim(5*od->nr, Dins[0], ms);
475 /* Calculate the order tensor S for each experiment via optimization */
476 for (ex = 0; ex < od->nex; ex++)
478 for (i = 0; i < 5; i++)
481 for (j = 0; j <= i; j++)
488 for (fa = 0; fa < nfa; fa += 3)
492 /* Here we update Dtav in t_fcdata using the data in history_t.
493 * Thus the results stay correct when this routine
494 * is called multiple times.
496 for (i = 0; i < 5; i++)
498 Dtav[d][i] = edt*hist->orire_Dtav[d*5+i] + edt_1*Dins[d][i];
502 type = forceatoms[fa];
503 ex = ip[type].orires.ex;
504 weight = ip[type].orires.kfac;
505 /* Calculate the vector rhs and half the matrix T for the 5 equations */
506 for (i = 0; i < 5; i++)
508 rhs[ex][i] += Dtav[d][i]*ip[type].orires.obs*weight;
509 for (j = 0; j <= i; j++)
511 T[ex][i][j] += Dtav[d][i]*Dtav[d][j]*weight;
516 /* Now we have all the data we can calculate S */
517 for (ex = 0; ex < od->nex; ex++)
519 /* Correct corrfac and copy one half of T to the other half */
520 for (i = 0; i < 5; i++)
522 rhs[ex][i] *= corrfac;
523 T[ex][i][i] *= sqr(corrfac);
524 for (j = 0; j < i; j++)
526 T[ex][i][j] *= sqr(corrfac);
527 T[ex][j][i] = T[ex][i][j];
530 m_inv_gen(T[ex], 5, T[ex]);
531 /* Calculate the orientation tensor S for this experiment */
537 for (i = 0; i < 5; i++)
539 S[ex][0][0] += 1.5*T[ex][0][i]*rhs[ex][i];
540 S[ex][0][1] += 1.5*T[ex][1][i]*rhs[ex][i];
541 S[ex][0][2] += 1.5*T[ex][2][i]*rhs[ex][i];
542 S[ex][1][1] += 1.5*T[ex][3][i]*rhs[ex][i];
543 S[ex][1][2] += 1.5*T[ex][4][i]*rhs[ex][i];
545 S[ex][1][0] = S[ex][0][1];
546 S[ex][2][0] = S[ex][0][2];
547 S[ex][2][1] = S[ex][1][2];
548 S[ex][2][2] = -S[ex][0][0] - S[ex][1][1];
555 for (fa = 0; fa < nfa; fa += 3)
557 type = forceatoms[fa];
558 ex = ip[type].orires.ex;
560 od->otav[d] = two_thr*
561 corrfac*(S[ex][0][0]*Dtav[d][0] + S[ex][0][1]*Dtav[d][1] +
562 S[ex][0][2]*Dtav[d][2] + S[ex][1][1]*Dtav[d][3] +
563 S[ex][1][2]*Dtav[d][4]);
566 od->oins[d] = two_thr*(S[ex][0][0]*Dins[d][0] + S[ex][0][1]*Dins[d][1] +
567 S[ex][0][2]*Dins[d][2] + S[ex][1][1]*Dins[d][3] +
568 S[ex][1][2]*Dins[d][4]);
572 /* When ensemble averaging is used recalculate the local orientation
573 * for output to the energy file.
575 od->oinsl[d] = two_thr*
576 (S[ex][0][0]*Dinsl[d][0] + S[ex][0][1]*Dinsl[d][1] +
577 S[ex][0][2]*Dinsl[d][2] + S[ex][1][1]*Dinsl[d][3] +
578 S[ex][1][2]*Dinsl[d][4]);
581 dev = od->otav[d] - ip[type].orires.obs;
583 wsv2 += ip[type].orires.kfac*sqr(dev);
584 sw += ip[type].orires.kfac;
588 od->rmsdev = sqrt(wsv2/sw);
590 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
591 for (ex = 0; ex < od->nex; ex++)
593 tmmul(R, S[ex], TMP);
599 /* Approx. 120*nfa/3 flops */
602 real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
603 const rvec x[], rvec f[], rvec fshift[],
604 const t_pbc *pbc, const t_graph *g,
605 real gmx_unused lambda, real gmx_unused *dvdlambda,
606 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
607 int gmx_unused *global_atom_index)
610 int fa, d, i, type, ex, power, ki = CENTRAL;
612 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
615 const t_oriresdata *od;
623 bTAV = (od->edt != 0);
628 /* Smoothly switch on the restraining when time averaging is used */
629 smooth_fc *= (1.0 - od->exp_min_t_tau);
633 for (fa = 0; fa < nfa; fa += 3)
635 type = forceatoms[fa];
636 ai = forceatoms[fa+1];
637 aj = forceatoms[fa+2];
640 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
644 rvec_sub(x[ai], x[aj], r);
647 invr = gmx_invsqrt(r2);
649 ex = ip[type].orires.ex;
650 power = ip[type].orires.power;
651 fc = smooth_fc*ip[type].orires.kfac;
652 dev = od->otav[d] - ip[type].orires.obs;
655 * there is no real potential when time averaging is applied
657 vtot += 0.5*fc*sqr(dev);
661 /* Calculate the force as the sqrt of tav times instantaneous */
662 devins = od->oins[d] - ip[type].orires.obs;
669 dev = sqrt(dev*devins);
677 pfac = fc*ip[type].orires.c*invr2;
678 for (i = 0; i < power; i++)
682 mvmul(od->S[ex], r, Sr);
683 for (i = 0; i < DIM; i++)
686 -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
691 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
695 for (i = 0; i < DIM; i++)
699 fshift[ki][i] += fij[i];
700 fshift[CENTRAL][i] -= fij[i];
708 /* Approx. 80*nfa/3 flops */
711 void update_orires_history(t_fcdata *fcd, history_t *hist)
719 /* Copy the new time averages that have been calculated
720 * in calc_orires_dev.
722 hist->orire_initf = od->exp_min_t_tau;
723 for (pair = 0; pair < od->nr; pair++)
725 for (i = 0; i < 5; i++)
727 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];