Moved copyrite.* to fileio from gmxlib and legacyheaders.
[gromacs.git] / src / gromacs / gmxlib / orires.cpp
bloba249021710966de953aa9f7b0ebfe6fe9753bda8
1 /*
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,2015, 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.
37 #include "gmxpre.h"
39 #include "orires.h"
41 #include <cmath>
43 #include "gromacs/fileio/copyrite.h"
44 #include "gromacs/gmxlib/main.h"
45 #include "gromacs/legacyheaders/network.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/legacyheaders/types/fcdata.h"
48 #include "gromacs/legacyheaders/types/ifunc.h"
49 #include "gromacs/legacyheaders/types/mdatom.h"
50 #include "gromacs/legacyheaders/types/state.h"
51 #include "gromacs/linearalgebra/nrjac.h"
52 #include "gromacs/math/do_fit.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/ishift.h"
55 #include "gromacs/pbcutil/mshift.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
62 void init_orires(FILE *fplog, const gmx_mtop_t *mtop,
63 rvec xref[],
64 const t_inputrec *ir,
65 const t_commrec *cr, t_oriresdata *od,
66 t_state *state)
68 int i, j, d, ex, nmol, *nr_ex;
69 double mtot;
70 rvec com;
71 gmx_mtop_ilistloop_t iloop;
72 t_ilist *il;
73 gmx_mtop_atomloop_all_t aloop;
74 t_atom *atom;
75 const gmx_multisim_t *ms;
77 od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
78 if (0 == od->nr)
80 /* Not doing orientation restraints */
81 return;
84 if (DOMAINDECOMP(cr))
86 gmx_fatal(FARGS, "Orientation restraints do not work with more than one domain (ie. MPI rank).");
88 /* Orientation restraints */
89 if (!MASTER(cr))
91 /* Nothing to do */
92 return;
94 ms = cr->ms;
96 od->fc = ir->orires_fc;
97 od->nex = 0;
98 od->S = NULL;
99 od->M = NULL;
100 od->eig = NULL;
101 od->v = NULL;
103 nr_ex = NULL;
105 iloop = gmx_mtop_ilistloop_init(mtop);
106 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
108 for (i = 0; i < il[F_ORIRES].nr; i += 3)
110 ex = mtop->ffparams.iparams[il[F_ORIRES].iatoms[i]].orires.ex;
111 if (ex >= od->nex)
113 srenew(nr_ex, ex+1);
114 for (j = od->nex; j < ex+1; j++)
116 nr_ex[j] = 0;
118 od->nex = ex+1;
120 nr_ex[ex]++;
123 snew(od->S, od->nex);
124 /* When not doing time averaging, the instaneous and time averaged data
125 * are indentical and the pointers can point to the same memory.
127 snew(od->Dinsl, od->nr);
128 if (ms)
130 snew(od->Dins, od->nr);
132 else
134 od->Dins = od->Dinsl;
137 if (ir->orires_tau == 0)
139 od->Dtav = od->Dins;
140 od->edt = 0.0;
141 od->edt_1 = 1.0;
143 else
145 snew(od->Dtav, od->nr);
146 od->edt = std::exp(-ir->delta_t/ir->orires_tau);
147 od->edt_1 = 1.0 - od->edt;
149 /* Extend the state with the orires history */
150 state->flags |= (1<<estORIRE_INITF);
151 state->hist.orire_initf = 1;
152 state->flags |= (1<<estORIRE_DTAV);
153 state->hist.norire_Dtav = od->nr*5;
154 snew(state->hist.orire_Dtav, state->hist.norire_Dtav);
157 snew(od->oinsl, od->nr);
158 if (ms)
160 snew(od->oins, od->nr);
162 else
164 od->oins = od->oinsl;
166 if (ir->orires_tau == 0)
168 od->otav = od->oins;
170 else
172 snew(od->otav, od->nr);
174 snew(od->tmp, od->nex);
175 snew(od->TMP, od->nex);
176 for (ex = 0; ex < od->nex; ex++)
178 snew(od->TMP[ex], 5);
179 for (i = 0; i < 5; i++)
181 snew(od->TMP[ex][i], 5);
185 od->nref = 0;
186 for (i = 0; i < mtop->natoms; i++)
188 if (ggrpnr(&mtop->groups, egcORFIT, i) == 0)
190 od->nref++;
193 snew(od->mref, od->nref);
194 snew(od->xref, od->nref);
195 snew(od->xtmp, od->nref);
197 snew(od->eig, od->nex*12);
199 /* Determine the reference structure on the master node.
200 * Copy it to the other nodes after checking multi compatibility,
201 * so we are sure the subsystems match before copying.
203 clear_rvec(com);
204 mtot = 0.0;
205 j = 0;
206 aloop = gmx_mtop_atomloop_all_init(mtop);
207 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
209 if (mtop->groups.grpnr[egcORFIT] == NULL ||
210 mtop->groups.grpnr[egcORFIT][i] == 0)
212 /* Not correct for free-energy with changing masses */
213 od->mref[j] = atom->m;
214 if (ms == NULL || MASTERSIM(ms))
216 copy_rvec(xref[i], od->xref[j]);
217 for (d = 0; d < DIM; d++)
219 com[d] += od->mref[j]*xref[i][d];
222 mtot += od->mref[j];
223 j++;
226 svmul(1.0/mtot, com, com);
227 if (ms == NULL || MASTERSIM(ms))
229 for (j = 0; j < od->nref; j++)
231 rvec_dec(od->xref[j], com);
235 fprintf(fplog, "Found %d orientation experiments\n", od->nex);
236 for (i = 0; i < od->nex; i++)
238 fprintf(fplog, " experiment %d has %d restraints\n", i+1, nr_ex[i]);
241 sfree(nr_ex);
243 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n",
244 od->nref, mtot);
246 if (ms)
248 fprintf(fplog, " the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
250 check_multi_int(fplog, ms, od->nr,
251 "the number of orientation restraints",
252 FALSE);
253 check_multi_int(fplog, ms, od->nref,
254 "the number of fit atoms for orientation restraining",
255 FALSE);
256 check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
257 /* Copy the reference coordinates from the master to the other nodes */
258 gmx_sum_sim(DIM*od->nref, od->xref[0], ms);
261 please_cite(fplog, "Hess2003");
264 void diagonalize_orires_tensors(t_oriresdata *od)
266 int ex, i, j, nrot, ord[DIM], t;
267 matrix S, TMP;
269 if (od->M == NULL)
271 snew(od->M, DIM);
272 for (i = 0; i < DIM; i++)
274 snew(od->M[i], DIM);
276 snew(od->eig_diag, DIM);
277 snew(od->v, DIM);
278 for (i = 0; i < DIM; i++)
280 snew(od->v[i], DIM);
284 for (ex = 0; ex < od->nex; ex++)
286 /* Rotate the S tensor back to the reference frame */
287 mmul(od->R, od->S[ex], TMP);
288 mtmul(TMP, od->R, S);
289 for (i = 0; i < DIM; i++)
291 for (j = 0; j < DIM; j++)
293 od->M[i][j] = S[i][j];
297 jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
299 for (i = 0; i < DIM; i++)
301 ord[i] = i;
303 for (i = 0; i < DIM; i++)
305 for (j = i+1; j < DIM; j++)
307 if (sqr(od->eig_diag[ord[j]]) > sqr(od->eig_diag[ord[i]]))
309 t = ord[i];
310 ord[i] = ord[j];
311 ord[j] = t;
316 for (i = 0; i < DIM; i++)
318 od->eig[ex*12 + i] = od->eig_diag[ord[i]];
320 for (i = 0; i < DIM; i++)
322 for (j = 0; j < DIM; j++)
324 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
330 void print_orires_log(FILE *log, t_oriresdata *od)
332 int ex, i;
333 real *eig;
335 diagonalize_orires_tensors(od);
337 for (ex = 0; ex < od->nex; ex++)
339 eig = od->eig + ex*12;
340 fprintf(log, " Orientation experiment %d:\n", ex+1);
341 fprintf(log, " order parameter: %g\n", eig[0]);
342 for (i = 0; i < DIM; i++)
344 fprintf(log, " eig: %6.3f %6.3f %6.3f %6.3f\n",
345 (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
346 eig[DIM+i*DIM+XX],
347 eig[DIM+i*DIM+YY],
348 eig[DIM+i*DIM+ZZ]);
350 fprintf(log, "\n");
354 real calc_orires_dev(const gmx_multisim_t *ms,
355 int nfa, const t_iatom forceatoms[], const t_iparams ip[],
356 const t_mdatoms *md, const rvec x[], const t_pbc *pbc,
357 t_fcdata *fcd, history_t *hist)
359 int fa, d, i, j, type, ex, nref;
360 real edt, edt_1, invn, pfac, r2, invr, corrfac, weight, wsv2, sw, dev;
361 tensor *S, R, TMP;
362 rvec5 *Dinsl, *Dins, *Dtav, *rhs;
363 real *mref, ***T;
364 double mtot;
365 rvec *xref, *xtmp, com, r_unrot, r;
366 t_oriresdata *od;
367 gmx_bool bTAV;
368 const real two_thr = 2.0/3.0;
370 od = &(fcd->orires);
372 if (od->nr == 0)
374 /* This means that this is not the master node */
375 gmx_fatal(FARGS, "Orientation restraints are only supported on the master rank, use fewer ranks");
378 bTAV = (od->edt != 0);
379 edt = od->edt;
380 edt_1 = od->edt_1;
381 S = od->S;
382 Dinsl = od->Dinsl;
383 Dins = od->Dins;
384 Dtav = od->Dtav;
385 T = od->TMP;
386 rhs = od->tmp;
387 nref = od->nref;
388 mref = od->mref;
389 xref = od->xref;
390 xtmp = od->xtmp;
392 if (bTAV)
394 od->exp_min_t_tau = hist->orire_initf*edt;
396 /* Correction factor to correct for the lack of history
397 * at short times.
399 corrfac = 1.0/(1.0 - od->exp_min_t_tau);
401 else
403 corrfac = 1.0;
406 if (ms)
408 invn = 1.0/ms->nsim;
410 else
412 invn = 1.0;
415 clear_rvec(com);
416 mtot = 0;
417 j = 0;
418 for (i = 0; i < md->nr; i++)
420 if (md->cORF[i] == 0)
422 copy_rvec(x[i], xtmp[j]);
423 mref[j] = md->massT[i];
424 for (d = 0; d < DIM; d++)
426 com[d] += mref[j]*xref[j][d];
428 mtot += mref[j];
429 j++;
432 svmul(1.0/mtot, com, com);
433 for (j = 0; j < nref; j++)
435 rvec_dec(xtmp[j], com);
437 /* Calculate the rotation matrix to rotate x to the reference orientation */
438 calc_fit_R(DIM, nref, mref, xref, xtmp, R);
439 copy_mat(R, od->R);
441 d = 0;
442 for (fa = 0; fa < nfa; fa += 3)
444 type = forceatoms[fa];
445 if (pbc)
447 pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
449 else
451 rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
453 mvmul(R, r_unrot, r);
454 r2 = norm2(r);
455 invr = gmx_invsqrt(r2);
456 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
457 pfac = ip[type].orires.c*invr*invr*3;
458 for (i = 0; i < ip[type].orires.power; i++)
460 pfac *= invr;
462 Dinsl[d][0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
463 Dinsl[d][1] = pfac*(2*r[0]*r[1]);
464 Dinsl[d][2] = pfac*(2*r[0]*r[2]);
465 Dinsl[d][3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
466 Dinsl[d][4] = pfac*(2*r[1]*r[2]);
468 if (ms)
470 for (i = 0; i < 5; i++)
472 Dins[d][i] = Dinsl[d][i]*invn;
476 d++;
479 if (ms)
481 gmx_sum_sim(5*od->nr, Dins[0], ms);
484 /* Calculate the order tensor S for each experiment via optimization */
485 for (ex = 0; ex < od->nex; ex++)
487 for (i = 0; i < 5; i++)
489 rhs[ex][i] = 0;
490 for (j = 0; j <= i; j++)
492 T[ex][i][j] = 0;
496 d = 0;
497 for (fa = 0; fa < nfa; fa += 3)
499 if (bTAV)
501 /* Here we update Dtav in t_fcdata using the data in history_t.
502 * Thus the results stay correct when this routine
503 * is called multiple times.
505 for (i = 0; i < 5; i++)
507 Dtav[d][i] = edt*hist->orire_Dtav[d*5+i] + edt_1*Dins[d][i];
511 type = forceatoms[fa];
512 ex = ip[type].orires.ex;
513 weight = ip[type].orires.kfac;
514 /* Calculate the vector rhs and half the matrix T for the 5 equations */
515 for (i = 0; i < 5; i++)
517 rhs[ex][i] += Dtav[d][i]*ip[type].orires.obs*weight;
518 for (j = 0; j <= i; j++)
520 T[ex][i][j] += Dtav[d][i]*Dtav[d][j]*weight;
523 d++;
525 /* Now we have all the data we can calculate S */
526 for (ex = 0; ex < od->nex; ex++)
528 /* Correct corrfac and copy one half of T to the other half */
529 for (i = 0; i < 5; i++)
531 rhs[ex][i] *= corrfac;
532 T[ex][i][i] *= sqr(corrfac);
533 for (j = 0; j < i; j++)
535 T[ex][i][j] *= sqr(corrfac);
536 T[ex][j][i] = T[ex][i][j];
539 m_inv_gen(T[ex], 5, T[ex]);
540 /* Calculate the orientation tensor S for this experiment */
541 S[ex][0][0] = 0;
542 S[ex][0][1] = 0;
543 S[ex][0][2] = 0;
544 S[ex][1][1] = 0;
545 S[ex][1][2] = 0;
546 for (i = 0; i < 5; i++)
548 S[ex][0][0] += 1.5*T[ex][0][i]*rhs[ex][i];
549 S[ex][0][1] += 1.5*T[ex][1][i]*rhs[ex][i];
550 S[ex][0][2] += 1.5*T[ex][2][i]*rhs[ex][i];
551 S[ex][1][1] += 1.5*T[ex][3][i]*rhs[ex][i];
552 S[ex][1][2] += 1.5*T[ex][4][i]*rhs[ex][i];
554 S[ex][1][0] = S[ex][0][1];
555 S[ex][2][0] = S[ex][0][2];
556 S[ex][2][1] = S[ex][1][2];
557 S[ex][2][2] = -S[ex][0][0] - S[ex][1][1];
560 wsv2 = 0;
561 sw = 0;
563 d = 0;
564 for (fa = 0; fa < nfa; fa += 3)
566 type = forceatoms[fa];
567 ex = ip[type].orires.ex;
569 od->otav[d] = two_thr*
570 corrfac*(S[ex][0][0]*Dtav[d][0] + S[ex][0][1]*Dtav[d][1] +
571 S[ex][0][2]*Dtav[d][2] + S[ex][1][1]*Dtav[d][3] +
572 S[ex][1][2]*Dtav[d][4]);
573 if (bTAV)
575 od->oins[d] = two_thr*(S[ex][0][0]*Dins[d][0] + S[ex][0][1]*Dins[d][1] +
576 S[ex][0][2]*Dins[d][2] + S[ex][1][1]*Dins[d][3] +
577 S[ex][1][2]*Dins[d][4]);
579 if (ms)
581 /* When ensemble averaging is used recalculate the local orientation
582 * for output to the energy file.
584 od->oinsl[d] = two_thr*
585 (S[ex][0][0]*Dinsl[d][0] + S[ex][0][1]*Dinsl[d][1] +
586 S[ex][0][2]*Dinsl[d][2] + S[ex][1][1]*Dinsl[d][3] +
587 S[ex][1][2]*Dinsl[d][4]);
590 dev = od->otav[d] - ip[type].orires.obs;
592 wsv2 += ip[type].orires.kfac*sqr(dev);
593 sw += ip[type].orires.kfac;
595 d++;
597 od->rmsdev = std::sqrt(wsv2/sw);
599 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
600 for (ex = 0; ex < od->nex; ex++)
602 tmmul(R, S[ex], TMP);
603 mmul(TMP, R, S[ex]);
606 return od->rmsdev;
608 /* Approx. 120*nfa/3 flops */
611 real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
612 const rvec x[], rvec f[], rvec fshift[],
613 const t_pbc *pbc, const t_graph *g,
614 real gmx_unused lambda, real gmx_unused *dvdlambda,
615 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
616 int gmx_unused *global_atom_index)
618 atom_id ai, aj;
619 int fa, d, i, type, ex, power, ki = CENTRAL;
620 ivec dt;
621 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
622 rvec r, Sr, fij;
623 real vtot;
624 const t_oriresdata *od;
625 gmx_bool bTAV;
627 vtot = 0;
628 od = &(fcd->orires);
630 if (od->fc != 0)
632 bTAV = (od->edt != 0);
634 smooth_fc = od->fc;
635 if (bTAV)
637 /* Smoothly switch on the restraining when time averaging is used */
638 smooth_fc *= (1.0 - od->exp_min_t_tau);
641 d = 0;
642 for (fa = 0; fa < nfa; fa += 3)
644 type = forceatoms[fa];
645 ai = forceatoms[fa+1];
646 aj = forceatoms[fa+2];
647 if (pbc)
649 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
651 else
653 rvec_sub(x[ai], x[aj], r);
655 r2 = norm2(r);
656 invr = gmx_invsqrt(r2);
657 invr2 = invr*invr;
658 ex = ip[type].orires.ex;
659 power = ip[type].orires.power;
660 fc = smooth_fc*ip[type].orires.kfac;
661 dev = od->otav[d] - ip[type].orires.obs;
663 /* NOTE:
664 * there is no real potential when time averaging is applied
666 vtot += 0.5*fc*sqr(dev);
668 if (bTAV)
670 /* Calculate the force as the sqrt of tav times instantaneous */
671 devins = od->oins[d] - ip[type].orires.obs;
672 if (dev*devins <= 0)
674 dev = 0;
676 else
678 dev = std::sqrt(dev*devins);
679 if (devins < 0)
681 dev = -dev;
686 pfac = fc*ip[type].orires.c*invr2;
687 for (i = 0; i < power; i++)
689 pfac *= invr;
691 mvmul(od->S[ex], r, Sr);
692 for (i = 0; i < DIM; i++)
694 fij[i] =
695 -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
698 if (g)
700 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
701 ki = IVEC2IS(dt);
704 for (i = 0; i < DIM; i++)
706 f[ai][i] += fij[i];
707 f[aj][i] -= fij[i];
708 fshift[ki][i] += fij[i];
709 fshift[CENTRAL][i] -= fij[i];
711 d++;
715 return vtot;
717 /* Approx. 80*nfa/3 flops */
720 void update_orires_history(t_fcdata *fcd, history_t *hist)
722 t_oriresdata *od;
723 int pair, i;
725 od = &(fcd->orires);
726 if (od->edt != 0)
728 /* Copy the new time averages that have been calculated
729 * in calc_orires_dev.
731 hist->orire_initf = od->exp_min_t_tau;
732 for (pair = 0; pair < od->nr; pair++)
734 for (i = 0; i < 5; i++)
736 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];