2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2020,2021, 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.
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.
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.
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.
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.
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.
36 * \brief Tests for the halo exchange
38 * The test sets up the rank topology and performs a coordinate halo
39 * exchange (for both CPU and GPU codepaths) for several 1D and 2D
40 * pulse configirations. Each pulse involves a few non-contiguous
41 * indices. The sending rank, atom number and spatial 3D index are
42 * encoded in the x values, to allow correctness checking following
47 * \author Alan Gray <alang@nvidia.com>
48 * \ingroup module_domdec
59 #include <gtest/gtest.h>
61 #include "gromacs/domdec/atomdistribution.h"
62 #include "gromacs/domdec/domdec_internal.h"
63 #include "gromacs/domdec/gpuhaloexchange.h"
65 # include "gromacs/gpu_utils/device_stream.h"
66 # include "gromacs/gpu_utils/devicebuffer.h"
68 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
69 #include "gromacs/gpu_utils/hostallocator.h"
70 #include "gromacs/mdtypes/inputrec.h"
72 #include "testutils/mpitest.h"
73 #include "testutils/test_hardware_environment.h"
82 /*! \brief Get encoded numerical value for sending rank, atom number and spatial 3D index
84 * \param [in] sendRank MPI rank of sender
85 * \param [in] atomNumber Atom number
86 * \param [in] spatial3dIndex Spatial 3D Index
88 * \returns Encoded value
90 float encodedValue(const int sendRank, const int atomNumber, const int spatial3dIndex)
92 return sendRank * 1000 + atomNumber * 100 + spatial3dIndex;
95 /*! \brief Initialize halo array
97 * \param [in] x Atom coordinate data array
98 * \param [in] numHomeAtoms Number of home atoms
99 * \param [in] numAtomsTotal Total number of atoms, including halo
101 void initHaloData(RVec* x, const int numHomeAtoms, const int numAtomsTotal)
104 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
106 for (int i = 0; i < numAtomsTotal; i++)
108 for (int j = 0; j < DIM; j++)
110 x[i][j] = i < numHomeAtoms ? encodedValue(rank, i, j) : -1;
115 /*! \brief Perform GPU halo exchange, including required setup and data transfers
117 * \param [in] dd Domain decomposition object
118 * \param [in] box Box matrix
119 * \param [in] h_x Atom coordinate data array on host
120 * \param [in] numAtomsTotal Total number of atoms, including halo
122 void gpuHalo(gmx_domdec_t* dd, matrix box, HostVector<RVec>* h_x, int numAtomsTotal)
124 #if (GMX_GPU_CUDA && GMX_THREAD_MPI)
125 // pin memory if possible
126 changePinningPolicy(h_x, PinningPolicy::PinnedIfSupported);
127 // Set up GPU hardware environment and assign this MPI rank to a device
129 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
130 int numDevices = getTestHardwareEnvironment()->getTestDeviceList().size();
131 const auto& testDevice = getTestHardwareEnvironment()->getTestDeviceList()[rank % numDevices];
132 const auto& deviceContext = testDevice->deviceContext();
133 setActiveDevice(testDevice->deviceInfo());
134 DeviceStream deviceStream(deviceContext, DeviceStreamPriority::Normal, false);
136 // Set up GPU buffer and copy input data from host
137 DeviceBuffer<RVec> d_x;
139 int d_x_size_alloc = -1;
140 reallocateDeviceBuffer(&d_x, numAtomsTotal, &d_x_size, &d_x_size_alloc, deviceContext);
142 copyToDeviceBuffer(&d_x, h_x->data(), 0, numAtomsTotal, deviceStream, GpuApiCallBehavior::Sync, nullptr);
144 const int numPulses = std::accumulate(
145 dd->comm->cd.begin(), dd->comm->cd.end(), 0, [](const int a, const auto& b) {
146 return a + b.numPulses();
148 const int numExtraConsumptions = GMX_THREAD_MPI ? 1 : 0;
149 // Will be consumed once for each pulse, and, with tMPI, once more for dim=0,pulse=0 case
150 GpuEventSynchronizer coordinatesReadyOnDeviceEvent(numPulses + numExtraConsumptions,
151 numPulses + numExtraConsumptions);
152 coordinatesReadyOnDeviceEvent.markEvent(deviceStream);
154 std::array<std::vector<GpuHaloExchange>, DIM> gpuHaloExchange;
156 // Create halo exchange objects
157 for (int d = 0; d < dd->ndim; d++)
159 for (int pulse = 0; pulse < dd->comm->cd[d].numPulses(); pulse++)
161 gpuHaloExchange[d].push_back(
162 GpuHaloExchange(dd, d, MPI_COMM_WORLD, deviceContext, pulse, nullptr));
166 // Perform GPU halo exchange
167 for (int d = 0; d < dd->ndim; d++)
169 for (int pulse = 0; pulse < dd->comm->cd[d].numPulses(); pulse++)
171 gpuHaloExchange[d][pulse].reinitHalo(d_x, nullptr);
172 gpuHaloExchange[d][pulse].communicateHaloCoordinates(box, &coordinatesReadyOnDeviceEvent);
175 // Barrier is needed to avoid other threads using events after its owner has exited and destroyed the context.
176 MPI_Barrier(MPI_COMM_WORLD);
178 deviceStream.synchronize();
180 // Copy results back to host
181 copyFromDeviceBuffer(
182 h_x->data(), &d_x, 0, numAtomsTotal, deviceStream, GpuApiCallBehavior::Sync, nullptr);
184 freeDeviceBuffer(d_x);
186 GMX_UNUSED_VALUE(dd);
187 GMX_UNUSED_VALUE(box);
188 GMX_UNUSED_VALUE(h_x);
189 GMX_UNUSED_VALUE(numAtomsTotal);
193 /*! \brief Define 1D rank topology with 4 MPI tasks
195 * \param [in] dd Domain decomposition object
197 void define1dRankTopology(gmx_domdec_t* dd)
200 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
202 const int numRanks = getNumberOfTestMpiRanks();
203 dd->neighbor[0][0] = (rank + 1) % numRanks;
204 dd->neighbor[0][1] = (rank == 0) ? (numRanks - 1) : rank - 1;
207 /*! \brief Define 2D rank topology with 4 MPI tasks
214 * \param [in] dd Domain decomposition object
216 void define2dRankTopology(gmx_domdec_t* dd)
220 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
225 dd->neighbor[0][0] = 1;
226 dd->neighbor[0][1] = 1;
227 dd->neighbor[1][0] = 2;
228 dd->neighbor[1][1] = 2;
231 dd->neighbor[0][0] = 0;
232 dd->neighbor[0][1] = 0;
233 dd->neighbor[1][0] = 3;
234 dd->neighbor[1][1] = 3;
237 dd->neighbor[0][0] = 3;
238 dd->neighbor[0][1] = 3;
239 dd->neighbor[1][0] = 0;
240 dd->neighbor[1][1] = 0;
243 dd->neighbor[0][0] = 2;
244 dd->neighbor[0][1] = 2;
245 dd->neighbor[1][0] = 1;
246 dd->neighbor[1][1] = 1;
251 /*! \brief Define a 1D halo with 1 pulses
253 * \param [in] dd Domain decomposition object
254 * \param [in] indvec Vector of index vectors
256 void define1dHaloWith1Pulse(gmx_domdec_t* dd, std::vector<gmx_domdec_ind_t>* indvec)
260 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
262 std::vector<int> indexvec;
263 gmx_domdec_ind_t ind;
269 // Set up indices involved in halo
273 dd->comm->cd[dimIndex].receiveInPlace = true;
274 dd->dim[dimIndex] = 0;
275 dd->ci[dimIndex] = rank;
277 // First pulse involves (arbitrary) indices 1 and 3
278 indexvec.push_back(1);
279 indexvec.push_back(3);
281 ind.index = indexvec;
282 ind.nsend[nzone + 1] = 2;
283 ind.nrecv[nzone + 1] = 2;
284 indvec->push_back(ind);
286 dd->comm->cd[dimIndex].ind = *indvec;
289 /*! \brief Define a 1D halo with 2 pulses
291 * \param [in] dd Domain decomposition object
292 * \param [in] indvec Vector of index vectors
294 void define1dHaloWith2Pulses(gmx_domdec_t* dd, std::vector<gmx_domdec_ind_t>* indvec)
298 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
300 std::vector<int> indexvec;
301 gmx_domdec_ind_t ind;
307 // Set up indices involved in halo
311 dd->comm->cd[dimIndex].receiveInPlace = true;
312 dd->dim[dimIndex] = 0;
313 dd->ci[dimIndex] = rank;
315 // First pulse involves (arbitrary) indices 1 and 3
316 indexvec.push_back(1);
317 indexvec.push_back(3);
319 ind.index = indexvec;
320 ind.nsend[nzone + 1] = 2;
321 ind.nrecv[nzone + 1] = 2;
322 indvec->push_back(ind);
324 // Add another pulse with (arbitrary) indices 4,5,7
327 indexvec.push_back(4);
328 indexvec.push_back(5);
329 indexvec.push_back(7);
331 ind.index = indexvec;
332 ind.nsend[nzone + 1] = 3;
333 ind.nrecv[nzone + 1] = 3;
334 indvec->push_back(ind);
336 dd->comm->cd[dimIndex].ind = *indvec;
339 /*! \brief Define a 2D halo with 1 pulse in each dimension
341 * \param [in] dd Domain decomposition object
342 * \param [in] indvec Vector of index vectors
344 void define2dHaloWith1PulseInEachDim(gmx_domdec_t* dd, std::vector<gmx_domdec_ind_t>* indvec)
348 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
350 std::vector<int> indexvec;
351 gmx_domdec_ind_t ind;
355 for (int dimIndex = 0; dimIndex < dd->ndim; dimIndex++)
358 // Set up indices involved in halo
362 dd->comm->cd[dimIndex].receiveInPlace = true;
363 dd->dim[dimIndex] = 0;
364 dd->ci[dimIndex] = rank;
366 // Single pulse involving (arbitrary) indices 1 and 3
367 indexvec.push_back(1);
368 indexvec.push_back(3);
370 ind.index = indexvec;
371 ind.nsend[nzone + 1] = 2;
372 ind.nrecv[nzone + 1] = 2;
373 indvec->push_back(ind);
375 dd->comm->cd[dimIndex].ind = *indvec;
381 /*! \brief Define a 2D halo with 2 pulses in the first dimension
383 * \param [in] dd Domain decomposition object
384 * \param [in] indvec Vector of index vectors
386 void define2dHaloWith2PulsesInDim1(gmx_domdec_t* dd, std::vector<gmx_domdec_ind_t>* indvec)
390 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
392 std::vector<int> indexvec;
393 gmx_domdec_ind_t ind;
397 for (int dimIndex = 0; dimIndex < dd->ndim; dimIndex++)
400 // Set up indices involved in halo
404 dd->comm->cd[dimIndex].receiveInPlace = true;
405 dd->dim[dimIndex] = 0;
406 dd->ci[dimIndex] = rank;
408 // First pulse involves (arbitrary) indices 1 and 3
409 indexvec.push_back(1);
410 indexvec.push_back(3);
412 ind.index = indexvec;
413 ind.nsend[nzone + 1] = 2;
414 ind.nrecv[nzone + 1] = 2;
415 indvec->push_back(ind);
417 if (dimIndex == 0) // Add another pulse with (arbitrary) indices 4,5,7
421 indexvec.push_back(4);
422 indexvec.push_back(5);
423 indexvec.push_back(7);
425 ind.index = indexvec;
426 ind.nsend[nzone + 1] = 3;
427 ind.nrecv[nzone + 1] = 3;
428 indvec->push_back(ind);
431 dd->comm->cd[dimIndex].ind = *indvec;
437 /*! \brief Check results for above-defined 1D halo with 1 pulse
439 * \param [in] x Atom coordinate data array
440 * \param [in] dd Domain decomposition object
441 * \param [in] numHomeAtoms Number of home atoms
443 void checkResults1dHaloWith1Pulse(const RVec* x, const gmx_domdec_t* dd, const int numHomeAtoms)
445 // Check results are expected from values encoded in x data
446 for (int j = 0; j < DIM; j++)
448 // First Pulse in first dim: atoms 1 and 3 from forward horizontal neighbour
449 EXPECT_EQ(x[numHomeAtoms][j], encodedValue(dd->neighbor[0][0], 1, j));
450 EXPECT_EQ(x[numHomeAtoms + 1][j], encodedValue(dd->neighbor[0][0], 3, j));
454 /*! \brief Check results for above-defined 1D halo with 2 pulses
456 * \param [in] x Atom coordinate data array
457 * \param [in] dd Domain decomposition object
458 * \param [in] numHomeAtoms Number of home atoms
460 void checkResults1dHaloWith2Pulses(const RVec* x, const gmx_domdec_t* dd, const int numHomeAtoms)
462 // Check results are expected from values encoded in x data
463 for (int j = 0; j < DIM; j++)
465 // First Pulse in first dim: atoms 1 and 3 from forward horizontal neighbour
466 EXPECT_EQ(x[numHomeAtoms][j], encodedValue(dd->neighbor[0][0], 1, j));
467 EXPECT_EQ(x[numHomeAtoms + 1][j], encodedValue(dd->neighbor[0][0], 3, j));
468 // Second Pulse in first dim: atoms 4,5,7 from forward horizontal neighbour
469 EXPECT_EQ(x[numHomeAtoms + 2][j], encodedValue(dd->neighbor[0][0], 4, j));
470 EXPECT_EQ(x[numHomeAtoms + 3][j], encodedValue(dd->neighbor[0][0], 5, j));
471 EXPECT_EQ(x[numHomeAtoms + 4][j], encodedValue(dd->neighbor[0][0], 7, j));
475 /*! \brief Check results for above-defined 2D halo with 1 pulse in each dimension
477 * \param [in] x Atom coordinate data array
478 * \param [in] dd Domain decomposition object
479 * \param [in] numHomeAtoms Number of home atoms
481 void checkResults2dHaloWith1PulseInEachDim(const RVec* x, const gmx_domdec_t* dd, const int numHomeAtoms)
483 // Check results are expected from values encoded in x data
484 for (int j = 0; j < DIM; j++)
486 // First Pulse in first dim: atoms 1 and 3 from forward horizontal neighbour
487 EXPECT_EQ(x[numHomeAtoms][j], encodedValue(dd->neighbor[0][0], 1, j));
488 EXPECT_EQ(x[numHomeAtoms + 1][j], encodedValue(dd->neighbor[0][0], 3, j));
489 // First Pulse in second dim: atoms 1 and 3 from forward vertical neighbour
490 EXPECT_EQ(x[numHomeAtoms + 2][j], encodedValue(dd->neighbor[1][0], 1, j));
491 EXPECT_EQ(x[numHomeAtoms + 3][j], encodedValue(dd->neighbor[1][0], 3, j));
495 /*! \brief Check results for above-defined 2D halo with 2 pulses in the first dimension
497 * \param [in] x Atom coordinate data array
498 * \param [in] dd Domain decomposition object
499 * \param [in] numHomeAtoms Number of home atoms
501 void checkResults2dHaloWith2PulsesInDim1(const RVec* x, const gmx_domdec_t* dd, const int numHomeAtoms)
503 // Check results are expected from values encoded in x data
504 for (int j = 0; j < DIM; j++)
506 // First Pulse in first dim: atoms 1 and 3 from forward horizontal neighbour
507 EXPECT_EQ(x[numHomeAtoms][j], encodedValue(dd->neighbor[0][0], 1, j));
508 EXPECT_EQ(x[numHomeAtoms + 1][j], encodedValue(dd->neighbor[0][0], 3, j));
509 // Second Pulse in first dim: atoms 4,5,7 from forward horizontal neighbour
510 EXPECT_EQ(x[numHomeAtoms + 2][j], encodedValue(dd->neighbor[0][0], 4, j));
511 EXPECT_EQ(x[numHomeAtoms + 3][j], encodedValue(dd->neighbor[0][0], 5, j));
512 EXPECT_EQ(x[numHomeAtoms + 4][j], encodedValue(dd->neighbor[0][0], 7, j));
513 // First Pulse in second dim: atoms 1 and 3 from forward vertical neighbour
514 EXPECT_EQ(x[numHomeAtoms + 5][j], encodedValue(dd->neighbor[1][0], 1, j));
515 EXPECT_EQ(x[numHomeAtoms + 6][j], encodedValue(dd->neighbor[1][0], 3, j));
519 TEST(HaloExchangeTest, Coordinates1dHaloWith1Pulse)
521 GMX_MPI_TEST(RequireRankCount<4>);
524 const int numHomeAtoms = 10;
525 const int numHaloAtoms = 2;
526 const int numAtomsTotal = numHomeAtoms + numHaloAtoms;
527 HostVector<RVec> h_x;
528 h_x.resize(numAtomsTotal);
530 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
535 dd.mpi_comm_all = MPI_COMM_WORLD;
536 gmx_domdec_comm_t comm;
538 dd.unitCellInfo.haveScrewPBC = false;
540 DDAtomRanges atomRanges;
541 atomRanges.setEnd(DDAtomRanges::Type::Home, numHomeAtoms);
542 dd.comm->atomRanges = atomRanges;
544 define1dRankTopology(&dd);
546 std::vector<gmx_domdec_ind_t> indvec;
547 define1dHaloWith1Pulse(&dd, &indvec);
549 // Perform halo exchange
550 matrix box = { { 0., 0., 0. } };
551 dd_move_x(&dd, box, static_cast<ArrayRef<RVec>>(h_x), nullptr);
554 checkResults1dHaloWith1Pulse(h_x.data(), &dd, numHomeAtoms);
556 if (GMX_GPU_CUDA && GMX_THREAD_MPI) // repeat with GPU halo codepath
558 // early return if no devices are available.
559 if (getTestHardwareEnvironment()->getTestDeviceList().empty())
564 // Re-initialize input
565 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
567 // Perform GPU halo exchange
568 gpuHalo(&dd, box, &h_x, numAtomsTotal);
571 checkResults1dHaloWith1Pulse(h_x.data(), &dd, numHomeAtoms);
575 TEST(HaloExchangeTest, Coordinates1dHaloWith2Pulses)
577 GMX_MPI_TEST(RequireRankCount<4>);
580 const int numHomeAtoms = 10;
581 const int numHaloAtoms = 5;
582 const int numAtomsTotal = numHomeAtoms + numHaloAtoms;
583 HostVector<RVec> h_x;
584 h_x.resize(numAtomsTotal);
586 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
591 dd.mpi_comm_all = MPI_COMM_WORLD;
592 gmx_domdec_comm_t comm;
594 dd.unitCellInfo.haveScrewPBC = false;
596 DDAtomRanges atomRanges;
597 atomRanges.setEnd(DDAtomRanges::Type::Home, numHomeAtoms);
598 dd.comm->atomRanges = atomRanges;
600 define1dRankTopology(&dd);
602 std::vector<gmx_domdec_ind_t> indvec;
603 define1dHaloWith2Pulses(&dd, &indvec);
605 // Perform halo exchange
606 matrix box = { { 0., 0., 0. } };
607 dd_move_x(&dd, box, static_cast<ArrayRef<RVec>>(h_x), nullptr);
610 checkResults1dHaloWith2Pulses(h_x.data(), &dd, numHomeAtoms);
612 if (GMX_GPU_CUDA && GMX_THREAD_MPI) // repeat with GPU halo codepath
614 // early return if no devices are available.
615 if (getTestHardwareEnvironment()->getTestDeviceList().empty())
620 // Re-initialize input
621 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
623 // Perform GPU halo exchange
624 gpuHalo(&dd, box, &h_x, numAtomsTotal);
627 checkResults1dHaloWith2Pulses(h_x.data(), &dd, numHomeAtoms);
632 TEST(HaloExchangeTest, Coordinates2dHaloWith1PulseInEachDim)
634 GMX_MPI_TEST(RequireRankCount<4>);
637 const int numHomeAtoms = 10;
638 const int numHaloAtoms = 4;
639 const int numAtomsTotal = numHomeAtoms + numHaloAtoms;
640 HostVector<RVec> h_x;
641 h_x.resize(numAtomsTotal);
643 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
648 dd.mpi_comm_all = MPI_COMM_WORLD;
649 gmx_domdec_comm_t comm;
651 dd.unitCellInfo.haveScrewPBC = false;
653 DDAtomRanges atomRanges;
654 atomRanges.setEnd(DDAtomRanges::Type::Home, numHomeAtoms);
655 dd.comm->atomRanges = atomRanges;
657 define2dRankTopology(&dd);
659 std::vector<gmx_domdec_ind_t> indvec;
660 define2dHaloWith1PulseInEachDim(&dd, &indvec);
662 // Perform halo exchange
663 matrix box = { { 0., 0., 0. } };
664 dd_move_x(&dd, box, static_cast<ArrayRef<RVec>>(h_x), nullptr);
667 checkResults2dHaloWith1PulseInEachDim(h_x.data(), &dd, numHomeAtoms);
669 if (GMX_GPU_CUDA && GMX_THREAD_MPI) // repeat with GPU halo codepath
671 // early return if no devices are available.
672 if (getTestHardwareEnvironment()->getTestDeviceList().empty())
677 // Re-initialize input
678 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
680 // Perform GPU halo exchange
681 gpuHalo(&dd, box, &h_x, numAtomsTotal);
684 checkResults2dHaloWith1PulseInEachDim(h_x.data(), &dd, numHomeAtoms);
688 TEST(HaloExchangeTest, Coordinates2dHaloWith2PulsesInDim1)
690 GMX_MPI_TEST(RequireRankCount<4>);
693 const int numHomeAtoms = 10;
694 const int numHaloAtoms = 7;
695 const int numAtomsTotal = numHomeAtoms + numHaloAtoms;
696 HostVector<RVec> h_x;
697 h_x.resize(numAtomsTotal);
699 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
704 dd.mpi_comm_all = MPI_COMM_WORLD;
705 gmx_domdec_comm_t comm;
707 dd.unitCellInfo.haveScrewPBC = false;
709 DDAtomRanges atomRanges;
710 atomRanges.setEnd(DDAtomRanges::Type::Home, numHomeAtoms);
711 dd.comm->atomRanges = atomRanges;
713 define2dRankTopology(&dd);
715 std::vector<gmx_domdec_ind_t> indvec;
716 define2dHaloWith2PulsesInDim1(&dd, &indvec);
718 // Perform halo exchange
719 matrix box = { { 0., 0., 0. } };
720 dd_move_x(&dd, box, static_cast<ArrayRef<RVec>>(h_x), nullptr);
723 checkResults2dHaloWith2PulsesInDim1(h_x.data(), &dd, numHomeAtoms);
725 if (GMX_GPU_CUDA && GMX_THREAD_MPI) // repeat with GPU halo codepath
727 // early return if no devices are available.
728 if (getTestHardwareEnvironment()->getTestDeviceList().empty())
733 // Re-initialize input
734 initHaloData(h_x.data(), numHomeAtoms, numAtomsTotal);
736 // Perform GPU halo exchange
737 gpuHalo(&dd, box, &h_x, numAtomsTotal);
740 checkResults2dHaloWith2PulsesInDim1(h_x.data(), &dd, numHomeAtoms);