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