2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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 This file implements functions to interact with the dynamic load
38 * balancing machinery.
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_domdec
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/utility/gmxassert.h"
52 #include "domdec_internal.h"
55 float dd_pme_f_ratio(const gmx_domdec_t
*dd
)
57 GMX_ASSERT(DDMASTER(dd
), "This function should only be called on the master rank");
59 if (dd
->comm
->load
[0].mdf
> 0 && dd
->comm
->cycl_n
[ddCyclPME
] > 0)
61 return dd
->comm
->load
[0].pme
/dd
->comm
->load
[0].mdf
;
69 void set_dlb_limits(gmx_domdec_t
*dd
)
72 for (int d
= 0; d
< dd
->ndim
; d
++)
74 /* Set the number of pulses to the value for DLB */
75 dd
->comm
->cd
[d
].ind
.resize(dd
->comm
->cd
[d
].np_dlb
);
77 dd
->comm
->cellsize_min
[dd
->dim
[d
]] =
78 dd
->comm
->cellsize_min_dlb
[dd
->dim
[d
]];
82 void dd_dlb_set_should_check_whether_to_turn_dlb_on(gmx_domdec_t
*dd
, gmx_bool bValue
)
84 if (dd
->comm
->dlbState
== DlbState::offCanTurnOn
)
86 dd
->comm
->bCheckWhetherToTurnDlbOn
= bValue
;
90 /* Store the DD partitioning count, so we can ignore cycle counts
91 * over the next nstlist steps, which are often slower.
93 dd
->comm
->ddPartioningCountFirstDlbOff
= dd
->ddp_count
;
98 gmx_bool
dd_dlb_get_should_check_whether_to_turn_dlb_on(gmx_domdec_t
*dd
)
100 if (dd
->comm
->dlbState
!= DlbState::offCanTurnOn
)
105 if (dd
->ddp_count
<= dd
->comm
->ddPartioningCountFirstDlbOff
)
107 /* We ignore the first nstlist steps at the start of the run
108 * or after PME load balancing or after turning DLB off, since
109 * these often have extra allocation or cache miss overhead.
114 if (dd
->comm
->cycl_n
[ddCyclStep
] == 0)
116 /* We can have zero timed steps when dd_partition_system is called
117 * more than once at the same step, e.g. with replica exchange.
118 * Turning on DLB would trigger an assertion failure later, but is
119 * also useless right after exchanging replicas.
124 /* We should check whether we should use DLB directly after
126 if (dd
->comm
->bCheckWhetherToTurnDlbOn
)
128 /* This flag was set when the PME load-balancing routines
129 unlocked DLB, and should now be cleared. */
130 dd_dlb_set_should_check_whether_to_turn_dlb_on(dd
, FALSE
);
133 /* We check whether we should use DLB every c_checkTurnDlbOnInterval
134 * partitionings (we do not do this every partioning, so that we
135 * avoid excessive communication). */
136 return dd
->comm
->n_load_have
% c_checkTurnDlbOnInterval
== c_checkTurnDlbOnInterval
- 1;
139 gmx_bool
dd_dlb_is_on(const gmx_domdec_t
*dd
)
141 return isDlbOn(dd
->comm
);
144 gmx_bool
dd_dlb_is_locked(const gmx_domdec_t
*dd
)
146 return (dd
->comm
->dlbState
== DlbState::offTemporarilyLocked
);
149 void dd_dlb_lock(gmx_domdec_t
*dd
)
151 /* We can only lock the DLB when it is set to auto, otherwise don't do anything */
152 if (dd
->comm
->dlbState
== DlbState::offCanTurnOn
)
154 dd
->comm
->dlbState
= DlbState::offTemporarilyLocked
;
158 void dd_dlb_unlock(gmx_domdec_t
*dd
)
160 /* We can only lock the DLB when it is set to auto, otherwise don't do anything */
161 if (dd
->comm
->dlbState
== DlbState::offTemporarilyLocked
)
163 dd
->comm
->dlbState
= DlbState::offCanTurnOn
;
164 dd_dlb_set_should_check_whether_to_turn_dlb_on(dd
, TRUE
);