Apply clang-format to source tree
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \internal \file
38  * \brief Defines the high-level constraint code.
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_mdlib
43  */
44 #include "gmxpre.h"
45
46 #include "constr.h"
47
48 #include <cassert>
49 #include <cmath>
50 #include <cstdlib>
51
52 #include <algorithm>
53
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/essentialdynamics/edsam.h"
57 #include "gromacs/fileio/confio.h"
58 #include "gromacs/fileio/gmxfio.h"
59 #include "gromacs/fileio/pdbio.h"
60 #include "gromacs/gmxlib/nrnb.h"
61 #include "gromacs/math/utilities.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdlib/lincs.h"
65 #include "gromacs/mdlib/settle.h"
66 #include "gromacs/mdlib/shake.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/pulling/pull.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/topology/block.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/pleasecite.h"
84 #include "gromacs/utility/smalloc.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_mdatoms&      md_p,
107          const t_commrec*      cr_p,
108          const gmx_multisim_t* ms,
109          t_nrnb*               nrnb,
110          gmx_wallcycle*        wcycle_p,
111          bool                  pbcHandlingRequired,
112          int                   numConstraints,
113          int                   numSettles);
114     ~Impl();
115     void setConstraints(const gmx_localtop_t& top, const t_mdatoms& md);
116     bool apply(bool               bLog,
117                bool               bEner,
118                int64_t            step,
119                int                delta_step,
120                real               step_scaling,
121                rvec*              x,
122                rvec*              xprime,
123                rvec*              min_proj,
124                const matrix       box,
125                real               lambda,
126                real*              dvdlambda,
127                rvec*              v,
128                tensor*            vir,
129                ConstraintVariable econq);
130     //! The total number of constraints.
131     int ncon_tot = 0;
132     //! The number of flexible constraints.
133     int nflexcon = 0;
134     //! A list of atoms to constraints for each moleculetype.
135     std::vector<t_blocka> at2con_mt;
136     //! A list of atoms to settles for each moleculetype
137     std::vector<std::vector<int>> at2settle_mt;
138     //! LINCS data.
139     Lincs* lincsd = nullptr; // TODO this should become a unique_ptr
140     //! SHAKE data.
141     shakedata* shaked = nullptr;
142     //! SETTLE data.
143     settledata* settled = nullptr;
144     //! The maximum number of warnings.
145     int maxwarn = 0;
146     //! The number of warnings for LINCS.
147     int warncount_lincs = 0;
148     //! The number of warnings for SETTLE.
149     int warncount_settle = 0;
150     //! The essential dynamics data.
151     gmx_edsam* ed = nullptr;
152
153     //! Thread-local virial contribution.
154     tensor* vir_r_m_dr_th = { nullptr };
155     //! Did a settle error occur?
156     bool* bSettleErrorHasOccurred = nullptr;
157
158     //! Pointer to the global topology - only used for printing warnings.
159     const gmx_mtop_t& mtop;
160     //! Parameters for the interactions in this domain.
161     const t_idef* idef = nullptr;
162     //! Data about atoms in this domain.
163     const t_mdatoms& md;
164     //! Whether we need to do pbc for handling bonds.
165     bool pbcHandlingRequired_ = false;
166
167     //! Logging support.
168     FILE* log = nullptr;
169     //! Communication support.
170     const t_commrec* cr = nullptr;
171     //! Multi-sim support.
172     const gmx_multisim_t* ms = nullptr;
173     //! Pulling code object, if any.
174     pull_t* pull_work = nullptr;
175     /*!\brief Input options.
176      *
177      * \todo Replace with IMdpOptions */
178     const t_inputrec& ir;
179     //! Flop counting support.
180     t_nrnb* nrnb = nullptr;
181     //! Tracks wallcycle usage.
182     gmx_wallcycle* wcycle;
183 };
184
185 Constraints::~Constraints() = default;
186
187 int Constraints::numFlexibleConstraints() const
188 {
189     return impl_->nflexcon;
190 }
191
192 bool Constraints::havePerturbedConstraints() const
193 {
194     const gmx_ffparams_t& ffparams = impl_->mtop.ffparams;
195
196     for (size_t i = 0; i < ffparams.functype.size(); i++)
197     {
198         if ((ffparams.functype[i] == F_CONSTR || ffparams.functype[i] == F_CONSTRNC)
199             && ffparams.iparams[i].constr.dA != ffparams.iparams[i].constr.dB)
200         {
201             return true;
202         }
203     }
204
205     return false;
206 }
207
208 //! Clears constraint quantities for atoms in nonlocal region.
209 static void clear_constraint_quantity_nonlocal(gmx_domdec_t* dd, rvec* q)
210 {
211     int nonlocal_at_start, nonlocal_at_end, at;
212
213     dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
214
215     for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
216     {
217         clear_rvec(q[at]);
218     }
219 }
220
221 void too_many_constraint_warnings(int eConstrAlg, int warncount)
222 {
223     gmx_fatal(
224             FARGS,
225             "Too many %s warnings (%d)\n"
226             "If you know what you are doing you can %s"
227             "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
228             "but normally it is better to fix the problem",
229             (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
230             (eConstrAlg == econtLINCS) ? "adjust the lincs warning threshold in your mdp file\nor " : "\n");
231 }
232
233 //! Writes out coordinates.
234 static void write_constr_pdb(const char*       fn,
235                              const char*       title,
236                              const gmx_mtop_t& mtop,
237                              int               start,
238                              int               homenr,
239                              const t_commrec*  cr,
240                              const rvec        x[],
241                              const matrix      box)
242 {
243     char          fname[STRLEN];
244     FILE*         out;
245     int           dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
246     gmx_domdec_t* dd;
247     const char *  anm, *resnm;
248
249     dd = nullptr;
250     if (DOMAINDECOMP(cr))
251     {
252         dd = cr->dd;
253         dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
254         start  = 0;
255         homenr = dd_ac1;
256     }
257
258     if (PAR(cr))
259     {
260         sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
261     }
262     else
263     {
264         sprintf(fname, "%s.pdb", fn);
265     }
266
267     out = gmx_fio_fopen(fname, "w");
268
269     fprintf(out, "TITLE     %s\n", title);
270     gmx_write_pdb_box(out, -1, box);
271     int molb = 0;
272     for (i = start; i < start + homenr; i++)
273     {
274         if (dd != nullptr)
275         {
276             if (i >= dd_numHomeAtoms(*dd) && i < dd_ac0)
277             {
278                 continue;
279             }
280             ii = dd->globalAtomIndices[i];
281         }
282         else
283         {
284             ii = i;
285         }
286         mtopGetAtomAndResidueName(mtop, ii, &molb, &anm, &resnr, &resnm, nullptr);
287         gmx_fprintf_pdb_atomline(out, epdbATOM, ii + 1, anm, ' ', resnm, ' ', resnr, ' ',
288                                  10 * x[i][XX], 10 * x[i][YY], 10 * x[i][ZZ], 1.0, 0.0, "");
289     }
290     fprintf(out, "TER\n");
291
292     gmx_fio_fclose(out);
293 }
294
295 //! Writes out domain contents to help diagnose crashes.
296 static void dump_confs(FILE*             log,
297                        int64_t           step,
298                        const gmx_mtop_t& mtop,
299                        int               start,
300                        int               homenr,
301                        const t_commrec*  cr,
302                        const rvec        x[],
303                        rvec              xprime[],
304                        const matrix      box)
305 {
306     char buf[STRLEN], buf2[22];
307
308     char* env = getenv("GMX_SUPPRESS_DUMP");
309     if (env)
310     {
311         return;
312     }
313
314     sprintf(buf, "step%sb", gmx_step_str(step, buf2));
315     write_constr_pdb(buf, "initial coordinates", mtop, start, homenr, cr, x, box);
316     sprintf(buf, "step%sc", gmx_step_str(step, buf2));
317     write_constr_pdb(buf, "coordinates after constraining", mtop, start, homenr, cr, xprime, box);
318     if (log)
319     {
320         fprintf(log, "Wrote pdb files with previous and current coordinates\n");
321     }
322     fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
323 }
324
325 bool Constraints::apply(bool               bLog,
326                         bool               bEner,
327                         int64_t            step,
328                         int                delta_step,
329                         real               step_scaling,
330                         rvec*              x,
331                         rvec*              xprime,
332                         rvec*              min_proj,
333                         const matrix       box,
334                         real               lambda,
335                         real*              dvdlambda,
336                         rvec*              v,
337                         tensor*            vir,
338                         ConstraintVariable econq)
339 {
340     return impl_->apply(bLog, bEner, step, delta_step, step_scaling, x, xprime, min_proj, box,
341                         lambda, dvdlambda, v, vir, econq);
342 }
343
344 bool Constraints::Impl::apply(bool               bLog,
345                               bool               bEner,
346                               int64_t            step,
347                               int                delta_step,
348                               real               step_scaling,
349                               rvec*              x,
350                               rvec*              xprime,
351                               rvec*              min_proj,
352                               const matrix       box,
353                               real               lambda,
354                               real*              dvdlambda,
355                               rvec*              v,
356                               tensor*            vir,
357                               ConstraintVariable econq)
358 {
359     bool   bOK, bDump;
360     int    start, homenr;
361     tensor vir_r_m_dr;
362     real   scaled_delta_t;
363     real   invdt, vir_fac = 0, t;
364     int    nsettle;
365     t_pbc  pbc, *pbc_null;
366     char   buf[22];
367     int    nth;
368
369     wallcycle_start(wcycle, ewcCONSTR);
370
371     if (econq == ConstraintVariable::ForceDispl && !EI_ENERGY_MINIMIZATION(ir.eI))
372     {
373         gmx_incons(
374                 "constrain called for forces displacements while not doing energy minimization, "
375                 "can not do this while the LINCS and SETTLE constraint connection matrices are "
376                 "mass weighted");
377     }
378
379     bOK   = TRUE;
380     bDump = FALSE;
381
382     start  = 0;
383     homenr = md.homenr;
384
385     scaled_delta_t = step_scaling * ir.delta_t;
386
387     /* Prepare time step for use in constraint implementations, and
388        avoid generating inf when ir.delta_t = 0. */
389     if (ir.delta_t == 0)
390     {
391         invdt = 0.0;
392     }
393     else
394     {
395         invdt = 1.0 / scaled_delta_t;
396     }
397
398     if (ir.efep != efepNO && EI_DYNAMICS(ir.eI))
399     {
400         /* Set the constraint lengths for the step at which this configuration
401          * is meant to be. The invmasses should not be changed.
402          */
403         lambda += delta_step * ir.fepvals->delta_lambda;
404     }
405
406     if (vir != nullptr)
407     {
408         clear_mat(vir_r_m_dr);
409     }
410     const t_ilist* settle = &idef->il[F_SETTLE];
411     nsettle               = settle->nr / (1 + NRAL(F_SETTLE));
412
413     if (nsettle > 0)
414     {
415         nth = gmx_omp_nthreads_get(emntSETTLE);
416     }
417     else
418     {
419         nth = 1;
420     }
421
422     /* We do not need full pbc when constraints do not cross update groups
423      * i.e. when dd->constraint_comm==NULL.
424      * Note that PBC for constraints is different from PBC for bondeds.
425      * For constraints there is both forward and backward communication.
426      */
427     if (ir.ePBC != epbcNONE && (cr->dd || pbcHandlingRequired_)
428         && !(cr->dd && cr->dd->constraint_comm == nullptr))
429     {
430         /* With pbc=screw the screw has been changed to a shift
431          * by the constraint coordinate communication routine,
432          * so that here we can use normal pbc.
433          */
434         pbc_null = set_pbc_dd(&pbc, ir.ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr, FALSE, box);
435     }
436     else
437     {
438         pbc_null = nullptr;
439     }
440
441     /* Communicate the coordinates required for the non-local constraints
442      * for LINCS and/or SETTLE.
443      */
444     if (cr->dd)
445     {
446         dd_move_x_constraints(cr->dd, box, x, xprime, econq == ConstraintVariable::Positions);
447
448         if (v != nullptr)
449         {
450             /* We need to initialize the non-local components of v.
451              * We never actually use these values, but we do increment them,
452              * so we should avoid uninitialized variables and overflows.
453              */
454             clear_constraint_quantity_nonlocal(cr->dd, v);
455         }
456     }
457
458     if (lincsd != nullptr)
459     {
460         bOK = constrain_lincs(bLog || bEner, ir, step, lincsd, md, cr, ms, x, xprime, min_proj, box,
461                               pbc_null, lambda, dvdlambda, invdt, v, vir != nullptr, vir_r_m_dr,
462                               econq, nrnb, maxwarn, &warncount_lincs);
463         if (!bOK && maxwarn < INT_MAX)
464         {
465             if (log != nullptr)
466             {
467                 fprintf(log, "Constraint error in algorithm %s at step %s\n",
468                         econstr_names[econtLINCS], gmx_step_str(step, buf));
469             }
470             bDump = TRUE;
471         }
472     }
473
474     if (shaked != nullptr)
475     {
476         bOK = constrain_shake(log, shaked, md.invmass, *idef, ir, x, xprime, min_proj, nrnb, lambda,
477                               dvdlambda, invdt, v, vir != nullptr, vir_r_m_dr, maxwarn < INT_MAX, econq);
478
479         if (!bOK && maxwarn < INT_MAX)
480         {
481             if (log != nullptr)
482             {
483                 fprintf(log, "Constraint error in algorithm %s at step %s\n",
484                         econstr_names[econtSHAKE], gmx_step_str(step, buf));
485             }
486             bDump = TRUE;
487         }
488     }
489
490     if (nsettle > 0)
491     {
492         bool bSettleErrorHasOccurred0 = false;
493
494         switch (econq)
495         {
496             case ConstraintVariable::Positions:
497 #pragma omp parallel for num_threads(nth) schedule(static)
498                 for (int th = 0; th < nth; th++)
499                 {
500                     try
501                     {
502                         if (th > 0)
503                         {
504                             clear_mat(vir_r_m_dr_th[th]);
505                         }
506
507                         csettle(settled, nth, th, pbc_null, x[0], xprime[0], invdt, v ? v[0] : nullptr,
508                                 vir != nullptr, th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th],
509                                 th == 0 ? &bSettleErrorHasOccurred0 : &bSettleErrorHasOccurred[th]);
510                     }
511                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
512                 }
513                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
514                 if (v != nullptr)
515                 {
516                     inc_nrnb(nrnb, eNR_CONSTR_V, nsettle * 3);
517                 }
518                 if (vir != nullptr)
519                 {
520                     inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle * 3);
521                 }
522                 break;
523             case ConstraintVariable::Velocities:
524             case ConstraintVariable::Derivative:
525             case ConstraintVariable::Force:
526             case ConstraintVariable::ForceDispl:
527 #pragma omp parallel for num_threads(nth) schedule(static)
528                 for (int th = 0; th < nth; th++)
529                 {
530                     try
531                     {
532                         int calcvir_atom_end;
533
534                         if (vir == nullptr)
535                         {
536                             calcvir_atom_end = 0;
537                         }
538                         else
539                         {
540                             calcvir_atom_end = md.homenr;
541                         }
542
543                         if (th > 0)
544                         {
545                             clear_mat(vir_r_m_dr_th[th]);
546                         }
547
548                         int start_th = (nsettle * th) / nth;
549                         int end_th   = (nsettle * (th + 1)) / nth;
550
551                         if (start_th >= 0 && end_th - start_th > 0)
552                         {
553                             settle_proj(settled, econq, end_th - start_th,
554                                         settle->iatoms + start_th * (1 + NRAL(F_SETTLE)), pbc_null,
555                                         x, xprime, min_proj, calcvir_atom_end,
556                                         th == 0 ? vir_r_m_dr : vir_r_m_dr_th[th]);
557                         }
558                     }
559                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
560                 }
561                 /* This is an overestimate */
562                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
563                 break;
564             case ConstraintVariable::Deriv_FlexCon:
565                 /* Nothing to do, since the are no flexible constraints in settles */
566                 break;
567             default: gmx_incons("Unknown constraint quantity for settle");
568         }
569
570         if (vir != nullptr)
571         {
572             /* Reduce the virial contributions over the threads */
573             for (int th = 1; th < nth; th++)
574             {
575                 m_add(vir_r_m_dr, vir_r_m_dr_th[th], vir_r_m_dr);
576             }
577         }
578
579         if (econq == ConstraintVariable::Positions)
580         {
581             for (int th = 1; th < nth; th++)
582             {
583                 bSettleErrorHasOccurred0 = bSettleErrorHasOccurred0 || bSettleErrorHasOccurred[th];
584             }
585
586             if (bSettleErrorHasOccurred0)
587             {
588                 char buf[STRLEN];
589                 sprintf(buf,
590                         "\nstep "
591                         "%" PRId64
592                         ": One or more water molecules can not be settled.\n"
593                         "Check for bad contacts and/or reduce the timestep if appropriate.\n",
594                         step);
595                 if (log)
596                 {
597                     fprintf(log, "%s", buf);
598                 }
599                 fprintf(stderr, "%s", buf);
600                 warncount_settle++;
601                 if (warncount_settle > maxwarn)
602                 {
603                     too_many_constraint_warnings(-1, warncount_settle);
604                 }
605                 bDump = TRUE;
606
607                 bOK = FALSE;
608             }
609         }
610     }
611
612     if (vir != nullptr)
613     {
614         /* The normal uses of constrain() pass step_scaling = 1.0.
615          * The call to constrain() for SD1 that passes step_scaling =
616          * 0.5 also passes vir = NULL, so cannot reach this
617          * assertion. This assertion should remain until someone knows
618          * that this path works for their intended purpose, and then
619          * they can use scaled_delta_t instead of ir.delta_t
620          * below. */
621         assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
622         switch (econq)
623         {
624             case ConstraintVariable::Positions: vir_fac = 0.5 / (ir.delta_t * ir.delta_t); break;
625             case ConstraintVariable::Velocities: vir_fac = 0.5 / ir.delta_t; break;
626             case ConstraintVariable::Force:
627             case ConstraintVariable::ForceDispl: vir_fac = 0.5; break;
628             default: gmx_incons("Unsupported constraint quantity for virial");
629         }
630
631         if (EI_VV(ir.eI))
632         {
633             vir_fac *= 2; /* only constraining over half the distance here */
634         }
635         for (int i = 0; i < DIM; i++)
636         {
637             for (int j = 0; j < DIM; j++)
638             {
639                 (*vir)[i][j] = vir_fac * vir_r_m_dr[i][j];
640             }
641         }
642     }
643
644     if (bDump)
645     {
646         dump_confs(log, step, mtop, start, homenr, cr, x, xprime, box);
647     }
648
649     if (econq == ConstraintVariable::Positions)
650     {
651         if (ir.bPull && pull_have_constraint(pull_work))
652         {
653             if (EI_DYNAMICS(ir.eI))
654             {
655                 t = ir.init_t + (step + delta_step) * ir.delta_t;
656             }
657             else
658             {
659                 t = ir.init_t;
660             }
661             set_pbc(&pbc, ir.ePBC, box);
662             pull_constraint(pull_work, &md, &pbc, cr, ir.delta_t, t, x, xprime, v, *vir);
663         }
664         if (ed && delta_step > 0)
665         {
666             /* apply the essential dynamics constraints here */
667             do_edsam(&ir, step, cr, xprime, v, box, ed);
668         }
669     }
670     wallcycle_stop(wcycle, ewcCONSTR);
671
672     if (v != nullptr && md.cFREEZE)
673     {
674         /* Set the velocities of frozen dimensions to zero */
675
676         int gmx_unused numThreads = gmx_omp_nthreads_get(emntUpdate);
677
678 #pragma omp parallel for num_threads(numThreads) schedule(static)
679         for (int i = 0; i < md.homenr; i++)
680         {
681             int freezeGroup = md.cFREEZE[i];
682
683             for (int d = 0; d < DIM; d++)
684             {
685                 if (ir.opts.nFreeze[freezeGroup][d])
686                 {
687                     v[i][d] = 0;
688                 }
689             }
690         }
691     }
692
693     return bOK;
694 }
695
696 ArrayRef<real> Constraints::rmsdData() const
697 {
698     if (impl_->lincsd)
699     {
700         return lincs_rmsdData(impl_->lincsd);
701     }
702     else
703     {
704         return {};
705     }
706 }
707
708 real Constraints::rmsd() const
709 {
710     if (impl_->lincsd)
711     {
712         return lincs_rmsd(impl_->lincsd);
713     }
714     else
715     {
716         return 0;
717     }
718 }
719
720 int Constraints::numConstraintsTotal()
721 {
722     return impl_->ncon_tot;
723 }
724
725 FlexibleConstraintTreatment flexibleConstraintTreatment(bool haveDynamicsIntegrator)
726 {
727     if (haveDynamicsIntegrator)
728     {
729         return FlexibleConstraintTreatment::Include;
730     }
731     else
732     {
733         return FlexibleConstraintTreatment::Exclude;
734     }
735 }
736
737 /*! \brief Returns a block struct to go from atoms to constraints
738  *
739  * The block struct will contain constraint indices with lower indices
740  * directly matching the order in F_CONSTR and higher indices matching
741  * the order in F_CONSTRNC offset by the number of constraints in F_CONSTR.
742  *
743  * \param[in]  numAtoms  The number of atoms to construct the list for
744  * \param[in]  ilists    The interaction lists, size F_NRE
745  * \param[in]  iparams   Interaction parameters, can be null when
746  * flexibleConstraintTreatment=Include \param[in]  flexibleConstraintTreatment  The flexible
747  * constraint treatment, see enum above \returns a block struct with all constraints for each atom
748  */
749 template<typename T>
750 static t_blocka makeAtomsToConstraintsList(int                         numAtoms,
751                                            const T*                    ilists,
752                                            const t_iparams*            iparams,
753                                            FlexibleConstraintTreatment flexibleConstraintTreatment)
754 {
755     GMX_ASSERT(flexibleConstraintTreatment == FlexibleConstraintTreatment::Include || iparams != nullptr,
756                "With flexible constraint detection we need valid iparams");
757
758     std::vector<int> count(numAtoms);
759
760     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
761     {
762         const T&  ilist  = ilists[ftype];
763         const int stride = 1 + NRAL(ftype);
764         for (int i = 0; i < ilist.size(); i += stride)
765         {
766             if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
767                 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
768             {
769                 for (int j = 1; j < 3; j++)
770                 {
771                     int a = ilist.iatoms[i + j];
772                     count[a]++;
773                 }
774             }
775         }
776     }
777
778     t_blocka at2con;
779     at2con.nr           = numAtoms;
780     at2con.nalloc_index = at2con.nr + 1;
781     snew(at2con.index, at2con.nalloc_index);
782     at2con.index[0] = 0;
783     for (int a = 0; a < numAtoms; a++)
784     {
785         at2con.index[a + 1] = at2con.index[a] + count[a];
786         count[a]            = 0;
787     }
788     at2con.nra      = at2con.index[at2con.nr];
789     at2con.nalloc_a = at2con.nra;
790     snew(at2con.a, at2con.nalloc_a);
791
792     /* The F_CONSTRNC constraints have constraint numbers
793      * that continue after the last F_CONSTR constraint.
794      */
795     int numConstraints = 0;
796     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
797     {
798         const T&  ilist  = ilists[ftype];
799         const int stride = 1 + NRAL(ftype);
800         for (int i = 0; i < ilist.size(); i += stride)
801         {
802             if (flexibleConstraintTreatment == FlexibleConstraintTreatment::Include
803                 || !isConstraintFlexible(iparams, ilist.iatoms[i]))
804             {
805                 for (int j = 1; j < 3; j++)
806                 {
807                     int a                                  = ilist.iatoms[i + j];
808                     at2con.a[at2con.index[a] + count[a]++] = numConstraints;
809                 }
810             }
811             numConstraints++;
812         }
813     }
814
815     return at2con;
816 }
817
818 t_blocka make_at2con(int                         numAtoms,
819                      const t_ilist*              ilist,
820                      const t_iparams*            iparams,
821                      FlexibleConstraintTreatment flexibleConstraintTreatment)
822 {
823     return makeAtomsToConstraintsList(numAtoms, ilist, iparams, flexibleConstraintTreatment);
824 }
825
826 t_blocka make_at2con(const gmx_moltype_t&           moltype,
827                      gmx::ArrayRef<const t_iparams> iparams,
828                      FlexibleConstraintTreatment    flexibleConstraintTreatment)
829 {
830     return makeAtomsToConstraintsList(moltype.atoms.nr, moltype.ilist.data(), iparams.data(),
831                                       flexibleConstraintTreatment);
832 }
833
834 //! Return the number of flexible constraints in the \c ilist and \c iparams.
835 template<typename T>
836 static int countFlexibleConstraintsTemplate(const T* ilist, const t_iparams* iparams)
837 {
838     int nflexcon = 0;
839     for (int ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
840     {
841         const int numIatomsPerConstraint = 3;
842         for (int i = 0; i < ilist[ftype].size(); i += numIatomsPerConstraint)
843         {
844             const int type = ilist[ftype].iatoms[i];
845             if (iparams[type].constr.dA == 0 && iparams[type].constr.dB == 0)
846             {
847                 nflexcon++;
848             }
849         }
850     }
851
852     return nflexcon;
853 }
854
855 int countFlexibleConstraints(const t_ilist* ilist, const t_iparams* iparams)
856 {
857     return countFlexibleConstraintsTemplate(ilist, iparams);
858 }
859
860 //! Returns the index of the settle to which each atom belongs.
861 static std::vector<int> make_at2settle(int natoms, const InteractionList& ilist)
862 {
863     /* Set all to no settle */
864     std::vector<int> at2s(natoms, -1);
865
866     const int stride = 1 + NRAL(F_SETTLE);
867
868     for (int s = 0; s < ilist.size(); s += stride)
869     {
870         at2s[ilist.iatoms[s + 1]] = s / stride;
871         at2s[ilist.iatoms[s + 2]] = s / stride;
872         at2s[ilist.iatoms[s + 3]] = s / stride;
873     }
874
875     return at2s;
876 }
877
878 void Constraints::Impl::setConstraints(const gmx_localtop_t& top, const t_mdatoms& md)
879 {
880     idef = &top.idef;
881
882     if (ncon_tot > 0)
883     {
884         /* With DD we might also need to call LINCS on a domain no constraints for
885          * communicating coordinates to other nodes that do have constraints.
886          */
887         if (ir.eConstrAlg == econtLINCS)
888         {
889             set_lincs(top.idef, md, EI_DYNAMICS(ir.eI), cr, lincsd);
890         }
891         if (ir.eConstrAlg == econtSHAKE)
892         {
893             if (cr->dd)
894             {
895                 // We are using the local topology, so there are only
896                 // F_CONSTR constraints.
897                 make_shake_sblock_dd(shaked, &idef->il[F_CONSTR], cr->dd);
898             }
899             else
900             {
901                 make_shake_sblock_serial(shaked, idef, md);
902             }
903         }
904     }
905
906     if (settled)
907     {
908         settle_set_constraints(settled, &idef->il[F_SETTLE], md);
909     }
910
911     /* Make a selection of the local atoms for essential dynamics */
912     if (ed && cr->dd)
913     {
914         dd_make_local_ed_indices(cr->dd, ed);
915     }
916 }
917
918 void Constraints::setConstraints(const gmx_localtop_t& top, const t_mdatoms& md)
919 {
920     impl_->setConstraints(top, md);
921 }
922
923 /*! \brief Makes a per-moleculetype container of mappings from atom
924  * indices to constraint indices.
925  *
926  * Note that flexible constraints are only enabled with a dynamical integrator. */
927 static std::vector<t_blocka> makeAtomToConstraintMappings(const gmx_mtop_t& mtop,
928                                                           FlexibleConstraintTreatment flexibleConstraintTreatment)
929 {
930     std::vector<t_blocka> mapping;
931     mapping.reserve(mtop.moltype.size());
932     for (const gmx_moltype_t& moltype : mtop.moltype)
933     {
934         mapping.push_back(make_at2con(moltype, mtop.ffparams.iparams, flexibleConstraintTreatment));
935     }
936     return mapping;
937 }
938
939 Constraints::Constraints(const gmx_mtop_t&     mtop,
940                          const t_inputrec&     ir,
941                          pull_t*               pull_work,
942                          FILE*                 log,
943                          const t_mdatoms&      md,
944                          const t_commrec*      cr,
945                          const gmx_multisim_t* ms,
946                          t_nrnb*               nrnb,
947                          gmx_wallcycle*        wcycle,
948                          bool                  pbcHandlingRequired,
949                          int                   numConstraints,
950                          int                   numSettles) :
951     impl_(new Impl(mtop, ir, pull_work, log, md, cr, ms, nrnb, wcycle, pbcHandlingRequired, numConstraints, numSettles))
952 {
953 }
954
955 Constraints::Impl::Impl(const gmx_mtop_t&     mtop_p,
956                         const t_inputrec&     ir_p,
957                         pull_t*               pull_work,
958                         FILE*                 log_p,
959                         const t_mdatoms&      md_p,
960                         const t_commrec*      cr_p,
961                         const gmx_multisim_t* ms_p,
962                         t_nrnb*               nrnb_p,
963                         gmx_wallcycle*        wcycle_p,
964                         bool                  pbcHandlingRequired,
965                         int                   numConstraints,
966                         int                   numSettles) :
967     ncon_tot(numConstraints),
968     mtop(mtop_p),
969     md(md_p),
970     pbcHandlingRequired_(pbcHandlingRequired),
971     log(log_p),
972     cr(cr_p),
973     ms(ms_p),
974     pull_work(pull_work),
975     ir(ir_p),
976     nrnb(nrnb_p),
977     wcycle(wcycle_p)
978 {
979     if (numConstraints + numSettles > 0 && ir.epc == epcMTTK)
980     {
981         gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
982     }
983
984     nflexcon = 0;
985     if (numConstraints > 0)
986     {
987         at2con_mt = makeAtomToConstraintMappings(mtop, flexibleConstraintTreatment(EI_DYNAMICS(ir.eI)));
988
989         for (const gmx_molblock_t& molblock : mtop.molblock)
990         {
991             int count = countFlexibleConstraintsTemplate(mtop.moltype[molblock.type].ilist.data(),
992                                                          mtop.ffparams.iparams.data());
993             nflexcon += molblock.nmol * count;
994         }
995
996         if (nflexcon > 0)
997         {
998             if (log)
999             {
1000                 fprintf(log, "There are %d flexible constraints\n", nflexcon);
1001                 if (ir.fc_stepsize == 0)
1002                 {
1003                     fprintf(log,
1004                             "\n"
1005                             "WARNING: step size for flexible constraining = 0\n"
1006                             "         All flexible constraints will be rigid.\n"
1007                             "         Will try to keep all flexible constraints at their original "
1008                             "length,\n"
1009                             "         but the lengths may exhibit some drift.\n\n");
1010                     nflexcon = 0;
1011                 }
1012             }
1013             if (nflexcon > 0)
1014             {
1015                 please_cite(log, "Hess2002");
1016             }
1017         }
1018
1019         if (ir.eConstrAlg == econtLINCS)
1020         {
1021             lincsd = init_lincs(log, mtop, nflexcon, at2con_mt,
1022                                 DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd), ir.nLincsIter,
1023                                 ir.nProjOrder);
1024         }
1025
1026         if (ir.eConstrAlg == econtSHAKE)
1027         {
1028             if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
1029             {
1030                 gmx_fatal(FARGS,
1031                           "SHAKE is not supported with domain decomposition and constraint that "
1032                           "cross domain boundaries, use LINCS");
1033             }
1034             if (nflexcon)
1035             {
1036                 gmx_fatal(FARGS,
1037                           "For this system also velocities and/or forces need to be constrained, "
1038                           "this can not be done with SHAKE, you should select LINCS");
1039             }
1040             please_cite(log, "Ryckaert77a");
1041             if (ir.bShakeSOR)
1042             {
1043                 please_cite(log, "Barth95a");
1044             }
1045
1046             shaked = shake_init();
1047         }
1048     }
1049
1050     if (numSettles > 0)
1051     {
1052         please_cite(log, "Miyamoto92a");
1053
1054         settled = settle_init(mtop);
1055
1056         /* Make an atom to settle index for use in domain decomposition */
1057         for (size_t mt = 0; mt < mtop.moltype.size(); mt++)
1058         {
1059             at2settle_mt.emplace_back(
1060                     make_at2settle(mtop.moltype[mt].atoms.nr, mtop.moltype[mt].ilist[F_SETTLE]));
1061         }
1062
1063         /* Allocate thread-local work arrays */
1064         int nthreads = gmx_omp_nthreads_get(emntSETTLE);
1065         if (nthreads > 1 && vir_r_m_dr_th == nullptr)
1066         {
1067             snew(vir_r_m_dr_th, nthreads);
1068             snew(bSettleErrorHasOccurred, nthreads);
1069         }
1070     }
1071
1072     maxwarn   = 999;
1073     char* env = getenv("GMX_MAXCONSTRWARN");
1074     if (env)
1075     {
1076         maxwarn = 0;
1077         sscanf(env, "%8d", &maxwarn);
1078         if (maxwarn < 0)
1079         {
1080             maxwarn = INT_MAX;
1081         }
1082         if (log)
1083         {
1084             fprintf(log, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1085         }
1086         if (MASTER(cr))
1087         {
1088             fprintf(stderr, "Setting the maximum number of constraint warnings to %d\n", maxwarn);
1089         }
1090     }
1091     warncount_lincs  = 0;
1092     warncount_settle = 0;
1093 }
1094
1095 Constraints::Impl::~Impl()
1096 {
1097     for (auto blocka : at2con_mt)
1098     {
1099         done_blocka(&blocka);
1100     }
1101     if (bSettleErrorHasOccurred != nullptr)
1102     {
1103         sfree(bSettleErrorHasOccurred);
1104     }
1105     if (vir_r_m_dr_th != nullptr)
1106     {
1107         sfree(vir_r_m_dr_th);
1108     }
1109     if (settled != nullptr)
1110     {
1111         settle_free(settled);
1112     }
1113     done_lincs(lincsd);
1114 }
1115
1116 void Constraints::saveEdsamPointer(gmx_edsam* ed)
1117 {
1118     impl_->ed = ed;
1119 }
1120
1121 ArrayRef<const t_blocka> Constraints::atom2constraints_moltype() const
1122 {
1123     return impl_->at2con_mt;
1124 }
1125
1126 ArrayRef<const std::vector<int>> Constraints::atom2settle_moltype() const
1127 {
1128     return impl_->at2settle_mt;
1129 }
1130
1131 void do_constrain_first(FILE*                     fplog,
1132                         gmx::Constraints*         constr,
1133                         const t_inputrec*         ir,
1134                         const t_mdatoms*          md,
1135                         int                       natoms,
1136                         ArrayRefWithPadding<RVec> x,
1137                         ArrayRefWithPadding<RVec> v,
1138                         const matrix              box,
1139                         real                      lambda)
1140 {
1141     int     i, m, start, end;
1142     int64_t step;
1143     real    dt = ir->delta_t;
1144     real    dvdl_dum;
1145     rvec*   savex;
1146
1147     auto xRvec = as_rvec_array(x.paddedArrayRef().data());
1148     auto vRvec = as_rvec_array(v.paddedArrayRef().data());
1149
1150     /* We need to allocate one element extra, since we might use
1151      * (unaligned) 4-wide SIMD loads to access rvec entries.
1152      */
1153     snew(savex, natoms + 1);
1154
1155     start = 0;
1156     end   = md->homenr;
1157
1158     if (debug)
1159     {
1160         fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n", start, md->homenr, end);
1161     }
1162     /* Do a first constrain to reset particles... */
1163     step = ir->init_step;
1164     if (fplog)
1165     {
1166         char buf[STEPSTRSIZE];
1167         fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n", gmx_step_str(step, buf));
1168     }
1169     dvdl_dum = 0;
1170
1171     /* constrain the current position */
1172     constr->apply(TRUE, FALSE, step, 0, 1.0, xRvec, xRvec, nullptr, box, lambda, &dvdl_dum, nullptr,
1173                   nullptr, gmx::ConstraintVariable::Positions);
1174     if (EI_VV(ir->eI))
1175     {
1176         /* constrain the inital velocity, and save it */
1177         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1178         constr->apply(TRUE, FALSE, step, 0, 1.0, xRvec, vRvec, vRvec, box, lambda, &dvdl_dum,
1179                       nullptr, nullptr, gmx::ConstraintVariable::Velocities);
1180     }
1181     /* constrain the inital velocities at t-dt/2 */
1182     if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
1183     {
1184         auto subX = x.paddedArrayRef().subArray(start, end);
1185         auto subV = v.paddedArrayRef().subArray(start, end);
1186         for (i = start; (i < end); i++)
1187         {
1188             for (m = 0; (m < DIM); m++)
1189             {
1190                 /* Reverse the velocity */
1191                 subV[i][m] = -subV[i][m];
1192                 /* Store the position at t-dt in buf */
1193                 savex[i][m] = subX[i][m] + dt * subV[i][m];
1194             }
1195         }
1196         /* Shake the positions at t=-dt with the positions at t=0
1197          * as reference coordinates.
1198          */
1199         if (fplog)
1200         {
1201             char buf[STEPSTRSIZE];
1202             fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n", gmx_step_str(step, buf));
1203         }
1204         dvdl_dum = 0;
1205         constr->apply(TRUE, FALSE, step, -1, 1.0, xRvec, savex, nullptr, box, lambda, &dvdl_dum,
1206                       vRvec, nullptr, gmx::ConstraintVariable::Positions);
1207
1208         for (i = start; i < end; i++)
1209         {
1210             for (m = 0; m < DIM; m++)
1211             {
1212                 /* Re-reverse the velocities */
1213                 subV[i][m] = -subV[i][m];
1214             }
1215         }
1216     }
1217     sfree(savex);
1218 }
1219
1220 } // namespace gmx