ecaaa2c85f5a0b1e92e2c8f8bd819e71afc5407d
[alexxy/gromacs.git] / src / gromacs / mdlib / constr.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 by 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 /*! \internal \file
39  * \brief Defines the high-level constraint code.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \author Mark Abraham <mark.j.abraham@gmail.com>
43  * \ingroup module_mdlib
44  */
45 #include "gmxpre.h"
46
47 #include "constr.h"
48
49 #include <cassert>
50 #include <cmath>
51 #include <cstdlib>
52
53 #include <algorithm>
54
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/essentialdynamics/edsam.h"
58 #include "gromacs/fileio/confio.h"
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/fileio/pdbio.h"
61 #include "gromacs/gmxlib/nrnb.h"
62 #include "gromacs/math/arrayrefwithpadding.h"
63 #include "gromacs/math/utilities.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/gmx_omp_nthreads.h"
66 #include "gromacs/mdlib/lincs.h"
67 #include "gromacs/mdlib/settle.h"
68 #include "gromacs/mdlib/shake.h"
69 #include "gromacs/mdtypes/commrec.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/md_enums.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/pulling/pull.h"
75 #include "gromacs/timing/wallcycle.h"
76 #include "gromacs/topology/ifunc.h"
77 #include "gromacs/topology/mtop_lookup.h"
78 #include "gromacs/topology/mtop_util.h"
79 #include "gromacs/utility/arrayref.h"
80 #include "gromacs/utility/exceptions.h"
81 #include "gromacs/utility/fatalerror.h"
82 #include "gromacs/utility/gmxassert.h"
83 #include "gromacs/utility/listoflists.h"
84 #include "gromacs/utility/pleasecite.h"
85 #include "gromacs/utility/txtdump.h"
86
87 namespace gmx
88 {
89
90 /* \brief Impl class for Constraints
91  *
92  * \todo Members like md, idef are valid only for the lifetime of a
93  * domain, which would be good to make clearer in the structure of the
94  * code. It should not be possible to call apply() if setConstraints()
95  * has not been called. For example, this could be achieved if
96  * setConstraints returned a valid object with such a method.  That
97  * still requires that we manage the lifetime of that object
98  * correctly, however. */
99 class Constraints::Impl
100 {
101 public:
102     Impl(const gmx_mtop_t&     mtop_p,
103          const t_inputrec&     ir_p,
104          pull_t*               pull_work,
105          FILE*                 log_p,
106          const t_commrec*      cr_p,
107          const gmx_multisim_t* ms,
108          t_nrnb*               nrnb,
109          gmx_wallcycle*        wcycle_p,
110          bool                  pbcHandlingRequired,
111          int                   numConstraints,
112          int                   numSettles);
113     ~Impl();
114     void setConstraints(gmx_localtop_t*                     top,
115                         int                                 numAtoms,
116                         int                                 numHomeAtoms,
117                         gmx::ArrayRef<const real>           masses,
118                         gmx::ArrayRef<const real>           inverseMasses,
119                         bool                                hasMassPerturbedAtoms,
120                         real                                lambda,
121                         gmx::ArrayRef<const unsigned short> cFREEZE);
122     bool apply(bool                      bLog,
123                bool                      bEner,
124                int64_t                   step,
125                int                       delta_step,
126                real                      step_scaling,
127                ArrayRefWithPadding<RVec> x,
128                ArrayRefWithPadding<RVec> xprime,
129                ArrayRef<RVec>            min_proj,
130                const matrix              box,
131                real                      lambda,
132                real*                     dvdlambda,
133                ArrayRefWithPadding<RVec> v,
134                bool                      computeVirial,
135                tensor                    constraintsVirial,
136                ConstraintVariable        econq);
137     //! The total number of constraints.
138     int ncon_tot = 0;
139     //! The number of flexible constraints.
140     int nflexcon = 0;
141     //! A list of atoms to constraints for each moleculetype.
142     std::vector<ListOfLists<int>> at2con_mt;
143     //! A list of atoms to settles for each moleculetype
144     std::vector<std::vector<int>> at2settle_mt;
145     //! LINCS data.
146     Lincs* lincsd = nullptr; // TODO this should become a unique_ptr
147     //! SHAKE data.
148     std::unique_ptr<shakedata> shaked;
149     //! SETTLE data.
150     std::unique_ptr<SettleData> settled;
151     //! The maximum number of warnings.
152     int maxwarn = 0;
153     //! The number of warnings for LINCS.
154     int warncount_lincs = 0;
155     //! The number of warnings for SETTLE.
156     int warncount_settle = 0;
157     //! The essential dynamics data.
158     gmx_edsam* ed = nullptr;
159
160     //! Thread-local virial contribution.
161     tensor* threadConstraintsVirial = { nullptr };
162     //! Did a settle error occur?
163     bool* bSettleErrorHasOccurred = nullptr;
164
165     //! Pointer to the global topology - only used for printing warnings.
166     const gmx_mtop_t& mtop;
167     //! Parameters for the interactions in this domain.
168     const InteractionDefinitions* idef = nullptr;
169     //! Total number of atoms.
170     int numAtoms_ = 0;
171     //! Number of local atoms.
172     int numHomeAtoms_ = 0;
173     //! Atoms masses.
174     gmx::ArrayRef<const real> masses_;
175     //! Inverse masses.
176     gmx::ArrayRef<const real> inverseMasses_;
177     //! If there are atoms with perturbed mass (for FEP).
178     bool hasMassPerturbedAtoms_;
179     //! FEP lambda value.
180     real lambda_;
181     //! Freeze groups data
182     gmx::ArrayRef<const unsigned short> cFREEZE_;
183     //! Whether we need to do pbc for handling bonds.
184     bool pbcHandlingRequired_ = false;
185
186     //! Logging support.
187     FILE* log = nullptr;
188     //! Communication support.
189     const t_commrec* cr = nullptr;
190     //! Multi-sim support.
191     const gmx_multisim_t* ms = nullptr;
192     //! Pulling code object, if any.
193     pull_t* pull_work = nullptr;
194     /*!\brief Input options.
195      *
196      * \todo Replace with IMdpOptions */
197     const t_inputrec& ir;
198     //! Flop counting support.
199     t_nrnb* nrnb = nullptr;
200     //! Tracks wallcycle usage.
201     gmx_wallcycle* wcycle;
202 };
203
204 Constraints::~Constraints() = default;
205
206 int Constraints::numFlexibleConstraints() const
207 {
208     return impl_->nflexcon;
209 }
210
211 bool Constraints::havePerturbedConstraints() const
212 {
213     const gmx_ffparams_t& ffparams = impl_->mtop.ffparams;
214
215     for (size_t i = 0; i < ffparams.functype.size(); i++)
216     {
217         if ((ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
218             && ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
219         {
220             return true;
221         }
222     }
223
224     return false;
225 }
226
227 //! Clears constraint quantities for atoms in nonlocal region.
228 static void clear_constraint_quantity_nonlocal(const gmx_domdec_t& dd, ArrayRef<RVec> q)
229 {
230     int nonlocal_at_start, nonlocal_at_end;
231
232     dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
233
234     for (int at = nonlocal_at_start; at < nonlocal_at_end; at++)
235     {
236         clear_rvec(q[at]);
237     }
238 }
239
240 void too_many_constraint_warnings(ConstraintAlgorithm eConstrAlg, int warncount)
241 {
242     gmx_fatal(FARGS,
243               "Too many %s warnings (%d)\n"
244               "If you know what you are doing you can %s"
245               "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
246               "but normally it is better to fix the problem",
247               (eConstrAlg == ConstraintAlgorithm::Lincs) ? "LINCS" : "SETTLE",
248               warncount,
249               (eConstrAlg == ConstraintAlgorithm::Lincs)
250                       ? "adjust the lincs warning threshold in your mdp file\nor "
251                       : "\n");
252 }
253
254 //! Writes out coordinates.
255 static void write_constr_pdb(const char*          fn,
256                              const char*          title,
257                              const gmx_mtop_t&    mtop,
258                              int                  start,
259                              int                  homenr,
260                              const t_commrec*     cr,
261                              ArrayRef<const RVec> x,
262                              const matrix         box)
263 {
264     char          fname[STRLEN];
265     FILE*         out;
266     int           dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
267     gmx_domdec_t* dd;
268     const char *  anm, *resnm;
269
270     dd = nullptr;
271     if (DOMAINDECOMP(cr))
272     {
273         dd = cr->dd;
274         dd_get_constraint_range(*dd, &dd_ac0, &dd_ac1);
275         start  = 0;
276         homenr = dd_ac1;
277     }
278
279     if (PAR(cr))
280     {
281         sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
282     }
283     else
284     {
285         sprintf(fname, "%s.pdb", fn);
286     }
287
288     out = gmx_fio_fopen(fname, "w");
289
290     fprintf(out, "TITLE     %s\n", title);
291     gmx_write_pdb_box(out, PbcType::Unset, box);
292     int molb = 0;
293     for (i = start; i < start + homenr; i++)
294     {
295         if (dd != nullptr)
296         {
297             if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
298             {
299                 continue;
300             }
301             ii = dd->globalAtomIndices[i];
302         }
303         else
304         {
305             ii = i;
306         }
307         mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
308         gmx_fprintf_pdb_atomline(out,
309                                  PdbRecordType::Atom,
310                                  ii + 1,
311                                  anm,
312                                  ' ',
313                                  resnm,
314                                  ' ',
315                                  resnr,
316                                  ' ',
317                                  10 * x[i][XX],
318                                  10 * x[i][YY],
319                                  10 * x[i][ZZ],
320                                  1.0,
321                                  0.0,
322                                  "");
323     }
324     fprintf(out, "TER\n");
325
326     gmx_fio_fclose(out);
327 }
328
329 //! Writes out domain contents to help diagnose crashes.
330 static void dump_confs(FILE*                log,
331                        int64_t              step,
332                        const gmx_mtop_t&    mtop,
333                        int                  start,
334                        int                  homenr,
335                        const t_commrec*     cr,
336                        ArrayRef<const RVec> x,
337                        ArrayRef<const RVec> xprime,
338                        const matrix         box)
339 {
340     char buf[STRLEN], buf2[22];
341
342     char* env = getenv("GMX_SUPPRESS_DUMP");
343     if (env)
344     {
345         return;
346     }
347
348     sprintf(buf, "step%sb", gmx_step_str(step, buf2));
349     write_constr_pdb(buf, "initial coordinates", mtop, start, homenr, cr, x, box);
350     sprintf(buf, "step%sc", gmx_step_str(step, buf2));
351     write_constr_pdb(buf, "coordinates after constraining", mtop, start, homenr, cr, xprime, box);
352     if (log)
353     {
354         fprintf(log, "Wrote pdb files with previous and current coordinates\n");
355     }
356     fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
357 }
358
359 bool Constraints::apply(bool                      bLog,
360                         bool                      bEner,
361                         int64_t                   step,
362                         int                       delta_step,
363                         real                      step_scaling,
364                         ArrayRefWithPadding<RVec> x,
365                         ArrayRefWithPadding<RVec> xprime,
366                         ArrayRef<RVec>            min_proj,
367                         const matrix              box,
368                         real                      lambda,
369                         real*                     dvdlambda,
370                         ArrayRefWithPadding<RVec> v,
371                         bool                      computeVirial,
372                         tensor                    constraintsVirial,
373                         ConstraintVariable        econq)
374 {
375     return impl_->apply(bLog,
376                         bEner,
377                         step,
378                         delta_step,
379                         step_scaling,
380                         std::move(x),
381                         std::move(xprime),
382                         min_proj,
383                         box,
384                         lambda,
385                         dvdlambda,
386                         std::move(v),
387                         computeVirial,
388                         constraintsVirial,
389                         econq);
390 }
391
392 bool Constraints::Impl::apply(bool                      bLog,
393                               bool                      bEner,
394                               int64_t                   step,
395                               int                       delta_step,
396                               real                      step_scaling,
397                               ArrayRefWithPadding<RVec> x,
398                               ArrayRefWithPadding<RVec> xprime,
399                               ArrayRef<RVec>            min_proj,
400                               const matrix              box,
401                               real                      lambda,
402                               real*                     dvdlambda,
403                               ArrayRefWithPadding<RVec> v,
404                               bool                      computeVirial,
405                               tensor                    constraintsVirial,
406                               ConstraintVariable        econq)
407 {
408     bool  bOK, bDump;
409     int   start;
410     real  scaled_delta_t;
411     real  invdt, vir_fac = 0, t;
412     int   nsettle;
413     t_pbc pbc, *pbc_null;
414     char  buf[22];
415     int   nth;
416
417     wallcycle_start(wcycle, WallCycleCounter::Constr);
418
419     if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
420     {
421         gmx_incons(
422                 "constrain called for forces displacements while not doing energy minimization, "
423                 "can not do this while the LINCS and SETTLE constraint connection matrices are "
424                 "mass weighted");
425     }
426
427     bOK   = TRUE;
428     bDump = FALSE;
429
430     start = 0;
431
432     scaled_delta_t = step_scaling * ir.delta_t;
433
434     /* Prepare time step for use in constraint implementations, and
435        avoid generating inf when ir.delta_t = 0. */
436     if (ir.delta_t == 0)
437     {
438         invdt = 0.0;
439     }
440     else
441     {
442         invdt = 1.0 / scaled_delta_t;
443     }
444
445     if (ir.efep != FreeEnergyPerturbationType::No && EI_DYNAMICS(ir.eI))
446     {
447         /* Set the constraint lengths for the step at which this configuration
448          * is meant to be. The invmasses should not be changed.
449          */
450         lambda += delta_step * ir.fepvals->delta_lambda;
451     }
452
453     if (computeVirial)
454     {
455         clear_mat(constraintsVirial);
456     }
457     const InteractionList& settle = idef->il[F_SETTLE];
458     nsettle                       = settle.size() / (1 + NRAL(F_SETTLE));
459
460     if (nsettle > 0)
461     {
462         nth = gmx_omp_nthreads_get(ModuleMultiThread::Settle);
463     }
464     else
465     {
466         nth = 1;
467     }
468
469     /* We do not need full pbc when constraints do not cross update groups
470      * i.e. when dd->constraint_comm==NULL.
471      * Note that PBC for constraints is different from PBC for bondeds.
472      * For constraints there is both forward and backward communication.
473      */
474     if (ir.pbcType != PbcType::No && (cr->dd || pbcHandlingRequired_)
475         && !(cr->dd && cr->dd->constraint_comm == nullptr))
476     {
477         /* With pbc=screw the screw has been changed to a shift
478          * by the constraint coordinate communication routine,
479          * so that here we can use normal pbc.
480          */
481         pbc_null = set_pbc_dd(&pbc, ir.pbcType, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, FALSE, box);
482     }
483     else
484     {
485         pbc_null = nullptr;
486     }
487
488     /* Communicate the coordinates required for the non-local constraints
489      * for LINCS and/or SETTLE.
490      */
491     if (cr->dd)
492     {
493         dd_move_x_constraints(cr->dd,
494                               box,
495                               x.unpaddedArrayRef(),
496                               xprime.unpaddedArrayRef(),
497                               econq == ConstraintVariable::Positions);
498
499         if (!v.empty())
500         {
501             /* We need to initialize the non-local components of v.
502              * We never actually use these values, but we do increment them,
503              * so we should avoid uninitialized variables and overflows.
504              */
505             clear_constraint_quantity_nonlocal(*cr->dd, v.unpaddedArrayRef());
506         }
507     }
508
509     if (lincsd != nullptr)
510     {
511         bOK = constrain_lincs(bLog || bEner,
512                               ir,
513                               step,
514                               lincsd,
515                               inverseMasses_,
516                               cr,
517                               ms,
518                               x,
519                               xprime,
520                               min_proj,
521                               box,
522                               pbc_null,
523                               hasMassPerturbedAtoms_,
524                               lambda,
525                               dvdlambda,
526                               invdt,
527                               v.unpaddedArrayRef(),
528                               computeVirial,
529                               constraintsVirial,
530                               econq,
531                               nrnb,
532                               maxwarn,
533                               &warncount_lincs);
534         if (!bOK && maxwarn < INT_MAX)
535         {
536             if (log != nullptr)
537             {
538                 fprintf(log,
539                         "Constraint error in algorithm %s at step %s\n",
540                         enumValueToString(ConstraintAlgorithm::Lincs),
541                         gmx_step_str(step, buf));
542             }
543             bDump = TRUE;
544         }
545     }
546
547     if (shaked != nullptr)
548     {
549         bOK = constrain_shake(log,
550                               shaked.get(),
551                               inverseMasses_,
552                               *idef,
553                               ir,
554                               x.unpaddedArrayRef(),
555                               xprime.unpaddedArrayRef(),
556                               min_proj,
557                               pbc_null,
558                               nrnb,
559                               lambda,
560                               dvdlambda,
561                               invdt,
562                               v.unpaddedArrayRef(),
563                               computeVirial,
564                               constraintsVirial,
565                               maxwarn < INT_MAX,
566                               econq);
567
568         if (!bOK && maxwarn < INT_MAX)
569         {
570             if (log != nullptr)
571             {
572                 fprintf(log,
573                         "Constraint error in algorithm %s at step %s\n",
574                         enumValueToString(ConstraintAlgorithm::Shake),
575                         gmx_step_str(step, buf));
576             }
577             bDump = TRUE;
578         }
579     }
580
581     if (nsettle > 0)
582     {
583         bool bSettleErrorHasOccurred0 = false;
584
585         switch (econq)
586         {
587             case ConstraintVariable::Positions:
588 #pragma omp parallel for num_threads(nth) schedule(static)
589                 for (int th = 0; th < nth; th++)
590                 {
591                     try
592                     {
593                         if (th > 0)
594                         {
595                             clear_mat(threadConstraintsVirial[th]);
596                         }
597
598                         csettle(*settled,
599                                 nth,
600                                 th,
601                                 pbc_null,
602                                 x,
603                                 xprime,
604                                 invdt,
605                                 v,
606                                 computeVirial,
607                                 th == 0 ? constraintsVirial : threadConstraintsVirial[th],
608                                 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
609                     }
610                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
611                 }
612                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
613                 if (!v.empty())
614                 {
615                     inc_nrnb(nrnb, eNR_CONSTR_V, nsettle * 3);
616                 }
617                 if (computeVirial)
618                 {
619                     inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle * 3);
620                 }
621                 break;
622             case ConstraintVariable::Velocities:
623             case ConstraintVariable::Derivative:
624             case ConstraintVariable::Force:
625             case ConstraintVariable::ForceDispl:
626 #pragma omp parallel for num_threads(nth) schedule(static)
627                 for (int th = 0; th < nth; th++)
628                 {
629                     try
630                     {
631                         int calcvir_atom_end;
632
633                         if (!computeVirial)
634                         {
635                             calcvir_atom_end = 0;
636                         }
637                         else
638                         {
639                             calcvir_atom_end = numHomeAtoms_;
640                         }
641
642                         if (th > 0)
643                         {
644                             clear_mat(threadConstraintsVirial[th]);
645                         }
646
647                         int start_th = (nsettle * th) / nth;
648                         int end_th   = (nsettle * (th + 1)) / nth;
649
650                         if (start_th >= 0 && end_th - start_th > 0)
651                         {
652                             settle_proj(*settled,
653                                         econq,
654                                         end_th - start_th,
655                                         settle.iatoms.data() + start_th * (1 + NRAL(F_SETTLE)),
656                                         pbc_null,
657                                         x.unpaddedArrayRef(),
658                                         xprime.unpaddedArrayRef(),
659                                         min_proj,
660                                         calcvir_atom_end,
661                                         th == 0 ? constraintsVirial : threadConstraintsVirial[th]);
662                         }
663                     }
664                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
665                 }
666                 /* This is an overestimate */
667                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
668                 break;
669             case ConstraintVariable::Deriv_FlexCon:
670                 /* Nothing to do, since the are no flexible constraints in settles */
671                 break;
672             default: gmx_incons("Unknown constraint quantity for settle");
673         }
674
675         if (computeVirial)
676         {
677             /* Reduce the virial contributions over the threads */
678             for (int th = 1; th < nth; th++)
679             {
680                 m_add(constraintsVirial, threadConstraintsVirial[th], constraintsVirial);
681             }
682         }
683
684         if (econq == ConstraintVariable::Positions)
685         {
686             for (int th = 1; th < nth; th++)
687             {
688                 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
689             }
690
691             if (bSettleErrorHasOccurred0)
692             {
693                 char buf[STRLEN];
694                 sprintf(buf,
695                         "\nstep "
696                         "%" PRId64
697                         ": One or more water molecules can not be settled.\n"
698                         "Check for bad contacts and/or reduce the timestep if appropriate.\n",
699                         step);
700                 if (log)
701                 {
702                     fprintf(log, "%s", buf);
703                 }
704                 fprintf(stderr, "%s", buf);
705                 warncount_settle++;
706                 if (warncount_settle > maxwarn)
707                 {
708                     too_many_constraint_warnings(ConstraintAlgorithm::Count, warncount_settle);
709                 }
710                 bDump = TRUE;
711
712                 bOK = FALSE;
713             }
714         }
715     }
716
717     if (computeVirial)
718     {
719         /* The normal uses of constrain() pass step_scaling = 1.0.
720          * The call to constrain() for SD1 that passes step_scaling =
721          * 0.5 also passes vir = NULL, so cannot reach this
722          * assertion. This assertion should remain until someone knows
723          * that this path works for their intended purpose, and then
724          * they can use scaled_delta_t instead of ir.delta_t
725          * below. */
726         assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
727         switch (econq)
728         {
729             case ConstraintVariable::Positions: vir_fac = 0.5 / (ir.delta_t * ir.delta_t); break;
730             case ConstraintVariable::Velocities: vir_fac = 0.5 / ir.delta_t; break;
731             case ConstraintVariable::Force:
732             case ConstraintVariable::ForceDispl: vir_fac = 0.5; break;
733             default: gmx_incons("Unsupported constraint quantity for virial");
734         }
735
736         if (EI_VV(ir.eI))
737         {
738             vir_fac *= 2; /* only constraining over half the distance here */
739         }
740         for (int i = 0; i < DIM; i++)
741         {
742             for (int j = 0; j < DIM; j++)
743             {
744                 constraintsVirial[i][j] *= vir_fac;
745             }
746         }
747     }
748
749     if (bDump)
750     {
751         dump_confs(log, step, mtop, start, numHomeAtoms_, cr, x.unpaddedArrayRef(), xprime.unpaddedArrayRef(), box);
752     }
753
754     if (econq == ConstraintVariable::Positions)
755     {
756         if (ir.bPull && pull_have_constraint(*pull_work))
757         {
758             if (EI_DYNAMICS(ir.eI))
759             {
760                 t = ir.init_t + (step + delta_step) * ir.delta_t;
761             }
762             else
763             {
764                 t = ir.init_t;
765             }
766             set_pbc(&pbc, ir.pbcType, box);
767             pull_constraint(pull_work,
768                             masses_,
769                             &pbc,
770                             cr,
771                             ir.delta_t,
772                             t,
773                             x.unpaddedArrayRef(),
774                             xprime.unpaddedArrayRef(),
775                             v.unpaddedArrayRef(),
776                             constraintsVirial);
777         }
778         if (ed && delta_step > 0)
779         {
780             /* apply the essential dynamics constraints here */
781             do_edsam(&ir, step, cr, xprime.unpaddedArrayRef(), v.unpaddedArrayRef(), box, ed);
782         }
783     }
784     wallcycle_stop(wcycle, WallCycleCounter::Constr);
785
786     const bool haveVelocities = (!v.empty() || econq == ConstraintVariable::Velocities);
787     if (haveVelocities && !cFREEZE_.empty())
788     {
789         /* Set the velocities of frozen dimensions to zero */
790         ArrayRef<RVec> vRef;
791         if (econq == ConstraintVariable::Velocities)
792         {
793             vRef = xprime.unpaddedArrayRef();
794         }
795         else
796         {
797             vRef = v.unpaddedArrayRef();
798         }
799
800         int gmx_unused numThreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
801
802 #pragma omp parallel for num_threads(numThreads) schedule(static)
803         for (int i = 0; i < numHomeAtoms_; i++)
804         {
805             int freezeGroup = cFREEZE_[i];
806
807             for (int d = 0; d < DIM; d++)
808             {
809                 if (ir.opts.nFreeze[freezeGroup][d])
810                 {
811                     vRef[i][d] = 0;
812                 }
813             }
814         }
815     }
816
817     return bOK;
818 } // namespace gmx
819
820 ArrayRef<real> Constraints::rmsdData() const
821 {
822     if (impl_->lincsd)
823     {
824         return lincs_rmsdData(impl_->lincsd);
825     }
826     else
827     {
828         return {};
829     }
830 }
831
832 real Constraints::rmsd() const
833 {
834     if (impl_->lincsd)
835     {
836         return lincs_rmsd(impl_->lincsd);
837     }
838     else
839     {
840         return 0;
841     }
842 }
843
844 int Constraints::numConstraintsTotal()
845 {
846     return impl_->ncon_tot;
847 }
848
849 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
850 {
851     if (haveDynamicsIntegrator)
852     {
853         return FlexibleConstraintTreatment::Include;
854     }
855     else
856     {
857         return FlexibleConstraintTreatment::Exclude;
858     }
859 }
860
861 /*! \brief Returns a block struct to go from atoms to constraints
862  *
863  * The block struct will contain constraint indices with lower indices
864  * directly matching the order in F_CONSTR and higher indices matching
865  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
866  *
867  * \param[in]  numAtoms  The number of atoms to construct the list for
868  * \param[in]  ilists    The interaction lists, size F_NRE
869  * \param[in]  iparams   Interaction parameters, can be null when
870  *                       \p flexibleConstraintTreatment==Include
871  * \param[in]  flexibleConstraintTreatment  The flexible constraint treatment,
872  *                                          see enum above
873  *
874  * \returns a block struct with all constraints for each atom
875  */
876 static ListOfLists<int> makeAtomsToConstraintsList(int                             numAtoms,
877                                                    ArrayRef<const InteractionList> ilists,
878                                                    ArrayRef<const t_iparams>       iparams,
879                                                    FlexibleConstraintTreatment flexibleConstraintTreatment)
880 {
881     GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || !iparams.empty(),
882                "With flexible constraint detection we need valid iparams");
883
884     std::vector<int> count(numAtoms);
885
886     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
887     {
888         const InteractionList& ilist  = ilists[ftype];
889         const int              stride = 1 + NRAL(ftype);
890         for (int i = 0; i < ilist.size(); i += stride)
891         {
892             if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
893                 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
894             {
895                 for (int j = 1; j < 3; j++)
896                 {
897                     int a = ilist.iatoms[i + j];
898                     count[a]++;
899                 }
900             }
901         }
902     }
903
904     std::vector<int> listRanges(numAtoms + 1);
905     for (int a = 0; a < numAtoms; a++)
906     {
907         listRanges[a + 1] = listRanges[a] + count[a];
908         count[a]          = 0;
909     }
910     std::vector<int> elements(listRanges[numAtoms]);
911
912     /* The F_CONSTRNC constraints have constraint numbers
913      * that continue after the last F_CONSTR constraint.
914      */
915     int numConstraints = 0;
916     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
917     {
918         const InteractionList& ilist  = ilists[ftype];
919         const int              stride = 1 + NRAL(ftype);
920         for (int i = 0; i < ilist.size(); i += stride)
921         {
922             if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
923                 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
924             {
925                 for (int j = 1; j < 3; j++)
926                 {
927                     const int a                          = ilist.iatoms[i + j];
928                     elements[listRanges[a] + count[a]++] = numConstraints;
929                 }
930             }
931             numConstraints++;
932         }
933     }
934
935     return ListOfLists<int>(std::move(listRanges), std::move(elements));
936 }
937
938 ListOfLists<int> make_at2con(int                             numAtoms,
939                              ArrayRef<const InteractionList> ilist,
940                              ArrayRef<const t_iparams>       iparams,
941                              FlexibleConstraintTreatment     flexibleConstraintTreatment)
942 {
943     return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
944 }
945
946 ListOfLists<int> make_at2con(const gmx_moltype_t&           moltype,
947                              gmx::ArrayRef<const t_iparams> iparams,
948                              FlexibleConstraintTreatment    flexibleConstraintTreatment)
949 {
950     return makeAtomsToConstraintsList(
951             moltype.atoms.nr, makeConstArrayRef(moltype.ilist), iparams, flexibleConstraintTreatment);
952 }
953
954 //! Return the number of flexible constraints in the \c ilist and \c iparams.
955 int countFlexibleConstraints(ArrayRef<const InteractionList> ilist, ArrayRef<const t_iparams> iparams)
956 {
957     int nflexcon = 0;
958     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
959     {
960         const int numIatomsPerConstraint = 3;
961         for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
962         {
963             const int type = ilist[ftype].iatoms[i];
964             if (iparams[type].constr.dA == 0 && iparams[type].constr.dB == 0)
965             {
966                 nflexcon++;
967             }
968         }
969     }
970
971     return nflexcon;
972 }
973
974 //! Returns the index of the settle to which each atom belongs.
975 static std::vector<int> make_at2settle(int natoms, const InteractionList& ilist)
976 {
977     /* Set all to no settle */
978     std::vector<int> at2s(natoms, -1);
979
980     const int stride = 1 + NRAL(F_SETTLE);
981
982     for (int s = 0; s < ilist.size(); s += stride)
983     {
984         at2s[ilist.iatoms[s + 1]] = s / stride;
985         at2s[ilist.iatoms[s + 2]] = s / stride;
986         at2s[ilist.iatoms[s + 3]] = s / stride;
987     }
988
989     return at2s;
990 }
991
992 void Constraints::Impl::setConstraints(gmx_localtop_t*                     top,
993                                        int                                 numAtoms,
994                                        int                                 numHomeAtoms,
995                                        gmx::ArrayRef<const real>           masses,
996                                        gmx::ArrayRef<const real>           inverseMasses,
997                                        const bool                          hasMassPerturbedAtoms,
998                                        const real                          lambda,
999                                        gmx::ArrayRef<const unsigned short> cFREEZE)
1000 {
1001     numAtoms_              = numAtoms;
1002     numHomeAtoms_          = numHomeAtoms;
1003     masses_                = masses;
1004     inverseMasses_         = inverseMasses;
1005     hasMassPerturbedAtoms_ = hasMassPerturbedAtoms;
1006     lambda_                = lambda;
1007     cFREEZE_               = cFREEZE;
1008
1009     idef = &top->idef;
1010
1011     if (ncon_tot > 0)
1012     {
1013         /* With DD we might also need to call LINCS on a domain no constraints for
1014          * communicating coordinates to other nodes that do have constraints.
1015          */
1016         if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1017         {
1018             set_lincs(*idef, numAtoms_, inverseMasses_, lambda_, EI_DYNAMICS(ir.eI), cr, lincsd);
1019         }
1020         if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1021         {
1022             if (cr->dd)
1023             {
1024                 // We are using the local topology, so there are only
1025                 // F_CONSTR constraints.
1026                 GMX_RELEASE_ASSERT(idef->il[F_CONSTRNC].empty(),
1027                                    "Here we should not have no-connect constraints");
1028                 make_shake_sblock_dd(shaked.get(), idef->il[F_CONSTR]);
1029             }
1030             else
1031             {
1032                 make_shake_sblock_serial(shaked.get(), &top->idef, numAtoms_);
1033             }
1034         }
1035     }
1036
1037     if (settled)
1038     {
1039         settled->setConstraints(idef->il[F_SETTLE], numHomeAtoms_, masses_, inverseMasses_);
1040     }
1041
1042     /* Make a selection of the local atoms for essential dynamics */
1043     if (ed && cr->dd)
1044     {
1045         dd_make_local_ed_indices(cr->dd, ed);
1046     }
1047 }
1048
1049 void Constraints::setConstraints(gmx_localtop_t*                     top,
1050                                  const int                           numAtoms,
1051                                  const int                           numHomeAtoms,
1052                                  gmx::ArrayRef<const real>           masses,
1053                                  gmx::ArrayRef<const real>           inverseMasses,
1054                                  const bool                          hasMassPerturbedAtoms,
1055                                  const real                          lambda,
1056                                  gmx::ArrayRef<const unsigned short> cFREEZE)
1057 {
1058     impl_->setConstraints(
1059             top, numAtoms, numHomeAtoms, masses, inverseMasses, hasMassPerturbedAtoms, lambda, cFREEZE);
1060 }
1061
1062 /*! \brief Makes a per-moleculetype container of mappings from atom
1063  * indices to constraint indices.
1064  *
1065  * Note that flexible constraints are only enabled with a dynamical integrator. */
1066 static std::vector<ListOfLists<int>> makeAtomToConstraintMappings(const gmx_mtop_t& mtop,
1067                                                                   FlexibleConstraintTreatment flexibleConstraintTreatment)
1068 {
1069     std::vector<ListOfLists<int>> mapping;
1070     mapping.reserve(mtop.moltype.size());
1071     for (const gmx_moltype_t& moltype : mtop.moltype)
1072     {
1073         mapping.push_back(make_at2con(moltype, mtop.ffparams.iparams, flexibleConstraintTreatment));
1074     }
1075     return mapping;
1076 }
1077
1078 Constraints::Constraints(const gmx_mtop_t&     mtop,
1079                          const t_inputrec&     ir,
1080                          pull_t*               pull_work,
1081                          FILE*                 log,
1082                          const t_commrec*      cr,
1083                          const gmx_multisim_t* ms,
1084                          t_nrnb*               nrnb,
1085                          gmx_wallcycle*        wcycle,
1086                          bool                  pbcHandlingRequired,
1087                          int                   numConstraints,
1088                          int                   numSettles) :
1089     impl_(new Impl(mtop, ir, pull_work, log, cr, ms, nrnb, wcycle, pbcHandlingRequired, numConstraints, numSettles))
1090 {
1091 }
1092
1093 Constraints::Impl::Impl(const gmx_mtop_t&     mtop_p,
1094                         const t_inputrec&     ir_p,
1095                         pull_t*               pull_work,
1096                         FILE*                 log_p,
1097                         const t_commrec*      cr_p,
1098                         const gmx_multisim_t* ms_p,
1099                         t_nrnb*               nrnb_p,
1100                         gmx_wallcycle*        wcycle_p,
1101                         bool                  pbcHandlingRequired,
1102                         int                   numConstraints,
1103                         int                   numSettles) :
1104     ncon_tot(numConstraints),
1105     mtop(mtop_p),
1106     pbcHandlingRequired_(pbcHandlingRequired),
1107     log(log_p),
1108     cr(cr_p),
1109     ms(ms_p),
1110     pull_work(pull_work),
1111     ir(ir_p),
1112     nrnb(nrnb_p),
1113     wcycle(wcycle_p)
1114 {
1115     if (numConstraints + numSettles > 0 && ir.epc == PressureCoupling::Mttk)
1116     {
1117         gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1118     }
1119
1120     nflexcon = 0;
1121     if (numConstraints > 0)
1122     {
1123         at2con_mt = makeAtomToConstraintMappings(mtop, flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
1124
1125         for (const gmx_molblock_t& molblock : mtop.molblock)
1126         {
1127             int count = countFlexibleConstraints(mtop.moltype[molblock.type].ilist, mtop.ffparams.iparams);
1128             nflexcon += molblock.nmol * count;
1129         }
1130
1131         if (nflexcon > 0)
1132         {
1133             if (log)
1134             {
1135                 fprintf(log, "There are %d flexible constraints\n", nflexcon);
1136                 if (ir.fc_stepsize == 0)
1137                 {
1138                     fprintf(log,
1139                             "\n"
1140                             "WARNING: step size for flexible constraining = 0\n"
1141                             "         All flexible constraints will be rigid.\n"
1142                             "         Will try to keep all flexible constraints at their original "
1143                             "length,\n"
1144                             "         but the lengths may exhibit some drift.\n\n");
1145                     nflexcon = 0;
1146                 }
1147             }
1148             if (nflexcon > 0)
1149             {
1150                 please_cite(log, "Hess2002");
1151             }
1152         }
1153
1154         if (ir.eConstrAlg == ConstraintAlgorithm::Lincs)
1155         {
1156             lincsd = init_lincs(log,
1157                                 mtop,
1158                                 nflexcon,
1159                                 at2con_mt,
1160                                 DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd),
1161                                 ir.nLincsIter,
1162                                 ir.nProjOrder);
1163         }
1164
1165         if (ir.eConstrAlg == ConstraintAlgorithm::Shake)
1166         {
1167             if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
1168             {
1169                 gmx_fatal(FARGS,
1170                           "SHAKE is not supported with domain decomposition and constraint that "
1171                           "cross domain boundaries, use LINCS");
1172             }
1173             if (nflexcon)
1174             {
1175                 gmx_fatal(FARGS,
1176                           "For this system also velocities and/or forces need to be constrained, "
1177                           "this can not be done with SHAKE, you should select LINCS");
1178             }
1179             please_cite(log, "Ryckaert77a");
1180             if (ir.bShakeSOR)
1181             {
1182                 please_cite(log, "Barth95a");
1183             }
1184
1185             shaked = std::make_unique<shakedata>();
1186         }
1187     }
1188
1189     if (numSettles > 0)
1190     {
1191         please_cite(log, "Miyamoto92a");
1192
1193         settled = std::make_unique<SettleData>(mtop);
1194
1195         /* Make an atom to settle index for use in domain decomposition */
1196         for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1197         {
1198             at2settle_mt.emplace_back(
1199                     make_at2settle(mtop.moltype[mt].atoms.nr, mtop.moltype[mt].ilist[F_SETTLE]));
1200         }
1201
1202         /* Allocate thread-local work arrays */
1203         int nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Settle);
1204         if (nthreads > 1 && threadConstraintsVirial == nullptr)
1205         {
1206             snew(threadConstraintsVirial, nthreads);
1207             snew(bSettleErrorHasOccurred, nthreads);
1208         }
1209     }
1210
1211     maxwarn   = 999;
1212     char* env = getenv("GMX_MAXCONSTRWARN");
1213     if (env)
1214     {
1215         maxwarn = 0;
1216         sscanf(env, "%8d", &maxwarn);
1217         if (maxwarn < 0)
1218         {
1219             maxwarn = INT_MAX;
1220         }
1221         if (log)
1222         {
1223             fprintf(log, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1224         }
1225         if (MASTER(cr))
1226         {
1227             fprintf(stderr, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1228         }
1229     }
1230     warncount_lincs  = 0;
1231     warncount_settle = 0;
1232 }
1233
1234 Constraints::Impl::~Impl()
1235 {
1236     if (bSettleErrorHasOccurred != nullptr)
1237     {
1238         sfree(bSettleErrorHasOccurred);
1239     }
1240     if (threadConstraintsVirial != nullptr)
1241     {
1242         sfree(threadConstraintsVirial);
1243     }
1244     done_lincs(lincsd);
1245 }
1246
1247 void Constraints::saveEdsamPointer(gmx_edsam* ed)
1248 {
1249     impl_->ed = ed;
1250 }
1251
1252 ArrayRef<const ListOfLists<int>> Constraints::atom2constraints_moltype() const
1253 {
1254     return impl_->at2con_mt;
1255 }
1256
1257 ArrayRef<const std::vector<int>> Constraints::atom2settle_moltype() const
1258 {
1259     return impl_->at2settle_mt;
1260 }
1261
1262 void do_constrain_first(FILE*                     fplog,
1263                         gmx::Constraints*         constr,
1264                         const t_inputrec*         ir,
1265                         int                       numAtoms,
1266                         int                       numHomeAtoms,
1267                         ArrayRefWithPadding<RVec> x,
1268                         ArrayRefWithPadding<RVec> v,
1269                         const matrix              box,
1270                         real                      lambda)
1271 {
1272     int     i, m, start, end;
1273     int64_t step;
1274     real    dt = ir->delta_t;
1275     real    dvdl_dum;
1276
1277     PaddedVector<RVec> savex(numAtoms);
1278
1279     start = 0;
1280     end   = numHomeAtoms;
1281
1282     if (debug)
1283     {
1284         fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n", start, numHomeAtoms, end);
1285     }
1286     /* Do a first constrain to reset particles... */
1287     step = ir->init_step;
1288     if (fplog)
1289     {
1290         char buf[STEPSTRSIZE];
1291         fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n", gmx_step_str(step, buf));
1292     }
1293     dvdl_dum = 0;
1294
1295     bool needsLogging  = true;
1296     bool computeEnergy = false;
1297     bool computeVirial = false;
1298     /* constrain the current position */
1299     constr->apply(needsLogging,
1300                   computeEnergy,
1301                   step,
1302                   0,
1303                   1.0,
1304                   x,
1305                   x,
1306                   {},
1307                   box,
1308                   lambda,
1309                   &dvdl_dum,
1310                   {},
1311                   computeVirial,
1312                   nullptr,
1313                   gmx::ConstraintVariable::Positions);
1314     if (EI_VV(ir->eI))
1315     {
1316         /* constrain the inital velocity, and save it */
1317         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1318         constr->apply(needsLogging,
1319                       computeEnergy,
1320                       step,
1321                       0,
1322                       1.0,
1323                       x,
1324                       v,
1325                       v.unpaddedArrayRef(),
1326                       box,
1327                       lambda,
1328                       &dvdl_dum,
1329                       {},
1330                       computeVirial,
1331                       nullptr,
1332                       gmx::ConstraintVariable::Velocities);
1333     }
1334     /* constrain the inital velocities at t-dt/2 */
1335     if (EI_STATE_VELOCITY(ir->eI) && ir->eI != IntegrationAlgorithm::VV)
1336     {
1337         auto subX = x.paddedArrayRef().subArray(start, end);
1338         auto subV = v.paddedArrayRef().subArray(start, end);
1339         for (i = start; (i < end); i++)
1340         {
1341             for (m = 0; (m < DIM); m++)
1342             {
1343                 /* Reverse the velocity */
1344                 subV[i][m] = -subV[i][m];
1345                 /* Store the position at t-dt in buf */
1346                 savex[i][m] = subX[i][m] + dt * subV[i][m];
1347             }
1348         }
1349         /* Shake the positions at t=-dt with the positions at t=0
1350          * as reference coordinates.
1351          */
1352         if (fplog)
1353         {
1354             char buf[STEPSTRSIZE];
1355             fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n", gmx_step_str(step, buf));
1356         }
1357         dvdl_dum = 0;
1358         constr->apply(needsLogging,
1359                       computeEnergy,
1360                       step,
1361                       -1,
1362                       1.0,
1363                       x,
1364                       savex.arrayRefWithPadding(),
1365                       {},
1366                       box,
1367                       lambda,
1368                       &dvdl_dum,
1369                       v,
1370                       computeVirial,
1371                       nullptr,
1372                       gmx::ConstraintVariable::Positions);
1373
1374         for (i = start; i < end; i++)
1375         {
1376             for (m = 0; m < DIM; m++)
1377             {
1378                 /* Re-reverse the velocities */
1379                 subV[i][m] = -subV[i][m];
1380             }
1381         }
1382     }
1383 }
1384
1385 void constrain_velocities(gmx::Constraints* constr,
1386                           bool              do_log,
1387                           bool              do_ene,
1388                           int64_t           step,
1389                           t_state*          state,
1390                           real*             dvdlambda,
1391                           gmx_bool          computeVirial,
1392                           tensor            constraintsVirial)
1393 {
1394     if (constr != nullptr)
1395     {
1396         constr->apply(do_log,
1397                       do_ene,
1398                       step,
1399                       1,
1400                       1.0,
1401                       state->x.arrayRefWithPadding(),
1402                       state->v.arrayRefWithPadding(),
1403                       state->v.arrayRefWithPadding().unpaddedArrayRef(),
1404                       state->box,
1405                       state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1406                       dvdlambda,
1407                       ArrayRefWithPadding<RVec>(),
1408                       computeVirial,
1409                       constraintsVirial,
1410                       ConstraintVariable::Velocities);
1411     }
1412 }
1413
1414 void constrain_coordinates(gmx::Constraints*         constr,
1415                            bool                      do_log,
1416                            bool                      do_ene,
1417                            int64_t                   step,
1418                            t_state*                  state,
1419                            ArrayRefWithPadding<RVec> xp,
1420                            real*                     dvdlambda,
1421                            gmx_bool                  computeVirial,
1422                            tensor                    constraintsVirial)
1423 {
1424     if (constr != nullptr)
1425     {
1426         constr->apply(do_log,
1427                       do_ene,
1428                       step,
1429                       1,
1430                       1.0,
1431                       state->x.arrayRefWithPadding(),
1432                       std::move(xp),
1433                       ArrayRef<RVec>(),
1434                       state->box,
1435                       state->lambda[FreeEnergyPerturbationCouplingType::Bonded],
1436                       dvdlambda,
1437                       state->v.arrayRefWithPadding(),
1438                       computeVirial,
1439                       constraintsVirial,
1440                       ConstraintVariable::Positions);
1441     }
1442 }
1443
1444 } // namespace gmx