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-2023 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 "ParticleRange.hpp"
40#include "PropagationMode.hpp"
41#include "accumulators/AutoUpdateAccumulators.hpp"
45#include "cells.hpp"
46#include "collision_detection/CollisionDetection.hpp"
47#include "communication.hpp"
48#include "errorhandling.hpp"
49#include "forces.hpp"
51#include "lb/utils.hpp"
55#include "npt.hpp"
56#include "rattle.hpp"
57#include "rotation.hpp"
58#include "signalhandling.hpp"
59#include "system/System.hpp"
60#include "thermostat.hpp"
62#include "virtual_sites/com.hpp"
65
66#include <boost/mpi/collectives/all_reduce.hpp>
67
68#ifdef ESPRESSO_CALIPER
69#include <caliper/cali.h>
70#endif
71
72#ifdef ESPRESSO_VALGRIND
73#include <callgrind.h>
74#endif
75
76#include <algorithm>
77#include <cassert>
78#include <cmath>
79#include <csignal>
80#include <functional>
81#include <limits>
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}
198
199void System::System::integrator_sanity_checks() const {
200 auto const thermo_switch = thermostat->thermo_switch;
201 if (time_step <= 0.) {
202 runtimeErrorMsg() << "time_step not set";
203 }
204 if (propagation->integ_switch == INTEG_METHOD_STEEPEST_DESCENT) {
205 if (thermo_switch != THERMO_OFF) {
207 << "The steepest descent integrator is incompatible with thermostats";
208 }
209 }
210 if (propagation->integ_switch == INTEG_METHOD_NVT) {
211 if (thermo_switch & (THERMO_NPT_ISO | THERMO_BROWNIAN | THERMO_SD)) {
212 runtimeErrorMsg() << "The VV integrator is incompatible with the "
213 "currently active combination of thermostats";
214 }
215 }
216#ifdef ESPRESSO_NPT
217 if (propagation->used_propagations & PropagationMode::TRANS_LANGEVIN_NPT) {
218 if (thermo_switch != THERMO_NPT_ISO) {
219 runtimeErrorMsg() << "The NpT integrator requires the NpT thermostat";
220 }
221 if (box_geo->type() == BoxType::LEES_EDWARDS) {
222 runtimeErrorMsg() << "The NpT integrator cannot use Lees-Edwards";
223 }
224 try {
225 nptiso->coulomb_dipole_sanity_checks(*this);
226 } catch (std::runtime_error const &err) {
227 runtimeErrorMsg() << err.what();
228 }
229 }
230#endif
231 if (propagation->used_propagations & PropagationMode::TRANS_BROWNIAN) {
232 if (thermo_switch != THERMO_BROWNIAN) {
233 runtimeErrorMsg() << "The BD integrator requires the BD thermostat";
234 }
235 }
236 if (propagation->used_propagations & PropagationMode::TRANS_STOKESIAN) {
237#ifdef ESPRESSO_STOKESIAN_DYNAMICS
238 if (thermo_switch != THERMO_SD) {
239 runtimeErrorMsg() << "The SD integrator requires the SD thermostat";
240 }
241#endif
242 }
243 if (lb.is_solver_set() and (propagation->used_propagations &
246 if (thermostat->lb == nullptr) {
247 runtimeErrorMsg() << "The LB integrator requires the LB thermostat";
248 }
249 }
250 if (bonded_ias->get_n_thermalized_bonds() >= 1 and
251 (thermostat->thermalized_bond == nullptr or
252 (thermo_switch & THERMO_BOND) == 0)) {
254 << "Thermalized bonds require the thermalized_bond thermostat";
255 }
256
257#ifdef ESPRESSO_ROTATION
258 for (auto const &p : cell_structure->local_particles()) {
259 using namespace PropagationMode;
260 if (p.can_rotate() and not p.is_virtual() and
261 (p.propagation() & (SYSTEM_DEFAULT | ROT_EULER | ROT_LANGEVIN |
262 ROT_BROWNIAN | ROT_STOKESIAN)) == 0) {
264 << "Rotating particles must have a rotation propagation mode enabled";
265 break;
266 }
267 }
268#endif // ESPRESSO_ROTATION
269
270#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
271 if ((thermo_switch & THERMO_LANGEVIN) == 0) {
272 for (auto const &p : cell_structure->local_particles()) {
273 if (p.stoner_wohlfarth_is_enabled()) {
274 runtimeErrorMsg() << "The thermal Stoner-Wohlfarth model requires the "
275 "Langevin thermostat";
276 break;
277 }
278 }
279 }
280#endif // ESPRESSO_THERMAL_STONER_WOHLFARTH
281}
282
283#ifdef ESPRESSO_WALBERLA
284void walberla_tau_sanity_checks(std::string method, double tau,
285 double time_step) {
286 if (time_step <= 0.) {
287 return;
288 }
289 // use float epsilon since tau may be a float
290 auto const eps = static_cast<double>(std::numeric_limits<float>::epsilon());
291 if ((tau - time_step) / (tau + time_step) < -eps)
292 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
293 ") must be >= MD time_step (" +
294 std::to_string(time_step) + ")");
295 auto const factor = tau / time_step;
296 if (std::fabs(std::round(factor) - factor) / factor > eps)
297 throw std::invalid_argument(method + " tau (" + std::to_string(tau) +
298 ") must be an integer multiple of the "
299 "MD time_step (" +
300 std::to_string(time_step) + "). Factor is " +
301 std::to_string(factor));
302}
303
304void walberla_agrid_sanity_checks(std::string method,
305 Utils::Vector3d const &geo_left,
306 Utils::Vector3d const &geo_right,
307 Utils::Vector3d const &lattice_left,
308 Utils::Vector3d const &lattice_right,
309 double agrid) {
310 // waLBerla and ESPResSo must agree on domain decomposition
311 auto const tol = agrid / 1E6;
312 if ((lattice_left - geo_left).norm2() > tol or
313 (lattice_right - geo_right).norm2() > tol) {
314 runtimeErrorMsg() << "\nMPI rank " << ::this_node << ": "
315 << "left ESPResSo: [" << geo_left << "], "
316 << "left waLBerla: [" << lattice_left << "]"
317 << "\nMPI rank " << ::this_node << ": "
318 << "right ESPResSo: [" << geo_right << "], "
319 << "right waLBerla: [" << lattice_right << "]"
320 << "\nfor method: " << method;
321 throw std::runtime_error(
322 "waLBerla and ESPResSo disagree about domain decomposition.");
323 }
324}
325#endif // ESPRESSO_WALBERLA
326
328#ifdef ESPRESSO_CALIPER
329 CALI_CXX_MARK_FUNCTION;
330#endif
331 auto &cell_structure = *system.cell_structure;
332 auto const offset = LeesEdwards::verlet_list_offset(
333 *system.box_geo, cell_structure.get_le_pos_offset_at_last_resort());
334 if (cell_structure.check_resort_required(offset)) {
335 cell_structure.set_resort_particles(Cells::RESORT_LOCAL);
336 }
337}
338
339/** @brief Calls the hook for propagation kernels before the force calculation
340 * @return whether or not to stop the integration loop early.
341 */
342static bool integrator_step_1(CellStructure &cell_structure,
343 Propagation const &propagation,
344 System::System &system, double time_step) {
345#ifdef ESPRESSO_CALIPER
346 CALI_CXX_MARK_FUNCTION;
347#endif
348 // steepest decent
350 return steepest_descent_step(cell_structure.local_particles());
351
352 auto const &thermostat = *system.thermostat;
353 auto const kT = thermostat.kT;
354 cell_structure.for_each_local_particle([&](Particle &p) {
355#ifdef ESPRESSO_VIRTUAL_SITES
356 // virtual sites are updated later in the integration loop
357 if (p.is_virtual())
358 return;
359#endif
360 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
361 if (propagation.should_propagate_with(
363 symplectic_euler_propagator_1(p, time_step);
365 symplectic_euler_propagator_1(p, time_step);
366#ifdef ESPRESSO_ROTATION
368 symplectic_euler_rotator_1(p, time_step);
369#endif
371 symplectic_euler_propagator_1(p, time_step);
372#ifdef ESPRESSO_ROTATION
374 symplectic_euler_rotator_1(p, time_step);
375#endif
376 } else {
377 if (propagation.should_propagate_with(
379 velocity_verlet_propagator_1(p, time_step);
381 velocity_verlet_propagator_1(p, time_step);
382#ifdef ESPRESSO_ROTATION
384 velocity_verlet_rotator_1(p, time_step);
385#endif
387 velocity_verlet_propagator_1(p, time_step);
388#ifdef ESPRESSO_ROTATION
390 velocity_verlet_rotator_1(p, time_step);
391#endif
392 }
394 brownian_dynamics_propagator(*thermostat.brownian, p, time_step, kT);
395#ifdef ESPRESSO_ROTATION
397 brownian_dynamics_rotator(*thermostat.brownian, p, time_step, kT);
398#endif
399 });
400
401#ifdef ESPRESSO_NPT
404 auto pred = PropagationPredicateNPT(propagation.default_propagation);
405 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
407 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
408 time_step, system);
409 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
411 cell_structure.local_particles().filter(pred), *thermostat.npt_iso,
412 time_step, system);
413 }
414 }
415#endif
416
417#ifdef ESPRESSO_STOKESIAN_DYNAMICS
420 auto pred = PropagationPredicateStokesian(propagation.default_propagation);
421 stokesian_dynamics_step_1(cell_structure.local_particles().filter(pred),
422 *thermostat.stokesian, time_step, kT);
423 }
424#endif // ESPRESSO_STOKESIAN_DYNAMICS
425
426 return false;
427}
428
429static void integrator_step_2(CellStructure &cell_structure,
430 Propagation const &propagation,
431 [[maybe_unused]] System::System &system,
432 double time_step) {
433#ifdef ESPRESSO_CALIPER
434 CALI_CXX_MARK_FUNCTION;
435#endif
437 return;
438
439 cell_structure.for_each_local_particle([&](Particle &p) {
440#ifdef ESPRESSO_VIRTUAL_SITES
441 // virtual sites are updated later in the integration loop
442 if (p.is_virtual())
443 return;
444#endif
445 if (propagation.integ_switch == INTEG_METHOD_SYMPLECTIC_EULER) {
446 if (propagation.should_propagate_with(
448 symplectic_euler_propagator_2(p, time_step);
450 symplectic_euler_propagator_2(p, time_step);
451#ifdef ESPRESSO_ROTATION
453 symplectic_euler_rotator_2(p, time_step);
454#endif
456 symplectic_euler_propagator_2(p, time_step);
457#ifdef ESPRESSO_ROTATION
459 symplectic_euler_rotator_2(p, time_step);
460#endif
461 } else {
462 if (propagation.should_propagate_with(
464 velocity_verlet_propagator_2(p, time_step);
466 velocity_verlet_propagator_2(p, time_step);
467#ifdef ESPRESSO_ROTATION
469 velocity_verlet_rotator_2(p, time_step);
470#endif
472 velocity_verlet_propagator_2(p, time_step);
473#ifdef ESPRESSO_ROTATION
475 velocity_verlet_rotator_2(p, time_step);
476#endif
477 }
478 });
479
480#ifdef ESPRESSO_NPT
483 auto pred = PropagationPredicateNPT(propagation.default_propagation);
484 if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_AND) {
486 cell_structure.local_particles().filter(pred), time_step, system);
487 } else if (propagation.integ_switch == INTEG_METHOD_NPT_ISO_MTK) {
489 cell_structure.local_particles().filter(pred), time_step, system);
490 }
491 }
492#endif
493}
494
495int System::System::integrate(int n_steps, int reuse_forces) {
496#ifdef ESPRESSO_CALIPER
497 CALI_CXX_MARK_FUNCTION;
498#endif
499 auto &propagation = *this->propagation;
500#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
501 auto const has_vs_rel = [&propagation]() {
502 return propagation.used_propagations &
506 };
507#endif
508#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
509 auto const has_vs_com = [&propagation]() {
510 return propagation.used_propagations &
512 };
513#endif
514#ifdef ESPRESSO_BOND_CONSTRAINT
515 auto const n_rigid_bonds = bonded_ias->get_n_rigid_bonds();
516#endif
517
518 // Prepare particle structure and run sanity checks of all active algorithms
519 propagation.update_default_propagation(thermostat->thermo_switch);
520 update_used_propagations();
521 on_integration_start();
522
523 // If any method vetoes (e.g. P3M not initialized), immediately bail out
525 return INTEG_ERROR_RUNTIME;
526
527 // Additional preparations for the first integration step
528 if (reuse_forces == INTEG_REUSE_FORCES_NEVER or
529 ((reuse_forces != INTEG_REUSE_FORCES_ALWAYS) and
530 propagation.recalc_forces)) {
531#ifdef ESPRESSO_CALIPER
532 CALI_MARK_BEGIN("Initial Force Calculation");
533#endif
534 thermostat->lb_coupling_deactivate();
535
536#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
537 if (has_vs_rel()) {
538 vs_relative_update_particles(*cell_structure, *box_geo);
539 }
540#endif
541#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
542 if (has_vs_com()) {
543 vs_com_update_particles(*cell_structure, *box_geo);
544 }
545#endif
546
547 // Communication step: distribute ghost positions
548 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
549
550 calculate_forces();
551
552 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
553#ifdef ESPRESSO_ROTATION
554 convert_initial_torques(cell_structure->local_particles());
555#endif
556 }
557
558#ifdef ESPRESSO_CALIPER
559 CALI_MARK_END("Initial Force Calculation");
560#endif
561 }
562
563 thermostat->lb_coupling_activate();
564
566 return INTEG_ERROR_RUNTIME;
567
568 // Keep track of the number of Verlet updates (i.e. particle resorts)
569 int n_verlet_updates = 0;
570
571 // Keep track of whether an interrupt signal was caught (only in singleton
572 // mode, since signal handlers are unreliable with more than 1 MPI rank)
573 auto const singleton_mode = comm_cart.size() == 1;
574 auto caught_sigint = false;
575 auto caught_error = false;
576
577 auto lb_active = false;
578 auto ek_active = false;
579 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
580 lb_active = lb.is_solver_set();
581 ek_active = ek.is_ready_for_propagation();
582 }
583 auto const calc_md_steps_per_tau = [this](double tau) {
584 return static_cast<int>(std::round(tau / time_step));
585 };
586
587#ifdef ESPRESSO_VALGRIND
588 CALLGRIND_START_INSTRUMENTATION;
589#endif
590 // Integration loop
591#ifdef ESPRESSO_CALIPER
592 CALI_CXX_MARK_LOOP_BEGIN(integration_loop, "Integration loop");
593#endif
594 int integrated_steps = 0;
595 for (int step = 0; step < n_steps; step++) {
596#ifdef ESPRESSO_CALIPER
597 CALI_CXX_MARK_LOOP_ITERATION(integration_loop, step);
598#endif
599
600#ifdef ESPRESSO_BOND_CONSTRAINT
601 if (n_rigid_bonds)
602 save_old_position(cell_structure->local_particles(),
603 cell_structure->ghost_particles());
604#endif
605
606 lees_edwards->update_box_params(*box_geo, sim_time);
607 bool early_exit =
608 integrator_step_1(*cell_structure, propagation, *this, time_step);
609 if (early_exit)
610 break;
611
612 sim_time += time_step;
613 if (box_geo->type() == BoxType::LEES_EDWARDS) {
614 auto const kernel = LeesEdwards::Push{*box_geo};
615 cell_structure->for_each_local_particle(
616 [&kernel](Particle &p) { kernel(p); });
617 }
618
619#ifdef ESPRESSO_NPT
620 if (not has_npt_enabled())
621#endif
622 {
624 }
625 // Propagate philox RNG counters
626 thermostat->philox_counter_increment();
627
628#ifdef ESPRESSO_BOND_CONSTRAINT
629 // Correct particle positions that participate in a rigid/constrained bond
630 if (n_rigid_bonds) {
631 correct_position_shake(*cell_structure, *box_geo, *bonded_ias);
632 }
633#endif
634
635#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
636 if (has_vs_rel()) {
637#ifdef ESPRESSO_NPT
638 if (has_npt_enabled()) {
639 cell_structure->update_ghosts_and_resort_particle(
641 }
642#endif // ESPRESSO_NPT
643 vs_relative_update_particles(*cell_structure, *box_geo);
644 }
645#endif // ESPRESSO_VIRTUAL_SITES_RELATIVE
646#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
647 if (has_vs_com()) {
648#ifdef ESPRESSO_NPT
649 if (has_npt_enabled()) {
650 cell_structure->update_ghosts_and_resort_particle(
652 }
653#endif // ESPRESSO_NPT
654 vs_com_update_particles(*cell_structure, *box_geo);
655 }
656#endif // ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
657
658 if (cell_structure->get_resort_particles() >= Cells::RESORT_LOCAL)
659 n_verlet_updates++;
660
661 // Communication step: distribute ghost positions
662 cell_structure->update_ghosts_and_resort_particle(get_global_ghost_flags());
663
664#ifdef ESPRESSO_THERMAL_STONER_WOHLFARTH
665 run_magnetodynamics(*cell_structure, *thermostat);
666#endif
667
668 calculate_forces();
669
670#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
671 if (thermostat->lb and
672 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
673 lb_tracers_add_particle_force_to_fluid(*cell_structure, *box_geo,
674 *local_geo, lb);
675 }
676#endif
677 integrator_step_2(*cell_structure, propagation, *this, time_step);
678 if (propagation.integ_switch == INTEG_METHOD_BD) {
680 }
681 if (box_geo->type() == BoxType::LEES_EDWARDS) {
682 auto const kernel = LeesEdwards::UpdateOffset{*box_geo};
683 cell_structure->for_each_local_particle(
684 [&kernel](Particle &p) { kernel(p); });
685 }
686#ifdef ESPRESSO_BOND_CONSTRAINT
687 if (n_rigid_bonds) {
688 correct_velocity_shake(*cell_structure, *box_geo, *bonded_ias);
689 }
690#endif
691
692 // propagate one-step functionalities
693 if (propagation.integ_switch != INTEG_METHOD_STEEPEST_DESCENT) {
694 if (lb_active and ek_active) {
695 // assume that they are coupled, which is not necessarily true
696 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
697 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
698
699 if (md_steps_per_lb_step != md_steps_per_ek_step) {
701 << "LB and EK are active but with different time steps.";
702 }
703
704 assert(lb.is_gpu() == ek.is_gpu());
705 assert(propagation.lb_skipped_md_steps ==
706 propagation.ek_skipped_md_steps);
707
708 propagation.lb_skipped_md_steps += 1;
709 propagation.ek_skipped_md_steps += 1;
710 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
711 propagation.lb_skipped_md_steps = 0;
712 propagation.ek_skipped_md_steps = 0;
713#ifdef ESPRESSO_CALIPER
714 CALI_MARK_BEGIN("lb_propagation");
715#endif
716 lb.propagate();
717 lb.ghost_communication_vel();
718#ifdef ESPRESSO_CALIPER
719 CALI_MARK_END("lb_propagation");
720 CALI_MARK_BEGIN("ek_propagation");
721#endif
722 ek.propagate();
723#ifdef ESPRESSO_CALIPER
724 CALI_MARK_END("ek_propagation");
725#endif
726 }
727 } else if (lb_active) {
728 auto const md_steps_per_lb_step = calc_md_steps_per_tau(lb.get_tau());
729 propagation.lb_skipped_md_steps += 1;
730 if (propagation.lb_skipped_md_steps >= md_steps_per_lb_step) {
731 propagation.lb_skipped_md_steps = 0;
732#ifdef ESPRESSO_CALIPER
733 CALI_MARK_BEGIN("lb_propagation");
734#endif
735 lb.propagate();
736#ifdef ESPRESSO_CALIPER
737 CALI_MARK_END("lb_propagation");
738#endif
739 }
740 } else if (ek_active) {
741 auto const md_steps_per_ek_step = calc_md_steps_per_tau(ek.get_tau());
742 propagation.ek_skipped_md_steps += 1;
743 if (propagation.ek_skipped_md_steps >= md_steps_per_ek_step) {
744 propagation.ek_skipped_md_steps = 0;
745#ifdef ESPRESSO_CALIPER
746 CALI_MARK_BEGIN("ek_propagation");
747#endif
748 ek.propagate();
749#ifdef ESPRESSO_CALIPER
750 CALI_MARK_END("ek_propagation");
751#endif
752 }
753 }
754 if (lb_active and (propagation.used_propagations &
756 thermostat->lb->rng_increment();
757 }
758
759#ifdef ESPRESSO_VIRTUAL_SITES_INERTIALESS_TRACERS
760 if (thermostat->lb and
761 (propagation.used_propagations & PropagationMode::TRANS_LB_TRACER)) {
762#ifdef ESPRESSO_CALIPER
763 CALI_MARK_BEGIN("lb_tracers_propagation");
764#endif
765 if (lb_active) {
766 lb.ghost_communication_vel();
767 }
768 lb_tracers_propagate(*cell_structure, lb, time_step);
769#ifdef ESPRESSO_CALIPER
770 CALI_MARK_END("lb_tracers_propagation");
771#endif
772 }
773#endif
774
775#ifdef ESPRESSO_COLLISION_DETECTION
776 collision_detection->handle_collisions();
777#endif
778 bond_breakage->process_queue(*this);
779 }
780
781 integrated_steps++;
782
784 caught_error = true;
785 break;
786 }
787
788 // Check if SIGINT has been caught.
789 if (singleton_mode and ctrl_C == 1) {
790 caught_sigint = true;
791 break;
792 }
793
794 } // for-loop over integration steps
795 if (lb_active) {
796 lb.ghost_communication();
797 }
798 lees_edwards->update_box_params(*box_geo, sim_time);
799#ifdef ESPRESSO_CALIPER
800 CALI_CXX_MARK_LOOP_END(integration_loop);
801#endif
802
803#ifdef ESPRESSO_VALGRIND
804 CALLGRIND_STOP_INSTRUMENTATION;
805#endif
806
807#ifdef ESPRESSO_VIRTUAL_SITES_RELATIVE
808 if (has_vs_rel()) {
809 vs_relative_update_particles(*cell_structure, *box_geo);
810 }
811#endif
812#ifdef ESPRESSO_VIRTUAL_SITES_CENTER_OF_MASS
813 if (has_vs_com()) {
814 vs_com_update_particles(*cell_structure, *box_geo);
815 }
816#endif
817
818 // Verlet list statistics
819 cell_structure->update_verlet_stats(n_steps, n_verlet_updates);
820
821#ifdef ESPRESSO_NPT
822 if (has_npt_enabled()) {
823 synchronize_npt_state();
824 }
825#endif
826 if (caught_sigint) {
827 ctrl_C = 0;
828 return INTEG_ERROR_SIGINT;
829 }
830 if (caught_error) {
831 return INTEG_ERROR_RUNTIME;
832 }
833#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
834 if (boost::mpi::all_reduce(::comm_cart, not cell_structure->use_verlet_list,
835 std::logical_or<>())) {
836 cell_structure->use_verlet_list = false;
837 }
838#endif
839 return integrated_steps;
840}
841
842int System::System::integrate_with_signal_handler(int n_steps, int reuse_forces,
843 bool update_accumulators) {
844 assert(n_steps >= 0);
845
846 // Override the signal handler so that the integrator obeys Ctrl+C
847 SignalHandler sa(SIGINT, [](int) { ctrl_C = 1; });
848
849 /* if skin wasn't set, do an educated guess now */
850 if (not cell_structure->is_verlet_skin_set()) {
851 try {
852 cell_structure->set_verlet_skin_heuristic();
853 } catch (...) {
854 if (comm_cart.rank() == 0) {
855 throw;
856 }
857 return INTEG_ERROR_RUNTIME;
858 }
859 }
860
861 if (not update_accumulators or n_steps == 0) {
862 return integrate(n_steps, reuse_forces);
863 }
864
865 for (int i = 0; i < n_steps;) {
866 /* Integrate to either the next accumulator update, or the
867 * end, depending on what comes first. */
868 auto const steps =
869 std::min((n_steps - i), auto_update_accumulators->next_update());
870
871 auto const local_retval = integrate(steps, reuse_forces);
872
873 // make sure all ranks exit when one rank fails
874 std::remove_const_t<decltype(local_retval)> global_retval;
875 boost::mpi::all_reduce(comm_cart, local_retval, global_retval,
876 std::plus<int>());
877 if (global_retval < 0) {
878 return global_retval; // propagate error code
879 }
880
881 reuse_forces = INTEG_REUSE_FORCES_ALWAYS;
882
883 (*auto_update_accumulators)(comm_cart, steps);
884
885 i += steps;
886 }
887
888 return 0;
889}
890
892 sim_time = value;
893 propagation->recalc_forces = true;
894 lees_edwards->update_box_params(*box_geo, sim_time);
895}
@ 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.
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:105
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()
Force calculation.
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:232
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:121
void convert_initial_torques(const ParticleRange &particles)
Convert torques to the body-fixed frame before the integration loop.
Definition rotation.cpp:188
This file contains all subroutines required to process rotational motion.
bool steepest_descent_step(const ParticleRange &particles)
Steepest descent integrator.
void stokesian_dynamics_step_1(ParticleRangeStokesian const &particles, StokesianThermostat const &stokesian, double time_step, double kT)
void run_magnetodynamics(CellStructure &cell_structure, Thermostat::Thermostat const &thermostat)
Run magnetodynamics update for local virtual particles.
Struct holding all information for one particle.
Definition Particle.hpp:450
auto is_virtual() const
Definition Particle.hpp:603
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 NpT isotropic for MTK approach.
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 for Andersen method.
void velocity_verlet_npt_Andersen_step_1(ParticleRangeNPT const &particles, IsotropicNptThermostat const &npt_iso, double time_step, System::System &system)
Special propagator for NpT isotropic for 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 for Andersen method.