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