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