2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
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.
37 * \brief Declares implementation functions and types for the domain
38 * decomposition module.
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_domdec
43 #ifndef GMX_DOMDEC_DOMDEC_INTERNAL_H
44 #define GMX_DOMDEC_DOMDEC_INTERNAL_H
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/mdlib/updategroupscog.h"
51 #include "gromacs/timing/cyclecounter.h"
52 #include "gromacs/topology/block.h"
58 #define DD_NLOAD_MAX 9
62 //! Indices to communicate in a dimension
63 struct gmx_domdec_ind_t
66 /*! \brief The numbers of charge groups to send and receive for each
67 * cell that requires communication, the last entry contains the total
68 * number of atoms that needs to be communicated.
70 int nsend[DD_MAXIZONE+2] = {};
71 int nrecv[DD_MAXIZONE+2] = {};
73 //! The charge groups to send
74 std::vector<int> index;
76 /* The atom range for non-in-place communication */
77 int cell2at0[DD_MAXIZONE] = {};
78 int cell2at1[DD_MAXIZONE] = {};
82 //! Things relating to index communication
83 struct gmx_domdec_comm_dim_t
85 /* Returns the number of grid pulses (the number of domains in the halo along this dimension) */
91 /**< For dlb, for use with edlbAUTO */
93 /**< The indices to communicate, size np */
94 std::vector<gmx_domdec_ind_t> ind;
95 /**< Can we receive data in place? */
96 bool receiveInPlace = false;
99 /*! \brief Load balancing data along a dim used on the master rank of that dim */
104 /**< State var.: max lower bound., incl. neighbors */
105 real cellFracLowerMax = 0;
106 /**< State var.: min upper bound., incl. neighbors */
107 real cellFracUpperMin = 0;
108 /**< Temp. var.: lower limit for cell boundary */
110 /**< Temp. var.: upper limit for cell boundary */
114 /**< Temp. var.: is this cell size at the limit */
115 std::vector<bool> isCellMin;
116 /**< State var.: cell boundaries, box relative */
117 std::vector<real> cellFrac;
118 /**< Temp. var.: old cell size */
119 std::vector<real> oldCellFrac;
121 std::vector<Bounds> bounds;
122 /**< State var.: is DLB limited in this row */
123 bool dlbIsLimited = false;
125 std::vector<real> buf_ncd;
128 /*! \brief Struct for managing cell sizes with DLB along a dimension */
129 struct DDCellsizesWithDlb
131 /**< Cell row root struct, only available on the first rank in a row */
132 std::unique_ptr<RowMaster> rowMaster;
133 /**< The cell sizes, in fractions, along a row, not available on the first rank in a row */
134 std::vector<real> fracRow;
135 /**< The lower corner, in fractions, in triclinic space */
137 /**< The upper corner, in fractions, in triclinic space */
139 /**< The maximum lower corner among all our neighbors */
140 real fracLowerMax = 0;
141 /**< The minimum upper corner among all our neighbors */
142 real fracUpperMin = 0;
145 /*! \brief Struct for compute load commuication
147 * Here floats are accurate enough, since these variables
148 * only influence the load balancing, not the actual MD results.
152 /**< The number of load recordings */
154 /**< Scan of the sum of load over dimensions */
155 float *load = nullptr;
156 /**< The sum of the load over the ranks up to our current dimension */
158 /**< The maximum over the ranks contributing to \p sum */
160 /**< Like \p sum, but takes the maximum when the load balancing is limited */
162 /**< Minimum cell volume, relative to the box */
164 /**< The PP time during which PME can overlap */
166 /**< The PME-only rank load */
168 /**< Bit flags that tell if DLB was limited, per dimension */
172 /*! \brief Data needed to sort an atom to the desired location in the local state */
175 /**< Neighborsearch grid cell index */
177 /**< Global atom/charge group index */
179 /**< Local atom/charge group index */
183 /*! \brief Temporary buffers for sorting atoms */
186 /**< Sorted array of indices */
187 std::vector<gmx_cgsort_t> sorted;
188 /**< Array of stationary atom/charge group indices */
189 std::vector<gmx_cgsort_t> stationary;
190 /**< Array of moved atom/charge group indices */
191 std::vector<gmx_cgsort_t> moved;
192 /**< Integer buffer for sorting */
193 std::vector<int> intBuffer;
196 /*! \brief Manages atom ranges and order for the local state atom vectors */
200 /*! \brief The local state atom order
202 * This enum determines the order of the atoms in the local state.
203 * ddnatHOME and ddnatZONE should be first and second,
204 * the others can be ordered as wanted.
206 enum class Type : int
208 Home, /**< The home atoms */
209 Zones, /**< All zones in the eighth shell */
210 Vsites, /**< Atoms communicated for virtual sites */
211 Constraints, /**< Atoms communicated for constraints */
212 Number /**< Not a count, only present for convenience */
215 /*! \brief Returns the start atom index for range \p rangeType */
216 int start(Type rangeType) const
218 if (rangeType == Type::Home)
224 return end_[static_cast<int>(rangeType) - 1];
228 /*! \brief Returns the end atom index for range \p rangeType */
229 int end(Type rangeType) const
231 return end_[static_cast<int>(rangeType)];
234 /*! \brief Returns the number of home atoms */
235 int numHomeAtoms() const
237 return end_[static_cast<int>(Type::Home)];
240 /*! \brief Returns the total number of atoms */
241 int numAtomsTotal() const
243 return end_[static_cast<int>(Type::Number) - 1];
246 /*! \brief Sets the end index of range \p rangeType to \p end
248 * This should be called either with Type::Home or with a type
249 * that is larger than that passed in the previous call to setEnd.
250 * A release assertion for these conditions are present.
252 void setEnd(Type rangeType,
255 GMX_RELEASE_ASSERT(rangeType == Type::Home || rangeType > lastTypeSet_, "Can only set either home or a larger type than the last one");
257 for (int i = static_cast<int>(rangeType); i < static_cast<int>(Type::Number); i++)
262 lastTypeSet_ = rangeType;
266 /*! \brief The list of end atom indices */
267 std::array<int, static_cast<int>(Type::Number)> end_;
268 Type lastTypeSet_ = Type::Number;
271 /*! \brief Enum of dynamic load balancing states
273 * Allowed DLB states and transitions
274 * - intialization at startup:
275 * -> offUser ("-dlb no")
276 * -> onUser ("-dlb yes")
277 * -> offCanTurnOn ("-dlb auto")
279 * - in automatic mode (i.e. initial state offCanTurnOn):
280 * offCanTurnOn -> onCanTurnOff
281 * offCanTurnOn -> offForever
282 * offCanTurnOn -> offTemporarilyLocked
283 * offTemporarilyLocked -> offCanTurnOn
284 * onCanTurnOff -> offCanTurnOn
288 offUser, /**< DLB is permanently off per user request */
289 offForever, /**< DLB is off due to a runtime condition (not supported or causes performance loss) and will never be turned on */
290 offCanTurnOn, /**< DLB is off and will turn on on imbalance */
291 offTemporarilyLocked, /**< DLB is off and temporarily can't turn on */
292 onCanTurnOff, /**< DLB is on and can turn off when slow */
293 onUser, /**< DLB is permanently on per user request */
294 nr /**< The number of DLB states */
297 /*! \brief The PME domain decomposition for one dimension */
300 /**< The dimension */
302 /**< Tells if DD and PME dims match */
303 gmx_bool dim_match = false;
304 /**< The number of PME ranks/domains in this dimension */
306 /**< Cell sizes for determining the PME comm. with SLB */
307 real *slb_dim_f = nullptr;
308 /**< The minimum pp node location, size nslab */
309 int *pp_min = nullptr;
310 /**< The maximum pp node location, size nslab */
311 int *pp_max = nullptr;
312 /**< The maximum shift for coordinate redistribution in PME */
318 /**< The minimum bottom of this zone */
320 /**< The maximum top of this zone */
322 /**< The minimum top of this zone */
324 /**< The maximum bottom communicaton height for this zone */
326 /**< The maximum top communicaton height for this zone */
328 /**< The bottom value of the first cell in this zone */
330 /**< The top value of the first cell in this zone */
332 /**< Bool disguised as a real, 1 when the above data has been set. 0 otherwise */
336 /*! \brief The number of reals in gmx_ddzone_t */
337 constexpr int c_ddzoneNumReals = 8;
339 template<typename T> class DDBufferAccess;
341 /*! \brief Temporary storage container that minimizes (re)allocation and clearing
343 * This is only the storage, actual access happens through DDBufferAccess.
344 * All methods check if the buffer is (not) in use.
350 /*! \brief Returns a buffer of size \p numElements, the elements are undefined */
351 gmx::ArrayRef<T> resize(size_t numElements)
353 GMX_ASSERT(isInUse_, "Should only operate on acquired buffers");
355 if (numElements > buffer_.size())
357 buffer_.resize(numElements);
360 return gmx::arrayRefFromArray(buffer_.data(), numElements);
363 /*! \brief Acquire the buffer for use with size set to \p numElements, the elements are undefined */
364 gmx::ArrayRef<T> acquire(size_t numElements)
366 GMX_RELEASE_ASSERT(!isInUse_, "Should only request free buffers");
369 return resize(numElements);
372 /*! \brief Releases the buffer, buffer_ should not be used after this */
375 GMX_RELEASE_ASSERT(isInUse_, "Should only release buffers in use");
379 std::vector<T> buffer_; /**< The actual memory buffer */
380 bool isInUse_ = false; /**< Flag that tells whether the buffer is in use */
382 friend class DDBufferAccess<T>;
385 /*! \brief Class that manages access to a temporary memory buffer */
390 /*! \brief Constructor, returns a buffer of size \p numElements, element values are undefined
392 * \note The actual memory buffer \p ddBuffer can not be used to
393 * create other DDBufferAccess objects until the one created
396 DDBufferAccess(DDBuffer<T> &ddBuffer,
397 size_t numElements) :
400 buffer = ddBuffer_.acquire(numElements);
408 /*! \brief Resizes the buffer to \p numElements, new elements are undefined
410 * \note The buffer arrayref is updated after this call.
412 void resize(size_t numElements)
414 buffer = ddBuffer_.resize(numElements);
418 DDBuffer<T> &ddBuffer_; /**< Reference to the storage class */
420 gmx::ArrayRef<T> buffer; /**< The access to the memory buffer */
423 /*! \brief Temporary buffer for setting up communiation over one pulse and all zones in the halo */
424 struct dd_comm_setup_work_t
426 /**< The local atom group indices to send */
427 std::vector<int> localAtomGroupBuffer;
428 /**< Buffer for collecting the global atom group indices to send */
429 std::vector<int> atomGroupBuffer;
430 /**< Buffer for collecting the atom group positions to send */
431 std::vector<gmx::RVec> positionBuffer;
432 /**< The number of atoms contained in the atom groups to send */
434 /**< The number of atom groups to send for the last zone */
438 /*! \brief Information about the simulated system */
441 //! True when update groups are used
442 bool useUpdateGroups = false;
443 //! Update atom grouping for each molecule type
444 std::vector<gmx::RangePartitioning> updateGroupingPerMoleculetype;
445 //! The maximum radius over all update groups
446 real maxUpdateGroupRadius;
448 //! Are there inter-domain bonded interactions?
449 bool haveInterDomainBondeds = false;
450 //! Are there inter-domain multi-body interactions?
451 bool haveInterDomainMultiBodyBondeds = false;
453 //! Cut-off for multi-body interactions
454 real minCutoffForMultiBody = 0;
455 //! Cut-off for non-bonded/2-body interactions
457 //! The lower limit for the DD cell size
458 real cellsizeLimit = 0;
460 //! Can atoms connected by constraints be assigned to different domains?
461 bool haveSplitConstraints = false;
462 //! Can atoms connected by settles be assigned to different domains?
463 bool haveSplitSettles = false;
464 //! Estimated communication range needed for constraints
465 real constraintCommunicationRange = 0;
467 //! Whether to only communicate atoms beyond the non-bonded cut-off when they are involved in bonded interactions with non-local atoms
468 bool filterBondedCommunication = false;
469 //! Whether to increase the multi-body cut-off beyond the minimum required
470 bool increaseMultiBodyCutoff = false;
473 /*! \brief Settings that affect the behavior of the domain decomposition
475 * These settings depend on options chosen by the user, set by enviroment
476 * variables, as well as hardware support. The initial DLB state also
477 * depends on the integrator.
479 * Note: Settings that depend on the simulated system are in DDSystemInfo.
483 //! Use MPI_Sendrecv communication instead of non-blocking calls
484 bool useSendRecv2 = false;
486 /* Information for managing the dynamic load balancing */
487 //! Maximum DLB scaling per load balancing step in percent
488 int dlb_scale_lim = 0;
489 //! Flop counter (0=no,1=yes,2=with (eFlop-1)*5% noise
492 //! Whether to order the DD dimensions from z to x
493 bool useDDOrderZYX = false;
495 //! Whether to use MPI Cartesian reordering of communicators, when supported (almost never)
496 bool useCartesianReorder = true;
498 //! Whether we should record the load
499 bool recordLoad = false;
502 //! Step interval for dumping the local+non-local atoms to pdb
504 //! Step interval for duming the DD grid to pdb
505 int nstDDDumpGrid = 0;
506 //! DD debug print level: 0, 1, 2
509 //! The DLB state at the start of the run
510 DlbState initialDlbState = DlbState::offCanTurnOn;
513 /*! \brief Information on how the DD ranks are set up */
516 /**< The number of particle-particle (non PME-only) ranks */
518 /**< The DD PP grid */
519 ivec numPPCells = { 0, 0, 0 };
521 /* PME and Cartesian communicator stuff */
522 /**< The number of decomposition dimensions for PME, 0: no PME */
523 int npmedecompdim = 0;
524 /**< The number of ranks doing PME (PP/PME or only PME) */
526 /**< The number of PME ranks/domains along x */
528 /**< The number of PME ranks/domains along y */
530 /**< Use Cartesian communication between PP and PME ranks */
531 gmx_bool bCartesianPP_PME = false;
532 /**< Cartesian grid for combinted PP+PME ranks */
534 /**< The number of dimensions for the PME setup that are Cartesian */
536 /**< The PME ranks, size npmenodes */
537 int *pmenodes = nullptr;
538 /**< The Cartesian index to sim rank conversion, used with bCartesianPP_PME */
539 int *ddindex2simnodeid = nullptr;
540 /**< The 1D or 2D PME domain decomposition setup */
541 gmx_ddpme_t ddpme[2];
543 /* The DD particle-particle nodes only */
544 /**< Use a Cartesian communicator for PP */
545 gmx_bool bCartesianPP = false;
546 /**< The Cartesian index to DD rank conversion, used with bCartesianPP */
547 int *ddindex2ddnodeid = nullptr;
550 /*! \brief Struct for domain decomposition communication
552 * This struct contains most information about domain decomposition
553 * communication setup, some communication buffers, some statistics
554 * and also the setup for the communication between particle-particle
555 * and PME only ranks.
557 * All arrays are indexed with 0 to dd->ndim (not Cartesian indexing),
558 * unless stated otherwise.
560 struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
562 /**< Constant parameters that control DD behavior */
563 DDSettings ddSettings;
565 /**< Information on how the DD ranks are set up */
566 DDRankSetup ddRankSetup;
568 /* The DLB state, used for reloading old states, during e.g. EM */
569 /**< The global charge groups, this defined the DD state (except for the DLB state) */
570 t_block cgs_gl = { };
572 /* Charge group / atom sorting */
573 /**< Data structure for cg/atom sorting */
574 std::unique_ptr<gmx_domdec_sort_t> sort;
576 //! Centers of mass of local update groups
577 std::unique_ptr<gmx::UpdateGroupsCog> updateGroupsCog;
579 /* Data for the optional filtering of communication of atoms for bonded interactions */
580 /**< Links between cg's through bonded interactions */
581 t_blocka *cglink = nullptr;
582 /**< Local cg availability, TODO: remove when group scheme is removed */
583 char *bLocalCG = nullptr;
585 /* The DLB state, possible values are defined above */
587 /* With dlbState=DlbState::offCanTurnOn, should we check if to DLB on at the next DD? */
588 gmx_bool bCheckWhetherToTurnDlbOn = false;
589 /* The first DD count since we are running without DLB */
590 int ddPartioningCountFirstDlbOff = 0;
592 /* Cell sizes for static load balancing, first index cartesian */
593 real **slb_frac = nullptr;
595 /**< Information about the simulated system */
596 DDSystemInfo systemInfo;
598 /* The width of the communicated boundaries */
599 /**< Cut-off for multi-body interactions, also 2-body bonded when \p cutoff_mody > \p cutoff */
600 real cutoff_mbody = 0;
601 /**< The minimum guaranteed cell-size, Cartesian indexing */
602 rvec cellsize_min = { };
603 /**< The minimum guaranteed cell-size with dlb=auto */
604 rvec cellsize_min_dlb = { };
605 /**< The lower limit for the DD cell size with DLB */
606 real cellsize_limit = 0;
607 /**< Effectively no NB cut-off limit with DLB for systems without PBC? */
608 gmx_bool bVacDLBNoLimit = false;
610 /** With PME load balancing we set limits on DLB */
611 gmx_bool bPMELoadBalDLBLimits = false;
612 /** DLB needs to take into account that we want to allow this maximum
613 * cut-off (for PME load balancing), this could limit cell boundaries.
615 real PMELoadBal_max_cutoff = 0;
617 /**< tric_dir from \p gmx_ddbox_t is only stored here because dd_get_ns_ranges needs it */
619 /**< box lower corner, required with dim's without pbc and -gcom */
621 /**< box size, required with dim's without pbc and -gcom */
624 /**< The DD cell lower corner, in triclinic space */
626 /**< The DD cell upper corner, in triclinic space */
629 /**< The old \p cell_x0, to check cg displacements */
630 rvec old_cell_x0 = { };
631 /**< The old \p cell_x1, to check cg displacements */
632 rvec old_cell_x1 = { };
634 /** The communication setup and charge group boundaries for the zones */
635 gmx_domdec_zones_t zones;
637 /* The zone limits for DD dimensions 1 and 2 (not 0), determined from
638 * cell boundaries of neighboring cells for staggered grids when using
639 * dynamic load balancing.
641 /**< Zone limits for dim 1 with staggered grids */
642 gmx_ddzone_t zone_d1[2];
643 /**< Zone limits for dim 2 with staggered grids */
644 gmx_ddzone_t zone_d2[2][2];
646 /** The coordinate/force communication setup and indices */
647 gmx_domdec_comm_dim_t cd[DIM];
648 /** The maximum number of cells to communicate with in one dimension */
651 /** Which cg distribution is stored on the master node,
652 * stored as DD partitioning call count.
654 int64_t master_cg_ddp_count = 0;
656 /** The number of cg's received from the direct neighbors */
657 int zone_ncg1[DD_MAXZONE] = {0};
659 /** The atom ranges in the local state */
660 DDAtomRanges atomRanges;
662 /** Array for signalling if atoms have moved to another domain */
663 std::vector<int> movedBuffer;
665 /** Communication int buffer for general use */
666 DDBuffer<int> intBuffer;
668 /** Communication rvec buffer for general use */
669 DDBuffer<gmx::RVec> rvecBuffer;
671 /* Temporary storage for thread parallel communication setup */
672 /**< Thread-local work data */
673 std::vector<dd_comm_setup_work_t> dth;
675 /* Communication buffer only used with multiple grid pulses */
676 /**< Another rvec comm. buffer */
677 DDBuffer<gmx::RVec> rvecBuffer2;
679 /* Communication buffers for local redistribution */
680 /**< Charge group flag comm. buffers */
681 std::array<std::vector<int>, DIM*2> cggl_flag;
682 /**< Charge group center comm. buffers */
683 std::array<std::vector<gmx::RVec>, DIM*2> cgcm_state;
685 /* Cell sizes for dynamic load balancing */
686 std::vector<DDCellsizesWithDlb> cellsizesWithDlb;
688 /* Stuff for load communication */
689 /**< The recorded load data */
690 domdec_load_t *load = nullptr;
691 /**< The number of MPI ranks sharing the GPU our rank is using */
692 int nrank_gpu_shared = 0;
694 /**< The MPI load communicator */
695 MPI_Comm *mpi_comm_load = nullptr;
696 /**< The MPI load communicator for ranks sharing a GPU */
697 MPI_Comm mpi_comm_gpu_shared;
700 /**< Struct for timing the force load balancing region */
701 BalanceRegion *balanceRegion = nullptr;
703 /* Cycle counters over nstlist steps */
704 /**< Total cycles counted */
705 float cycl[ddCyclNr] = { };
706 /**< The number of cycle recordings */
707 int cycl_n[ddCyclNr] = { };
708 /**< The maximum cycle count */
709 float cycl_max[ddCyclNr] = { };
710 /**< Total flops counted */
712 /**< The number of flop recordings */
714 /** How many times did we have load measurements */
716 /** How many times have we collected the load measurements */
717 int n_load_collect = 0;
719 /* Cycle count history for DLB checks */
720 /**< The averaged cycles per step over the last nstlist step before turning on DLB */
721 float cyclesPerStepBeforeDLB = 0;
722 /**< The running average of the cycles per step during DLB */
723 float cyclesPerStepDlbExpAverage = 0;
724 /**< Have we turned off DLB (after turning DLB on)? */
725 bool haveTurnedOffDlb = false;
726 /**< The DD step at which we last measured that DLB off was faster than DLB on, 0 if there was no such step */
727 int64_t dlbSlowerPartitioningCount = 0;
729 /* Statistics for atoms */
730 /**< The atoms per range, summed over the steps */
731 double sum_nat[static_cast<int>(DDAtomRanges::Type::Number)] = { };
733 /* Statistics for calls and times */
734 /**< The number of partioning calls */
736 /**< The number of load recordings */
738 /**< Total MD step time */
739 double load_step = 0.0;
740 /**< Total PP force time */
741 double load_sum = 0.0;
742 /**< Max \p load_sum over the ranks */
743 double load_max = 0.0;
744 /**< Was load balancing limited, per DD dim */
746 /**< Total time on PP done during PME overlap time */
747 double load_mdf = 0.0;
748 /**< Total time on our PME-only rank */
749 double load_pme = 0.0;
751 /** The last partition step */
752 int64_t partition_step = INT_MIN;
755 /*! \brief DD zone permutation
757 * Zone permutation from the Cartesian x-major/z-minor order to an order
758 * that leads to consecutive charge groups for neighbor searching.
759 * TODO: remove when the group scheme is removed
761 static const int zone_perm[3][4] = { {0, 0, 0, 0}, {1, 0, 0, 0}, {3, 0, 1, 2} };
763 /*! \brief DD zone reordering to Cartesian order
765 * Index to reorder the zone such that the end up in Cartesian order
766 * with dimension index 0 major and dimension index 2 minor.
768 static const int zone_reorder_cartesian[DD_MAXZONE] = { 0, 1, 3, 2, 5, 4, 6, 7 };
770 /* dd_zo and dd_zp3 is set up such that i zones with non-zero
771 * components see only j zones with that component 0.
774 /*! \brief Returns the DD cut-off distance for multi-body interactions */
775 real dd_cutoff_multibody(const gmx_domdec_t *dd);
777 /*! \brief Returns the DD cut-off distance for two-body interactions */
778 real dd_cutoff_twobody(const gmx_domdec_t *dd);
780 /*! \brief Returns the domain index given the number of domains and the domain coordinates
782 * This order is required to minimize the coordinate communication in PME
783 * which uses decomposition in the x direction.
785 static inline int dd_index(const ivec numDomains,
786 const ivec domainCoordinates)
788 return ((domainCoordinates[XX]*numDomains[YY] + domainCoordinates[YY])*numDomains[ZZ]) + domainCoordinates[ZZ];
791 /*! Returns the size of the buffer to hold fractional cell boundaries for DD dimension index dimIndex */
792 static inline int ddCellFractionBufferSize(const gmx_domdec_t *dd,
795 return dd->nc[dd->dim[dimIndex ]] + 1 + dimIndex*2 + 1 + dimIndex;
798 /*! \brief Maximum number of ranks for using send/recv for state scattering and gathering
800 * Use separate MPI send and receive commands
801 * when #ranks <= c_maxNumRanksUseSendRecvForScatterAndGather
802 * This saves memory (and some copying for small #ranks).
803 * For high parallelization scatter and gather calls are used.
805 static constexpr int c_maxNumRanksUseSendRecvForScatterAndGather = 4;
807 /*! \brief Make DD cells larger by this factor than the limit to avoid rounding issues */
808 static constexpr double DD_CELL_MARGIN = 1.0001;
810 /*! \brief Factor for checking DD cell size limitation during DLB, should be in between 1 and DD_CELL_MARGIN */
811 static constexpr double DD_CELL_MARGIN2 = 1.00005;
813 /*! \brief With pressure scaling, keep cell sizes 2% above the limit to allow for some scaling */
814 static constexpr double DD_PRES_SCALE_MARGIN = 1.02;