Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
tuning.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2022 The ESPResSo project
3 * Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010
4 * Max-Planck-Institute for Polymer Research, Theory Group
5 *
6 * This file is part of ESPResSo.
7 *
8 * ESPResSo is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * ESPResSo is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the GNU General Public License
19 * along with this program. If not, see <http://www.gnu.org/licenses/>.
20 */
21
22#include "tuning.hpp"
23
25#include "communication.hpp"
26#include "errorhandling.hpp"
27#include "integrate.hpp"
29#include "system/System.hpp"
30
32
33#include <boost/mpi/collectives/all_reduce.hpp>
34#include <boost/mpi/collectives/broadcast.hpp>
35
36#include <mpi.h>
37
38#include <algorithm>
39#include <functional>
40#include <string>
41
42std::string TuningFailed::get_first_error() const {
43 using namespace ErrorHandling;
45 auto message = std::string("tuning failed: an exception was thrown while "
46 "benchmarking the integration loop");
47 for (auto const &warning : queued_warnings) {
48 if (warning.level() == RuntimeError::ErrorLevel::ERROR) {
49 message += " (" + warning.what() + ")";
50 break;
51 }
52 }
53 return message;
54}
55
57 if (acc.avg() <= 5 * MPI_Wtick()) {
59 << "Clock resolution is too low to reliably time integration.";
60 }
61 if (acc.sig() >= 0.1 * acc.avg()) {
62 runtimeWarningMsg() << "Statistics of tuning samples is very bad.";
63 }
64}
65
67 auto const error_code = system.integrate(0, reuse_forces);
68 if (error_code == INTEG_ERROR_RUNTIME) {
69 throw TuningFailed{};
70 }
71}
72
75
76 // check if the system can be integrated with the current parameters
78
79 // measure force calculation time
80 for (int i = 0; i < int_steps; i++) {
81 auto const tick = MPI_Wtime();
83 auto const tock = MPI_Wtime();
84 running_average.add_sample((tock - tick));
85 }
86
87 if (this_node == 0) {
89 }
90
91 /* MPI returns in seconds, returned value should be in ms. */
92 auto retval = 1000. * running_average.avg();
93 boost::mpi::broadcast(::comm_cart, retval, 0);
94 return retval;
95}
96
97/**
98 * \brief Time the integration.
99 * This times the integration and
100 * propagates the system.
101 *
102 * @param system The system to tune.
103 * @param int_steps Number of steps to integrate.
104 * @return Time per integration in ms.
105 */
107 auto const error_code_init =
110 return -1;
111 }
112
113 /* perform force calculation test */
114 auto const tick = MPI_Wtime();
115 auto const error_code = system.integrate(int_steps, INTEG_REUSE_FORCES_NEVER);
116 auto const tock = MPI_Wtime();
117 if (error_code == INTEG_ERROR_RUNTIME) {
118 return -1;
119 }
120
121 /* MPI returns in seconds, returned value should be in ms. */
122 return 1000. * (tock - tick) / int_steps;
123}
124
125namespace System {
126void System::tune_verlet_skin(double min_skin, double max_skin, double tol,
127 int int_steps, bool adjust_max_skin) {
128
129 double a = min_skin;
130 double b = max_skin;
131
132 /* The maximal skin is the remainder from the required cutoff to
133 * the maximal range that can be supported by the cell system, but
134 * never larger than half the box size. */
135 auto const max_permissible_skin = std::min(
136 std::ranges::min(cell_structure->max_cutoff()) - maximal_cutoff(),
137 0.5 * std::ranges::max(box_geo->length()));
138
141
142 while (fabs(a - b) > tol) {
143 cell_structure->set_verlet_skin(a);
144 auto const time_a = time_calc(*this, int_steps);
145
146 cell_structure->set_verlet_skin(b);
147 auto const time_b = time_calc(*this, int_steps);
148
149 if (time_a > time_b) {
150 a = 0.5 * (a + b);
151 } else {
152 b = 0.5 * (a + b);
153 }
154 }
155 cell_structure->set_verlet_skin(0.5 * (a + b));
156}
157} // namespace System
Main system class.
void tune_verlet_skin(double min_skin, double max_skin, double tol, int int_steps, bool adjust_max_skin)
Tune the Verlet skin.
Definition tuning.cpp:126
Keep running average and variance.
Scalar avg() const
Average of the samples.
Scalar sig() const
Standard deviation of the samples.
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
boost::mpi::communicator comm_cart
The communicator.
int this_node
The number of this node.
This file contains the errorhandling code for severe errors, like a broken bond or illegal parameter ...
#define runtimeWarningMsg()
Molecular dynamics integrator.
#define INTEG_ERROR_RUNTIME
Definition integrate.hpp:42
#define INTEG_REUSE_FORCES_CONDITIONALLY
recalculate forces only if Propagation::recalc_forces is set
Definition integrate.hpp:51
#define INTEG_REUSE_FORCES_NEVER
recalculate forces unconditionally (mostly used for timing)
Definition integrate.hpp:49
std::vector< RuntimeError > mpi_gather_runtime_errors_all(bool is_head_node)
Gather messages on main rank.
Various procedures concerning interactions between particles.
static void run_full_force_calc(System::System &system, int reuse_forces)
Definition tuning.cpp:66
static void check_statistics(Utils::Statistics::RunningAverage< double > &acc)
Definition tuning.cpp:56
static double time_calc(System::System &system, int int_steps)
Time the integration.
Definition tuning.cpp:106
double benchmark_integration_step(System::System &system, int int_steps)
Benchmark the integration loop.
Definition tuning.cpp:73