Change gmx_moltype_t.ilist to InteractionLists
[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/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/fcdata.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/pbcutil/ishift.h"
55 #include "gromacs/pbcutil/mshift.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/topology/mtop_util.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/pleasecite.h"
62 #include "gromacs/utility/smalloc.h"
63
64 // TODO This implementation of ensemble orientation restraints is nasty because
65 // a user can't just do multi-sim with single-sim orientation restraints.
66
67 void init_orires(FILE                 *fplog,
68                  const gmx_mtop_t     *mtop,
69                  const t_inputrec     *ir,
70                  const t_commrec      *cr,
71                  const gmx_multisim_t *ms,
72                  t_state              *globalState,
73                  t_oriresdata         *od)
74 {
75     od->nr = gmx_mtop_ftype_count(mtop, F_ORIRES);
76     if (0 == od->nr)
77     {
78         /* Not doing orientation restraints */
79         return;
80     }
81
82     const int numFitParams = 5;
83     if (od->nr <= numFitParams)
84     {
85         gmx_fatal(FARGS, "The system has %d orientation restraints, but at least %d are required, since there are %d fitting parameters.",
86                   od->nr, numFitParams + 1, numFitParams);
87     }
88
89     if (ir->bPeriodicMols)
90     {
91         /* Since we apply fitting, we need to make molecules whole and this
92          * can not be done when periodic molecules are present.
93          */
94         gmx_fatal(FARGS, "Orientation restraints can not be applied when periodic molecules are present in the system");
95     }
96
97     if (PAR(cr))
98     {
99         gmx_fatal(FARGS, "Orientation restraints do not work with MPI parallelization. Choose 1 MPI rank, if possible.");
100     }
101
102     GMX_RELEASE_ASSERT(globalState != nullptr, "We need a valid global state in init_orires");
103
104     od->fc  = ir->orires_fc;
105     od->nex = 0;
106     od->S   = nullptr;
107     od->M   = nullptr;
108     od->eig = nullptr;
109     od->v   = nullptr;
110
111     int                   *nr_ex   = nullptr;
112     int                    typeMin = INT_MAX;
113     int                    typeMax = 0;
114     gmx_mtop_ilistloop_t   iloop   = gmx_mtop_ilistloop_init(mtop);
115     int                    nmol;
116     while (const InteractionLists *il = gmx_mtop_ilistloop_next(iloop, &nmol))
117     {
118         if (nmol > 1)
119         {
120             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);
121         }
122
123         for (int i = 0; i < (*il)[F_ORIRES].size(); i += 3)
124         {
125             int type = (*il)[F_ORIRES].iatoms[i];
126             int ex   = mtop->ffparams.iparams[type].orires.ex;
127             if (ex >= od->nex)
128             {
129                 srenew(nr_ex, ex+1);
130                 for (int j = od->nex; j < ex+1; j++)
131                 {
132                     nr_ex[j] = 0;
133                 }
134                 od->nex = ex+1;
135             }
136             GMX_ASSERT(nr_ex, "Check for allocated nr_ex to keep the static analyzer happy");
137             nr_ex[ex]++;
138
139             typeMin = std::min(typeMin, type);
140             typeMax = std::max(typeMax, type);
141         }
142     }
143     /* With domain decomposition we use the type index for indexing in global arrays */
144     GMX_RELEASE_ASSERT(typeMax - typeMin + 1 == od->nr, "All orientation restraint parameter entries in the topology should be consecutive");
145     /* Store typeMin so we can index array with the type offset */
146     od->typeMin = typeMin;
147
148     snew(od->S, od->nex);
149     /* When not doing time averaging, the instaneous and time averaged data
150      * are indentical and the pointers can point to the same memory.
151      */
152     snew(od->Dinsl, od->nr);
153
154     if (ms)
155     {
156         snew(od->Dins, od->nr);
157     }
158     else
159     {
160         od->Dins = od->Dinsl;
161     }
162
163     if (ir->orires_tau == 0)
164     {
165         od->Dtav  = od->Dins;
166         od->edt   = 0.0;
167         od->edt_1 = 1.0;
168     }
169     else
170     {
171         snew(od->Dtav, od->nr);
172         od->edt   = std::exp(-ir->delta_t/ir->orires_tau);
173         od->edt_1 = 1.0 - od->edt;
174
175         /* Extend the state with the orires history */
176         globalState->flags           |= (1<<estORIRE_INITF);
177         globalState->hist.orire_initf = 1;
178         globalState->flags           |= (1<<estORIRE_DTAV);
179         globalState->hist.norire_Dtav = od->nr*5;
180         snew(globalState->hist.orire_Dtav, globalState->hist.norire_Dtav);
181     }
182
183     snew(od->oinsl, od->nr);
184     if (ms)
185     {
186         snew(od->oins, od->nr);
187     }
188     else
189     {
190         od->oins = od->oinsl;
191     }
192     if (ir->orires_tau == 0)
193     {
194         od->otav = od->oins;
195     }
196     else
197     {
198         snew(od->otav, od->nr);
199     }
200     snew(od->tmpEq, od->nex);
201
202     od->nref = 0;
203     for (int i = 0; i < mtop->natoms; i++)
204     {
205         if (getGroupType(&mtop->groups, egcORFIT, i) == 0)
206         {
207             od->nref++;
208         }
209     }
210     snew(od->mref, od->nref);
211     snew(od->xref, od->nref);
212     snew(od->xtmp, od->nref);
213
214     snew(od->eig, od->nex*12);
215
216     /* Determine the reference structure on the master node.
217      * Copy it to the other nodes after checking multi compatibility,
218      * so we are sure the subsystems match before copying.
219      */
220     rvec                     com   = { 0, 0, 0 };
221     double                   mtot  = 0.0;
222     int                      j     = 0;
223     gmx_mtop_atomloop_all_t  aloop = gmx_mtop_atomloop_all_init(mtop);
224     int                      i     = -1;
225     const t_atom            *atom;
226     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
227     {
228         if (mtop->groups.grpnr[egcORFIT] == nullptr ||
229             mtop->groups.grpnr[egcORFIT][i] == 0)
230         {
231             /* Not correct for free-energy with changing masses */
232             od->mref[j] = atom->m;
233             // Note that only one rank per sim is supported.
234             if (isMasterSim(ms))
235             {
236                 copy_rvec(globalState->x[i], od->xref[j]);
237                 for (int d = 0; d < DIM; d++)
238                 {
239                     com[d] += od->mref[j]*globalState->x[i][d];
240                 }
241             }
242             mtot += od->mref[j];
243             j++;
244         }
245     }
246     svmul(1.0/mtot, com, com);
247     if (isMasterSim(ms))
248     {
249         for (int j = 0; j < od->nref; j++)
250         {
251             rvec_dec(od->xref[j], com);
252         }
253     }
254
255     fprintf(fplog, "Found %d orientation experiments\n", od->nex);
256     for (int i = 0; i < od->nex; i++)
257     {
258         fprintf(fplog, "  experiment %d has %d restraints\n", i+1, nr_ex[i]);
259     }
260
261     sfree(nr_ex);
262
263     fprintf(fplog, "  the fit group consists of %d atoms and has total mass %g\n",
264             od->nref, mtot);
265
266     if (ms)
267     {
268         fprintf(fplog, "  the orientation restraints are ensemble averaged over %d systems\n", ms->nsim);
269
270         check_multi_int(fplog, ms, od->nr,
271                         "the number of orientation restraints",
272                         FALSE);
273         check_multi_int(fplog, ms, od->nref,
274                         "the number of fit atoms for orientation restraining",
275                         FALSE);
276         check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
277         /* Copy the reference coordinates from the master to the other nodes */
278         gmx_sum_sim(DIM*od->nref, od->xref[0], ms);
279     }
280
281     please_cite(fplog, "Hess2003");
282 }
283
284 void diagonalize_orires_tensors(t_oriresdata *od)
285 {
286     if (od->M == nullptr)
287     {
288         snew(od->M, DIM);
289         for (int i = 0; i < DIM; i++)
290         {
291             snew(od->M[i], DIM);
292         }
293         snew(od->eig_diag, DIM);
294         snew(od->v, DIM);
295         for (int i = 0; i < DIM; i++)
296         {
297             snew(od->v[i], DIM);
298         }
299     }
300
301     for (int ex = 0; ex < od->nex; ex++)
302     {
303         /* Rotate the S tensor back to the reference frame */
304         matrix S, TMP;
305         mmul(od->R, od->S[ex], TMP);
306         mtmul(TMP, od->R, S);
307         for (int i = 0; i < DIM; i++)
308         {
309             for (int j = 0; j < DIM; j++)
310             {
311                 od->M[i][j] = S[i][j];
312             }
313         }
314
315         int nrot;
316         jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
317
318         int ord[DIM];
319         for (int i = 0; i < DIM; i++)
320         {
321             ord[i] = i;
322         }
323         for (int i = 0; i < DIM; i++)
324         {
325             for (int j = i+1; j < DIM; j++)
326             {
327                 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
328                 {
329                     int t  = ord[i];
330                     ord[i] = ord[j];
331                     ord[j] = t;
332                 }
333             }
334         }
335
336         for (int i = 0; i < DIM; i++)
337         {
338             od->eig[ex*12 + i] = od->eig_diag[ord[i]];
339         }
340         for (int i = 0; i < DIM; i++)
341         {
342             for (int j = 0; j < DIM; j++)
343             {
344                 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
345             }
346         }
347     }
348 }
349
350 void print_orires_log(FILE *log, t_oriresdata *od)
351 {
352     real *eig;
353
354     diagonalize_orires_tensors(od);
355
356     for (int ex = 0; ex < od->nex; ex++)
357     {
358         eig = od->eig + ex*12;
359         fprintf(log, "  Orientation experiment %d:\n", ex+1);
360         fprintf(log, "    order parameter: %g\n", eig[0]);
361         for (int i = 0; i < DIM; i++)
362         {
363             fprintf(log, "    eig: %6.3f   %6.3f %6.3f %6.3f\n",
364                     (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
365                     eig[DIM+i*DIM+XX],
366                     eig[DIM+i*DIM+YY],
367                     eig[DIM+i*DIM+ZZ]);
368         }
369         fprintf(log, "\n");
370     }
371 }
372
373 real calc_orires_dev(const gmx_multisim_t *ms,
374                      int nfa, const t_iatom forceatoms[], const t_iparams ip[],
375                      const t_mdatoms *md, const rvec x[], const t_pbc *pbc,
376                      t_fcdata *fcd, history_t *hist)
377 {
378     int              nref;
379     real             edt, edt_1, invn, pfac, r2, invr, corrfac, wsv2, sw, dev;
380     OriresMatEq     *matEq;
381     real            *mref;
382     double           mtot;
383     rvec            *xref, *xtmp, com, r_unrot, r;
384     t_oriresdata    *od;
385     gmx_bool         bTAV;
386     const real       two_thr = 2.0/3.0;
387
388     od = &(fcd->orires);
389
390     if (od->nr == 0)
391     {
392         /* This means that this is not the master node */
393         gmx_fatal(FARGS, "Orientation restraints are only supported on the master rank, use fewer ranks");
394     }
395
396     bTAV  = (od->edt != 0);
397     edt   = od->edt;
398     edt_1 = od->edt_1;
399     matEq = od->tmpEq;
400     nref  = od->nref;
401     mref  = od->mref;
402     xref  = od->xref;
403     xtmp  = od->xtmp;
404
405     if (bTAV)
406     {
407         od->exp_min_t_tau = hist->orire_initf*edt;
408
409         /* Correction factor to correct for the lack of history
410          * at short times.
411          */
412         corrfac = 1.0/(1.0 - od->exp_min_t_tau);
413     }
414     else
415     {
416         corrfac = 1.0;
417     }
418
419     if (ms)
420     {
421         invn = 1.0/ms->nsim;
422     }
423     else
424     {
425         invn = 1.0;
426     }
427
428     clear_rvec(com);
429     mtot  = 0;
430     int j = 0;
431     for (int i = 0; i < md->nr; i++)
432     {
433         if (md->cORF[i] == 0)
434         {
435             copy_rvec(x[i], xtmp[j]);
436             mref[j] = md->massT[i];
437             for (int d = 0; d < DIM; d++)
438             {
439                 com[d] += mref[j]*xtmp[j][d];
440             }
441             mtot += mref[j];
442             j++;
443         }
444     }
445     svmul(1.0/mtot, com, com);
446     for (int j = 0; j < nref; j++)
447     {
448         rvec_dec(xtmp[j], com);
449     }
450     /* Calculate the rotation matrix to rotate x to the reference orientation */
451     calc_fit_R(DIM, nref, mref, xref, xtmp, od->R);
452
453     for (int fa = 0; fa < nfa; fa += 3)
454     {
455         const int type           = forceatoms[fa];
456         const int restraintIndex = type - od->typeMin;
457         if (pbc)
458         {
459             pbc_dx_aiuc(pbc, x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
460         }
461         else
462         {
463             rvec_sub(x[forceatoms[fa+1]], x[forceatoms[fa+2]], r_unrot);
464         }
465         mvmul(od->R, r_unrot, r);
466         r2   = norm2(r);
467         invr = gmx::invsqrt(r2);
468         /* Calculate the prefactor for the D tensor, this includes the factor 3! */
469         pfac = ip[type].orires.c*invr*invr*3;
470         for (int i = 0; i < ip[type].orires.power; i++)
471         {
472             pfac *= invr;
473         }
474         rvec5 &Dinsl = od->Dinsl[restraintIndex];
475         Dinsl[0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
476         Dinsl[1] = pfac*(2*r[0]*r[1]);
477         Dinsl[2] = pfac*(2*r[0]*r[2]);
478         Dinsl[3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
479         Dinsl[4] = pfac*(2*r[1]*r[2]);
480
481         if (ms)
482         {
483             for (int i = 0; i < 5; i++)
484             {
485                 od->Dins[restraintIndex][i] = Dinsl[i]*invn;
486             }
487         }
488     }
489
490     if (ms)
491     {
492         gmx_sum_sim(5*od->nr, od->Dins[0], ms);
493     }
494
495     /* Calculate the order tensor S for each experiment via optimization */
496     for (int ex = 0; ex < od->nex; ex++)
497     {
498         for (int i = 0; i < 5; i++)
499         {
500             matEq[ex].rhs[i] = 0;
501             for (int j = 0; j <= i; j++)
502             {
503                 matEq[ex].mat[i][j] = 0;
504             }
505         }
506     }
507
508     for (int fa = 0; fa < nfa; fa += 3)
509     {
510         const int  type           = forceatoms[fa];
511         const int  restraintIndex = type - od->typeMin;
512         rvec5     &Dtav           = od->Dtav[restraintIndex];
513         if (bTAV)
514         {
515             /* Here we update Dtav in t_fcdata using the data in history_t.
516              * Thus the results stay correct when this routine
517              * is called multiple times.
518              */
519             for (int i = 0; i < 5; i++)
520             {
521                 Dtav[i] =
522                     edt*hist->orire_Dtav[restraintIndex*5 + i] +
523                     edt_1*od->Dins[restraintIndex][i];
524             }
525         }
526
527         int  ex     = ip[type].orires.ex;
528         real weight = ip[type].orires.kfac;
529         /* Calculate the vector rhs and half the matrix T for the 5 equations */
530         for (int i = 0; i < 5; i++)
531         {
532             matEq[ex].rhs[i] += Dtav[i]*ip[type].orires.obs*weight;
533             for (int j = 0; j <= i; j++)
534             {
535                 matEq[ex].mat[i][j] += Dtav[i]*Dtav[j]*weight;
536             }
537         }
538     }
539     /* Now we have all the data we can calculate S */
540     for (int ex = 0; ex < od->nex; ex++)
541     {
542         OriresMatEq &eq = matEq[ex];
543         /* Correct corrfac and copy one half of T to the other half */
544         for (int i = 0; i < 5; i++)
545         {
546             eq.rhs[i]    *= corrfac;
547             eq.mat[i][i] *= gmx::square(corrfac);
548             for (int j = 0; j < i; j++)
549             {
550                 eq.mat[i][j] *= gmx::square(corrfac);
551                 eq.mat[j][i]  = eq.mat[i][j];
552             }
553         }
554         m_inv_gen(&eq.mat[0][0], 5, &eq.mat[0][0]);
555         /* Calculate the orientation tensor S for this experiment */
556         matrix &S = od->S[ex];
557         S[0][0] = 0;
558         S[0][1] = 0;
559         S[0][2] = 0;
560         S[1][1] = 0;
561         S[1][2] = 0;
562         for (int i = 0; i < 5; i++)
563         {
564             S[0][0] += 1.5*eq.mat[0][i]*eq.rhs[i];
565             S[0][1] += 1.5*eq.mat[1][i]*eq.rhs[i];
566             S[0][2] += 1.5*eq.mat[2][i]*eq.rhs[i];
567             S[1][1] += 1.5*eq.mat[3][i]*eq.rhs[i];
568             S[1][2] += 1.5*eq.mat[4][i]*eq.rhs[i];
569         }
570         S[1][0] = S[0][1];
571         S[2][0] = S[0][2];
572         S[2][1] = S[1][2];
573         S[2][2] = -S[0][0] - S[1][1];
574     }
575
576     const matrix *S = od->S;
577
578     wsv2            = 0;
579     sw              = 0;
580
581     for (int fa = 0; fa < nfa; fa += 3)
582     {
583         const int    type           = forceatoms[fa];
584         const int    restraintIndex = type - od->typeMin;
585         const int    ex             = ip[type].orires.ex;
586
587         const rvec5 &Dtav           = od->Dtav[restraintIndex];
588         od->otav[restraintIndex]    = two_thr*
589             corrfac*(S[ex][0][0]*Dtav[0] + S[ex][0][1]*Dtav[1] +
590                      S[ex][0][2]*Dtav[2] + S[ex][1][1]*Dtav[3] +
591                      S[ex][1][2]*Dtav[4]);
592         if (bTAV)
593         {
594             const rvec5 &Dins = od->Dins[restraintIndex];
595             od->oins[restraintIndex] = two_thr*
596                 (S[ex][0][0]*Dins[0] + S[ex][0][1]*Dins[1] +
597                  S[ex][0][2]*Dins[2] + S[ex][1][1]*Dins[3] +
598                  S[ex][1][2]*Dins[4]);
599         }
600         if (ms)
601         {
602             /* When ensemble averaging is used recalculate the local orientation
603              * for output to the energy file.
604              */
605             const rvec5 &Dinsl = od->Dinsl[restraintIndex];
606             od->oinsl[restraintIndex] = two_thr*
607                 (S[ex][0][0]*Dinsl[0] + S[ex][0][1]*Dinsl[1] +
608                  S[ex][0][2]*Dinsl[2] + S[ex][1][1]*Dinsl[3] +
609                  S[ex][1][2]*Dinsl[4]);
610         }
611
612         dev = od->otav[restraintIndex] - ip[type].orires.obs;
613
614         wsv2 += ip[type].orires.kfac*gmx::square(dev);
615         sw   += ip[type].orires.kfac;
616     }
617     od->rmsdev = std::sqrt(wsv2/sw);
618
619     /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
620     for (int ex = 0; ex < od->nex; ex++)
621     {
622         matrix RS;
623         tmmul(od->R, od->S[ex], RS);
624         mmul(RS, od->R, od->S[ex]);
625     }
626
627     return od->rmsdev;
628
629     /* Approx. 120*nfa/3 flops */
630 }
631
632 real orires(int nfa, const t_iatom forceatoms[], const t_iparams ip[],
633             const rvec x[], rvec4 f[], rvec fshift[],
634             const t_pbc *pbc, const t_graph *g,
635             real gmx_unused lambda, real gmx_unused *dvdlambda,
636             const t_mdatoms gmx_unused *md, t_fcdata *fcd,
637             int gmx_unused *global_atom_index)
638 {
639     int                 ex, power, ki = CENTRAL;
640     ivec                dt;
641     real                r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
642     rvec                r, Sr, fij;
643     real                vtot;
644     const t_oriresdata *od;
645     gmx_bool            bTAV;
646
647     vtot = 0;
648     od   = &(fcd->orires);
649
650     if (od->fc != 0)
651     {
652         bTAV = (od->edt != 0);
653
654         smooth_fc = od->fc;
655         if (bTAV)
656         {
657             /* Smoothly switch on the restraining when time averaging is used */
658             smooth_fc *= (1.0 - od->exp_min_t_tau);
659         }
660
661         for (int fa = 0; fa < nfa; fa += 3)
662         {
663             const int type           = forceatoms[fa];
664             const int ai             = forceatoms[fa + 1];
665             const int aj             = forceatoms[fa + 2];
666             const int restraintIndex = type - od->typeMin;
667             if (pbc)
668             {
669                 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
670             }
671             else
672             {
673                 rvec_sub(x[ai], x[aj], r);
674             }
675             r2    = norm2(r);
676             invr  = gmx::invsqrt(r2);
677             invr2 = invr*invr;
678             ex    = ip[type].orires.ex;
679             power = ip[type].orires.power;
680             fc    = smooth_fc*ip[type].orires.kfac;
681             dev   = od->otav[restraintIndex] - ip[type].orires.obs;
682
683             /* NOTE:
684              * there is no real potential when time averaging is applied
685              */
686             vtot += 0.5*fc*gmx::square(dev);
687
688             if (bTAV)
689             {
690                 /* Calculate the force as the sqrt of tav times instantaneous */
691                 devins = od->oins[restraintIndex] - ip[type].orires.obs;
692                 if (dev*devins <= 0)
693                 {
694                     dev = 0;
695                 }
696                 else
697                 {
698                     dev = std::sqrt(dev*devins);
699                     if (devins < 0)
700                     {
701                         dev = -dev;
702                     }
703                 }
704             }
705
706             pfac  = fc*ip[type].orires.c*invr2;
707             for (int i = 0; i < power; i++)
708             {
709                 pfac *= invr;
710             }
711             mvmul(od->S[ex], r, Sr);
712             for (int i = 0; i < DIM; i++)
713             {
714                 fij[i] =
715                     -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr, r)*r[i]);
716             }
717
718             if (g)
719             {
720                 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
721                 ki = IVEC2IS(dt);
722             }
723
724             for (int i = 0; i < DIM; i++)
725             {
726                 f[ai][i]           += fij[i];
727                 f[aj][i]           -= fij[i];
728                 fshift[ki][i]      += fij[i];
729                 fshift[CENTRAL][i] -= fij[i];
730             }
731         }
732     }
733
734     return vtot;
735
736     /* Approx. 80*nfa/3 flops */
737 }
738
739 void update_orires_history(const t_fcdata *fcd, history_t *hist)
740 {
741     const t_oriresdata *od = &(fcd->orires);
742
743     if (od->edt != 0)
744     {
745         /* Copy the new time averages that have been calculated
746          *  in calc_orires_dev.
747          */
748         hist->orire_initf = od->exp_min_t_tau;
749         for (int pair = 0; pair < od->nr; pair++)
750         {
751             for (int i = 0; i < 5; i++)
752             {
753                 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];
754             }
755         }
756     }
757 }