ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
integrate.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2026 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/** \file
23 * Molecular dynamics integrator.
24 *
25 * For more information about the integrator
26 * see \ref integrate.hpp "integrate.hpp".
27 */
28
29#include "integrate.hpp"
37
38#include "BoxGeometry.hpp"
39#include "PropagationMode.hpp"
40#include "accumulators/AutoUpdateAccumulators.hpp"
44#include "cells.hpp"
45#include "collision_detection/CollisionDetection.hpp"
46#include "communication.hpp"
47#include "errorhandling.hpp"
49#include "lb/utils.hpp"
52#include "npt.hpp"
53#include "rattle.hpp"
54#include "rotation.hpp"
55#include "signalhandling.hpp"
57#include "system/System.hpp"
59#include "thermostat.hpp"
61#include "virtual_sites/com.hpp"
64
66
67#include <boost/mpi/collectives/all_reduce.hpp>
68
69#ifdef ESPRESSO_CALIPER
70#include <caliper/cali.h>
71#endif
72
73#ifdef ESPRESSO_VALGRIND
74#include <callgrind.h>
75#endif
76
77#include <algorithm>
78#include <cassert>
79#include <cmath>
80#include <csignal>
81#include <functional>
82#include <limits>
83#include <sstream>
84#include <stdexcept>
85#include <string>
86#include <utility>
87
88#ifdef ESPRESSO_WALBERLA
89#ifdef ESPRESSO_WALBERLA_STATIC_ASSERT
90#error "waLberla headers should not be visible to the ESPResSo core"
91#endif
92#endif
93
94namespace {
95volatile std::sig_atomic_t ctrl_C = 0;
96} // namespace
97
98namespace LeesEdwards {
99
100/**
101 * @brief Update the Lees-Edwards parameters of the box geometry
102 * for the current simulation time.
103 */
104void LeesEdwards::update_box_params(BoxGeometry &box_geo, double sim_time) {
105 if (box_geo.type() == BoxType::LEES_EDWARDS) {
106 assert(m_protocol != nullptr);
107 box_geo.lees_edwards_update(get_pos_offset(sim_time, *m_protocol),
108 get_shear_velocity(sim_time, *m_protocol));
109 }
110}
111
112void LeesEdwards::set_protocol(std::shared_ptr<ActiveProtocol> protocol) {
113 auto &system = get_system();
114 auto &cell_structure = *system.cell_structure;
115 auto &box_geo = *system.box_geo;
116 box_geo.set_type(BoxType::LEES_EDWARDS);
117 m_protocol = std::move(protocol);
118 update_box_params(box_geo, system.get_sim_time());
119 system.propagation->recalc_forces = true;
120 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
121}
122
124 auto &system = get_system();
125 auto &cell_structure = *system.cell_structure;
126 auto &box_geo = *system.box_geo;
127 m_protocol = nullptr;
128 box_geo.set_type(BoxType::CUBOID);
129 system.propagation->recalc_forces = true;
130 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
131}
132
133} // namespace LeesEdwards
134
136 switch (integ_switch) {
139 break;
140 case INTEG_METHOD_NVT:
142 // NOLINTNEXTLINE(bugprone-branch-clone)
143 if ((thermo_switch & THERMO_LB) and (thermo_switch & THERMO_LANGEVIN)) {
145#ifdef ESPRESSO_ROTATION
147#endif
148 } else if (thermo_switch & THERMO_LB) {
150#ifdef ESPRESSO_ROTATION
152#endif
153 } else if (thermo_switch & THERMO_LANGEVIN) {
155#ifdef ESPRESSO_ROTATION
157#endif
158 } else {
160#ifdef ESPRESSO_ROTATION
162#endif
163 }
164 break;
165 }
166#ifdef ESPRESSO_NPT
170 break;
171#endif
172 case INTEG_METHOD_BD:
174#ifdef ESPRESSO_ROTATION
176#endif
177 break;
178#ifdef ESPRESSO_STOKESIAN_DYNAMICS
179 case INTEG_METHOD_SD:
181 break;
182#endif // ESPRESSO_STOKESIAN_DYNAMICS
183 default:
184 throw std::runtime_error("Unknown value for integ_switch");
185 }
186}
187
189 int used_propagations = PropagationMode::NONE;
190 for (auto &p : cell_structure->local_particles()) {
191 used_propagations |= p.propagation();
192 }
193 if (used_propagations & PropagationMode::SYSTEM_DEFAULT) {
194 used_propagations |= propagation->default_propagation;
195 }
196 used_propagations = boost::mpi::all_reduce(::comm_cart, used_propagations,
197 std::bit_or<int>());
198 propagation->used_propagations = used_propagations;
199 propagation->recalc_used_propagations = false;
200}
201
202void System::System::integrator_sanity_checks() const {
203 auto const thermo_switch = thermostat->thermo_switch;
204 if (time_step <= 0.) {
205 runtimeErrorMsg() << "time_step not set";
206 }
207 if (propagation->integ_switch == INTEG_METHOD_STEEPEST_DESCENT) {
208 if (thermo_switch != THERMO_OFF) {
210 << "The steepest descent integrator is incompatible with thermostats";
211 }
212 }
213 if (propagation->integ_switch == INTEG_METHOD_NVT) {
214 if (thermo_switch & (THERMO_NPT_ISO | THERMO_BROWNIAN | THERMO_SD)) {
215 runtimeErrorMsg() << "The VV integrator is incompatible with the "
216 "currently active combination of thermostats";
217 }
218 }
219#ifdef ESPRESSO_NPT
220 if (propagation->used_propagations & PropagationMode::TRANS_LANGEVIN_NPT) {
221 if (thermo_switch != THERMO_NPT_ISO) {
222 runtimeErrorMsg() << "The NpT integrator requires the NpT thermostat";
223 }
224 if (box_geo->type() == BoxType::LEES_EDWARDS) {
225 runtimeErrorMsg() << "The NpT integrator cannot use Lees-Edwards";
226 }
227 try {
228 nptiso->coulomb_dipole_sanity_checks(*this);
229 } catch (std::runtime_error const &err) {
230 runtimeErrorMsg() << err.what();
231 }
232 }
233#endif
234 if (propagation->used_propagations & PropagationMode::TRANS_BROWNIAN) {
235 if (thermo_switch != THERMO_BROWNIAN) {
236 runtimeErrorMsg() << "The BD integrator requires the BD thermostat";
237 }
238 }
239 if (propagation->used_propagations & PropagationMode::TRANS_STOKESIAN) {
240#ifdef ESPRESSO_STOKESIAN_DYNAMICS
241 if (thermo_switch != THERMO_SD) {
242 runtimeErrorMsg() << "The SD integrator requires the SD thermostat";
243 }
244#endif
245 }
246 if (lb.is_solver_set() and (propagation->used_propagations &
249 if (thermostat->lb == nullptr) {
250 runtimeErrorMsg() << "The LB integrator requires the LB thermostat";
251 }
252 }
253 if (bonded_ias->get_n_thermalized_bonds() >= 1 and
254 (thermostat->thermalized_bond == nullptr or
255 (thermo_switch & THERMO_BOND) == 0)) {
257 << "Thermalized bonds require the thermalized_bond thermostat";
258 }
259
260#ifdef ESPRESSO_ROTATION
261 for (auto const &p : cell_structure->local_particles()) {
262 using namespace PropagationMode;
263 if (p.can_rotate() and not p.is_virtual() and
264 (p.propagation() & (SYSTEM_DEFAULT | ROT_EULER | ROT_LANGEVIN |
265 ROT_BROWNIAN | ROT_STOKESIAN)) == 0) {
267 << "Rotating particles must have a rotation propagation mode enabled";
268 break;
269 }
270 }
271#endif // ESPRESSO_ROTATION
272
273#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
274#ifdef ESPRESSO_EXTERNAL_FORCES
275 if (propagation->used_propagations &
277 for (auto const &p : cell_structure->local_particles()) {
278 using namespace PropagationMode;
279 if ((p.propagation() & TRANS_VS_CENTER_OF_MASS) and
280 p.has_fixed_coordinates()) {
281 runtimeErrorMsg() << "VS COM particles cannot be fixed in space";
282 break;
283 }
284 }
285 }
286#endif // ESPRESSO_EXTERNAL_FORCES
287#ifdef ESPRESSO_BOND_CONSTRAINT
288 if (bonded_ias->get_n_rigid_bonds()) {
289 using namespace PropagationMode;
290 for (auto const &p : cell_structure->local_particles()) {
291 if (p.propagation() & TRANS_VS_CENTER_OF_MASS) {
292 for (auto const bond : p.bonds()) {
293 if (std::holds_alternative<RigidBond>(
294 *bonded_ias->at(bond.bond_id()))) {
295 runtimeErrorMsg() << "VS COM particles cannot use rigid bonds";
296 break;
297 }
298 }
299 }
300 }
301 }
302#endif // ESPRESSO_BOND_CONSTRAINT
303#endif // ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
304
305#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
306 if ((thermo_switch & THERMO_LANGEVIN) == 0) {
307 for (auto const &p : cell_structure->local_particles()) {
308 if (p.stoner_wohlfarth_is_enabled()) {
309 runtimeErrorMsg() << "The thermal Stoner-Wohlfarth model requires the "
310 "Langevin thermostat";
311 break;
312 }
313 }
314 }
315#endif // ESPRESSO_THERMAL_STONER_WOHLFARTH
316}
317
318#ifdef ESPRESSO_WALBERLA
319void walberla_tau_sanity_checks(std::string method, double tau,
320 double time_step) {
321 if (time_step <= 0.) {
322 return;
323 }
324 // use float epsilon since tau may be a float
325 auto const eps = static_cast<double>(std::numeric_limits<float>::epsilon());
326 if ((tau - time_step) / (tau + time_step) < -eps)
327 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
328 ") must be >= MD time_step (" +
329 std::to_string(time_step) + ")");
330 auto const factor = tau / time_step;
331 if (std::fabs(std::round(factor) - factor) / factor > eps)
332 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
333 ") must be an integer multiple of the "
334 "MD time_step (" +
335 std::to_string(time_step) + "). Factor is " +
336 std::to_string(factor));
337}
338
339void walberla_agrid_sanity_checks(std::string method,
340 Utils::Vector3d const &geo_left,
341 Utils::Vector3d const &geo_right,
342 Utils::Vector3d const &lattice_left,
343 Utils::Vector3d const &lattice_right,
344 double agrid) {
345 // waLBerla and ESPResSo must agree on domain decomposition
346 auto const tol = agrid / 1E6;
347 if ((lattice_left - geo_left).norm2() > tol or
348 (lattice_right - geo_right).norm2() > tol) {
349 std::stringstream error_msg;
350 error_msg << "waLBerla and ESPResSo disagree about domain decomposition"
351 << "\nMPI rank " << ::this_node << ": "
352 << "left ESPResSo: [" << geo_left << "], "
353 << "left waLBerla: [" << lattice_left << "]"
354 << "\nMPI rank " << ::this_node << ": "
355 << "right ESPResSo: [" << geo_right << "], "
356 << "right waLBerla: [" << lattice_right << "]"
357 << "\nfor method: " << method;
358 throw std::runtime_error(error_msg.str());
359 }
360}
361#endif // ESPRESSO_WALBERLA
362
364#ifdef ESPRESSO_CALIPER
365 CALI_CXX_MARK_FUNCTION;
366#endif
367 auto &cell_structure = *system.cell_structure;
368 auto const offset = LeesEdwards::verlet_list_offset(
369 *system.box_geo, cell_structure.get_le_pos_offset_at_last_resort());
370 if (cell_structure.check_resort_required(offset)) {
371 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
372 }
373}
374
375/** @brief Calls the hook for propagation kernels before the force calculation
376 * @return whether or not to stop the integration loop early.
377 */
378static bool integrator_step_1(CellStructure &cell_structure,
379 Propagation const &propagation,
380 System::System &system, double time_step) {
381#ifdef ESPRESSO_CALIPER
382 CALI_CXX_MARK_FUNCTION;
383#endif
384 // steepest decent
386 return system.steepest_descent->propagate(cell_structure);
387
388 auto const &thermostat = *system.thermostat;
389 auto const kT = thermostat.kT;
390 cell_structure.for_each_local_particle([&](Particle &p) {
391#ifdef ESPRESSO_VIRTUAL_SITES
392 // virtual sites are updated later in the integration loop
393 if (p.is_virtual())
394 return;
395#endif
396 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
397 if (propagation.should_propagate_with(
399 symplectic_euler_propagator_1(p, time_step);
401 symplectic_euler_propagator_1(p, time_step);
402#ifdef ESPRESSO_ROTATION
404 symplectic_euler_rotator_1(p, time_step);
405#endif
407 symplectic_euler_propagator_1(p, time_step);
408#ifdef ESPRESSO_ROTATION
410 symplectic_euler_rotator_1(p, time_step);
411#endif
412 } else {
413 if (propagation.should_propagate_with(
415 velocity_verlet_propagator_1(p, time_step);
417 velocity_verlet_propagator_1(p, time_step);
418#ifdef ESPRESSO_ROTATION
420 velocity_verlet_rotator_1(p, time_step);
421#endif
423 velocity_verlet_propagator_1(p, time_step);
424#ifdef ESPRESSO_ROTATION
426 velocity_verlet_rotator_1(p, time_step);
427#endif
428 }
430 brownian_dynamics_propagator(*thermostat.brownian, p, time_step, kT);
431#ifdef ESPRESSO_ROTATION
433 brownian_dynamics_rotator(*thermostat.brownian, p, time_step, kT);
434#endif
435 });
436
437#ifdef ESPRESSO_NPT
440 auto pred = PropagationPredicateNPT(propagation.default_propagation);
441 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
443 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
444 time_step, system);
445 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
447 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
448 time_step, system);
449 }
450 }
451#endif
452
453#ifdef ESPRESSO_STOKESIAN_DYNAMICS
456 auto pred = PropagationPredicateStokesian(propagation.default_propagation);
457 stokesian_dynamics_step_1(cell_structure.local_particles().filter(pred),
458 *system.stokesian_dynamics, *thermostat.stokesian,
459 time_step, kT);
460 }
461#endif // ESPRESSO_STOKESIAN_DYNAMICS
462
463 return false;
464}
465
466static void integrator_step_2(CellStructure &cell_structure,
467 Propagation const &propagation,
468 [[maybe_unused]] System::System &system,
469 double time_step) {
470#ifdef ESPRESSO_CALIPER
471 CALI_CXX_MARK_FUNCTION;
472#endif
474 return;
475
476 cell_structure.for_each_local_particle([&](Particle &p) {
477#ifdef ESPRESSO_VIRTUAL_SITES
478 // virtual sites are updated later in the integration loop
479 if (p.is_virtual())
480 return;
481#endif
482 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
483 if (propagation.should_propagate_with(
485 symplectic_euler_propagator_2(p, time_step);
487 symplectic_euler_propagator_2(p, time_step);
488#ifdef ESPRESSO_ROTATION
490 symplectic_euler_rotator_2(p, time_step);
491#endif
493 symplectic_euler_propagator_2(p, time_step);
494#ifdef ESPRESSO_ROTATION
496 symplectic_euler_rotator_2(p, time_step);
497#endif
498 } else {
499 if (propagation.should_propagate_with(
501 velocity_verlet_propagator_2(p, time_step);
503 velocity_verlet_propagator_2(p, time_step);
504#ifdef ESPRESSO_ROTATION
506 velocity_verlet_rotator_2(p, time_step);
507#endif
509 velocity_verlet_propagator_2(p, time_step);
510#ifdef ESPRESSO_ROTATION
512 velocity_verlet_rotator_2(p, time_step);
513#endif
514 }
515 });
516
517#ifdef ESPRESSO_NPT
520 auto pred = PropagationPredicateNPT(propagation.default_propagation);
521 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
523 cell_structure.local_particles().filter(pred), time_step, system);
524 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
526 cell_structure.local_particles().filter(pred), time_step, system);
527 }
528 }
529#endif
530}
531
532int System::System::integrate(int n_steps, int reuse_forces) {
533#ifdef ESPRESSO_CALIPER
534 CALI_CXX_MARK_FUNCTION;
535#endif
536 auto &propagation = *this->propagation;
537#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
538 auto const has_vs_rel = [&propagation]() {
539 return propagation.used_propagations &
543 };
544#endif
545#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
546 auto const has_vs_com = [&propagation]() {
547 return propagation.used_propagations &
549 };
550#endif
551#ifdef ESPRESSO_BOND_CONSTRAINT
552 auto const n_rigid_bonds = bonded_ias->get_n_rigid_bonds();
553#endif
554
555 // Prepare particle structure and run sanity checks of all active algorithms
556 propagation.update_default_propagation(thermostat->thermo_switch);
557 update_used_propagations();
558 on_integration_start();
559
560 // If any method vetoes (e.g. P3M not initialized), immediately bail out
562 return INTEG_ERROR_RUNTIME;
563
564 // Additional preparations for the first integration step
565 if (reuse_forces == INTEG_REUSE_FORCES_NEVER or
566 ((reuse_forces != INTEG_REUSE_FORCES_ALWAYS) and
567 propagation.recalc_forces)) {
568#ifdef ESPRESSO_CALIPER
569 CALI_MARK_BEGIN("Initial Force Calculation");
570#endif
571 thermostat->lb_coupling_deactivate();
572
573#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
574 if (has_vs_rel()) {
575 vs_relative_update_particles(*cell_structure, *box_geo);
576 }
577#endif
578#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
579 if (has_vs_com()) {
580 vs_com_update_particles(*cell_structure, *box_geo);
581 }
582#endif
583
584 // Communication step: distribute ghost positions
585 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
586
587 calculate_forces();
588
589 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
590#ifdef ESPRESSO_ROTATION
591 convert_initial_torques(cell_structure->local_particles());
592#endif
593 }
594
595#ifdef ESPRESSO_CALIPER
596 CALI_MARK_END("Initial Force Calculation");
597#endif
598 }
599
600 thermostat->lb_coupling_activate();
601
603 return INTEG_ERROR_RUNTIME;
604
605 // Keep track of the number of Verlet updates (i.e. particle resorts)
606 int n_verlet_updates = 0;
607
608 // Keep track of whether an interrupt signal was caught (only in singleton
609 // mode, since signal handlers are unreliable with more than 1 MPI rank)
610 auto const singleton_mode = comm_cart.size() == 1;
611 auto caught_sigint = false;
612 auto caught_error = false;
613
614 auto lb_active = false;
615 auto ek_active = false;
616 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
617 lb_active = lb.is_solver_set();
618 ek_active = ek.is_ready_for_propagation();
619 }
620 auto const calc_md_steps_per_tau = [this](double tau) {
621 return static_cast<int>(std::round(tau / time_step));
622 };
623
624#ifdef ESPRESSO_VALGRIND
625 CALLGRIND_START_INSTRUMENTATION;
626#endif
627 // Integration loop
628#ifdef ESPRESSO_CALIPER
629 CALI_CXX_MARK_LOOP_BEGIN(integration_loop, "Integration loop");
630#endif
631 int integrated_steps = 0;
632 for (int step = 0; step < n_steps; step++) {
633#ifdef ESPRESSO_CALIPER
634 CALI_CXX_MARK_LOOP_ITERATION(integration_loop, step);
635#endif
636
637#ifdef ESPRESSO_BOND_CONSTRAINT
638 if (n_rigid_bonds)
639 save_old_position(cell_structure->local_particles(),
640 cell_structure->ghost_particles());
641#endif
642
643 lees_edwards->update_box_params(*box_geo, sim_time);
644 bool early_exit =
645 integrator_step_1(*cell_structure, propagation, *this, time_step);
646 if (early_exit)
647 break;
648
649 sim_time += time_step;
650 if (box_geo->type() == BoxType::LEES_EDWARDS) {
651 auto const kernel = LeesEdwards::Push{*box_geo};
652 cell_structure->for_each_local_particle(
653 [&kernel](Particle &p) { kernel(p); });
654 }
655
656#ifdef ESPRESSO_NPT
657 if (not has_npt_enabled())
658#endif
659 {
661 }
662 // Propagate philox RNG counters
663 thermostat->philox_counter_increment();
664
665#ifdef ESPRESSO_BOND_CONSTRAINT
666 // Correct particle positions that participate in a rigid/constrained bond
667 if (n_rigid_bonds) {
668 correct_position_shake(*cell_structure, *box_geo, *bonded_ias);
669 }
670#endif
671
672#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
673 if (has_vs_rel()) {
674#ifdef ESPRESSO_NPT
675 if (has_npt_enabled()) {
676 cell_structure->update_ghosts_and_resort_particle(
678 }
679#endif // ESPRESSO_NPT
680 vs_relative_update_particles(*cell_structure, *box_geo);
681 }
682#endif // ESPRESSO_VIRTUAL_SITES_RELATIVE
683#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
684 if (has_vs_com()) {
685#ifdef ESPRESSO_NPT
686 if (has_npt_enabled()) {
687 cell_structure->update_ghosts_and_resort_particle(
689 }
690#endif // ESPRESSO_NPT
691 vs_com_update_particles(*cell_structure, *box_geo);
692 }
693#endif // ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
694
695 if (cell_structure->get_resort_particles() >= Cells::RESORT_LOCAL)
696 n_verlet_updates++;
697
698 // Communication step: distribute ghost positions
699 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
700
701#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
702 integrate_magnetodynamics();
703#endif
704
705 calculate_forces();
706
707#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
708 if (thermostat->lb and
709 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
710 lb_tracers_add_particle_force_to_fluid(*cell_structure, *box_geo,
711 *local_geo, lb);
712 }
713#endif
714 integrator_step_2(*cell_structure, propagation, *this, time_step);
715 if (propagation.integ_switch == INTEG_METHOD_BD) {
717 }
718 if (box_geo->type() == BoxType::LEES_EDWARDS) {
719 auto const kernel = LeesEdwards::UpdateOffset{*box_geo};
720 cell_structure->for_each_local_particle(
721 [&kernel](Particle &p) { kernel(p); });
722 }
723#ifdef ESPRESSO_BOND_CONSTRAINT
724 if (n_rigid_bonds) {
725 correct_velocity_shake(*cell_structure, *box_geo, *bonded_ias);
726 }
727#endif
728
729 // propagate one-step functionalities
730 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
731 if (lb_active and ek_active) {
732 // assume that they are coupled, which is not necessarily true
733 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
734 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
735
736 if (md_steps_per_lb_step != md_steps_per_ek_step) {
738 << "LB and EK are active but with different time steps.";
739 }
740
741 assert(lb.is_gpu() == ek.is_gpu());
742 assert(propagation.lb_skipped_md_steps ==
743 propagation.ek_skipped_md_steps);
744
745 propagation.lb_skipped_md_steps += 1;
746 propagation.ek_skipped_md_steps += 1;
747 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
748 propagation.lb_skipped_md_steps = 0;
749 propagation.ek_skipped_md_steps = 0;
750#ifdef ESPRESSO_CALIPER
751 CALI_MARK_BEGIN("lb_propagation");
752#endif
753 lb.propagate();
754 lb.ghost_communication_vel();
755#ifdef ESPRESSO_CALIPER
756 CALI_MARK_END("lb_propagation");
757#endif
758#ifdef ESPRESSO_CALIPER
759 CALI_MARK_BEGIN("ek_propagation");
760#endif
761 ek.propagate();
762#ifdef ESPRESSO_CALIPER
763 CALI_MARK_END("ek_propagation");
764#endif
765 }
766 } else if (lb_active) {
767 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
768 propagation.lb_skipped_md_steps += 1;
769 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
770 propagation.lb_skipped_md_steps = 0;
771#ifdef ESPRESSO_CALIPER
772 CALI_MARK_BEGIN("lb_propagation");
773#endif
774 lb.propagate();
775#ifdef ESPRESSO_CALIPER
776 CALI_MARK_END("lb_propagation");
777#endif
778 }
779 } else if (ek_active) {
780 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
781 propagation.ek_skipped_md_steps += 1;
782 if (propagation.ek_skipped_md_steps >= md_steps_per_ek_step) {
783 propagation.ek_skipped_md_steps = 0;
784#ifdef ESPRESSO_CALIPER
785 CALI_MARK_BEGIN("ek_propagation");
786#endif
787 ek.propagate();
788#ifdef ESPRESSO_CALIPER
789 CALI_MARK_END("ek_propagation");
790#endif
791 }
792 }
793 if (lb_active and (propagation.used_propagations &
795 thermostat->lb->rng_increment();
796 }
797
798#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
799 if (thermostat->lb and
800 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
801#ifdef ESPRESSO_CALIPER
802 CALI_MARK_BEGIN("lb_tracers_propagation");
803#endif
804 if (lb_active) {
805 lb.ghost_communication_vel();
806 }
807 lb_tracers_propagate(*cell_structure, lb, time_step);
808#ifdef ESPRESSO_CALIPER
809 CALI_MARK_END("lb_tracers_propagation");
810#endif
811 }
812#endif
813
814#ifdef ESPRESSO_COLLISION_DETECTION
815 cell_structure->clear_new_bonds();
816 collision_detection->handle_collisions();
817 cell_structure->rebuild_bond_list();
818#endif
819 bond_breakage->process_queue(*this);
820 }
821
822 integrated_steps++;
823
825 caught_error = true;
826 break;
827 }
828
829 // Check if SIGINT has been caught.
830 if (singleton_mode and ctrl_C == 1) {
831 caught_sigint = true;
832 break;
833 }
834
835 } // for-loop over integration steps
836 if (lb_active) {
837 lb.ghost_communication();
838 }
839 lees_edwards->update_box_params(*box_geo, sim_time);
840#ifdef ESPRESSO_CALIPER
841 CALI_CXX_MARK_LOOP_END(integration_loop);
842#endif
843
844#ifdef ESPRESSO_VALGRIND
845 CALLGRIND_STOP_INSTRUMENTATION;
846#endif
847
848#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
849 if (has_vs_rel()) {
850 vs_relative_update_particles(*cell_structure, *box_geo);
851 }
852#endif
853#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
854 if (has_vs_com()) {
855 vs_com_update_particles(*cell_structure, *box_geo);
856 }
857#endif
858
859 // Verlet list statistics
860 cell_structure->update_verlet_stats(n_steps, n_verlet_updates);
861
862#ifdef ESPRESSO_NPT
863 if (has_npt_enabled()) {
864 synchronize_npt_state();
865 }
866#endif
867 if (caught_sigint) {
868 ctrl_C = 0;
869 return INTEG_ERROR_SIGINT;
870 }
871 if (caught_error) {
872 return INTEG_ERROR_RUNTIME;
873 }
874 if (boost::mpi::all_reduce(::comm_cart, not cell_structure->use_verlet_list,
875 std::logical_or<>())) {
876 cell_structure->use_verlet_list = false;
877 }
878 return integrated_steps;
879}
880
881int System::System::integrate_with_signal_handler(int n_steps, int reuse_forces,
882 bool update_accumulators) {
883 assert(n_steps >= 0);
884
885 // Override the signal handler so that the integrator obeys Ctrl+C
886 SignalHandler sa(SIGINT, [](int) { ctrl_C = 1; });
887
888 /* if skin wasn't set, do an educated guess now */
889 if (not cell_structure->is_verlet_skin_set()) {
890 try {
891 cell_structure->set_verlet_skin_heuristic();
892 } catch (...) {
893 if (comm_cart.rank() == 0) {
894 throw;
895 }
896 return INTEG_ERROR_RUNTIME;
897 }
898 }
899
900 if (not update_accumulators or n_steps == 0) {
901 return integrate(n_steps, reuse_forces);
902 }
903
904 for (int i = 0; i < n_steps;) {
905 /* Integrate to either the next accumulator update, or the
906 * end, depending on what comes first. */
907 auto const steps =
908 std::min((n_steps - i), auto_update_accumulators->next_update());
909
910 auto const local_retval = integrate(steps, reuse_forces);
911
912 // make sure all ranks exit when one rank fails
913 std::remove_const_t<decltype(local_retval)> global_retval;
914 boost::mpi::all_reduce(comm_cart, local_retval, global_retval,
915 std::plus<int>());
916 if (global_retval < 0) {
917 return global_retval; // propagate error code
918 }
919
920 reuse_forces = INTEG_REUSE_FORCES_ALWAYS;
921
922 (*auto_update_accumulators)(comm_cart, steps);
923
924 i += steps;
925 }
926
927 return 0;
928}
929
931 sim_time = value;
932 propagation->recalc_forces = true;
933 lees_edwards->update_box_params(*box_geo, sim_time);
934}
@ LEES_EDWARDS
@ INTEG_METHOD_NPT_ISO_AND
@ INTEG_METHOD_SD
@ INTEG_METHOD_STEEPEST_DESCENT
@ INTEG_METHOD_NVT
@ INTEG_METHOD_SYMPLECTIC_EULER
@ INTEG_METHOD_BD
@ INTEG_METHOD_NPT_ISO_MTK
@ THERMO_SD
@ THERMO_BROWNIAN
@ THERMO_BOND
@ THERMO_LB
@ THERMO_LANGEVIN
@ THERMO_NPT_ISO
@ THERMO_OFF
Data structures for bonded interactions.
This file contains everything related to the global cell structure / cell system.
void lees_edwards_update(double pos_offset, double shear_velocity)
Update the Lees-Edwards parameters of the box geometry for the current simulation time.
BoxType type() const
Describes a cell structure / cell system.
void for_each_local_particle(ParticleUnaryOp &&f, bool parallel=true) const
Run a kernel on all local particles.
ParticleRange local_particles() const
void update_box_params(BoxGeometry &box_geo, double sim_time)
Update the Lees-Edwards parameters of the box geometry for the current simulation time.
void set_protocol(std::shared_ptr< ActiveProtocol > protocol)
Set a new Lees-Edwards protocol.
void unset_protocol()
Delete the currently active Lees-Edwards protocol.
ParticleRangeFiltered< Predicate > filter(Predicate pred) const
int default_propagation
void update_default_propagation(int thermo_switch)
bool should_propagate_with(Particle const &p, int mode) const
int used_propagations
RAII guard for signal handling.
Main system class.
std::shared_ptr< StokesianDynamics > stokesian_dynamics
void update_used_propagations()
Update the global propagation bitmask.
void set_sim_time(double value)
Set sim_time.
int integrate_with_signal_handler(int n_steps, int reuse_forces, bool update_accumulators)
std::shared_ptr< SteepestDescent > steepest_descent
int integrate(int n_steps, int reuse_forces)
Integrate equations of motion.
std::shared_ptr< Thermostat::Thermostat > thermostat
std::shared_ptr< CellStructure > cell_structure
std::shared_ptr< BoxGeometry > box_geo
void vs_com_update_particles(CellStructure &cell_structure, BoxGeometry const &box_geo)
Definition com.cpp:130
boost::mpi::communicator comm_cart
The communicator.
int this_node
The number of this node.
int check_runtime_errors(boost::mpi::communicator const &comm)
Count runtime errors on all nodes.
This file contains the errorhandling code for severe errors, like a broken bond or illegal parameter ...
#define runtimeErrorMsg()
static bool integrator_step_1(CellStructure &cell_structure, Propagation const &propagation, System::System &system, double time_step)
Calls the hook for propagation kernels before the force calculation.
static void resort_particles_if_needed(System::System &system)
static void integrator_step_2(CellStructure &cell_structure, Propagation const &propagation, System::System &system, double time_step)
void walberla_agrid_sanity_checks(std::string method, Utils::Vector3d const &geo_left, Utils::Vector3d const &geo_right, Utils::Vector3d const &lattice_left, Utils::Vector3d const &lattice_right, double agrid)
void walberla_tau_sanity_checks(std::string method, double tau, double time_step)
Molecular dynamics integrator.
#define INTEG_ERROR_RUNTIME
Definition integrate.hpp:42
#define INTEG_ERROR_SIGINT
Definition integrate.hpp:43
#define INTEG_REUSE_FORCES_NEVER
recalculate forces unconditionally (mostly used for timing)
Definition integrate.hpp:49
#define INTEG_REUSE_FORCES_ALWAYS
do not recalculate forces (mostly when reading checkpoints with forces)
Definition integrate.hpp:53
void brownian_dynamics_rotator(BrownianThermostat const &brownian, Particle &p, double time_step, double kT)
void brownian_dynamics_propagator(BrownianThermostat const &brownian, Particle &p, double time_step, double kT)
void lb_tracers_propagate(CellStructure &cell_structure, LB::Solver const &lb, double time_step)
void lb_tracers_add_particle_force_to_fluid(CellStructure &cell_structure, BoxGeometry const &box_geo, LocalBox const &local_box, LB::Solver &lb)
@ DATA_PART_PROPERTIES
Particle::p.
Utils::Vector3d verlet_list_offset(BoxGeometry const &box, double pos_offset_at_last_resort)
double get_shear_velocity(double time, ActiveProtocol const &protocol)
Calculation of current velocity.
double get_pos_offset(double time, ActiveProtocol const &protocol)
Definition protocols.hpp:94
volatile std::sig_atomic_t ctrl_C
Definition integrate.cpp:95
Various procedures concerning interactions between particles.
Exports for the NpT code.
void correct_velocity_shake(CellStructure &cs, BoxGeometry const &box_geo, BondedInteractionsMap const &bonded_ias)
Correction of current velocities using RATTLE algorithm.
Definition rattle.cpp:229
void save_old_position(const ParticleRange &particles, const ParticleRange &ghost_particles)
copy current position
Definition rattle.cpp:59
void correct_position_shake(CellStructure &cs, BoxGeometry const &box_geo, BondedInteractionsMap const &bonded_ias)
Propagate velocity and position while using SHAKE algorithm for bond constraint.
Definition rattle.cpp:157
RATTLE algorithm ().
void vs_relative_update_particles(CellStructure &cell_structure, BoxGeometry const &box_geo)
Definition relative.cpp:120
void convert_initial_torques(const ParticleRange &particles)
Convert torques to the body-fixed frame before the integration loop.
Definition rotation.cpp:192
This file contains all subroutines required to process rotational motion.
See for the Stokesian dynamics method used here.
void stokesian_dynamics_step_1(ParticleRangeStokesian const &particles, StokesianDynamics const &integrator, StokesianThermostat const &stokesian, double time_step, double kT)
Struct holding all information for one particle.
Definition Particle.hpp:435
auto is_virtual() const
Definition Particle.hpp:588
void symplectic_euler_rotator_2(Particle &, double)
void symplectic_euler_rotator_1(Particle &p, double time_step)
void symplectic_euler_propagator_2(Particle &, double)
Final integration step of the Symplectic Euler integrator For symplectic Euler, there is no second st...
void symplectic_euler_propagator_1(Particle &p, double time_step)
Propagate the velocities and positions.
void velocity_verlet_rotator_1(Particle &p, double time_step)
void velocity_verlet_propagator_2(Particle &p, double time_step)
Final integration step of the Velocity Verlet integrator.
void velocity_verlet_propagator_1(Particle &p, double time_step)
Propagate the velocities and positions.
void velocity_verlet_rotator_2(Particle &p, double time_step)
void velocity_verlet_npt_MTK_step_1(ParticleRangeNPT const &particles, IsotropicNptThermostat const &npt_iso, double time_step, System::System &system)
Special propagator for velocity Verlet NpT with the Andersen method.
void velocity_verlet_npt_MTK_step_2(ParticleRangeNPT const &particles, double time_step, System::System &system)
Final integration step of the velocity Verlet NpT integrator with the MTK method.
void velocity_verlet_npt_Andersen_step_1(ParticleRangeNPT const &particles, IsotropicNptThermostat const &npt_iso, double time_step, System::System &system)
Special propagator for velocity Verlet NpT with the Andersen method.
void velocity_verlet_npt_Andersen_step_2(ParticleRangeNPT const &particles, double time_step, System::System &system)
Final integration step of the velocity Verlet NpT integrator with the Andersen method.