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