04faddad435aa3791c6a9662dc3e28ec5b3ec5e4
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_specatomcomm.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,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_specatomcomm.h"
48
49 #include <cassert>
50
51 #include <algorithm>
52
53 #include "gromacs/domdec/dlb.h"
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_network.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/domdec/ga2la.h"
58 #include "gromacs/domdec/hashedmap.h"
59 #include "gromacs/domdec/partition.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65
66 void dd_move_f_specat(gmx_domdec_t* dd, gmx_domdec_specat_comm_t* spac, rvec* f, rvec* fshift)
67 {
68     gmx_specatsend_t* spas;
69     rvec*             vbuf;
70     int               n, n0, n1, dim, dir;
71     ivec              vis;
72     int               is;
73     gmx_bool          bPBC, bScrew;
74
75     n = spac->at_end;
76     for (int d = dd->ndim - 1; d >= 0; d--)
77     {
78         dim = dd->dim[d];
79         if (dd->nc[dim] > 2)
80         {
81             /* Pulse the grid forward and backward */
82             spas = spac->spas[d];
83             n0   = spas[0].nrecv;
84             n1   = spas[1].nrecv;
85             n -= n1 + n0;
86             vbuf = as_rvec_array(spac->vbuf.data());
87             /* Send and receive the coordinates */
88             dd_sendrecv2_rvec(dd, d, f + n + n1, n0, vbuf, spas[0].a.size(), f + n, n1,
89                               vbuf + spas[0].a.size(), spas[1].a.size());
90             for (dir = 0; dir < 2; dir++)
91             {
92                 bPBC = ((dir == 0 && dd->ci[dim] == 0) || (dir == 1 && dd->ci[dim] == dd->nc[dim] - 1));
93                 bScrew = (bPBC && dd->unitCellInfo.haveScrewPBC && dim == XX);
94
95                 spas = &spac->spas[d][dir];
96                 /* Sum the buffer into the required forces */
97                 if (!bPBC || (!bScrew && fshift == nullptr))
98                 {
99                     for (int a : spas->a)
100                     {
101                         rvec_inc(f[a], *vbuf);
102                         vbuf++;
103                     }
104                 }
105                 else
106                 {
107                     clear_ivec(vis);
108                     vis[dim] = (dir == 0 ? 1 : -1);
109                     is       = IVEC2IS(vis);
110                     if (!bScrew)
111                     {
112                         /* Sum and add to shift forces */
113                         for (int a : spas->a)
114                         {
115                             rvec_inc(f[a], *vbuf);
116                             rvec_inc(fshift[is], *vbuf);
117                             vbuf++;
118                         }
119                     }
120                     else
121                     {
122                         /* Rotate the forces */
123                         for (int a : spas->a)
124                         {
125                             f[a][XX] += (*vbuf)[XX];
126                             f[a][YY] -= (*vbuf)[YY];
127                             f[a][ZZ] -= (*vbuf)[ZZ];
128                             if (fshift)
129                             {
130                                 rvec_inc(fshift[is], *vbuf);
131                             }
132                             vbuf++;
133                         }
134                     }
135                 }
136             }
137         }
138         else
139         {
140             /* Two cells, so we only need to communicate one way */
141             spas = &spac->spas[d][0];
142             n -= spas->nrecv;
143             /* Send and receive the coordinates */
144             ddSendrecv(dd, d, dddirForward, f + n, spas->nrecv, as_rvec_array(spac->vbuf.data()),
145                        spas->a.size());
146             /* Sum the buffer into the required forces */
147             if (dd->unitCellInfo.haveScrewPBC && dim == XX
148                 && (dd->ci[dim] == 0 || dd->ci[dim] == dd->nc[dim] - 1))
149             {
150                 int i = 0;
151                 for (int a : spas->a)
152                 {
153                     /* Rotate the force */
154                     f[a][XX] += spac->vbuf[i][XX];
155                     f[a][YY] -= spac->vbuf[i][YY];
156                     f[a][ZZ] -= spac->vbuf[i][ZZ];
157                     i++;
158                 }
159             }
160             else
161             {
162                 int i = 0;
163                 for (int a : spas->a)
164                 {
165                     rvec_inc(f[a], spac->vbuf[i]);
166                     i++;
167                 }
168             }
169         }
170     }
171 }
172
173 void dd_move_x_specat(gmx_domdec_t* dd, gmx_domdec_specat_comm_t* spac, const matrix box, rvec* x0, rvec* x1, gmx_bool bX1IsCoord)
174 {
175     gmx_specatsend_t* spas;
176     int               nvec, v, n, nn, ns0, ns1, nr0, nr1, nr, d, dim, dir, i;
177     gmx_bool          bPBC, bScrew = FALSE;
178     rvec              shift = { 0, 0, 0 };
179
180     nvec = 1;
181     if (x1 != nullptr)
182     {
183         nvec++;
184     }
185
186     n = spac->at_start;
187     for (d = 0; d < dd->ndim; d++)
188     {
189         dim = dd->dim[d];
190         if (dd->nc[dim] > 2)
191         {
192             /* Pulse the grid forward and backward */
193             rvec* vbuf = as_rvec_array(spac->vbuf.data());
194             for (dir = 0; dir < 2; dir++)
195             {
196                 if (dir == 0 && dd->ci[dim] == 0)
197                 {
198                     bPBC   = TRUE;
199                     bScrew = (dd->unitCellInfo.haveScrewPBC && dim == XX);
200                     copy_rvec(box[dim], shift);
201                 }
202                 else if (dir == 1 && dd->ci[dim] == dd->nc[dim] - 1)
203                 {
204                     bPBC   = TRUE;
205                     bScrew = (dd->unitCellInfo.haveScrewPBC && dim == XX);
206                     for (i = 0; i < DIM; i++)
207                     {
208                         shift[i] = -box[dim][i];
209                     }
210                 }
211                 else
212                 {
213                     bPBC   = FALSE;
214                     bScrew = FALSE;
215                 }
216                 spas = &spac->spas[d][dir];
217                 for (v = 0; v < nvec; v++)
218                 {
219                     rvec* x = (v == 0 ? x0 : x1);
220                     /* Copy the required coordinates to the send buffer */
221                     if (!bPBC || (v == 1 && !bX1IsCoord))
222                     {
223                         /* Only copy */
224                         for (int a : spas->a)
225                         {
226                             copy_rvec(x[a], *vbuf);
227                             vbuf++;
228                         }
229                     }
230                     else if (!bScrew)
231                     {
232                         /* Shift coordinates */
233                         for (int a : spas->a)
234                         {
235                             rvec_add(x[a], shift, *vbuf);
236                             vbuf++;
237                         }
238                     }
239                     else
240                     {
241                         /* Shift and rotate coordinates */
242                         for (int a : spas->a)
243                         {
244                             (*vbuf)[XX] = x[a][XX] + shift[XX];
245                             (*vbuf)[YY] = box[YY][YY] - x[a][YY] + shift[YY];
246                             (*vbuf)[ZZ] = box[ZZ][ZZ] - x[a][ZZ] + shift[ZZ];
247                             vbuf++;
248                         }
249                     }
250                 }
251             }
252             /* Send and receive the coordinates */
253             spas = spac->spas[d];
254             ns0  = spas[0].a.size();
255             nr0  = spas[0].nrecv;
256             ns1  = spas[1].a.size();
257             nr1  = spas[1].nrecv;
258             if (nvec == 1)
259             {
260                 rvec* vbuf = as_rvec_array(spac->vbuf.data());
261                 dd_sendrecv2_rvec(dd, d, vbuf + ns0, ns1, x0 + n, nr1, vbuf, ns0, x0 + n + nr1, nr0);
262             }
263             else
264             {
265                 rvec* vbuf = as_rvec_array(spac->vbuf.data());
266                 /* Communicate both vectors in one buffer */
267                 rvec* rbuf = as_rvec_array(spac->vbuf2.data());
268                 dd_sendrecv2_rvec(dd, d, vbuf + 2 * ns0, 2 * ns1, rbuf, 2 * nr1, vbuf, 2 * ns0,
269                                   rbuf + 2 * nr1, 2 * nr0);
270                 /* Split the buffer into the two vectors */
271                 nn = n;
272                 for (dir = 1; dir >= 0; dir--)
273                 {
274                     nr = spas[dir].nrecv;
275                     for (v = 0; v < 2; v++)
276                     {
277                         rvec* x = (v == 0 ? x0 : x1);
278                         for (i = 0; i < nr; i++)
279                         {
280                             copy_rvec(*rbuf, x[nn + i]);
281                             rbuf++;
282                         }
283                     }
284                     nn += nr;
285                 }
286             }
287             n += nr0 + nr1;
288         }
289         else
290         {
291             spas = &spac->spas[d][0];
292             /* Copy the required coordinates to the send buffer */
293             rvec* vbuf = as_rvec_array(spac->vbuf.data());
294             for (v = 0; v < nvec; v++)
295             {
296                 rvec* x = (v == 0 ? x0 : x1);
297                 if (dd->unitCellInfo.haveScrewPBC && dim == XX
298                     && (dd->ci[XX] == 0 || dd->ci[XX] == dd->nc[XX] - 1))
299                 {
300                     /* Here we only perform the rotation, the rest of the pbc
301                      * is handled in the constraint or viste routines.
302                      */
303                     for (int a : spas->a)
304                     {
305                         (*vbuf)[XX] = x[a][XX];
306                         (*vbuf)[YY] = box[YY][YY] - x[a][YY];
307                         (*vbuf)[ZZ] = box[ZZ][ZZ] - x[a][ZZ];
308                         vbuf++;
309                     }
310                 }
311                 else
312                 {
313                     for (int a : spas->a)
314                     {
315                         copy_rvec(x[a], *vbuf);
316                         vbuf++;
317                     }
318                 }
319             }
320             /* Send and receive the coordinates */
321             if (nvec == 1)
322             {
323                 rvec* vbuf = as_rvec_array(spac->vbuf.data());
324                 ddSendrecv(dd, d, dddirBackward, vbuf, spas->a.size(), x0 + n, spas->nrecv);
325             }
326             else
327             {
328                 rvec* vbuf = as_rvec_array(spac->vbuf.data());
329                 /* Communicate both vectors in one buffer */
330                 rvec* rbuf = as_rvec_array(spac->vbuf2.data());
331                 ddSendrecv(dd, d, dddirBackward, vbuf, 2 * spas->a.size(), rbuf, 2 * spas->nrecv);
332                 /* Split the buffer into the two vectors */
333                 nr = spas[0].nrecv;
334                 for (v = 0; v < 2; v++)
335                 {
336                     rvec* x = (v == 0 ? x0 : x1);
337                     for (i = 0; i < nr; i++)
338                     {
339                         copy_rvec(*rbuf, x[n + i]);
340                         rbuf++;
341                     }
342                 }
343             }
344             n += spas->nrecv;
345         }
346     }
347 }
348
349 int setup_specat_communication(gmx_domdec_t*             dd,
350                                std::vector<int>*         ireq,
351                                gmx_domdec_specat_comm_t* spac,
352                                gmx::HashedMap<int>*      ga2la_specat,
353                                int                       at_start,
354                                int                       vbuf_fac,
355                                const char*               specat_type,
356                                const char*               add_err)
357 {
358     int               nsend[2], nlast, nsend_zero[2] = { 0, 0 }, *nsend_ptr;
359     int               dim, ndir, nr, ns, nrecv_local, n0, start, buf[2];
360     int               nat_tot_specat, nat_tot_prev;
361     gmx_bool          bPBC;
362     gmx_specatsend_t* spas;
363
364     if (debug)
365     {
366         fprintf(debug, "Begin setup_specat_communication for %s\n", specat_type);
367     }
368
369     /* nsend[0]: the number of atoms requested by this node only,
370      *           we communicate this for more efficients checks
371      * nsend[1]: the total number of requested atoms
372      */
373     const int numRequested = ireq->size();
374     nsend[0]               = ireq->size();
375     nsend[1]               = nsend[0];
376     nlast                  = nsend[1];
377     for (int d = dd->ndim - 1; d >= 0; d--)
378     {
379         /* Pulse the grid forward and backward */
380         dim  = dd->dim[d];
381         bPBC = (dim < dd->unitCellInfo.npbcdim);
382         if (dd->nc[dim] == 2)
383         {
384             /* Only 2 cells, so we only need to communicate once */
385             ndir = 1;
386         }
387         else
388         {
389             ndir = 2;
390         }
391         for (int dir = 0; dir < ndir; dir++)
392         {
393             if (!bPBC && dd->nc[dim] > 2
394                 && ((dir == 0 && dd->ci[dim] == dd->nc[dim] - 1) || (dir == 1 && dd->ci[dim] == 0)))
395             {
396                 /* No pbc: the fist/last cell should not request atoms */
397                 nsend_ptr = nsend_zero;
398             }
399             else
400             {
401                 nsend_ptr = nsend;
402             }
403             /* Communicate the number of indices */
404             ddSendrecv(dd, d, dir == 0 ? dddirForward : dddirBackward, nsend_ptr, 2,
405                        spac->nreq[d][dir], 2);
406             nr = spac->nreq[d][dir][1];
407             ireq->resize(nlast + nr);
408             /* Communicate the indices */
409             ddSendrecv(dd, d, dir == 0 ? dddirForward : dddirBackward, ireq->data(), nsend_ptr[1],
410                        ireq->data() + nlast, nr);
411             nlast += nr;
412         }
413         nsend[1] = nlast;
414     }
415     if (debug)
416     {
417         fprintf(debug, "Communicated the counts\n");
418     }
419
420     /* Search for the requested atoms and communicate the indices we have */
421     nat_tot_specat = at_start;
422     nrecv_local    = 0;
423     for (int d = 0; d < dd->ndim; d++)
424     {
425         /* Pulse the grid forward and backward */
426         if (dd->dim[d] >= dd->unitCellInfo.npbcdim || dd->nc[dd->dim[d]] > 2)
427         {
428             ndir = 2;
429         }
430         else
431         {
432             ndir = 1;
433         }
434         nat_tot_prev = nat_tot_specat;
435         for (int dir = ndir - 1; dir >= 0; dir--)
436         {
437             /* To avoid cost of clearing by resize(), we only increase size */
438             if (static_cast<size_t>(nat_tot_specat) > spac->sendAtom.size())
439             {
440                 /* Note: resize initializes new elements to false, which is actually needed here */
441                 spac->sendAtom.resize(nat_tot_specat);
442             }
443             spas = &spac->spas[d][dir];
444             n0   = spac->nreq[d][dir][0];
445             nr   = spac->nreq[d][dir][1];
446             if (debug)
447             {
448                 fprintf(debug, "dim=%d, dir=%d, searching for %d atoms\n", d, dir, nr);
449             }
450             start = nlast - nr;
451             spas->a.clear();
452             spac->ibuf.clear();
453             nsend[0] = 0;
454             for (int i = 0; i < nr; i++)
455             {
456                 const int indr = (*ireq)[start + i];
457                 int       ind;
458                 /* Check if this is a home atom and if so ind will be set */
459                 if (const int* homeIndex = dd->ga2la->findHome(indr))
460                 {
461                     ind = *homeIndex;
462                 }
463                 else
464                 {
465                     /* Search in the communicated atoms */
466                     if (const int* a = ga2la_specat->find(indr))
467                     {
468                         ind = *a;
469                     }
470                     else
471                     {
472                         ind = -1;
473                     }
474                 }
475                 if (ind >= 0)
476                 {
477                     if (i < n0 || !spac->sendAtom[ind])
478                     {
479                         /* Store the local index so we know which coordinates
480                          * to send out later.
481                          */
482                         spas->a.push_back(ind);
483                         spac->sendAtom[ind] = true;
484                         /* Store the global index so we can send it now */
485                         spac->ibuf.push_back(indr);
486                         if (i < n0)
487                         {
488                             nsend[0]++;
489                         }
490                     }
491                 }
492             }
493             nlast = start;
494             /* Clear the local flags */
495             for (int a : spas->a)
496             {
497                 spac->sendAtom[a] = false;
498             }
499             /* Send and receive the number of indices to communicate */
500             nsend[1] = spas->a.size();
501             ddSendrecv(dd, d, dir == 0 ? dddirBackward : dddirForward, nsend, 2, buf, 2);
502             if (debug)
503             {
504                 fprintf(debug,
505                         "Send to rank %d, %d (%d) indices, "
506                         "receive from rank %d, %d (%d) indices\n",
507                         dd->neighbor[d][1 - dir], nsend[1], nsend[0], dd->neighbor[d][dir], buf[1],
508                         buf[0]);
509                 if (gmx_debug_at)
510                 {
511                     for (int i : spac->ibuf)
512                     {
513                         fprintf(debug, " %d", i + 1);
514                     }
515                     fprintf(debug, "\n");
516                 }
517             }
518             nrecv_local += buf[0];
519             spas->nrecv = buf[1];
520             dd->globalAtomIndices.resize(nat_tot_specat + spas->nrecv);
521             /* Send and receive the indices */
522             ddSendrecv(dd, d, dir == 0 ? dddirBackward : dddirForward, spac->ibuf.data(),
523                        spac->ibuf.size(), dd->globalAtomIndices.data() + nat_tot_specat, spas->nrecv);
524             nat_tot_specat += spas->nrecv;
525         }
526
527         /* Increase the x/f communication buffer sizes, when necessary */
528         ns = spac->spas[d][0].a.size();
529         nr = spac->spas[d][0].nrecv;
530         if (ndir == 2)
531         {
532             ns += spac->spas[d][1].a.size();
533             nr += spac->spas[d][1].nrecv;
534         }
535         if (vbuf_fac * ns > gmx::index(spac->vbuf.size()))
536         {
537             spac->vbuf.resize(vbuf_fac * ns);
538         }
539         if (vbuf_fac == 2 && vbuf_fac * nr > gmx::index(spac->vbuf2.size()))
540         {
541             spac->vbuf2.resize(vbuf_fac * nr);
542         }
543
544         /* Make a global to local index for the communication atoms */
545         for (int i = nat_tot_prev; i < nat_tot_specat; i++)
546         {
547             ga2la_specat->insert_or_assign(dd->globalAtomIndices[i], i);
548         }
549     }
550
551     /* Check that in the end we got the number of atoms we asked for */
552     if (nrecv_local != numRequested)
553     {
554         if (debug)
555         {
556             fprintf(debug, "Requested %d, received %d (tot recv %d)\n", numRequested, nrecv_local,
557                     nat_tot_specat - at_start);
558             if (gmx_debug_at)
559             {
560                 for (int i = 0; i < numRequested; i++)
561                 {
562                     const int* ind = ga2la_specat->find((*ireq)[i]);
563                     fprintf(debug, " %s%d", ind ? "" : "!", (*ireq)[i] + 1);
564                 }
565                 fprintf(debug, "\n");
566             }
567         }
568         fprintf(stderr, "\nDD cell %d %d %d: Neighboring cells do not have atoms:", dd->ci[XX],
569                 dd->ci[YY], dd->ci[ZZ]);
570         for (int i = 0; i < numRequested; i++)
571         {
572             if (!ga2la_specat->find((*ireq)[i]))
573             {
574                 fprintf(stderr, " %d", (*ireq)[i] + 1);
575             }
576         }
577         fprintf(stderr, "\n");
578         gmx_fatal(FARGS,
579                   "DD cell %d %d %d could only obtain %d of the %d atoms that are connected via "
580                   "%ss from the neighboring cells. This probably means your %s lengths are too "
581                   "long compared to the domain decomposition cell size. Decrease the number of "
582                   "domain decomposition grid cells%s%s.",
583                   dd->ci[XX], dd->ci[YY], dd->ci[ZZ], nrecv_local, numRequested, specat_type,
584                   specat_type, add_err, dd_dlb_is_on(dd) ? " or use the -rcon option of mdrun" : "");
585     }
586
587     spac->at_start = at_start;
588     spac->at_end   = nat_tot_specat;
589
590     if (debug)
591     {
592         fprintf(debug, "Done setup_specat_communication\n");
593     }
594
595     return nat_tot_specat;
596 }