1ec9cef352005a8a2a7b866571c4482ab9a89f8d
[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.orire_Dtav.resize(od->nr * 5);
198     }
199
200     snew(od->oinsl, od->nr);
201     if (ms)
202     {
203         snew(od->oins, od->nr);
204     }
205     else
206     {
207         od->oins = od->oinsl;
208     }
209     if (ir->orires_tau == 0)
210     {
211         od->otav = od->oins;
212     }
213     else
214     {
215         snew(od->otav, od->nr);
216     }
217     snew(od->tmpEq, od->nex);
218
219     od->nref = 0;
220     for (int i = 0; i < mtop.natoms; i++)
221     {
222         if (getGroupType(mtop.groups, SimulationAtomGroupType::OrientationRestraintsFit, i) == 0)
223         {
224             od->nref++;
225         }
226     }
227     snew(od->mref, od->nref);
228     snew(od->xref, od->nref);
229     snew(od->xtmp, od->nref);
230
231     snew(od->eig, od->nex * 12);
232
233     /* Determine the reference structure on the master node.
234      * Copy it to the other nodes after checking multi compatibility,
235      * so we are sure the subsystems match before copying.
236      */
237     auto   x    = makeArrayRef(globalState->x);
238     rvec   com  = { 0, 0, 0 };
239     double mtot = 0.0;
240     int    j    = 0;
241     for (const AtomProxy atomP : AtomRange(mtop))
242     {
243         const t_atom& local = atomP.atom();
244         int           i     = atomP.globalAtomNumber();
245         if (mtop.groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit].empty()
246             || mtop.groups.groupNumbers[SimulationAtomGroupType::OrientationRestraintsFit][i] == 0)
247         {
248             /* Not correct for free-energy with changing masses */
249             od->mref[j] = local.m;
250             // Note that only one rank per sim is supported.
251             if (isMasterSim(ms))
252             {
253                 copy_rvec(x[i], od->xref[j]);
254                 for (int d = 0; d < DIM; d++)
255                 {
256                     com[d] += od->mref[j] * x[i][d];
257                 }
258             }
259             mtot += od->mref[j];
260             j++;
261         }
262     }
263     svmul(1.0 / mtot, com, com);
264     if (isMasterSim(ms))
265     {
266         for (int j = 0; j < od->nref; j++)
267         {
268             rvec_dec(od->xref[j], com);
269         }
270     }
271
272     fprintf(fplog, "Found %d orientation experiments\n", od->nex);
273     for (int i = 0; i < od->nex; i++)
274     {
275         fprintf(fplog, "  experiment %d has %d restraints\n", i + 1, nr_ex[i]);
276     }
277
278     sfree(nr_ex);
279
280     fprintf(fplog, "  the fit group consists of %d atoms and has total mass %g\n", od->nref, mtot);
281
282     if (ms)
283     {
284         fprintf(fplog, "  the orientation restraints are ensemble averaged over %d systems\n", ms->numSimulations_);
285
286         check_multi_int(fplog, ms, od->nr, "the number of orientation restraints", FALSE);
287         check_multi_int(fplog, ms, od->nref, "the number of fit atoms for orientation restraining", FALSE);
288         check_multi_int(fplog, ms, ir->nsteps, "nsteps", FALSE);
289         /* Copy the reference coordinates from the master to the other nodes */
290         gmx_sum_sim(DIM * od->nref, od->xref[0], ms);
291     }
292
293     please_cite(fplog, "Hess2003");
294 }
295
296 void diagonalize_orires_tensors(t_oriresdata* od)
297 {
298     if (od->M == nullptr)
299     {
300         snew(od->M, DIM);
301         for (int i = 0; i < DIM; i++)
302         {
303             snew(od->M[i], DIM);
304         }
305         snew(od->eig_diag, DIM);
306         snew(od->v, DIM);
307         for (int i = 0; i < DIM; i++)
308         {
309             snew(od->v[i], DIM);
310         }
311     }
312
313     for (int ex = 0; ex < od->nex; ex++)
314     {
315         /* Rotate the S tensor back to the reference frame */
316         matrix S, TMP;
317         mmul(od->R, od->S[ex], TMP);
318         mtmul(TMP, od->R, S);
319         for (int i = 0; i < DIM; i++)
320         {
321             for (int j = 0; j < DIM; j++)
322             {
323                 od->M[i][j] = S[i][j];
324             }
325         }
326
327         int nrot;
328         jacobi(od->M, DIM, od->eig_diag, od->v, &nrot);
329
330         int ord[DIM];
331         for (int i = 0; i < DIM; i++)
332         {
333             ord[i] = i;
334         }
335         for (int i = 0; i < DIM; i++)
336         {
337             for (int j = i + 1; j < DIM; j++)
338             {
339                 if (gmx::square(od->eig_diag[ord[j]]) > gmx::square(od->eig_diag[ord[i]]))
340                 {
341                     int t  = ord[i];
342                     ord[i] = ord[j];
343                     ord[j] = t;
344                 }
345             }
346         }
347
348         for (int i = 0; i < DIM; i++)
349         {
350             od->eig[ex * 12 + i] = od->eig_diag[ord[i]];
351         }
352         for (int i = 0; i < DIM; i++)
353         {
354             for (int j = 0; j < DIM; j++)
355             {
356                 od->eig[ex * 12 + 3 + 3 * i + j] = od->v[j][ord[i]];
357             }
358         }
359     }
360 }
361
362 void print_orires_log(FILE* log, t_oriresdata* od)
363 {
364     real* eig;
365
366     diagonalize_orires_tensors(od);
367
368     for (int ex = 0; ex < od->nex; ex++)
369     {
370         eig = od->eig + ex * 12;
371         fprintf(log, "  Orientation experiment %d:\n", ex + 1);
372         fprintf(log, "    order parameter: %g\n", eig[0]);
373         for (int i = 0; i < DIM; i++)
374         {
375             fprintf(log,
376                     "    eig: %6.3f   %6.3f %6.3f %6.3f\n",
377                     (eig[0] != 0) ? eig[i] / eig[0] : eig[i],
378                     eig[DIM + i * DIM + XX],
379                     eig[DIM + i * DIM + YY],
380                     eig[DIM + i * DIM + ZZ]);
381         }
382         fprintf(log, "\n");
383     }
384 }
385
386 real calc_orires_dev(const gmx_multisim_t* ms,
387                      int                   nfa,
388                      const t_iatom         forceatoms[],
389                      const t_iparams       ip[],
390                      const t_mdatoms*      md,
391                      ArrayRef<const RVec>  xWholeMolecules,
392                      const rvec            x[],
393                      const t_pbc*          pbc,
394                      t_oriresdata*         od,
395                      const history_t*      hist)
396 {
397     int          nref;
398     real         edt, edt_1, invn, pfac, r2, invr, corrfac, wsv2, sw, dev;
399     OriresMatEq* matEq;
400     real*        mref;
401     double       mtot;
402     rvec *       xref, *xtmp, com, r_unrot, r;
403     gmx_bool     bTAV;
404     const real   two_thr = 2.0 / 3.0;
405
406     if (od->nr == 0)
407     {
408         /* This means that this is not the master node */
409         gmx_fatal(FARGS,
410                   "Orientation restraints are only supported on the master rank, use fewer ranks");
411     }
412
413     bTAV  = (od->edt != 0);
414     edt   = od->edt;
415     edt_1 = od->edt_1;
416     matEq = od->tmpEq;
417     nref  = od->nref;
418     mref  = od->mref;
419     xref  = od->xref;
420     xtmp  = od->xtmp;
421
422     if (bTAV)
423     {
424         od->exp_min_t_tau = hist->orire_initf * edt;
425
426         /* Correction factor to correct for the lack of history
427          * at short times.
428          */
429         corrfac = 1.0 / (1.0 - od->exp_min_t_tau);
430     }
431     else
432     {
433         corrfac = 1.0;
434     }
435
436     if (ms)
437     {
438         invn = 1.0 / ms->numSimulations_;
439     }
440     else
441     {
442         invn = 1.0;
443     }
444
445     clear_rvec(com);
446     mtot        = 0;
447     int   j     = 0;
448     auto* massT = md->massT;
449     auto* cORF  = md->cORF;
450     for (int i = 0; i < md->nr; i++)
451     {
452         if (cORF[i] == 0)
453         {
454             copy_rvec(xWholeMolecules[i], xtmp[j]);
455             mref[j] = massT[i];
456             for (int d = 0; d < DIM; d++)
457             {
458                 com[d] += mref[j] * xtmp[j][d];
459             }
460             mtot += mref[j];
461             j++;
462         }
463     }
464     svmul(1.0 / mtot, com, com);
465     for (int j = 0; j < nref; j++)
466     {
467         rvec_dec(xtmp[j], com);
468     }
469     /* Calculate the rotation matrix to rotate x to the reference orientation */
470     calc_fit_R(DIM, nref, mref, xref, xtmp, od->R);
471
472     for (int fa = 0; fa < nfa; fa += 3)
473     {
474         const int type           = forceatoms[fa];
475         const int restraintIndex = type - od->typeMin;
476         if (pbc)
477         {
478             pbc_dx_aiuc(pbc, x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
479         }
480         else
481         {
482             rvec_sub(x[forceatoms[fa + 1]], x[forceatoms[fa + 2]], r_unrot);
483         }
484         mvmul(od->R, r_unrot, r);
485         r2   = norm2(r);
486         invr = gmx::invsqrt(r2);
487         /* Calculate the prefactor for the D tensor, this includes the factor 3! */
488         pfac = ip[type].orires.c * invr * invr * 3;
489         for (int i = 0; i < ip[type].orires.power; i++)
490         {
491             pfac *= invr;
492         }
493         rvec5& Dinsl = od->Dinsl[restraintIndex];
494         Dinsl[0]     = pfac * (2 * r[0] * r[0] + r[1] * r[1] - r2);
495         Dinsl[1]     = pfac * (2 * r[0] * r[1]);
496         Dinsl[2]     = pfac * (2 * r[0] * r[2]);
497         Dinsl[3]     = pfac * (2 * r[1] * r[1] + r[0] * r[0] - r2);
498         Dinsl[4]     = pfac * (2 * r[1] * r[2]);
499
500         if (ms)
501         {
502             for (int i = 0; i < 5; i++)
503             {
504                 od->Dins[restraintIndex][i] = Dinsl[i] * invn;
505             }
506         }
507     }
508
509     if (ms)
510     {
511         gmx_sum_sim(5 * od->nr, od->Dins[0], ms);
512     }
513
514     /* Calculate the order tensor S for each experiment via optimization */
515     for (int ex = 0; ex < od->nex; ex++)
516     {
517         for (int i = 0; i < 5; i++)
518         {
519             matEq[ex].rhs[i] = 0;
520             for (int j = 0; j <= i; j++)
521             {
522                 matEq[ex].mat[i][j] = 0;
523             }
524         }
525     }
526
527     for (int fa = 0; fa < nfa; fa += 3)
528     {
529         const int type           = forceatoms[fa];
530         const int restraintIndex = type - od->typeMin;
531         rvec5&    Dtav           = od->Dtav[restraintIndex];
532         if (bTAV)
533         {
534             /* Here we update Dtav in t_fcdata using the data in history_t.
535              * Thus the results stay correct when this routine
536              * is called multiple times.
537              */
538             for (int i = 0; i < 5; i++)
539             {
540                 Dtav[i] = edt * hist->orire_Dtav[restraintIndex * 5 + i]
541                           + edt_1 * od->Dins[restraintIndex][i];
542             }
543         }
544
545         int  ex     = ip[type].orires.ex;
546         real weight = ip[type].orires.kfac;
547         /* Calculate the vector rhs and half the matrix T for the 5 equations */
548         for (int i = 0; i < 5; i++)
549         {
550             matEq[ex].rhs[i] += Dtav[i] * ip[type].orires.obs * weight;
551             for (int j = 0; j <= i; j++)
552             {
553                 matEq[ex].mat[i][j] += Dtav[i] * Dtav[j] * weight;
554             }
555         }
556     }
557     /* Now we have all the data we can calculate S */
558     for (int ex = 0; ex < od->nex; ex++)
559     {
560         OriresMatEq& eq = matEq[ex];
561         /* Correct corrfac and copy one half of T to the other half */
562         for (int i = 0; i < 5; i++)
563         {
564             eq.rhs[i] *= corrfac;
565             eq.mat[i][i] *= gmx::square(corrfac);
566             for (int j = 0; j < i; j++)
567             {
568                 eq.mat[i][j] *= gmx::square(corrfac);
569                 eq.mat[j][i] = eq.mat[i][j];
570             }
571         }
572         m_inv_gen(&eq.mat[0][0], 5, &eq.mat[0][0]);
573         /* Calculate the orientation tensor S for this experiment */
574         matrix& S = od->S[ex];
575         S[0][0]   = 0;
576         S[0][1]   = 0;
577         S[0][2]   = 0;
578         S[1][1]   = 0;
579         S[1][2]   = 0;
580         for (int i = 0; i < 5; i++)
581         {
582             S[0][0] += 1.5 * eq.mat[0][i] * eq.rhs[i];
583             S[0][1] += 1.5 * eq.mat[1][i] * eq.rhs[i];
584             S[0][2] += 1.5 * eq.mat[2][i] * eq.rhs[i];
585             S[1][1] += 1.5 * eq.mat[3][i] * eq.rhs[i];
586             S[1][2] += 1.5 * eq.mat[4][i] * eq.rhs[i];
587         }
588         S[1][0] = S[0][1];
589         S[2][0] = S[0][2];
590         S[2][1] = S[1][2];
591         S[2][2] = -S[0][0] - S[1][1];
592     }
593
594     const matrix* S = od->S;
595
596     wsv2 = 0;
597     sw   = 0;
598
599     for (int fa = 0; fa < nfa; fa += 3)
600     {
601         const int type           = forceatoms[fa];
602         const int restraintIndex = type - od->typeMin;
603         const int ex             = ip[type].orires.ex;
604
605         const rvec5& Dtav = od->Dtav[restraintIndex];
606         od->otav[restraintIndex] =
607                 two_thr * corrfac
608                 * (S[ex][0][0] * Dtav[0] + S[ex][0][1] * Dtav[1] + S[ex][0][2] * Dtav[2]
609                    + S[ex][1][1] * Dtav[3] + S[ex][1][2] * Dtav[4]);
610         if (bTAV)
611         {
612             const rvec5& Dins = od->Dins[restraintIndex];
613             od->oins[restraintIndex] =
614                     two_thr
615                     * (S[ex][0][0] * Dins[0] + S[ex][0][1] * Dins[1] + S[ex][0][2] * Dins[2]
616                        + S[ex][1][1] * Dins[3] + S[ex][1][2] * Dins[4]);
617         }
618         if (ms)
619         {
620             /* When ensemble averaging is used recalculate the local orientation
621              * for output to the energy file.
622              */
623             const rvec5& Dinsl = od->Dinsl[restraintIndex];
624             od->oinsl[restraintIndex] =
625                     two_thr
626                     * (S[ex][0][0] * Dinsl[0] + S[ex][0][1] * Dinsl[1] + S[ex][0][2] * Dinsl[2]
627                        + S[ex][1][1] * Dinsl[3] + S[ex][1][2] * Dinsl[4]);
628         }
629
630         dev = od->otav[restraintIndex] - ip[type].orires.obs;
631
632         wsv2 += ip[type].orires.kfac * gmx::square(dev);
633         sw += ip[type].orires.kfac;
634     }
635     od->rmsdev = std::sqrt(wsv2 / sw);
636
637     /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
638     for (int ex = 0; ex < od->nex; ex++)
639     {
640         matrix RS;
641         tmmul(od->R, od->S[ex], RS);
642         mmul(RS, od->R, od->S[ex]);
643     }
644
645     return od->rmsdev;
646
647     /* Approx. 120*nfa/3 flops */
648 }
649
650 real orires(int             nfa,
651             const t_iatom   forceatoms[],
652             const t_iparams ip[],
653             const rvec      x[],
654             rvec4           f[],
655             rvec            fshift[],
656             const t_pbc*    pbc,
657             real gmx_unused lambda,
658             real gmx_unused* dvdlambda,
659             gmx::ArrayRef<const real> /*charge*/,
660             t_fcdata gmx_unused* fcd,
661             t_disresdata gmx_unused* disresdata,
662             t_oriresdata*            oriresdata,
663             int gmx_unused* global_atom_index)
664 {
665     int      ex, power, ki = gmx::c_centralShiftIndex;
666     real     r2, invr, invr2, fc, smooth_fc, dev, devins, pfac;
667     rvec     r, Sr, fij;
668     real     vtot;
669     gmx_bool bTAV;
670
671     vtot = 0;
672
673     if (oriresdata->fc != 0)
674     {
675         bTAV = (oriresdata->edt != 0);
676
677         smooth_fc = oriresdata->fc;
678         if (bTAV)
679         {
680             /* Smoothly switch on the restraining when time averaging is used */
681             smooth_fc *= (1.0 - oriresdata->exp_min_t_tau);
682         }
683
684         for (int fa = 0; fa < nfa; fa += 3)
685         {
686             const int type           = forceatoms[fa];
687             const int ai             = forceatoms[fa + 1];
688             const int aj             = forceatoms[fa + 2];
689             const int restraintIndex = type - oriresdata->typeMin;
690             if (pbc)
691             {
692                 ki = pbc_dx_aiuc(pbc, x[ai], x[aj], r);
693             }
694             else
695             {
696                 rvec_sub(x[ai], x[aj], r);
697             }
698             r2    = norm2(r);
699             invr  = gmx::invsqrt(r2);
700             invr2 = invr * invr;
701             ex    = ip[type].orires.ex;
702             power = ip[type].orires.power;
703             fc    = smooth_fc * ip[type].orires.kfac;
704             dev   = oriresdata->otav[restraintIndex] - ip[type].orires.obs;
705
706             /* NOTE:
707              * there is no real potential when time averaging is applied
708              */
709             vtot += 0.5 * fc * gmx::square(dev);
710
711             if (bTAV)
712             {
713                 /* Calculate the force as the sqrt of tav times instantaneous */
714                 devins = oriresdata->oins[restraintIndex] - ip[type].orires.obs;
715                 if (dev * devins <= 0)
716                 {
717                     dev = 0;
718                 }
719                 else
720                 {
721                     dev = std::sqrt(dev * devins);
722                     if (devins < 0)
723                     {
724                         dev = -dev;
725                     }
726                 }
727             }
728
729             pfac = fc * ip[type].orires.c * invr2;
730             for (int i = 0; i < power; i++)
731             {
732                 pfac *= invr;
733             }
734             mvmul(oriresdata->S[ex], r, Sr);
735             for (int i = 0; i < DIM; i++)
736             {
737                 fij[i] = -pfac * dev * (4 * Sr[i] - 2 * (2 + power) * invr2 * iprod(Sr, r) * r[i]);
738             }
739
740             for (int i = 0; i < DIM; i++)
741             {
742                 f[ai][i] += fij[i];
743                 f[aj][i] -= fij[i];
744                 if (fshift)
745                 {
746                     fshift[ki][i] += fij[i];
747                     fshift[gmx::c_centralShiftIndex][i] -= fij[i];
748                 }
749             }
750         }
751     }
752
753     return vtot;
754
755     /* Approx. 80*nfa/3 flops */
756 }
757
758 void update_orires_history(const t_oriresdata& od, history_t* hist)
759 {
760     if (od.edt != 0)
761     {
762         /* Copy the new time averages that have been calculated
763          *  in calc_orires_dev.
764          */
765         hist->orire_initf = od.exp_min_t_tau;
766         for (int pair = 0; pair < od.nr; pair++)
767         {
768             for (int i = 0; i < 5; i++)
769             {
770                 hist->orire_Dtav[pair * 5 + i] = od.Dtav[pair][i];
771             }
772         }
773     }
774 }