2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Declares implementation functions and types for the domain
39 * decomposition module.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
44 #ifndef GMX_DOMDEC_DOMDEC_INTERNAL_H
45 #define GMX_DOMDEC_DOMDEC_INTERNAL_H
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/domdec/domdec_struct.h"
51 #include "gromacs/mdlib/updategroupscog.h"
52 #include "gromacs/timing/cyclecounter.h"
53 #include "gromacs/topology/block.h"
59 #define DD_NLOAD_MAX 9
63 //! Indices to communicate in a dimension
64 struct gmx_domdec_ind_t
67 /*! \brief The numbers of charge groups to send and receive for each
68 * cell that requires communication, the last entry contains the total
69 * number of atoms that needs to be communicated.
71 int nsend[DD_MAXIZONE + 2] = {};
72 int nrecv[DD_MAXIZONE + 2] = {};
74 //! The charge groups to send
75 std::vector<int> index;
77 /* The atom range for non-in-place communication */
78 int cell2at0[DD_MAXIZONE] = {};
79 int cell2at1[DD_MAXIZONE] = {};
83 //! Things relating to index communication
84 struct gmx_domdec_comm_dim_t
86 /* Returns the number of grid pulses (the number of domains in the halo along this dimension) */
87 int numPulses() const { return ind.size(); }
89 /**< For dlb, for use with edlbAUTO */
91 /**< The indices to communicate, size np */
92 std::vector<gmx_domdec_ind_t> ind;
93 /**< Can we receive data in place? */
94 bool receiveInPlace = false;
97 /*! \brief Load balancing data along a dim used on the master rank of that dim */
102 /**< State var.: max lower bound., incl. neighbors */
103 real cellFracLowerMax = 0;
104 /**< State var.: min upper bound., incl. neighbors */
105 real cellFracUpperMin = 0;
106 /**< Temp. var.: lower limit for cell boundary */
108 /**< Temp. var.: upper limit for cell boundary */
112 /**< Temp. var.: is this cell size at the limit */
113 std::vector<bool> isCellMin;
114 /**< State var.: cell boundaries, box relative */
115 std::vector<real> cellFrac;
116 /**< Temp. var.: old cell size */
117 std::vector<real> oldCellFrac;
119 std::vector<Bounds> bounds;
120 /**< State var.: is DLB limited in this row */
121 bool dlbIsLimited = false;
123 std::vector<real> buf_ncd;
126 /*! \brief Struct for managing cell sizes with DLB along a dimension */
127 struct DDCellsizesWithDlb
129 /**< Cell row root struct, only available on the first rank in a row */
130 std::unique_ptr<RowMaster> rowMaster;
131 /**< The cell sizes, in fractions, along a row, not available on the first rank in a row */
132 std::vector<real> fracRow;
133 /**< The lower corner, in fractions, in triclinic space */
135 /**< The upper corner, in fractions, in triclinic space */
137 /**< The maximum lower corner among all our neighbors */
138 real fracLowerMax = 0;
139 /**< The minimum upper corner among all our neighbors */
140 real fracUpperMin = 0;
143 /*! \brief Struct for compute load commuication
145 * Here floats are accurate enough, since these variables
146 * only influence the load balancing, not the actual MD results.
148 typedef struct domdec_load
150 /**< The number of load recordings */
152 /**< Scan of the sum of load over dimensions */
153 float* load = nullptr;
154 /**< The sum of the load over the ranks up to our current dimension */
156 /**< The maximum over the ranks contributing to \p sum */
158 /**< Like \p sum, but takes the maximum when the load balancing is limited */
160 /**< Minimum cell volume, relative to the box */
162 /**< The PP time during which PME can overlap */
164 /**< The PME-only rank load */
166 /**< Bit flags that tell if DLB was limited, per dimension */
170 /*! \brief Data needed to sort an atom to the desired location in the local state */
171 typedef struct gmx_cgsort
173 /**< Neighborsearch grid cell index */
175 /**< Global atom/charge group index */
177 /**< Local atom/charge group index */
181 /*! \brief Temporary buffers for sorting atoms */
182 typedef struct gmx_domdec_sort
184 /**< Sorted array of indices */
185 std::vector<gmx_cgsort_t> sorted;
186 /**< Array of stationary atom/charge group indices */
187 std::vector<gmx_cgsort_t> stationary;
188 /**< Array of moved atom/charge group indices */
189 std::vector<gmx_cgsort_t> moved;
190 /**< Integer buffer for sorting */
191 std::vector<int> intBuffer;
194 /*! \brief Manages atom ranges and order for the local state atom vectors */
198 /*! \brief The local state atom order
200 * This enum determines the order of the atoms in the local state.
201 * ddnatHOME and ddnatZONE should be first and second,
202 * the others can be ordered as wanted.
204 enum class Type : int
206 Home, /**< The home atoms */
207 Zones, /**< All zones in the eighth shell */
208 Vsites, /**< Atoms communicated for virtual sites */
209 Constraints, /**< Atoms communicated for constraints */
210 Number /**< Not a count, only present for convenience */
213 /*! \brief Returns the start atom index for range \p rangeType */
214 int start(Type rangeType) const
216 if (rangeType == Type::Home)
222 return end_[static_cast<int>(rangeType) - 1];
226 /*! \brief Returns the end atom index for range \p rangeType */
227 int end(Type rangeType) const { return end_[static_cast<int>(rangeType)]; }
229 /*! \brief Returns the number of home atoms */
230 int numHomeAtoms() const { return end_[static_cast<int>(Type::Home)]; }
232 /*! \brief Returns the total number of atoms */
233 int numAtomsTotal() const { return end_[static_cast<int>(Type::Number) - 1]; }
235 /*! \brief Sets the end index of range \p rangeType to \p end
237 * This should be called either with Type::Home or with a type
238 * that is larger than that passed in the previous call to setEnd.
239 * A release assertion for these conditions are present.
241 void setEnd(Type rangeType, int end)
243 GMX_RELEASE_ASSERT(rangeType == Type::Home || rangeType > lastTypeSet_,
244 "Can only set either home or a larger type than the last one");
246 for (int i = static_cast<int>(rangeType); i < static_cast<int>(Type::Number); i++)
251 lastTypeSet_ = rangeType;
255 /*! \brief The list of end atom indices */
256 std::array<int, static_cast<int>(Type::Number)> end_;
257 Type lastTypeSet_ = Type::Number;
260 /*! \brief Enum of dynamic load balancing states
262 * Allowed DLB states and transitions
263 * - intialization at startup:
264 * -> offUser ("-dlb no")
265 * -> onUser ("-dlb yes")
266 * -> offCanTurnOn ("-dlb auto")
268 * - in automatic mode (i.e. initial state offCanTurnOn):
269 * offCanTurnOn -> onCanTurnOff
270 * offCanTurnOn -> offForever
271 * offCanTurnOn -> offTemporarilyLocked
272 * offTemporarilyLocked -> offCanTurnOn
273 * onCanTurnOff -> offCanTurnOn
277 offUser, /**< DLB is permanently off per user request */
278 offForever, /**< DLB is off due to a runtime condition (not supported or causes performance loss) and will never be turned on */
279 offCanTurnOn, /**< DLB is off and will turn on on imbalance */
280 offTemporarilyLocked, /**< DLB is off and temporarily can't turn on */
281 onCanTurnOff, /**< DLB is on and can turn off when slow */
282 onUser, /**< DLB is permanently on per user request */
283 nr /**< The number of DLB states */
286 /*! \brief The PME domain decomposition for one dimension */
287 typedef struct gmx_ddpme
289 /**< The dimension */
291 /**< Tells if DD and PME dims match */
292 gmx_bool dim_match = false;
293 /**< The number of PME ranks/domains in this dimension */
295 /**< Cell sizes for determining the PME comm. with SLB */
296 real* slb_dim_f = nullptr;
297 /**< The minimum pp node location, size nslab */
298 int* pp_min = nullptr;
299 /**< The maximum pp node location, size nslab */
300 int* pp_max = nullptr;
301 /**< The maximum shift for coordinate redistribution in PME */
307 /**< The minimum bottom of this zone */
309 /**< The maximum top of this zone */
311 /**< The minimum top of this zone */
313 /**< The maximum bottom communicaton height for this zone */
315 /**< The maximum top communicaton height for this zone */
317 /**< The bottom value of the first cell in this zone */
319 /**< The top value of the first cell in this zone */
321 /**< Bool disguised as a real, 1 when the above data has been set. 0 otherwise */
325 /*! \brief The number of reals in gmx_ddzone_t */
326 constexpr int c_ddzoneNumReals = 8;
329 class DDBufferAccess;
331 /*! \brief Temporary storage container that minimizes (re)allocation and clearing
333 * This is only the storage, actual access happens through DDBufferAccess.
334 * All methods check if the buffer is (not) in use.
340 /*! \brief Returns a buffer of size \p numElements, the elements are undefined */
341 gmx::ArrayRef<T> resize(size_t numElements)
343 GMX_ASSERT(isInUse_, "Should only operate on acquired buffers");
345 if (numElements > buffer_.size())
347 buffer_.resize(numElements);
350 return gmx::arrayRefFromArray(buffer_.data(), numElements);
353 /*! \brief Acquire the buffer for use with size set to \p numElements, the elements are undefined */
354 gmx::ArrayRef<T> acquire(size_t numElements)
356 GMX_RELEASE_ASSERT(!isInUse_, "Should only request free buffers");
359 return resize(numElements);
362 /*! \brief Releases the buffer, buffer_ should not be used after this */
365 GMX_RELEASE_ASSERT(isInUse_, "Should only release buffers in use");
369 std::vector<T> buffer_; /**< The actual memory buffer */
370 bool isInUse_ = false; /**< Flag that tells whether the buffer is in use */
372 friend class DDBufferAccess<T>;
375 /*! \brief Class that manages access to a temporary memory buffer */
380 /*! \brief Constructor, returns a buffer of size \p numElements, element values are undefined
382 * \note The actual memory buffer \p ddBuffer can not be used to
383 * create other DDBufferAccess objects until the one created
386 DDBufferAccess(DDBuffer<T>& ddBuffer, size_t numElements) : ddBuffer_(ddBuffer)
388 buffer = ddBuffer_.acquire(numElements);
391 ~DDBufferAccess() { ddBuffer_.release(); }
393 /*! \brief Resizes the buffer to \p numElements, new elements are undefined
395 * \note The buffer arrayref is updated after this call.
397 void resize(size_t numElements) { buffer = ddBuffer_.resize(numElements); }
400 DDBuffer<T>& ddBuffer_; /**< Reference to the storage class */
402 gmx::ArrayRef<T> buffer; /**< The access to the memory buffer */
405 /*! \brief Temporary buffer for setting up communiation over one pulse and all zones in the halo */
406 struct dd_comm_setup_work_t
408 /**< The local atom group indices to send */
409 std::vector<int> localAtomGroupBuffer;
410 /**< Buffer for collecting the global atom group indices to send */
411 std::vector<int> atomGroupBuffer;
412 /**< Buffer for collecting the atom group positions to send */
413 std::vector<gmx::RVec> positionBuffer;
414 /**< The number of atoms contained in the atom groups to send */
416 /**< The number of atom groups to send for the last zone */
420 /*! \brief Information about the simulated system */
423 //! True when update groups are used
424 bool useUpdateGroups = false;
425 //! Update atom grouping for each molecule type
426 std::vector<gmx::RangePartitioning> updateGroupingPerMoleculetype;
427 //! The maximum radius over all update groups
428 real maxUpdateGroupRadius;
430 //! Are molecules always whole, i.e. not broken by PBC?
431 bool moleculesAreAlwaysWhole = false;
432 //! Are there inter-domain bonded interactions?
433 bool haveInterDomainBondeds = false;
434 //! Are there inter-domain multi-body interactions?
435 bool haveInterDomainMultiBodyBondeds = false;
437 //! Cut-off for multi-body interactions
438 real minCutoffForMultiBody = 0;
439 //! Cut-off for non-bonded/2-body interactions
441 //! The lower limit for the DD cell size
442 real cellsizeLimit = 0;
444 //! Can atoms connected by constraints be assigned to different domains?
445 bool haveSplitConstraints = false;
446 //! Can atoms connected by settles be assigned to different domains?
447 bool haveSplitSettles = false;
448 //! Estimated communication range needed for constraints
449 real constraintCommunicationRange = 0;
451 //! Whether to only communicate atoms beyond the non-bonded cut-off when they are involved in bonded interactions with non-local atoms
452 bool filterBondedCommunication = false;
453 //! Whether to increase the multi-body cut-off beyond the minimum required
454 bool increaseMultiBodyCutoff = false;
457 /*! \brief Settings that affect the behavior of the domain decomposition
459 * These settings depend on options chosen by the user, set by enviroment
460 * variables, as well as hardware support. The initial DLB state also
461 * depends on the integrator.
463 * Note: Settings that depend on the simulated system are in DDSystemInfo.
467 //! Use MPI_Sendrecv communication instead of non-blocking calls
468 bool useSendRecv2 = false;
470 /* Information for managing the dynamic load balancing */
471 //! Maximum DLB scaling per load balancing step in percent
472 int dlb_scale_lim = 0;
473 //! Flop counter (0=no,1=yes,2=with (eFlop-1)*5% noise
476 //! Request 1D domain decomposition
479 //! Whether to order the DD dimensions from z to x
480 bool useDDOrderZYX = false;
482 //! Whether to use MPI Cartesian reordering of communicators, when supported (almost never)
483 bool useCartesianReorder = true;
485 //! Whether we should record the load
486 bool recordLoad = false;
489 //! Step interval for dumping the local+non-local atoms to pdb
491 //! Step interval for duming the DD grid to pdb
492 int nstDDDumpGrid = 0;
493 //! DD debug print level: 0, 1, 2
496 //! The DLB state at the start of the run
497 DlbState initialDlbState = DlbState::offCanTurnOn;
500 /*! \brief Information on how the DD ranks are set up */
503 /**< The number of particle-particle (non PME-only) ranks */
505 /**< The DD PP grid */
506 ivec numPPCells = { 0, 0, 0 };
508 /* PME and Cartesian communicator stuff */
509 bool usePmeOnlyRanks = false;
510 /**< The number of decomposition dimensions for PME, 0: no PME */
511 int npmedecompdim = 0;
512 /**< The number of ranks doing PME (PP/PME or only PME) */
513 int numRanksDoingPme = 0;
514 /**< The number of PME ranks/domains along x */
516 /**< The number of PME ranks/domains along y */
518 /**< The 1D or 2D PME domain decomposition setup */
519 gmx_ddpme_t ddpme[2];
522 /*! \brief Information on Cartesian MPI setup of the DD ranks */
523 struct CartesianRankSetup
525 /**< Use Cartesian communication between PP and PME ranks */
526 bool bCartesianPP_PME = false;
527 /**< Cartesian grid for combinted PP+PME ranks */
529 /**< The number of dimensions for the PME setup that are Cartesian */
531 /**< The Cartesian index to sim rank conversion, used with bCartesianPP_PME */
532 std::vector<int> ddindex2simnodeid;
534 /* The DD particle-particle nodes only */
535 /**< Use a Cartesian communicator for PP */
536 bool bCartesianPP = false;
537 /**< The Cartesian index to DD rank conversion, used with bCartesianPP */
538 std::vector<int> ddindex2ddnodeid;
541 /*! \brief Struct for domain decomposition communication
543 * This struct contains most information about domain decomposition
544 * communication setup, some communication buffers, some statistics
545 * and also the setup for the communication between particle-particle
546 * and PME only ranks.
548 * All arrays are indexed with 0 to dd->ndim (not Cartesian indexing),
549 * unless stated otherwise.
551 struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
553 /**< Constant parameters that control DD behavior */
554 DDSettings ddSettings;
556 /**< Information on how the DD ranks are set up */
557 DDRankSetup ddRankSetup;
558 /**< Information on the Cartesian part of the DD rank setup */
559 CartesianRankSetup cartesianRankSetup;
561 /* Charge group / atom sorting */
562 /**< Data structure for cg/atom sorting */
563 std::unique_ptr<gmx_domdec_sort_t> sort;
565 //! Centers of mass of local update groups
566 std::unique_ptr<gmx::UpdateGroupsCog> updateGroupsCog;
568 /* Data for the optional filtering of communication of atoms for bonded interactions */
569 /**< Links between atoms through bonded interactions */
570 t_blocka* bondedLinks = nullptr;
572 /* The DLB state, possible values are defined above */
574 /* With dlbState=DlbState::offCanTurnOn, should we check if to DLB on at the next DD? */
575 gmx_bool bCheckWhetherToTurnDlbOn = false;
576 /* The first DD count since we are running without DLB */
577 int ddPartioningCountFirstDlbOff = 0;
579 /* Cell sizes for static load balancing, first index cartesian */
580 real** slb_frac = nullptr;
582 /**< Information about the simulated system */
583 DDSystemInfo systemInfo;
585 /* The width of the communicated boundaries */
586 /**< Cut-off for multi-body interactions, also 2-body bonded when \p cutoff_mody > \p cutoff */
587 real cutoff_mbody = 0;
588 /**< The minimum guaranteed cell-size, Cartesian indexing */
589 rvec cellsize_min = {};
590 /**< The minimum guaranteed cell-size with dlb=auto */
591 rvec cellsize_min_dlb = {};
592 /**< The lower limit for the DD cell size with DLB */
593 real cellsize_limit = 0;
594 /**< Effectively no NB cut-off limit with DLB for systems without PBC? */
595 gmx_bool bVacDLBNoLimit = false;
597 /** With PME load balancing we set limits on DLB */
598 gmx_bool bPMELoadBalDLBLimits = false;
599 /** DLB needs to take into account that we want to allow this maximum
600 * cut-off (for PME load balancing), this could limit cell boundaries.
602 real PMELoadBal_max_cutoff = 0;
604 /**< box lower corner, required with dim's without pbc and -gcom */
606 /**< box size, required with dim's without pbc and -gcom */
609 /**< The DD cell lower corner, in triclinic space */
611 /**< The DD cell upper corner, in triclinic space */
614 /**< The old \p cell_x0, to check cg displacements */
615 rvec old_cell_x0 = {};
616 /**< The old \p cell_x1, to check cg displacements */
617 rvec old_cell_x1 = {};
619 /** The communication setup and charge group boundaries for the zones */
620 gmx_domdec_zones_t zones;
622 /* The zone limits for DD dimensions 1 and 2 (not 0), determined from
623 * cell boundaries of neighboring cells for staggered grids when using
624 * dynamic load balancing.
626 /**< Zone limits for dim 1 with staggered grids */
627 gmx_ddzone_t zone_d1[2];
628 /**< Zone limits for dim 2 with staggered grids */
629 gmx_ddzone_t zone_d2[2][2];
631 /** The coordinate/force communication setup and indices */
632 gmx_domdec_comm_dim_t cd[DIM];
633 /** Restricts the maximum number of cells to communicate with in one dimension
635 * Dynamic load balancing is not permitted to change sizes if it
636 * would violate this restriction. */
639 /** Which cg distribution is stored on the master node,
640 * stored as DD partitioning call count.
642 int64_t master_cg_ddp_count = 0;
644 /** The number of cg's received from the direct neighbors */
645 int zone_ncg1[DD_MAXZONE] = { 0 };
647 /** The atom ranges in the local state */
648 DDAtomRanges atomRanges;
650 /** Array for signalling if atoms have moved to another domain */
651 std::vector<int> movedBuffer;
653 /** Communication int buffer for general use */
654 DDBuffer<int> intBuffer;
656 /** Communication rvec buffer for general use */
657 DDBuffer<gmx::RVec> rvecBuffer;
659 /* Temporary storage for thread parallel communication setup */
660 /**< Thread-local work data */
661 std::vector<dd_comm_setup_work_t> dth;
663 /* Communication buffer only used with multiple grid pulses */
664 /**< Another rvec comm. buffer */
665 DDBuffer<gmx::RVec> rvecBuffer2;
667 /* Communication buffers for local redistribution */
668 /**< Charge group flag comm. buffers */
669 std::array<std::vector<int>, DIM * 2> cggl_flag;
670 /**< Charge group center comm. buffers */
671 std::array<std::vector<gmx::RVec>, DIM * 2> cgcm_state;
673 /* Cell sizes for dynamic load balancing */
674 std::vector<DDCellsizesWithDlb> cellsizesWithDlb;
676 /* Stuff for load communication */
677 /**< The recorded load data */
678 domdec_load_t* load = nullptr;
679 /**< The number of MPI ranks sharing the GPU our rank is using */
680 int nrank_gpu_shared = 0;
682 /**< The MPI load communicator */
683 MPI_Comm* mpi_comm_load = nullptr;
684 /**< The MPI load communicator for ranks sharing a GPU */
685 MPI_Comm mpi_comm_gpu_shared;
688 /**< Struct for timing the force load balancing region */
689 BalanceRegion* balanceRegion = nullptr;
691 /* Cycle counters over nstlist steps */
692 /**< Total cycles counted */
693 float cycl[ddCyclNr] = {};
694 /**< The number of cycle recordings */
695 int cycl_n[ddCyclNr] = {};
696 /**< The maximum cycle count */
697 float cycl_max[ddCyclNr] = {};
698 /**< Total flops counted */
700 /**< The number of flop recordings */
702 /** How many times did we have load measurements */
704 /** How many times have we collected the load measurements */
705 int n_load_collect = 0;
707 /* Cycle count history for DLB checks */
708 /**< The averaged cycles per step over the last nstlist step before turning on DLB */
709 float cyclesPerStepBeforeDLB = 0;
710 /**< The running average of the cycles per step during DLB */
711 float cyclesPerStepDlbExpAverage = 0;
712 /**< Have we turned off DLB (after turning DLB on)? */
713 bool haveTurnedOffDlb = false;
714 /**< The DD step at which we last measured that DLB off was faster than DLB on, 0 if there was no such step */
715 int64_t dlbSlowerPartitioningCount = 0;
717 /* Statistics for atoms */
718 /**< The atoms per range, summed over the steps */
719 double sum_nat[static_cast<int>(DDAtomRanges::Type::Number)] = {};
721 /* Statistics for calls and times */
722 /**< The number of partioning calls */
724 /**< The number of load recordings */
726 /**< Total MD step time */
727 double load_step = 0.0;
728 /**< Total PP force time */
729 double load_sum = 0.0;
730 /**< Max \p load_sum over the ranks */
731 double load_max = 0.0;
732 /**< Was load balancing limited, per DD dim */
734 /**< Total time on PP done during PME overlap time */
735 double load_mdf = 0.0;
736 /**< Total time on our PME-only rank */
737 double load_pme = 0.0;
739 /** The last partition step */
740 int64_t partition_step = INT_MIN;
743 /*! \brief DD zone permutation
745 * Zone permutation from the Cartesian x-major/z-minor order to an order
746 * that leads to consecutive charge groups for neighbor searching.
747 * TODO: It should be possible to remove this now that the group scheme is removed
749 static const int zone_perm[3][4] = { { 0, 0, 0, 0 }, { 1, 0, 0, 0 }, { 3, 0, 1, 2 } };
751 /*! \brief DD zone reordering to Cartesian order
753 * Index to reorder the zone such that the end up in Cartesian order
754 * with dimension index 0 major and dimension index 2 minor.
756 static const int zone_reorder_cartesian[DD_MAXZONE] = { 0, 1, 3, 2, 5, 4, 6, 7 };
758 /* dd_zo and dd_zp3 is set up such that i zones with non-zero
759 * components see only j zones with that component 0.
762 /*! \brief Returns the DD cut-off distance for multi-body interactions */
763 real dd_cutoff_multibody(const gmx_domdec_t* dd);
765 /*! \brief Returns the DD cut-off distance for two-body interactions */
766 real dd_cutoff_twobody(const gmx_domdec_t* dd);
768 /*! \brief Returns the domain index given the number of domains and the domain coordinates
770 * This order is required to minimize the coordinate communication in PME
771 * which uses decomposition in the x direction.
773 static inline int dd_index(const ivec numDomains, const ivec domainCoordinates)
775 return ((domainCoordinates[XX] * numDomains[YY] + domainCoordinates[YY]) * numDomains[ZZ])
776 + domainCoordinates[ZZ];
779 /*! Returns the size of the buffer to hold fractional cell boundaries for DD dimension index dimIndex */
780 static inline int ddCellFractionBufferSize(const gmx_domdec_t* dd, int dimIndex)
782 return dd->numCells[dd->dim[dimIndex]] + 1 + dimIndex * 2 + 1 + dimIndex;
785 /*! \brief Maximum number of ranks for using send/recv for state scattering and gathering
787 * Use separate MPI send and receive commands
788 * when #ranks <= c_maxNumRanksUseSendRecvForScatterAndGather
789 * This saves memory (and some copying for small #ranks).
790 * For high parallelization scatter and gather calls are used.
792 static constexpr int c_maxNumRanksUseSendRecvForScatterAndGather = 4;
794 /*! \brief Make DD cells larger by this factor than the limit to avoid rounding issues */
795 static constexpr double DD_CELL_MARGIN = 1.0001;
797 /*! \brief Factor for checking DD cell size limitation during DLB, should be in between 1 and DD_CELL_MARGIN */
798 static constexpr double DD_CELL_MARGIN2 = 1.00005;
800 /*! \brief With pressure scaling, keep cell sizes 2% above the limit to allow for some scaling */
801 static constexpr double DD_PRES_SCALE_MARGIN = 1.02;