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
65#include <boost/mpi/collectives/all_reduce.hpp>
66
67#ifdef ESPRESSO_CALIPER
68#include <caliper/cali.h>
69#endif
70
71#ifdef ESPRESSO_VALGRIND
72#include <callgrind.h>
73#endif
74
75#include <algorithm>
76#include <cassert>
77#include <cmath>
78#include <csignal>
79#include <functional>
80#include <limits>
81#include <sstream>
82#include <stdexcept>
83#include <string>
84#include <utility>
85
86#ifdef ESPRESSO_WALBERLA
87#ifdef ESPRESSO_WALBERLA_STATIC_ASSERT
88#error "waLberla headers should not be visible to the ESPResSo core"
89#endif
90#endif
91
92namespace {
93volatile std::sig_atomic_t ctrl_C = 0;
94} // namespace
95
96namespace LeesEdwards {
97
98/**
99 * @brief Update the Lees-Edwards parameters of the box geometry
100 * for the current simulation time.
101 */
102void LeesEdwards::update_box_params(BoxGeometry &box_geo, double sim_time) {
103 if (box_geo.type() == BoxType::LEES_EDWARDS) {
104 assert(m_protocol != nullptr);
105 box_geo.lees_edwards_update(get_pos_offset(sim_time, *m_protocol),
106 get_shear_velocity(sim_time, *m_protocol));
107 }
108}
109
110void LeesEdwards::set_protocol(std::shared_ptr<ActiveProtocol> protocol) {
111 auto &system = get_system();
112 auto &cell_structure = *system.cell_structure;
113 auto &box_geo = *system.box_geo;
114 box_geo.set_type(BoxType::LEES_EDWARDS);
115 m_protocol = std::move(protocol);
116 update_box_params(box_geo, system.get_sim_time());
117 system.propagation->recalc_forces = true;
118 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
119}
120
122 auto &system = get_system();
123 auto &cell_structure = *system.cell_structure;
124 auto &box_geo = *system.box_geo;
125 m_protocol = nullptr;
126 box_geo.set_type(BoxType::CUBOID);
127 system.propagation->recalc_forces = true;
128 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
129}
130
131} // namespace LeesEdwards
132
134 switch (integ_switch) {
137 break;
138 case INTEG_METHOD_NVT:
140 // NOLINTNEXTLINE(bugprone-branch-clone)
141 if ((thermo_switch & THERMO_LB) and (thermo_switch & THERMO_LANGEVIN)) {
143#ifdef ESPRESSO_ROTATION
145#endif
146 } else if (thermo_switch & THERMO_LB) {
148#ifdef ESPRESSO_ROTATION
150#endif
151 } else if (thermo_switch & THERMO_LANGEVIN) {
153#ifdef ESPRESSO_ROTATION
155#endif
156 } else {
158#ifdef ESPRESSO_ROTATION
160#endif
161 }
162 break;
163 }
164#ifdef ESPRESSO_NPT
168 break;
169#endif
170 case INTEG_METHOD_BD:
172#ifdef ESPRESSO_ROTATION
174#endif
175 break;
176#ifdef ESPRESSO_STOKESIAN_DYNAMICS
177 case INTEG_METHOD_SD:
179 break;
180#endif // ESPRESSO_STOKESIAN_DYNAMICS
181 default:
182 throw std::runtime_error("Unknown value for integ_switch");
183 }
184}
185
187 int used_propagations = PropagationMode::NONE;
188 for (auto &p : cell_structure->local_particles()) {
189 used_propagations |= p.propagation();
190 }
191 if (used_propagations & PropagationMode::SYSTEM_DEFAULT) {
192 used_propagations |= propagation->default_propagation;
193 }
194 used_propagations = boost::mpi::all_reduce(::comm_cart, used_propagations,
195 std::bit_or<int>());
196 propagation->used_propagations = used_propagations;
197 propagation->recalc_used_propagations = false;
198}
199
200void System::System::integrator_sanity_checks() const {
201 auto const thermo_switch = thermostat->thermo_switch;
202 if (time_step <= 0.) {
203 runtimeErrorMsg() << "time_step not set";
204 }
205 if (propagation->integ_switch == INTEG_METHOD_STEEPEST_DESCENT) {
206 if (thermo_switch != THERMO_OFF) {
208 << "The steepest descent integrator is incompatible with thermostats";
209 }
210 }
211 if (propagation->integ_switch == INTEG_METHOD_NVT) {
212 if (thermo_switch & (THERMO_NPT_ISO | THERMO_BROWNIAN | THERMO_SD)) {
213 runtimeErrorMsg() << "The VV integrator is incompatible with the "
214 "currently active combination of thermostats";
215 }
216 }
217#ifdef ESPRESSO_NPT
218 if (propagation->used_propagations & PropagationMode::TRANS_LANGEVIN_NPT) {
219 if (thermo_switch != THERMO_NPT_ISO) {
220 runtimeErrorMsg() << "The NpT integrator requires the NpT thermostat";
221 }
222 if (box_geo->type() == BoxType::LEES_EDWARDS) {
223 runtimeErrorMsg() << "The NpT integrator cannot use Lees-Edwards";
224 }
225 try {
226 nptiso->coulomb_dipole_sanity_checks(*this);
227 } catch (std::runtime_error const &err) {
228 runtimeErrorMsg() << err.what();
229 }
230 }
231#endif
232 if (propagation->used_propagations & PropagationMode::TRANS_BROWNIAN) {
233 if (thermo_switch != THERMO_BROWNIAN) {
234 runtimeErrorMsg() << "The BD integrator requires the BD thermostat";
235 }
236 }
237 if (propagation->used_propagations & PropagationMode::TRANS_STOKESIAN) {
238#ifdef ESPRESSO_STOKESIAN_DYNAMICS
239 if (thermo_switch != THERMO_SD) {
240 runtimeErrorMsg() << "The SD integrator requires the SD thermostat";
241 }
242#endif
243 }
244 if (lb.is_solver_set() and (propagation->used_propagations &
247 if (thermostat->lb == nullptr) {
248 runtimeErrorMsg() << "The LB integrator requires the LB thermostat";
249 }
250 }
251 if (bonded_ias->get_n_thermalized_bonds() >= 1 and
252 (thermostat->thermalized_bond == nullptr or
253 (thermo_switch & THERMO_BOND) == 0)) {
255 << "Thermalized bonds require the thermalized_bond thermostat";
256 }
257
258#ifdef ESPRESSO_ROTATION
259 for (auto const &p : cell_structure->local_particles()) {
260 using namespace PropagationMode;
261 if (p.can_rotate() and not p.is_virtual() and
262 (p.propagation() & (SYSTEM_DEFAULT | ROT_EULER | ROT_LANGEVIN |
263 ROT_BROWNIAN | ROT_STOKESIAN)) == 0) {
265 << "Rotating particles must have a rotation propagation mode enabled";
266 break;
267 }
268 }
269#endif // ESPRESSO_ROTATION
270
271#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
272#ifdef ESPRESSO_EXTERNAL_FORCES
273 if (propagation->used_propagations &
275 for (auto const &p : cell_structure->local_particles()) {
276 using namespace PropagationMode;
277 if ((p.propagation() & TRANS_VS_CENTER_OF_MASS) and
278 p.has_fixed_coordinates()) {
279 runtimeErrorMsg() << "VS COM particles cannot be fixed in space";
280 break;
281 }
282 }
283 }
284#endif // ESPRESSO_EXTERNAL_FORCES
285#ifdef ESPRESSO_BOND_CONSTRAINT
286 if (bonded_ias->get_n_rigid_bonds()) {
287 using namespace PropagationMode;
288 for (auto const &p : cell_structure->local_particles()) {
289 if (p.propagation() & TRANS_VS_CENTER_OF_MASS) {
290 for (auto const bond : p.bonds()) {
291 if (std::holds_alternative<RigidBond>(
292 *bonded_ias->at(bond.bond_id()))) {
293 runtimeErrorMsg() << "VS COM particles cannot use rigid bonds";
294 break;
295 }
296 }
297 }
298 }
299 }
300#endif // ESPRESSO_BOND_CONSTRAINT
301#endif // ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
302
303#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
304 if ((thermo_switch & THERMO_LANGEVIN) == 0) {
305 for (auto const &p : cell_structure->local_particles()) {
306 if (p.stoner_wohlfarth_is_enabled()) {
307 runtimeErrorMsg() << "The thermal Stoner-Wohlfarth model requires the "
308 "Langevin thermostat";
309 break;
310 }
311 }
312 }
313#endif // ESPRESSO_THERMAL_STONER_WOHLFARTH
314}
315
316#ifdef ESPRESSO_WALBERLA
317void walberla_tau_sanity_checks(std::string method, double tau,
318 double time_step) {
319 if (time_step <= 0.) {
320 return;
321 }
322 // use float epsilon since tau may be a float
323 auto const eps = static_cast<double>(std::numeric_limits<float>::epsilon());
324 if ((tau - time_step) / (tau + time_step) < -eps)
325 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
326 ") must be >= MD time_step (" +
327 std::to_string(time_step) + ")");
328 auto const factor = tau / time_step;
329 if (std::fabs(std::round(factor) - factor) / factor > eps)
330 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
331 ") must be an integer multiple of the "
332 "MD time_step (" +
333 std::to_string(time_step) + "). Factor is " +
334 std::to_string(factor));
335}
336
342 double agrid) {
343 // waLBerla and ESPResSo must agree on domain decomposition
344 auto const tol = agrid / 1E6;
345 if ((lattice_left - geo_left).norm2() > tol or
346 (lattice_right - geo_right).norm2() > tol) {
347 std::stringstream error_msg;
348 error_msg << "waLBerla and ESPResSo disagree about domain decomposition"
349 << "\nMPI rank " << ::this_node << ": "
350 << "left ESPResSo: [" << geo_left << "], "
351 << "left waLBerla: [" << lattice_left << "]"
352 << "\nMPI rank " << ::this_node << ": "
353 << "right ESPResSo: [" << geo_right << "], "
354 << "right waLBerla: [" << lattice_right << "]"
355 << "\nfor method: " << method;
356 throw std::runtime_error(error_msg.str());
357 }
358}
359#endif // ESPRESSO_WALBERLA
360
362#ifdef ESPRESSO_CALIPER
364#endif
365 auto &cell_structure = *system.cell_structure;
366 auto const offset = LeesEdwards::verlet_list_offset(
367 *system.box_geo, cell_structure.get_le_pos_offset_at_last_resort());
368 if (cell_structure.check_resort_required(offset)) {
369 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
370 }
371}
372
373/** @brief Calls the hook for propagation kernels before the force calculation
374 * @return whether or not to stop the integration loop early.
375 */
376static bool integrator_step_1(CellStructure &cell_structure,
377 Propagation const &propagation,
378 System::System &system, double time_step) {
379#ifdef ESPRESSO_CALIPER
381#endif
382 // steepest decent
384 return system.steepest_descent->propagate(cell_structure);
385
386 auto const &thermostat = *system.thermostat;
387 auto const kT = thermostat.kT;
388 cell_structure.for_each_local_particle([&](Particle &p) {
389#ifdef ESPRESSO_VIRTUAL_SITES
390 // virtual sites are updated later in the integration loop
391 if (p.is_virtual())
392 return;
393#endif
394 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
395 if (propagation.should_propagate_with(
397 symplectic_euler_propagator_1(p, time_step);
399 symplectic_euler_propagator_1(p, time_step);
400#ifdef ESPRESSO_ROTATION
402 symplectic_euler_rotator_1(p, time_step);
403#endif
405 symplectic_euler_propagator_1(p, time_step);
406#ifdef ESPRESSO_ROTATION
408 symplectic_euler_rotator_1(p, time_step);
409#endif
410 } else {
411 if (propagation.should_propagate_with(
413 velocity_verlet_propagator_1(p, time_step);
415 velocity_verlet_propagator_1(p, time_step);
416#ifdef ESPRESSO_ROTATION
418 velocity_verlet_rotator_1(p, time_step);
419#endif
421 velocity_verlet_propagator_1(p, time_step);
422#ifdef ESPRESSO_ROTATION
424 velocity_verlet_rotator_1(p, time_step);
425#endif
426 }
428 brownian_dynamics_propagator(*thermostat.brownian, p, time_step, kT);
429#ifdef ESPRESSO_ROTATION
431 brownian_dynamics_rotator(*thermostat.brownian, p, time_step, kT);
432#endif
433 });
434
435#ifdef ESPRESSO_NPT
439 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
441 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
442 time_step, system);
443 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
445 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
446 time_step, system);
447 }
448 }
449#endif
450
451#ifdef ESPRESSO_STOKESIAN_DYNAMICS
456 *system.stokesian_dynamics, *thermostat.stokesian,
457 time_step, kT);
458 }
459#endif // ESPRESSO_STOKESIAN_DYNAMICS
460
461 return false;
462}
463
464static void integrator_step_2(CellStructure &cell_structure,
465 Propagation const &propagation,
467 double time_step) {
468#ifdef ESPRESSO_CALIPER
470#endif
472 return;
473
474 cell_structure.for_each_local_particle([&](Particle &p) {
475#ifdef ESPRESSO_VIRTUAL_SITES
476 // virtual sites are updated later in the integration loop
477 if (p.is_virtual())
478 return;
479#endif
480 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
481 if (propagation.should_propagate_with(
483 symplectic_euler_propagator_2(p, time_step);
485 symplectic_euler_propagator_2(p, time_step);
486#ifdef ESPRESSO_ROTATION
488 symplectic_euler_rotator_2(p, time_step);
489#endif
491 symplectic_euler_propagator_2(p, time_step);
492#ifdef ESPRESSO_ROTATION
494 symplectic_euler_rotator_2(p, time_step);
495#endif
496 } else {
497 if (propagation.should_propagate_with(
499 velocity_verlet_propagator_2(p, time_step);
501 velocity_verlet_propagator_2(p, time_step);
502#ifdef ESPRESSO_ROTATION
504 velocity_verlet_rotator_2(p, time_step);
505#endif
507 velocity_verlet_propagator_2(p, time_step);
508#ifdef ESPRESSO_ROTATION
510 velocity_verlet_rotator_2(p, time_step);
511#endif
512 }
513 });
514
515#ifdef ESPRESSO_NPT
519 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
521 cell_structure.local_particles().filter(pred), time_step, system);
522 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
524 cell_structure.local_particles().filter(pred), time_step, system);
525 }
526 }
527#endif
528}
529
531#ifdef ESPRESSO_CALIPER
533#endif
534 auto &propagation = *this->propagation;
535#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
536 auto const has_vs_rel = [&propagation]() {
537 return propagation.used_propagations &
541 };
542#endif
543#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
544 auto const has_vs_com = [&propagation]() {
545 return propagation.used_propagations &
547 };
548#endif
549#ifdef ESPRESSO_BOND_CONSTRAINT
550 auto const n_rigid_bonds = bonded_ias->get_n_rigid_bonds();
551#endif
552
553 // Prepare particle structure and run sanity checks of all active algorithms
554 propagation.update_default_propagation(thermostat->thermo_switch);
555 update_used_propagations();
556 on_integration_start();
557
558 // If any method vetoes (e.g. P3M not initialized), immediately bail out
560 return INTEG_ERROR_RUNTIME;
561
562 // Additional preparations for the first integration step
565 propagation.recalc_forces)) {
566#ifdef ESPRESSO_CALIPER
567 CALI_MARK_BEGIN("Initial Force Calculation");
568#endif
569 thermostat->lb_coupling_deactivate();
570
571#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
572 if (has_vs_rel()) {
573 vs_relative_update_particles(*cell_structure, *box_geo);
574 }
575#endif
576#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
577 if (has_vs_com()) {
578 vs_com_update_particles(*cell_structure, *box_geo);
579 }
580#endif
581
582 // Communication step: distribute ghost positions
583 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
584
585 calculate_forces();
586
587 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
588#ifdef ESPRESSO_ROTATION
589 convert_initial_torques(cell_structure->local_particles());
590#endif
591 }
592
593#ifdef ESPRESSO_CALIPER
594 CALI_MARK_END("Initial Force Calculation");
595#endif
596 }
597
598 thermostat->lb_coupling_activate();
599
601 return INTEG_ERROR_RUNTIME;
602
603 // Keep track of the number of Verlet updates (i.e. particle resorts)
604 int n_verlet_updates = 0;
605
606 // Keep track of whether an interrupt signal was caught (only in singleton
607 // mode, since signal handlers are unreliable with more than 1 MPI rank)
608 auto const singleton_mode = comm_cart.size() == 1;
609 auto caught_sigint = false;
610 auto caught_error = false;
611
612 auto lb_active = false;
613 auto ek_active = false;
614 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
615 lb_active = lb.is_solver_set();
616 ek_active = ek.is_ready_for_propagation();
617 }
618 auto const calc_md_steps_per_tau = [this](double tau) {
619 return static_cast<int>(std::round(tau / time_step));
620 };
621
622#ifdef ESPRESSO_VALGRIND
624#endif
625 // Integration loop
626#ifdef ESPRESSO_CALIPER
627 CALI_CXX_MARK_LOOP_BEGIN(integration_loop, "Integration loop");
628#endif
629 int integrated_steps = 0;
630 for (int step = 0; step < n_steps; step++) {
631#ifdef ESPRESSO_CALIPER
633#endif
634
635#ifdef ESPRESSO_BOND_CONSTRAINT
636 if (n_rigid_bonds)
637 save_old_position(cell_structure->local_particles(),
638 cell_structure->ghost_particles());
639#endif
640
641 lees_edwards->update_box_params(*box_geo, sim_time);
642 bool early_exit =
643 integrator_step_1(*cell_structure, propagation, *this, time_step);
644 if (early_exit)
645 break;
646
647 sim_time += time_step;
648 if (box_geo->type() == BoxType::LEES_EDWARDS) {
649 auto const kernel = LeesEdwards::Push{*box_geo};
650 cell_structure->for_each_local_particle(
651 [&kernel](Particle &p) { kernel(p); });
652 }
653
654#ifdef ESPRESSO_NPT
655 if (not has_npt_enabled())
656#endif
657 {
659 }
660 // Propagate philox RNG counters
661 thermostat->philox_counter_increment();
662
663#ifdef ESPRESSO_BOND_CONSTRAINT
664 // Correct particle positions that participate in a rigid/constrained bond
665 if (n_rigid_bonds) {
666 correct_position_shake(*cell_structure, *box_geo, *bonded_ias);
667 }
668#endif
669
670#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
671 if (has_vs_rel()) {
672#ifdef ESPRESSO_NPT
673 if (has_npt_enabled()) {
674 cell_structure->update_ghosts_and_resort_particle(
676 }
677#endif // ESPRESSO_NPT
678 vs_relative_update_particles(*cell_structure, *box_geo);
679 }
680#endif // ESPRESSO_VIRTUAL_SITES_RELATIVE
681#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
682 if (has_vs_com()) {
683#ifdef ESPRESSO_NPT
684 if (has_npt_enabled()) {
685 cell_structure->update_ghosts_and_resort_particle(
687 }
688#endif // ESPRESSO_NPT
689 vs_com_update_particles(*cell_structure, *box_geo);
690 }
691#endif // ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
692
693 if (cell_structure->get_resort_particles() >= Cells::RESORT_LOCAL)
695
696 // Communication step: distribute ghost positions
697 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
698
699#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
700 integrate_magnetodynamics();
701#endif
702
703 calculate_forces();
704
705#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
706 if (thermostat->lb and
707 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
708 lb_tracers_add_particle_force_to_fluid(*cell_structure, *box_geo,
709 *local_geo, lb);
710 }
711#endif
712 integrator_step_2(*cell_structure, propagation, *this, time_step);
713 if (propagation.integ_switch == INTEG_METHOD_BD) {
715 }
716 if (box_geo->type() == BoxType::LEES_EDWARDS) {
717 auto const kernel = LeesEdwards::UpdateOffset{*box_geo};
718 cell_structure->for_each_local_particle(
719 [&kernel](Particle &p) { kernel(p); });
720 }
721#ifdef ESPRESSO_BOND_CONSTRAINT
722 if (n_rigid_bonds) {
723 correct_velocity_shake(*cell_structure, *box_geo, *bonded_ias);
724 }
725#endif
726
727 // propagate one-step functionalities
728 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
729 if (lb_active and ek_active) {
730 // assume that they are coupled, which is not necessarily true
731 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
732 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
733
736 << "LB and EK are active but with different time steps.";
737 }
738
739 assert(lb.is_gpu() == ek.is_gpu());
740 assert(propagation.lb_skipped_md_steps ==
741 propagation.ek_skipped_md_steps);
742
743 propagation.lb_skipped_md_steps += 1;
744 propagation.ek_skipped_md_steps += 1;
745 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
746 propagation.lb_skipped_md_steps = 0;
747 propagation.ek_skipped_md_steps = 0;
748#ifdef ESPRESSO_CALIPER
749 CALI_MARK_BEGIN("lb_propagation");
750#endif
751 lb.propagate();
752 lb.ghost_communication_vel();
753#ifdef ESPRESSO_CALIPER
754 CALI_MARK_END("lb_propagation");
755 CALI_MARK_BEGIN("ek_propagation");
756#endif
757 ek.propagate();
758#ifdef ESPRESSO_CALIPER
759 CALI_MARK_END("ek_propagation");
760#endif
761 }
762 } else if (lb_active) {
763 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
764 propagation.lb_skipped_md_steps += 1;
765 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
766 propagation.lb_skipped_md_steps = 0;
767#ifdef ESPRESSO_CALIPER
768 CALI_MARK_BEGIN("lb_propagation");
769#endif
770 lb.propagate();
771#ifdef ESPRESSO_CALIPER
772 CALI_MARK_END("lb_propagation");
773#endif
774 }
775 } else if (ek_active) {
776 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
777 propagation.ek_skipped_md_steps += 1;
778 if (propagation.ek_skipped_md_steps >= md_steps_per_ek_step) {
779 propagation.ek_skipped_md_steps = 0;
780#ifdef ESPRESSO_CALIPER
781 CALI_MARK_BEGIN("ek_propagation");
782#endif
783 ek.propagate();
784#ifdef ESPRESSO_CALIPER
785 CALI_MARK_END("ek_propagation");
786#endif
787 }
788 }
789 if (lb_active and (propagation.used_propagations &
791 thermostat->lb->rng_increment();
792 }
793
794#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
795 if (thermostat->lb and
796 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
797#ifdef ESPRESSO_CALIPER
798 CALI_MARK_BEGIN("lb_tracers_propagation");
799#endif
800 if (lb_active) {
801 lb.ghost_communication_vel();
802 }
803 lb_tracers_propagate(*cell_structure, lb, time_step);
804#ifdef ESPRESSO_CALIPER
805 CALI_MARK_END("lb_tracers_propagation");
806#endif
807 }
808#endif
809
810#ifdef ESPRESSO_COLLISION_DETECTION
811 collision_detection->handle_collisions();
812#endif
813 bond_breakage->process_queue(*this);
814 }
815
817
819 caught_error = true;
820 break;
821 }
822
823 // Check if SIGINT has been caught.
824 if (singleton_mode and ctrl_C == 1) {
825 caught_sigint = true;
826 break;
827 }
828
829 } // for-loop over integration steps
830 if (lb_active) {
831 lb.ghost_communication();
832 }
833 lees_edwards->update_box_params(*box_geo, sim_time);
834#ifdef ESPRESSO_CALIPER
836#endif
837
838#ifdef ESPRESSO_VALGRIND
840#endif
841
842#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
843 if (has_vs_rel()) {
844 vs_relative_update_particles(*cell_structure, *box_geo);
845 }
846#endif
847#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
848 if (has_vs_com()) {
849 vs_com_update_particles(*cell_structure, *box_geo);
850 }
851#endif
852
853 // Verlet list statistics
854 cell_structure->update_verlet_stats(n_steps, n_verlet_updates);
855
856#ifdef ESPRESSO_NPT
857 if (has_npt_enabled()) {
858 synchronize_npt_state();
859 }
860#endif
861 if (caught_sigint) {
862 ctrl_C = 0;
863 return INTEG_ERROR_SIGINT;
864 }
865 if (caught_error) {
866 return INTEG_ERROR_RUNTIME;
867 }
868#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
869 if (boost::mpi::all_reduce(::comm_cart, not cell_structure->use_verlet_list,
870 std::logical_or<>())) {
871 cell_structure->use_verlet_list = false;
872 }
873#endif
874 return integrated_steps;
875}
876
878 bool update_accumulators) {
879 assert(n_steps >= 0);
880
881 // Override the signal handler so that the integrator obeys Ctrl+C
882 SignalHandler sa(SIGINT, [](int) { ctrl_C = 1; });
883
884 /* if skin wasn't set, do an educated guess now */
885 if (not cell_structure->is_verlet_skin_set()) {
886 try {
887 cell_structure->set_verlet_skin_heuristic();
888 } catch (...) {
889 if (comm_cart.rank() == 0) {
890 throw;
891 }
892 return INTEG_ERROR_RUNTIME;
893 }
894 }
895
897 return integrate(n_steps, reuse_forces);
898 }
899
900 for (int i = 0; i < n_steps;) {
901 /* Integrate to either the next accumulator update, or the
902 * end, depending on what comes first. */
903 auto const steps =
904 std::min((n_steps - i), auto_update_accumulators->next_update());
905
906 auto const local_retval = integrate(steps, reuse_forces);
907
908 // make sure all ranks exit when one rank fails
909 std::remove_const_t<decltype(local_retval)> global_retval;
910 boost::mpi::all_reduce(comm_cart, local_retval, global_retval,
911 std::plus<int>());
912 if (global_retval < 0) {
913 return global_retval; // propagate error code
914 }
915
917
918 (*auto_update_accumulators)(comm_cart, steps);
919
920 i += steps;
921 }
922
923 return 0;
924}
925
927 sim_time = value;
928 propagation->recalc_forces = true;
929 lees_edwards->update_box_params(*box_geo, sim_time);
930}
@ 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.
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)
int integrate(int n_steps, int reuse_forces)
Integrate equations of motion.
void vs_com_update_particles(CellStructure &cell_structure, BoxGeometry const &box_geo)
Definition com.cpp:130
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.
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:88
volatile std::sig_atomic_t ctrl_C
Definition integrate.cpp:93
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.