Convert gmx_mtop_t to C++
[gromacs.git] / src / gromacs / listed-forces / orires.cpp
blobb91d5d001a4e5321dc8a7d2ee47d9d6c1fe9eca5
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,2016,2017,2018, 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 <climits>
42 #include <cmath>
44 #include "gromacs/gmxlib/network.h"
45 #include "gromacs/linearalgebra/nrjac.h"
46 #include "gromacs/math/do_fit.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/main.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/mdtypes/fcdata.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/mdtypes/state.h"
55 #include "gromacs/pbcutil/ishift.h"
56 #include "gromacs/pbcutil/mshift.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/topology/topology.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/pleasecite.h"
63 #include "gromacs/utility/smalloc.h"
65 // TODO This implementation of ensemble orientation restraints is nasty because
66 // a user can't just do multi-sim with single-sim orientation restraints.
68 void init_orires(FILE *fplog,
69 const gmx_mtop_t *mtop,
70 const t_inputrec *ir,
71 const t_commrec *cr,
72 const gmx_multisim_t *ms,
73 t_state *globalState,
74 t_oriresdata *od)
76 od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
77 if (0 == od->nr)
79 /* Not doing orientation restraints */
80 return;
83 const int numFitParams = 5;
84 if (od->nr <= numFitParams)
86 gmx_fatal(FARGS, "The system has %d orientation restraints, but at least %d are required, since there are %d fitting parameters.",
87 od->nr, numFitParams + 1, numFitParams);
90 if (ir->bPeriodicMols)
92 /* Since we apply fitting, we need to make molecules whole and this
93 * can not be done when periodic molecules are present.
95 gmx_fatal(FARGS, "Orientation restraints can not be applied when periodic molecules are present in the system");
98 if (PAR(cr))
100 gmx_fatal(FARGS, "Orientation restraints do not work with MPI parallelization. Choose 1 MPI rank, if possible.");
103 GMX_RELEASE_ASSERT(globalState != nullptr, "We need a valid global state in init_orires");
105 od->fc = ir->orires_fc;
106 od->nex = 0;
107 od->S = nullptr;
108 od->M = nullptr;
109 od->eig = nullptr;
110 od->v = nullptr;
112 int *nr_ex = nullptr;
113 int typeMin = INT_MAX;
114 int typeMax = 0;
115 gmx_mtop_ilistloop_t iloop = gmx_mtop_ilistloop_init(mtop);
116 const t_ilist *il;
117 int nmol;
118 while (gmx_mtop_ilistloop_next(iloop, &il, &nmol))
120 if (nmol > 1)
122 gmx_fatal(FARGS, "Found %d copies of a molecule with orientation restrains while the current code only supports a single copy. If you want to ensemble average, run multiple copies of the system using the multi-sim feature of mdrun.", nmol);
125 for (int i = 0; i < il[F_ORIRES].nr; i += 3)
127 int type = il[F_ORIRES].iatoms[i];
128 int ex = mtop->ffparams.iparams[type].orires.ex;
129 if (ex >= od->nex)
131 srenew(nr_ex, ex+1);
132 for (int j = od->nex; j < ex+1; j++)
134 nr_ex[j] = 0;
136 od->nex = ex+1;
138 GMX_ASSERT(nr_ex, "Check for allocated nr_ex to keep the static analyzer happy");
139 nr_ex[ex]++;
141 typeMin = std::min(typeMin, type);
142 typeMax = std::max(typeMax, type);
145 /* With domain decomposition we use the type index for indexing in global arrays */
146 GMX_RELEASE_ASSERT(typeMax - typeMin + 1 == od->nr, "All orientation restraint parameter entries in the topology should be consecutive");
147 /* Store typeMin so we can index array with the type offset */
148 od->typeMin = typeMin;
150 snew(od->S, od->nex);
151 /* When not doing time averaging, the instaneous and time averaged data
152 * are indentical and the pointers can point to the same memory.
154 snew(od->Dinsl, od->nr);
156 if (ms)
158 snew(od->Dins, od->nr);
160 else
162 od->Dins = od->Dinsl;
165 if (ir->orires_tau == 0)
167 od->Dtav = od->Dins;
168 od->edt = 0.0;
169 od->edt_1 = 1.0;
171 else
173 snew(od->Dtav, od->nr);
174 od->edt = std::exp(-ir->delta_t/ir->orires_tau);
175 od->edt_1 = 1.0 - od->edt;
177 /* Extend the state with the orires history */
178 globalState->flags |= (1<<estORIRE_INITF);
179 globalState->hist.orire_initf = 1;
180 globalState->flags |= (1<<estORIRE_DTAV);
181 globalState->hist.norire_Dtav = od->nr*5;
182 snew(globalState->hist.orire_Dtav, globalState->hist.norire_Dtav);
185 snew(od->oinsl, od->nr);
186 if (ms)
188 snew(od->oins, od->nr);
190 else
192 od->oins = od->oinsl;
194 if (ir->orires_tau == 0)
196 od->otav = od->oins;
198 else
200 snew(od->otav, od->nr);
202 snew(od->tmpEq, od->nex);
204 od->nref = 0;
205 for (int i = 0; i < mtop->natoms; i++)
207 if (ggrpnr(&mtop->groups, egcORFIT, i) == 0)
209 od->nref++;
212 snew(od->mref, od->nref);
213 snew(od->xref, od->nref);
214 snew(od->xtmp, od->nref);
216 snew(od->eig, od->nex*12);
218 /* Determine the reference structure on the master node.
219 * Copy it to the other nodes after checking multi compatibility,
220 * so we are sure the subsystems match before copying.
222 rvec com = { 0, 0, 0 };
223 double mtot = 0.0;
224 int j = 0;
225 gmx_mtop_atomloop_all_t aloop = gmx_mtop_atomloop_all_init(mtop);
226 int i = -1;
227 const t_atom *atom;
228 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
230 if (mtop->groups.grpnr[egcORFIT] == nullptr ||
231 mtop->groups.grpnr[egcORFIT][i] == 0)
233 /* Not correct for free-energy with changing masses */
234 od->mref[j] = atom->m;
235 // Note that only one rank per sim is supported.
236 if (isMasterSim(ms))
238 copy_rvec(globalState->x[i], od->xref[j]);
239 for (int d = 0; d < DIM; d++)
241 com[d] += od->mref[j]*globalState->x[i][d];
244 mtot += od->mref[j];
245 j++;
248 svmul(1.0/mtot, com, com);
249 if (isMasterSim(ms))
251 for (int j = 0; j < od->nref; j++)
253 rvec_dec(od->xref[j], com);
257 fprintf(fplog, "Found %d orientation experiments\n", od->nex);
258 for (int i = 0; i < od->nex; i++)
260 fprintf(fplog, " experiment %d has %d restraints\n", i+1, nr_ex[i]);
263 sfree(nr_ex);
265 fprintf(fplog, " the fit group consists of %d atoms and has total mass %g\n",
266 od->nref, mtot);
268 if (ms)
270 fprintf(fplog, " the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
272 check_multi_int(fplog, ms, od->nr,
273 "the number of orientation restraints",
274 FALSE);
275 check_multi_int(fplog, ms, od->nref,
276 "the number of fit atoms for orientation restraining",
277 FALSE);
278 check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
279 /* Copy the reference coordinates from the master to the other nodes */
280 gmx_sum_sim(DIM*od->nref, od->xref[0], ms);
283 please_cite(fplog, "Hess2003");
286 void diagonalize_orires_tensors(t_oriresdata *od)
288 if (od->M == nullptr)
290 snew(od->M, DIM);
291 for (int i = 0; i < DIM; i++)
293 snew(od->M[i], DIM);
295 snew(od->eig_diag, DIM);
296 snew(od->v, DIM);
297 for (int i = 0; i < DIM; i++)
299 snew(od->v[i], DIM);
303 for (int ex = 0; ex < od->nex; ex++)
305 /* Rotate the S tensor back to the reference frame */
306 matrix S, TMP;
307 mmul(od->R, od->S[ex], TMP);
308 mtmul(TMP, od->R, S);
309 for (int i = 0; i < DIM; i++)
311 for (int j = 0; j < DIM; j++)
313 od->M[i][j] = S[i][j];
317 int nrot;
318 jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
320 int ord[DIM];
321 for (int i = 0; i < DIM; i++)
323 ord[i] = i;
325 for (int i = 0; i < DIM; i++)
327 for (int j = i+1; j < DIM; j++)
329 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
331 int t = ord[i];
332 ord[i] = ord[j];
333 ord[j] = t;
338 for (int i = 0; i < DIM; i++)
340 od->eig[ex*12 + i] = od->eig_diag[ord[i]];
342 for (int i = 0; i < DIM; i++)
344 for (int j = 0; j < DIM; j++)
346 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
352 void print_orires_log(FILE *log, t_oriresdata *od)
354 real *eig;
356 diagonalize_orires_tensors(od);
358 for (int ex = 0; ex < od->nex; ex++)
360 eig = od->eig + ex*12;
361 fprintf(log, " Orientation experiment %d:\n", ex+1);
362 fprintf(log, " order parameter: %g\n", eig[0]);
363 for (int i = 0; i < DIM; i++)
365 fprintf(log, " eig: %6.3f %6.3f %6.3f %6.3f\n",
366 (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
367 eig[DIM+i*DIM+XX],
368 eig[DIM+i*DIM+YY],
369 eig[DIM+i*DIM+ZZ]);
371 fprintf(log, "\n");
375 real calc_orires_dev(const gmx_multisim_t *ms,
376 int nfa, const t_iatom forceatoms[], const t_iparams ip[],
377 const t_mdatoms *md, const rvec x[], const t_pbc *pbc,
378 t_fcdata *fcd, history_t *hist)
380 int nref;
381 real edt, edt_1, invn, pfac, r2, invr, corrfac, wsv2, sw, dev;
382 OriresMatEq *matEq;
383 real *mref;
384 double mtot;
385 rvec *xref, *xtmp, com, r_unrot, r;
386 t_oriresdata *od;
387 gmx_bool bTAV;
388 const real two_thr = 2.0/3.0;
390 od = &(fcd->orires);
392 if (od->nr == 0)
394 /* This means that this is not the master node */
395 gmx_fatal(FARGS, "Orientation restraints are only supported on the master rank, use fewer ranks");
398 bTAV = (od->edt != 0);
399 edt = od->edt;
400 edt_1 = od->edt_1;
401 matEq = od->tmpEq;
402 nref = od->nref;
403 mref = od->mref;
404 xref = od->xref;
405 xtmp = od->xtmp;
407 if (bTAV)
409 od->exp_min_t_tau = hist->orire_initf*edt;
411 /* Correction factor to correct for the lack of history
412 * at short times.
414 corrfac = 1.0/(1.0 - od->exp_min_t_tau);
416 else
418 corrfac = 1.0;
421 if (ms)
423 invn = 1.0/ms->nsim;
425 else
427 invn = 1.0;
430 clear_rvec(com);
431 mtot = 0;
432 int j = 0;
433 for (int i = 0; i < md->nr; i++)
435 if (md->cORF[i] == 0)
437 copy_rvec(x[i], xtmp[j]);
438 mref[j] = md->massT[i];
439 for (int d = 0; d < DIM; d++)
441 com[d] += mref[j]*xtmp[j][d];
443 mtot += mref[j];
444 j++;
447 svmul(1.0/mtot, com, com);
448 for (int j = 0; j < nref; j++)
450 rvec_dec(xtmp[j], com);
452 /* Calculate the rotation matrix to rotate x to the reference orientation */
453 calc_fit_R(DIM, nref, mref, xref, xtmp, od->R);
455 for (int fa = 0; fa < nfa; fa += 3)
457 const int type = forceatoms[fa];
458 const int restraintIndex = type - od->typeMin;
459 if (pbc)
461 pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
463 else
465 rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
467 mvmul(od->R, r_unrot, r);
468 r2 = norm2(r);
469 invr = gmx::invsqrt(r2);
470 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
471 pfac = ip[type].orires.c*invr*invr*3;
472 for (int i = 0; i < ip[type].orires.power; i++)
474 pfac *= invr;
476 rvec5 &Dinsl = od->Dinsl[restraintIndex];
477 Dinsl[0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
478 Dinsl[1] = pfac*(2*r[0]*r[1]);
479 Dinsl[2] = pfac*(2*r[0]*r[2]);
480 Dinsl[3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
481 Dinsl[4] = pfac*(2*r[1]*r[2]);
483 if (ms)
485 for (int i = 0; i < 5; i++)
487 od->Dins[restraintIndex][i] = Dinsl[i]*invn;
492 if (ms)
494 gmx_sum_sim(5*od->nr, od->Dins[0], ms);
497 /* Calculate the order tensor S for each experiment via optimization */
498 for (int ex = 0; ex < od->nex; ex++)
500 for (int i = 0; i < 5; i++)
502 matEq[ex].rhs[i] = 0;
503 for (int j = 0; j <= i; j++)
505 matEq[ex].mat[i][j] = 0;
510 for (int fa = 0; fa < nfa; fa += 3)
512 const int type = forceatoms[fa];
513 const int restraintIndex = type - od->typeMin;
514 rvec5 &Dtav = od->Dtav[restraintIndex];
515 if (bTAV)
517 /* Here we update Dtav in t_fcdata using the data in history_t.
518 * Thus the results stay correct when this routine
519 * is called multiple times.
521 for (int i = 0; i < 5; i++)
523 Dtav[i] =
524 edt*hist->orire_Dtav[restraintIndex*5 + i] +
525 edt_1*od->Dins[restraintIndex][i];
529 int ex = ip[type].orires.ex;
530 real weight = ip[type].orires.kfac;
531 /* Calculate the vector rhs and half the matrix T for the 5 equations */
532 for (int i = 0; i < 5; i++)
534 matEq[ex].rhs[i] += Dtav[i]*ip[type].orires.obs*weight;
535 for (int j = 0; j <= i; j++)
537 matEq[ex].mat[i][j] += Dtav[i]*Dtav[j]*weight;
541 /* Now we have all the data we can calculate S */
542 for (int ex = 0; ex < od->nex; ex++)
544 OriresMatEq &eq = matEq[ex];
545 /* Correct corrfac and copy one half of T to the other half */
546 for (int i = 0; i < 5; i++)
548 eq.rhs[i] *= corrfac;
549 eq.mat[i][i] *= gmx::square(corrfac);
550 for (int j = 0; j < i; j++)
552 eq.mat[i][j] *= gmx::square(corrfac);
553 eq.mat[j][i] = eq.mat[i][j];
556 m_inv_gen(&eq.mat[0][0], 5, &eq.mat[0][0]);
557 /* Calculate the orientation tensor S for this experiment */
558 matrix &S = od->S[ex];
559 S[0][0] = 0;
560 S[0][1] = 0;
561 S[0][2] = 0;
562 S[1][1] = 0;
563 S[1][2] = 0;
564 for (int i = 0; i < 5; i++)
566 S[0][0] += 1.5*eq.mat[0][i]*eq.rhs[i];
567 S[0][1] += 1.5*eq.mat[1][i]*eq.rhs[i];
568 S[0][2] += 1.5*eq.mat[2][i]*eq.rhs[i];
569 S[1][1] += 1.5*eq.mat[3][i]*eq.rhs[i];
570 S[1][2] += 1.5*eq.mat[4][i]*eq.rhs[i];
572 S[1][0] = S[0][1];
573 S[2][0] = S[0][2];
574 S[2][1] = S[1][2];
575 S[2][2] = -S[0][0] - S[1][1];
578 const matrix *S = od->S;
580 wsv2 = 0;
581 sw = 0;
583 for (int fa = 0; fa < nfa; fa += 3)
585 const int type = forceatoms[fa];
586 const int restraintIndex = type - od->typeMin;
587 const int ex = ip[type].orires.ex;
589 const rvec5 &Dtav = od->Dtav[restraintIndex];
590 od->otav[restraintIndex] = two_thr*
591 corrfac*(S[ex][0][0]*Dtav[0] + S[ex][0][1]*Dtav[1] +
592 S[ex][0][2]*Dtav[2] + S[ex][1][1]*Dtav[3] +
593 S[ex][1][2]*Dtav[4]);
594 if (bTAV)
596 const rvec5 &Dins = od->Dins[restraintIndex];
597 od->oins[restraintIndex] = two_thr*
598 (S[ex][0][0]*Dins[0] + S[ex][0][1]*Dins[1] +
599 S[ex][0][2]*Dins[2] + S[ex][1][1]*Dins[3] +
600 S[ex][1][2]*Dins[4]);
602 if (ms)
604 /* When ensemble averaging is used recalculate the local orientation
605 * for output to the energy file.
607 const rvec5 &Dinsl = od->Dinsl[restraintIndex];
608 od->oinsl[restraintIndex] = two_thr*
609 (S[ex][0][0]*Dinsl[0] + S[ex][0][1]*Dinsl[1] +
610 S[ex][0][2]*Dinsl[2] + S[ex][1][1]*Dinsl[3] +
611 S[ex][1][2]*Dinsl[4]);
614 dev = od->otav[restraintIndex] - ip[type].orires.obs;
616 wsv2 += ip[type].orires.kfac*gmx::square(dev);
617 sw += ip[type].orires.kfac;
619 od->rmsdev = std::sqrt(wsv2/sw);
621 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
622 for (int ex = 0; ex < od->nex; ex++)
624 matrix RS;
625 tmmul(od->R, od->S[ex], RS);
626 mmul(RS, od->R, od->S[ex]);
629 return od->rmsdev;
631 /* Approx. 120*nfa/3 flops */
634 real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
635 const rvec x[], rvec4 f[], rvec fshift[],
636 const t_pbc *pbc, const t_graph *g,
637 real gmx_unused lambda, real gmx_unused *dvdlambda,
638 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
639 int gmx_unused *global_atom_index)
641 int ex, power, ki = CENTRAL;
642 ivec dt;
643 real r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
644 rvec r, Sr, fij;
645 real vtot;
646 const t_oriresdata *od;
647 gmx_bool bTAV;
649 vtot = 0;
650 od = &(fcd->orires);
652 if (od->fc != 0)
654 bTAV = (od->edt != 0);
656 smooth_fc = od->fc;
657 if (bTAV)
659 /* Smoothly switch on the restraining when time averaging is used */
660 smooth_fc *= (1.0 - od->exp_min_t_tau);
663 for (int fa = 0; fa < nfa; fa += 3)
665 const int type = forceatoms[fa];
666 const int ai = forceatoms[fa + 1];
667 const int aj = forceatoms[fa + 2];
668 const int restraintIndex = type - od->typeMin;
669 if (pbc)
671 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
673 else
675 rvec_sub(x[ai], x[aj], r);
677 r2 = norm2(r);
678 invr = gmx::invsqrt(r2);
679 invr2 = invr*invr;
680 ex = ip[type].orires.ex;
681 power = ip[type].orires.power;
682 fc = smooth_fc*ip[type].orires.kfac;
683 dev = od->otav[restraintIndex] - ip[type].orires.obs;
685 /* NOTE:
686 * there is no real potential when time averaging is applied
688 vtot += 0.5*fc*gmx::square(dev);
690 if (bTAV)
692 /* Calculate the force as the sqrt of tav times instantaneous */
693 devins = od->oins[restraintIndex] - ip[type].orires.obs;
694 if (dev*devins <= 0)
696 dev = 0;
698 else
700 dev = std::sqrt(dev*devins);
701 if (devins < 0)
703 dev = -dev;
708 pfac = fc*ip[type].orires.c*invr2;
709 for (int i = 0; i < power; i++)
711 pfac *= invr;
713 mvmul(od->S[ex], r, Sr);
714 for (int i = 0; i < DIM; i++)
716 fij[i] =
717 -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
720 if (g)
722 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
723 ki = IVEC2IS(dt);
726 for (int i = 0; i < DIM; i++)
728 f[ai][i] += fij[i];
729 f[aj][i] -= fij[i];
730 fshift[ki][i] += fij[i];
731 fshift[CENTRAL][i] -= fij[i];
736 return vtot;
738 /* Approx. 80*nfa/3 flops */
741 void update_orires_history(t_fcdata *fcd, history_t *hist)
743 t_oriresdata *od = &(fcd->orires);
745 if (od->edt != 0)
747 /* Copy the new time averages that have been calculated
748 * in calc_orires_dev.
750 hist->orire_initf = od->exp_min_t_tau;
751 for (int pair = 0; pair < od->nr; pair++)
753 for (int i = 0; i < 5; i++)
755 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];