1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
50 #include "mtop_util.h"
52 void init_orires(FILE *fplog, const gmx_mtop_t *mtop,
55 const gmx_multisim_t *ms, t_oriresdata *od,
58 int i, j, d, ex, nmol, nr, *nr_ex;
61 gmx_mtop_ilistloop_t iloop;
63 gmx_mtop_atomloop_all_t aloop;
66 od->fc = ir->orires_fc;
74 od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
82 iloop = gmx_mtop_ilistloop_init(mtop);
83 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
85 for (i = 0; i < il[F_ORIRES].nr; i += 3)
87 ex = mtop->ffparams.iparams[il[F_ORIRES].iatoms[i]].orires.ex;
91 for (j = od->nex; j < ex+1; j++)
100 snew(od->S, od->nex);
101 /* When not doing time averaging, the instaneous and time averaged data
102 * are indentical and the pointers can point to the same memory.
104 snew(od->Dinsl, od->nr);
107 snew(od->Dins, od->nr);
111 od->Dins = od->Dinsl;
114 if (ir->orires_tau == 0)
122 snew(od->Dtav, od->nr);
123 od->edt = exp(-ir->delta_t/ir->orires_tau);
124 od->edt_1 = 1.0 - od->edt;
126 /* Extend the state with the orires history */
127 state->flags |= (1<<estORIRE_INITF);
128 state->hist.orire_initf = 1;
129 state->flags |= (1<<estORIRE_DTAV);
130 state->hist.norire_Dtav = od->nr*5;
131 snew(state->hist.orire_Dtav, state->hist.norire_Dtav);
134 snew(od->oinsl, od->nr);
137 snew(od->oins, od->nr);
141 od->oins = od->oinsl;
143 if (ir->orires_tau == 0)
149 snew(od->otav, od->nr);
151 snew(od->tmp, od->nex);
152 snew(od->TMP, od->nex);
153 for (ex = 0; ex < od->nex; ex++)
155 snew(od->TMP[ex], 5);
156 for (i = 0; i < 5; i++)
158 snew(od->TMP[ex][i], 5);
163 for (i = 0; i < mtop->natoms; i++)
165 if (ggrpnr(&mtop->groups, egcORFIT, i) == 0)
170 snew(od->mref, od->nref);
171 snew(od->xref, od->nref);
172 snew(od->xtmp, od->nref);
174 snew(od->eig, od->nex*12);
176 /* Determine the reference structure on the master node.
177 * Copy it to the other nodes after checking multi compatibility,
178 * so we are sure the subsystems match before copying.
183 aloop = gmx_mtop_atomloop_all_init(mtop);
184 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
186 if (mtop->groups.grpnr[egcORFIT] == NULL ||
187 mtop->groups.grpnr[egcORFIT][i] == 0)
189 /* Not correct for free-energy with changing masses */
190 od->mref[j] = atom->m;
191 if (ms == NULL || MASTERSIM(ms))
193 copy_rvec(xref[i], od->xref[j]);
194 for (d = 0; d < DIM; d++)
196 com[d] += od->mref[j]*xref[i][d];
203 svmul(1.0/mtot, com, com);
204 if (ms == NULL || MASTERSIM(ms))
206 for (j = 0; j < od->nref; j++)
208 rvec_dec(od->xref[j], com);
212 fprintf(fplog, "Found %d orientation experiments\n", od->nex);
213 for (i = 0; i < od->nex; i++)
215 fprintf(fplog, " experiment %d has %d restraints\n", i+1, nr_ex[i]);
220 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n",
225 fprintf(fplog, " the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
227 check_multi_int(fplog, ms, od->nr,
228 "the number of orientation restraints",
230 check_multi_int(fplog, ms, od->nref,
231 "the number of fit atoms for orientation restraining",
233 check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
234 /* Copy the reference coordinates from the master to the other nodes */
235 gmx_sum_sim(DIM*od->nref, od->xref[0], ms);
238 please_cite(fplog, "Hess2003");
241 void diagonalize_orires_tensors(t_oriresdata *od)
243 int ex, i, j, nrot, ord[DIM], t;
249 for (i = 0; i < DIM; i++)
253 snew(od->eig_diag, DIM);
255 for (i = 0; i < DIM; i++)
261 for (ex = 0; ex < od->nex; ex++)
263 /* Rotate the S tensor back to the reference frame */
264 mmul(od->R, od->S[ex], TMP);
265 mtmul(TMP, od->R, S);
266 for (i = 0; i < DIM; i++)
268 for (j = 0; j < DIM; j++)
270 od->M[i][j] = S[i][j];
274 jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
276 for (i = 0; i < DIM; i++)
280 for (i = 0; i < DIM; i++)
282 for (j = i+1; j < DIM; j++)
284 if (sqr(od->eig_diag[ord[j]]) > sqr(od->eig_diag[ord[i]]))
293 for (i = 0; i < DIM; i++)
295 od->eig[ex*12 + i] = od->eig_diag[ord[i]];
297 for (i = 0; i < DIM; i++)
299 for (j = 0; j < DIM; j++)
301 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
307 void print_orires_log(FILE *log, t_oriresdata *od)
312 diagonalize_orires_tensors(od);
314 for (ex = 0; ex < od->nex; ex++)
316 eig = od->eig + ex*12;
317 fprintf(log, " Orientation experiment %d:\n", ex+1);
318 fprintf(log, " order parameter: %g\n", eig[0]);
319 for (i = 0; i < DIM; i++)
321 fprintf(log, " eig: %6.3f %6.3f %6.3f %6.3f\n",
322 (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
331 real calc_orires_dev(const gmx_multisim_t *ms,
332 int nfa, const t_iatom forceatoms[], const t_iparams ip[],
333 const t_mdatoms *md, const rvec x[], const t_pbc *pbc,
334 t_fcdata *fcd, history_t *hist)
336 int fa, d, i, j, type, ex, nref;
337 real edt, edt_1, invn, pfac, r2, invr, corrfac, weight, wsv2, sw, dev;
339 rvec5 *Dinsl, *Dins, *Dtav, *rhs;
342 rvec *xref, *xtmp, com, r_unrot, r;
345 const real two_thr = 2.0/3.0;
351 /* This means that this is not the master node */
352 gmx_fatal(FARGS, "Orientation restraints are only supported on the master node, use less processors");
355 bTAV = (od->edt != 0);
371 od->exp_min_t_tau = hist->orire_initf*edt;
373 /* Correction factor to correct for the lack of history
376 corrfac = 1.0/(1.0 - od->exp_min_t_tau);
395 for (i = 0; i < md->nr; i++)
397 if (md->cORF[i] == 0)
399 copy_rvec(x[i], xtmp[j]);
400 mref[j] = md->massT[i];
401 for (d = 0; d < DIM; d++)
403 com[d] += mref[j]*xref[j][d];
409 svmul(1.0/mtot, com, com);
410 for (j = 0; j < nref; j++)
412 rvec_dec(xtmp[j], com);
414 /* Calculate the rotation matrix to rotate x to the reference orientation */
415 calc_fit_R(DIM, nref, mref, xref, xtmp, R);
419 for (fa = 0; fa < nfa; fa += 3)
421 type = forceatoms[fa];
424 pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
428 rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
430 mvmul(R, r_unrot, r);
432 invr = gmx_invsqrt(r2);
433 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
434 pfac = ip[type].orires.c*invr*invr*3;
435 for (i = 0; i < ip[type].orires.power; i++)
439 Dinsl[d][0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
440 Dinsl[d][1] = pfac*(2*r[0]*r[1]);
441 Dinsl[d][2] = pfac*(2*r[0]*r[2]);
442 Dinsl[d][3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
443 Dinsl[d][4] = pfac*(2*r[1]*r[2]);
447 for (i = 0; i < 5; i++)
449 Dins[d][i] = Dinsl[d][i]*invn;
458 gmx_sum_sim(5*od->nr, Dins[0], ms);
461 /* Calculate the order tensor S for each experiment via optimization */
462 for (ex = 0; ex < od->nex; ex++)
464 for (i = 0; i < 5; i++)
467 for (j = 0; j <= i; j++)
474 for (fa = 0; fa < nfa; fa += 3)
478 /* Here we update Dtav in t_fcdata using the data in history_t.
479 * Thus the results stay correct when this routine
480 * is called multiple times.
482 for (i = 0; i < 5; i++)
484 Dtav[d][i] = edt*hist->orire_Dtav[d*5+i] + edt_1*Dins[d][i];
488 type = forceatoms[fa];
489 ex = ip[type].orires.ex;
490 weight = ip[type].orires.kfac;
491 /* Calculate the vector rhs and half the matrix T for the 5 equations */
492 for (i = 0; i < 5; i++)
494 rhs[ex][i] += Dtav[d][i]*ip[type].orires.obs*weight;
495 for (j = 0; j <= i; j++)
497 T[ex][i][j] += Dtav[d][i]*Dtav[d][j]*weight;
502 /* Now we have all the data we can calculate S */
503 for (ex = 0; ex < od->nex; ex++)
505 /* Correct corrfac and copy one half of T to the other half */
506 for (i = 0; i < 5; i++)
508 rhs[ex][i] *= corrfac;
509 T[ex][i][i] *= sqr(corrfac);
510 for (j = 0; j < i; j++)
512 T[ex][i][j] *= sqr(corrfac);
513 T[ex][j][i] = T[ex][i][j];
516 m_inv_gen(T[ex], 5, T[ex]);
517 /* Calculate the orientation tensor S for this experiment */
523 for (i = 0; i < 5; i++)
525 S[ex][0][0] += 1.5*T[ex][0][i]*rhs[ex][i];
526 S[ex][0][1] += 1.5*T[ex][1][i]*rhs[ex][i];
527 S[ex][0][2] += 1.5*T[ex][2][i]*rhs[ex][i];
528 S[ex][1][1] += 1.5*T[ex][3][i]*rhs[ex][i];
529 S[ex][1][2] += 1.5*T[ex][4][i]*rhs[ex][i];
531 S[ex][1][0] = S[ex][0][1];
532 S[ex][2][0] = S[ex][0][2];
533 S[ex][2][1] = S[ex][1][2];
534 S[ex][2][2] = -S[ex][0][0] - S[ex][1][1];
541 for (fa = 0; fa < nfa; fa += 3)
543 type = forceatoms[fa];
544 ex = ip[type].orires.ex;
546 od->otav[d] = two_thr*
547 corrfac*(S[ex][0][0]*Dtav[d][0] + S[ex][0][1]*Dtav[d][1] +
548 S[ex][0][2]*Dtav[d][2] + S[ex][1][1]*Dtav[d][3] +
549 S[ex][1][2]*Dtav[d][4]);
552 od->oins[d] = two_thr*(S[ex][0][0]*Dins[d][0] + S[ex][0][1]*Dins[d][1] +
553 S[ex][0][2]*Dins[d][2] + S[ex][1][1]*Dins[d][3] +
554 S[ex][1][2]*Dins[d][4]);
558 /* When ensemble averaging is used recalculate the local orientation
559 * for output to the energy file.
561 od->oinsl[d] = two_thr*
562 (S[ex][0][0]*Dinsl[d][0] + S[ex][0][1]*Dinsl[d][1] +
563 S[ex][0][2]*Dinsl[d][2] + S[ex][1][1]*Dinsl[d][3] +
564 S[ex][1][2]*Dinsl[d][4]);
567 dev = od->otav[d] - ip[type].orires.obs;
569 wsv2 += ip[type].orires.kfac*sqr(dev);
570 sw += ip[type].orires.kfac;
574 od->rmsdev = sqrt(wsv2/sw);
576 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
577 for (ex = 0; ex < od->nex; ex++)
579 tmmul(R, S[ex], TMP);
585 /* Approx. 120*nfa/3 flops */
588 real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
589 const rvec x[], rvec f[], rvec fshift[],
590 const t_pbc *pbc, const t_graph *g,
591 real gmx_unused lambda, real gmx_unused *dvdlambda,
592 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
593 int gmx_unused *global_atom_index)
596 int fa, d, i, type, ex, power, ki = CENTRAL;
598 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
601 const t_oriresdata *od;
609 bTAV = (od->edt != 0);
614 /* Smoothly switch on the restraining when time averaging is used */
615 smooth_fc *= (1.0 - od->exp_min_t_tau);
619 for (fa = 0; fa < nfa; fa += 3)
621 type = forceatoms[fa];
622 ai = forceatoms[fa+1];
623 aj = forceatoms[fa+2];
626 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
630 rvec_sub(x[ai], x[aj], r);
633 invr = gmx_invsqrt(r2);
635 ex = ip[type].orires.ex;
636 power = ip[type].orires.power;
637 fc = smooth_fc*ip[type].orires.kfac;
638 dev = od->otav[d] - ip[type].orires.obs;
641 * there is no real potential when time averaging is applied
643 vtot += 0.5*fc*sqr(dev);
647 /* Calculate the force as the sqrt of tav times instantaneous */
648 devins = od->oins[d] - ip[type].orires.obs;
655 dev = sqrt(dev*devins);
663 pfac = fc*ip[type].orires.c*invr2;
664 for (i = 0; i < power; i++)
668 mvmul(od->S[ex], r, Sr);
669 for (i = 0; i < DIM; i++)
672 -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
677 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
681 for (i = 0; i < DIM; i++)
685 fshift[ki][i] += fij[i];
686 fshift[CENTRAL][i] -= fij[i];
694 /* Approx. 80*nfa/3 flops */
697 void update_orires_history(t_fcdata *fcd, history_t *hist)
705 /* Copy the new time averages that have been calculated
706 * in calc_orires_dev.
708 hist->orire_initf = od->exp_min_t_tau;
709 for (pair = 0; pair < od->nr; pair++)
711 for (i = 0; i < 5; i++)
713 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];