Convert gmx_mtop_t to C++
[alexxy/gromacs.git] / src / gromacs / listed-forces / orires.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "orires.h"
40
41 #include <climits>
42 #include <cmath>
43
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"
64
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.
67
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)
75 {
76     od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
77     if (0 == od->nr)
78     {
79         /* Not doing orientation restraints */
80         return;
81     }
82
83     const int numFitParams = 5;
84     if (od->nr <= numFitParams)
85     {
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);
88     }
89
90     if (ir->bPeriodicMols)
91     {
92         /* Since we apply fitting, we need to make molecules whole and this
93          * can not be done when periodic molecules are present.
94          */
95         gmx_fatal(FARGS, "Orientation restraints can not be applied when periodic molecules are present in the system");
96     }
97
98     if (PAR(cr))
99     {
100         gmx_fatal(FARGS, "Orientation restraints do not work with MPI parallelization. Choose 1 MPI rank, if possible.");
101     }
102
103     GMX_RELEASE_ASSERT(globalState != nullptr, "We need a valid global state in init_orires");
104
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;
111
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))
119     {
120         if (nmol > 1)
121         {
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);
123         }
124
125         for (int i = 0; i < il[F_ORIRES].nr; i += 3)
126         {
127             int type = il[F_ORIRES].iatoms[i];
128             int ex   = mtop->ffparams.iparams[type].orires.ex;
129             if (ex >= od->nex)
130             {
131                 srenew(nr_ex, ex+1);
132                 for (int j = od->nex; j < ex+1; j++)
133                 {
134                     nr_ex[j] = 0;
135                 }
136                 od->nex = ex+1;
137             }
138             GMX_ASSERT(nr_ex, "Check for allocated nr_ex to keep the static analyzer happy");
139             nr_ex[ex]++;
140
141             typeMin = std::min(typeMin, type);
142             typeMax = std::max(typeMax, type);
143         }
144     }
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;
149
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.
153      */
154     snew(od->Dinsl, od->nr);
155
156     if (ms)
157     {
158         snew(od->Dins, od->nr);
159     }
160     else
161     {
162         od->Dins = od->Dinsl;
163     }
164
165     if (ir->orires_tau == 0)
166     {
167         od->Dtav  = od->Dins;
168         od->edt   = 0.0;
169         od->edt_1 = 1.0;
170     }
171     else
172     {
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;
176
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);
183     }
184
185     snew(od->oinsl, od->nr);
186     if (ms)
187     {
188         snew(od->oins, od->nr);
189     }
190     else
191     {
192         od->oins = od->oinsl;
193     }
194     if (ir->orires_tau == 0)
195     {
196         od->otav = od->oins;
197     }
198     else
199     {
200         snew(od->otav, od->nr);
201     }
202     snew(od->tmpEq, od->nex);
203
204     od->nref = 0;
205     for (int i = 0; i < mtop->natoms; i++)
206     {
207         if (ggrpnr(&mtop->groups, egcORFIT, i) == 0)
208         {
209             od->nref++;
210         }
211     }
212     snew(od->mref, od->nref);
213     snew(od->xref, od->nref);
214     snew(od->xtmp, od->nref);
215
216     snew(od->eig, od->nex*12);
217
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.
221      */
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))
229     {
230         if (mtop->groups.grpnr[egcORFIT] == nullptr ||
231             mtop->groups.grpnr[egcORFIT][i] == 0)
232         {
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))
237             {
238                 copy_rvec(globalState->x[i], od->xref[j]);
239                 for (int d = 0; d < DIM; d++)
240                 {
241                     com[d] += od->mref[j]*globalState->x[i][d];
242                 }
243             }
244             mtot += od->mref[j];
245             j++;
246         }
247     }
248     svmul(1.0/mtot, com, com);
249     if (isMasterSim(ms))
250     {
251         for (int j = 0; j < od->nref; j++)
252         {
253             rvec_dec(od->xref[j], com);
254         }
255     }
256
257     fprintf(fplog, "Found %d orientation experiments\n", od->nex);
258     for (int i = 0; i < od->nex; i++)
259     {
260         fprintf(fplog, "  experiment %d has %d restraints\n", i+1, nr_ex[i]);
261     }
262
263     sfree(nr_ex);
264
265     fprintf(fplog, "  the fit group consists of %d atoms and has total mass %g\n",
266             od->nref, mtot);
267
268     if (ms)
269     {
270         fprintf(fplog, "  the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
271
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);
281     }
282
283     please_cite(fplog, "Hess2003");
284 }
285
286 void diagonalize_orires_tensors(t_oriresdata *od)
287 {
288     if (od->M == nullptr)
289     {
290         snew(od->M, DIM);
291         for (int i = 0; i < DIM; i++)
292         {
293             snew(od->M[i], DIM);
294         }
295         snew(od->eig_diag, DIM);
296         snew(od->v, DIM);
297         for (int i = 0; i < DIM; i++)
298         {
299             snew(od->v[i], DIM);
300         }
301     }
302
303     for (int ex = 0; ex < od->nex; ex++)
304     {
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++)
310         {
311             for (int j = 0; j < DIM; j++)
312             {
313                 od->M[i][j] = S[i][j];
314             }
315         }
316
317         int nrot;
318         jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
319
320         int ord[DIM];
321         for (int i = 0; i < DIM; i++)
322         {
323             ord[i] = i;
324         }
325         for (int i = 0; i < DIM; i++)
326         {
327             for (int j = i+1; j < DIM; j++)
328             {
329                 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
330                 {
331                     int t  = ord[i];
332                     ord[i] = ord[j];
333                     ord[j] = t;
334                 }
335             }
336         }
337
338         for (int i = 0; i < DIM; i++)
339         {
340             od->eig[ex*12 + i] = od->eig_diag[ord[i]];
341         }
342         for (int i = 0; i < DIM; i++)
343         {
344             for (int j = 0; j < DIM; j++)
345             {
346                 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
347             }
348         }
349     }
350 }
351
352 void print_orires_log(FILE *log, t_oriresdata *od)
353 {
354     real *eig;
355
356     diagonalize_orires_tensors(od);
357
358     for (int ex = 0; ex < od->nex; ex++)
359     {
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++)
364         {
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]);
370         }
371         fprintf(log, "\n");
372     }
373 }
374
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)
379 {
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;
389
390     od = &(fcd->orires);
391
392     if (od->nr == 0)
393     {
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");
396     }
397
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;
406
407     if (bTAV)
408     {
409         od->exp_min_t_tau = hist->orire_initf*edt;
410
411         /* Correction factor to correct for the lack of history
412          * at short times.
413          */
414         corrfac = 1.0/(1.0 - od->exp_min_t_tau);
415     }
416     else
417     {
418         corrfac = 1.0;
419     }
420
421     if (ms)
422     {
423         invn = 1.0/ms->nsim;
424     }
425     else
426     {
427         invn = 1.0;
428     }
429
430     clear_rvec(com);
431     mtot  = 0;
432     int j = 0;
433     for (int i = 0; i < md->nr; i++)
434     {
435         if (md->cORF[i] == 0)
436         {
437             copy_rvec(x[i], xtmp[j]);
438             mref[j] = md->massT[i];
439             for (int d = 0; d < DIM; d++)
440             {
441                 com[d] += mref[j]*xtmp[j][d];
442             }
443             mtot += mref[j];
444             j++;
445         }
446     }
447     svmul(1.0/mtot, com, com);
448     for (int j = 0; j < nref; j++)
449     {
450         rvec_dec(xtmp[j], com);
451     }
452     /* Calculate the rotation matrix to rotate x to the reference orientation */
453     calc_fit_R(DIM, nref, mref, xref, xtmp, od->R);
454
455     for (int fa = 0; fa < nfa; fa += 3)
456     {
457         const int type           = forceatoms[fa];
458         const int restraintIndex = type - od->typeMin;
459         if (pbc)
460         {
461             pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
462         }
463         else
464         {
465             rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
466         }
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++)
473         {
474             pfac *= invr;
475         }
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]);
482
483         if (ms)
484         {
485             for (int i = 0; i < 5; i++)
486             {
487                 od->Dins[restraintIndex][i] = Dinsl[i]*invn;
488             }
489         }
490     }
491
492     if (ms)
493     {
494         gmx_sum_sim(5*od->nr, od->Dins[0], ms);
495     }
496
497     /* Calculate the order tensor S for each experiment via optimization */
498     for (int ex = 0; ex < od->nex; ex++)
499     {
500         for (int i = 0; i < 5; i++)
501         {
502             matEq[ex].rhs[i] = 0;
503             for (int j = 0; j <= i; j++)
504             {
505                 matEq[ex].mat[i][j] = 0;
506             }
507         }
508     }
509
510     for (int fa = 0; fa < nfa; fa += 3)
511     {
512         const int  type           = forceatoms[fa];
513         const int  restraintIndex = type - od->typeMin;
514         rvec5     &Dtav           = od->Dtav[restraintIndex];
515         if (bTAV)
516         {
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.
520              */
521             for (int i = 0; i < 5; i++)
522             {
523                 Dtav[i] =
524                     edt*hist->orire_Dtav[restraintIndex*5 + i] +
525                     edt_1*od->Dins[restraintIndex][i];
526             }
527         }
528
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++)
533         {
534             matEq[ex].rhs[i] += Dtav[i]*ip[type].orires.obs*weight;
535             for (int j = 0; j <= i; j++)
536             {
537                 matEq[ex].mat[i][j] += Dtav[i]*Dtav[j]*weight;
538             }
539         }
540     }
541     /* Now we have all the data we can calculate S */
542     for (int ex = 0; ex < od->nex; ex++)
543     {
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++)
547         {
548             eq.rhs[i]    *= corrfac;
549             eq.mat[i][i] *= gmx::square(corrfac);
550             for (int j = 0; j < i; j++)
551             {
552                 eq.mat[i][j] *= gmx::square(corrfac);
553                 eq.mat[j][i]  = eq.mat[i][j];
554             }
555         }
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++)
565         {
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];
571         }
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];
576     }
577
578     const matrix *S = od->S;
579
580     wsv2            = 0;
581     sw              = 0;
582
583     for (int fa = 0; fa < nfa; fa += 3)
584     {
585         const int    type           = forceatoms[fa];
586         const int    restraintIndex = type - od->typeMin;
587         const int    ex             = ip[type].orires.ex;
588
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)
595         {
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]);
601         }
602         if (ms)
603         {
604             /* When ensemble averaging is used recalculate the local orientation
605              * for output to the energy file.
606              */
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]);
612         }
613
614         dev = od->otav[restraintIndex] - ip[type].orires.obs;
615
616         wsv2 += ip[type].orires.kfac*gmx::square(dev);
617         sw   += ip[type].orires.kfac;
618     }
619     od->rmsdev = std::sqrt(wsv2/sw);
620
621     /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
622     for (int ex = 0; ex < od->nex; ex++)
623     {
624         matrix RS;
625         tmmul(od->R, od->S[ex], RS);
626         mmul(RS, od->R, od->S[ex]);
627     }
628
629     return od->rmsdev;
630
631     /* Approx. 120*nfa/3 flops */
632 }
633
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)
640 {
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;
648
649     vtot = 0;
650     od   = &(fcd->orires);
651
652     if (od->fc != 0)
653     {
654         bTAV = (od->edt != 0);
655
656         smooth_fc = od->fc;
657         if (bTAV)
658         {
659             /* Smoothly switch on the restraining when time averaging is used */
660             smooth_fc *= (1.0 - od->exp_min_t_tau);
661         }
662
663         for (int fa = 0; fa < nfa; fa += 3)
664         {
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)
670             {
671                 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
672             }
673             else
674             {
675                 rvec_sub(x[ai], x[aj], r);
676             }
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;
684
685             /* NOTE:
686              * there is no real potential when time averaging is applied
687              */
688             vtot += 0.5*fc*gmx::square(dev);
689
690             if (bTAV)
691             {
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)
695                 {
696                     dev = 0;
697                 }
698                 else
699                 {
700                     dev = std::sqrt(dev*devins);
701                     if (devins < 0)
702                     {
703                         dev = -dev;
704                     }
705                 }
706             }
707
708             pfac  = fc*ip[type].orires.c*invr2;
709             for (int i = 0; i < power; i++)
710             {
711                 pfac *= invr;
712             }
713             mvmul(od->S[ex], r, Sr);
714             for (int i = 0; i < DIM; i++)
715             {
716                 fij[i] =
717                     -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
718             }
719
720             if (g)
721             {
722                 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
723                 ki = IVEC2IS(dt);
724             }
725
726             for (int i = 0; i < DIM; i++)
727             {
728                 f[ai][i]           += fij[i];
729                 f[aj][i]           -= fij[i];
730                 fshift[ki][i]      += fij[i];
731                 fshift[CENTRAL][i] -= fij[i];
732             }
733         }
734     }
735
736     return vtot;
737
738     /* Approx. 80*nfa/3 flops */
739 }
740
741 void update_orires_history(t_fcdata *fcd, history_t *hist)
742 {
743     t_oriresdata *od = &(fcd->orires);
744
745     if (od->edt != 0)
746     {
747         /* Copy the new time averages that have been calculated
748          *  in calc_orires_dev.
749          */
750         hist->orire_initf = od->exp_min_t_tau;
751         for (int pair = 0; pair < od->nr; pair++)
752         {
753             for (int i = 0; i < 5; i++)
754             {
755                 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];
756             }
757         }
758     }
759 }