Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_constraints.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2006,2007,2008,2009,2010 by the GROMACS development team.
5  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
6  * Copyright (c) 2017,2018,2019,2020, 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
38 /*! \internal \file
39  *
40  * \brief This file implements functions for domdec to use
41  * while managing inter-atomic constraints.
42  *
43  * \author Berk Hess <hess@kth.se>
44  * \ingroup module_domdec
45  */
46
47 #include "gmxpre.h"
48
49 #include "domdec_constraints.h"
50
51 #include <cassert>
52
53 #include <algorithm>
54 #include <memory>
55
56 #include "gromacs/domdec/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/domdec/hashedmap.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdlib/constr.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/topology/ifunc.h"
68 #include "gromacs/topology/mtop_lookup.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/listoflists.h"
73
74 #include "domdec_internal.h"
75 #include "domdec_specatomcomm.h"
76
77 using gmx::ListOfLists;
78
79 /*! \brief Struct used during constraint setup with domain decomposition */
80 struct gmx_domdec_constraints_t
81 {
82     //! @cond Doxygen_Suppress
83     std::vector<int> molb_con_offset; /**< Offset in the constraint array for each molblock */
84     std::vector<int> molb_ncon_mol; /**< The number of constraints per molecule for each molblock */
85
86     int ncon; /**< The fully local and conneced constraints */
87     /* The global constraint number, only required for clearing gc_req */
88     std::vector<int> con_gl;     /**< Global constraint indices for local constraints */
89     std::vector<int> con_nlocat; /**< Number of local atoms (2/1/0) for each constraint */
90
91     std::vector<bool> gc_req; /**< Boolean that tells if a global constraint index has been requested; note: size global #constraints */
92
93     /* Hash table for keeping track of requests */
94     std::unique_ptr<gmx::HashedMap<int>> ga2la; /**< Global to local communicated constraint atom only index */
95
96     /* Multi-threading stuff */
97     int                          nthread; /**< Number of threads used for DD constraint setup */
98     std::vector<InteractionList> ils;     /**< Constraint ilist working arrays, size \p nthread */
99
100     /* Buffers for requesting atoms */
101     std::vector<std::vector<int>> requestedGlobalAtomIndices; /**< Buffers for requesting global atom indices, one per thread */
102
103     //! @endcond
104 };
105
106 void dd_move_x_constraints(gmx_domdec_t*            dd,
107                            const matrix             box,
108                            gmx::ArrayRef<gmx::RVec> x0,
109                            gmx::ArrayRef<gmx::RVec> x1,
110                            gmx_bool                 bX1IsCoord)
111 {
112     if (dd->constraint_comm)
113     {
114         dd_move_x_specat(
115                 dd, dd->constraint_comm, box, as_rvec_array(x0.data()), as_rvec_array(x1.data()), bX1IsCoord);
116
117         ddReopenBalanceRegionCpu(dd);
118     }
119 }
120
121 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd)
122 {
123     if (dd && dd->constraints)
124     {
125         return dd->constraints->con_nlocat;
126     }
127     else
128     {
129         return {};
130     }
131 }
132
133 void dd_clear_local_constraint_indices(gmx_domdec_t* dd)
134 {
135     gmx_domdec_constraints_t* dc = dd->constraints;
136
137     std::fill(dc->gc_req.begin(), dc->gc_req.end(), false);
138
139     if (dd->constraint_comm)
140     {
141         dc->ga2la->clear();
142     }
143 }
144
145 /*! \brief Walks over the constraints out from the local atoms into the non-local atoms and adds them to a list */
146 static void walk_out(int                       con,
147                      int                       con_offset,
148                      int                       a,
149                      int                       offset,
150                      int                       nrec,
151                      gmx::ArrayRef<const int>  ia1,
152                      gmx::ArrayRef<const int>  ia2,
153                      const ListOfLists<int>&   at2con,
154                      const gmx_ga2la_t&        ga2la,
155                      gmx_bool                  bHomeConnect,
156                      gmx_domdec_constraints_t* dc,
157                      gmx_domdec_specat_comm_t* dcc,
158                      InteractionList*          il_local,
159                      std::vector<int>*         ireq)
160 {
161     if (!dc->gc_req[con_offset + con])
162     {
163         /* Add this non-home constraint to the list */
164         dc->con_gl.push_back(con_offset + con);
165         dc->con_nlocat.push_back(bHomeConnect ? 1 : 0);
166         dc->gc_req[con_offset + con]     = true;
167         const int*         iap           = constr_iatomptr(ia1, ia2, con);
168         const int          parameterType = iap[0];
169         const int          a1_gl         = offset + iap[1];
170         const int          a2_gl         = offset + iap[2];
171         std::array<int, 2> atoms;
172         /* The following indexing code can probably be optizimed */
173         if (const int* a_loc = ga2la.findHome(a1_gl))
174         {
175             atoms[0] = *a_loc;
176         }
177         else
178         {
179             /* We set this index later */
180             atoms[0] = -a1_gl - 1;
181         }
182         if (const int* a_loc = ga2la.findHome(a2_gl))
183         {
184             atoms[1] = *a_loc;
185         }
186         else
187         {
188             /* We set this index later */
189             atoms[1] = -a2_gl - 1;
190         }
191         il_local->push_back(parameterType, atoms);
192         dc->ncon++;
193     }
194     /* Check to not ask for the same atom more than once */
195     if (!dc->ga2la->find(offset + a))
196     {
197         assert(dcc);
198         /* Add this non-home atom to the list */
199         ireq->push_back(offset + a);
200         /* Temporarily mark with -2, we get the index later */
201         dc->ga2la->insert(offset + a, -2);
202     }
203
204     if (nrec > 0)
205     {
206         /* Loop over the constraint connected to atom a */
207         for (const int coni : at2con[a])
208         {
209             if (coni != con)
210             {
211                 /* Walk further */
212                 const int* iap = constr_iatomptr(ia1, ia2, coni);
213                 int        b;
214                 if (a == iap[1])
215                 {
216                     b = iap[2];
217                 }
218                 else
219                 {
220                     b = iap[1];
221                 }
222                 if (!ga2la.findHome(offset + b))
223                 {
224                     walk_out(coni, con_offset, b, offset, nrec - 1, ia1, ia2, at2con, ga2la, FALSE, dc, dcc, il_local, ireq);
225                 }
226             }
227         }
228     }
229 }
230
231 /*! \brief Looks up SETTLE constraints for a range of charge-groups */
232 static void atoms_to_settles(gmx_domdec_t*                         dd,
233                              const gmx_mtop_t*                     mtop,
234                              const int*                            cginfo,
235                              gmx::ArrayRef<const std::vector<int>> at2settle_mt,
236                              int                                   cg_start,
237                              int                                   cg_end,
238                              InteractionList*                      ils_local,
239                              std::vector<int>*                     ireq)
240 {
241     const gmx_ga2la_t& ga2la = *dd->ga2la;
242     int                nral  = NRAL(F_SETTLE);
243
244     int mb = 0;
245     for (int a = cg_start; a < cg_end; a++)
246     {
247         if (GET_CGINFO_SETTLE(cginfo[a]))
248         {
249             int a_gl = dd->globalAtomIndices[a];
250             int a_mol;
251             mtopGetMolblockIndex(mtop, a_gl, &mb, nullptr, &a_mol);
252
253             const gmx_molblock_t* molb   = &mtop->molblock[mb];
254             int                   settle = at2settle_mt[molb->type][a_mol];
255
256             if (settle >= 0)
257             {
258                 int offset = a_gl - a_mol;
259
260                 const int* ia1 = mtop->moltype[molb->type].ilist[F_SETTLE].iatoms.data();
261
262                 int      a_gls[3];
263                 gmx_bool bAssign = FALSE;
264                 int      nlocal  = 0;
265                 for (int sa = 0; sa < nral; sa++)
266                 {
267                     int a_glsa = offset + ia1[settle * (1 + nral) + 1 + sa];
268                     a_gls[sa]  = a_glsa;
269                     if (ga2la.findHome(a_glsa))
270                     {
271                         if (nlocal == 0 && a_gl == a_glsa)
272                         {
273                             bAssign = TRUE;
274                         }
275                         nlocal++;
276                     }
277                 }
278
279                 if (bAssign)
280                 {
281                     const int          parameterType = ia1[settle * 4];
282                     std::array<int, 3> atoms;
283                     for (int sa = 0; sa < nral; sa++)
284                     {
285                         if (const int* a_loc = ga2la.findHome(a_gls[sa]))
286                         {
287                             atoms[sa] = *a_loc;
288                         }
289                         else
290                         {
291                             atoms[sa] = -a_gls[sa] - 1;
292                             /* Add this non-home atom to the list */
293                             ireq->push_back(a_gls[sa]);
294                             /* A check on double atom requests is
295                              * not required for settle.
296                              */
297                         }
298                     }
299                     ils_local->push_back(parameterType, atoms);
300                 }
301             }
302         }
303     }
304 }
305
306 /*! \brief Looks up constraint for the local atoms */
307 static void atoms_to_constraints(gmx_domdec_t*                         dd,
308                                  const gmx_mtop_t*                     mtop,
309                                  const int*                            cginfo,
310                                  gmx::ArrayRef<const ListOfLists<int>> at2con_mt,
311                                  int                                   nrec,
312                                  InteractionList*                      ilc_local,
313                                  std::vector<int>*                     ireq)
314 {
315     gmx_domdec_constraints_t* dc  = dd->constraints;
316     gmx_domdec_specat_comm_t* dcc = dd->constraint_comm;
317
318     const gmx_ga2la_t& ga2la = *dd->ga2la;
319
320     dc->con_gl.clear();
321     dc->con_nlocat.clear();
322
323     int mb    = 0;
324     int nhome = 0;
325     for (int a = 0; a < dd->ncg_home; a++)
326     {
327         if (GET_CGINFO_CONSTR(cginfo[a]))
328         {
329             int a_gl = dd->globalAtomIndices[a];
330             int molnr, a_mol;
331             mtopGetMolblockIndex(mtop, a_gl, &mb, &molnr, &a_mol);
332
333             const gmx_molblock_t& molb = mtop->molblock[mb];
334
335             gmx::ArrayRef<const int> ia1 = mtop->moltype[molb.type].ilist[F_CONSTR].iatoms;
336             gmx::ArrayRef<const int> ia2 = mtop->moltype[molb.type].ilist[F_CONSTRNC].iatoms;
337
338             /* Calculate the global constraint number offset for the molecule.
339              * This is only required for the global index to make sure
340              * that we use each constraint only once.
341              */
342             const int con_offset = dc->molb_con_offset[mb] + molnr * dc->molb_ncon_mol[mb];
343
344             /* The global atom number offset for this molecule */
345             const int offset = a_gl - a_mol;
346             /* Loop over the constraints connected to atom a_mol in the molecule */
347             const auto& at2con = at2con_mt[molb.type];
348             for (const int con : at2con[a_mol])
349             {
350                 const int* iap = constr_iatomptr(ia1, ia2, con);
351                 int        b_mol;
352                 if (a_mol == iap[1])
353                 {
354                     b_mol = iap[2];
355                 }
356                 else
357                 {
358                     b_mol = iap[1];
359                 }
360                 if (const int* a_loc = ga2la.findHome(offset + b_mol))
361                 {
362                     /* Add this fully home constraint at the first atom */
363                     if (a_mol < b_mol)
364                     {
365                         dc->con_gl.push_back(con_offset + con);
366                         dc->con_nlocat.push_back(2);
367                         const int          b_lo          = *a_loc;
368                         const int          parameterType = iap[0];
369                         std::array<int, 2> atoms;
370                         atoms[0] = (a_gl == iap[1] ? a : b_lo);
371                         atoms[1] = (a_gl == iap[1] ? b_lo : a);
372                         ilc_local->push_back(parameterType, atoms);
373                         dc->ncon++;
374                         nhome++;
375                     }
376                 }
377                 else
378                 {
379                     /* We need the nrec constraints coupled to this constraint,
380                      * so we need to walk out of the home cell by nrec+1 atoms,
381                      * since already atom bg is not locally present.
382                      * Therefore we call walk_out with nrec recursions to go
383                      * after this first call.
384                      */
385                     walk_out(con, con_offset, b_mol, offset, nrec, ia1, ia2, at2con, ga2la, TRUE, dc, dcc, ilc_local, ireq);
386                 }
387             }
388         }
389     }
390
391     GMX_ASSERT(dc->con_gl.size() == static_cast<size_t>(dc->ncon),
392                "con_gl size should match the number of constraints");
393     GMX_ASSERT(dc->con_nlocat.size() == static_cast<size_t>(dc->ncon),
394                "con_nlocat size should match the number of constraints");
395
396     if (debug)
397     {
398         fprintf(debug,
399                 "Constraints: home %3d border %3d atoms: %3zu\n",
400                 nhome,
401                 dc->ncon - nhome,
402                 dd->constraint_comm ? ireq->size() : 0);
403     }
404 }
405
406 int dd_make_local_constraints(gmx_domdec_t*                  dd,
407                               int                            at_start,
408                               const struct gmx_mtop_t*       mtop,
409                               const int*                     cginfo,
410                               gmx::Constraints*              constr,
411                               int                            nrec,
412                               gmx::ArrayRef<InteractionList> il_local)
413 {
414     gmx_domdec_constraints_t* dc;
415     InteractionList *         ilc_local, *ils_local;
416     gmx::HashedMap<int>*      ga2la_specat;
417     int                       at_end, i, j;
418
419     // This code should not be called unless this condition is true,
420     // because that's the only time init_domdec_constraints is
421     // called...
422     GMX_RELEASE_ASSERT(dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles,
423                        "dd_make_local_constraints called when there are no local constraints");
424     // ... and init_domdec_constraints always sets
425     // dd->constraint_comm...
426     GMX_RELEASE_ASSERT(
427             dd->constraint_comm,
428             "Invalid use of dd_make_local_constraints before construction of constraint_comm");
429     // ... which static analysis needs to be reassured about, because
430     // otherwise, when dd->splitSettles is
431     // true. dd->constraint_comm is unilaterally dereferenced before
432     // the call to atoms_to_settles.
433
434     dc = dd->constraints;
435
436     ilc_local = &il_local[F_CONSTR];
437     ils_local = &il_local[F_SETTLE];
438
439     dc->ncon = 0;
440     gmx::ArrayRef<const ListOfLists<int>> at2con_mt;
441     std::vector<int>*                     ireq = nullptr;
442     ilc_local->clear();
443     if (dd->constraint_comm)
444     {
445         // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
446         GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
447         at2con_mt = constr->atom2constraints_moltype();
448         ireq      = &dc->requestedGlobalAtomIndices[0];
449         ireq->clear();
450     }
451
452     gmx::ArrayRef<const std::vector<int>> at2settle_mt;
453     /* When settle works inside charge groups, we assigned them already */
454     if (dd->comm->systemInfo.haveSplitSettles)
455     {
456         // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
457         GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
458         at2settle_mt = constr->atom2settle_moltype();
459         ils_local->clear();
460     }
461
462     if (at2settle_mt.empty())
463     {
464         atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
465     }
466     else
467     {
468         int t0_set;
469
470         /* Do the constraints, if present, on the first thread.
471          * Do the settles on all other threads.
472          */
473         t0_set = ((!at2con_mt.empty() && dc->nthread > 1) ? 1 : 0);
474
475 #pragma omp parallel for num_threads(dc->nthread) schedule(static)
476         for (int thread = 0; thread < dc->nthread; thread++)
477         {
478             try
479             {
480                 if (!at2con_mt.empty() && thread == 0)
481                 {
482                     atoms_to_constraints(dd, mtop, cginfo, at2con_mt, nrec, ilc_local, ireq);
483                 }
484
485                 if (thread >= t0_set)
486                 {
487                     int              cg0, cg1;
488                     InteractionList* ilst;
489
490                     /* Distribute the settle check+assignments over
491                      * dc->nthread or dc->nthread-1 threads.
492                      */
493                     cg0 = (dd->ncg_home * (thread - t0_set)) / (dc->nthread - t0_set);
494                     cg1 = (dd->ncg_home * (thread - t0_set + 1)) / (dc->nthread - t0_set);
495
496                     if (thread == t0_set)
497                     {
498                         ilst = ils_local;
499                     }
500                     else
501                     {
502                         ilst = &dc->ils[thread];
503                     }
504                     ilst->clear();
505
506                     std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
507                     if (thread > 0)
508                     {
509                         ireqt.clear();
510                     }
511
512                     atoms_to_settles(dd, mtop, cginfo, at2settle_mt, cg0, cg1, ilst, &ireqt);
513                 }
514             }
515             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
516         }
517
518         /* Combine the generate settles and requested indices */
519         for (int thread = 1; thread < dc->nthread; thread++)
520         {
521             if (thread > t0_set)
522             {
523                 ils_local->append(dc->ils[thread]);
524             }
525
526             const std::vector<int>& ireqt = dc->requestedGlobalAtomIndices[thread];
527             ireq->insert(ireq->end(), ireqt.begin(), ireqt.end());
528         }
529
530         if (debug)
531         {
532             fprintf(debug, "Settles: total %3d\n", ils_local->size() / 4);
533         }
534     }
535
536     if (dd->constraint_comm)
537     {
538         int nral1;
539
540         at_end = setup_specat_communication(dd,
541                                             ireq,
542                                             dd->constraint_comm,
543                                             dd->constraints->ga2la.get(),
544                                             at_start,
545                                             2,
546                                             "constraint",
547                                             " or lincs-order");
548
549         /* Fill in the missing indices */
550         ga2la_specat = dd->constraints->ga2la.get();
551
552         nral1 = 1 + NRAL(F_CONSTR);
553         for (i = 0; i < ilc_local->size(); i += nral1)
554         {
555             int* iap = ilc_local->iatoms.data() + i;
556             for (j = 1; j < nral1; j++)
557             {
558                 if (iap[j] < 0)
559                 {
560                     const int* a = ga2la_specat->find(-iap[j] - 1);
561                     GMX_ASSERT(a, "We have checked before that this atom index has been set");
562                     iap[j] = *a;
563                 }
564             }
565         }
566
567         nral1 = 1 + NRAL(F_SETTLE);
568         for (i = 0; i < ils_local->size(); i += nral1)
569         {
570             int* iap = ils_local->iatoms.data() + i;
571             for (j = 1; j < nral1; j++)
572             {
573                 if (iap[j] < 0)
574                 {
575                     const int* a = ga2la_specat->find(-iap[j] - 1);
576                     GMX_ASSERT(a, "We have checked before that this atom index has been set");
577                     iap[j] = *a;
578                 }
579             }
580         }
581     }
582     else
583     {
584         // Currently unreachable
585         at_end = at_start;
586     }
587
588     return at_end;
589 }
590
591 void init_domdec_constraints(gmx_domdec_t* dd, const gmx_mtop_t* mtop)
592 {
593     gmx_domdec_constraints_t* dc;
594     const gmx_molblock_t*     molb;
595
596     if (debug)
597     {
598         fprintf(debug, "Begin init_domdec_constraints\n");
599     }
600
601     dd->constraints = new gmx_domdec_constraints_t;
602     dc              = dd->constraints;
603
604     dc->molb_con_offset.resize(mtop->molblock.size());
605     dc->molb_ncon_mol.resize(mtop->molblock.size());
606
607     int ncon = 0;
608     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
609     {
610         molb                    = &mtop->molblock[mb];
611         dc->molb_con_offset[mb] = ncon;
612         dc->molb_ncon_mol[mb]   = mtop->moltype[molb->type].ilist[F_CONSTR].size() / 3
613                                 + mtop->moltype[molb->type].ilist[F_CONSTRNC].size() / 3;
614         ncon += molb->nmol * dc->molb_ncon_mol[mb];
615     }
616
617     if (ncon > 0)
618     {
619         dc->gc_req.resize(ncon);
620     }
621
622     /* Use a hash table for the global to local index.
623      * The number of keys is a rough estimate, it will be optimized later.
624      */
625     int numKeysEstimate = std::min(mtop->natoms / 20, mtop->natoms / (2 * dd->nnodes));
626     dc->ga2la           = std::make_unique<gmx::HashedMap<int>>(numKeysEstimate);
627
628     dc->nthread = gmx_omp_nthreads_get(emntDomdec);
629     dc->ils.resize(dc->nthread);
630
631     dd->constraint_comm = new gmx_domdec_specat_comm_t;
632
633     dc->requestedGlobalAtomIndices.resize(dc->nthread);
634 }