pure-cpp 1.0.0
A C++ physics simulation benchmark comparing performance with Python implementations
space.hpp
Go to the documentation of this file.
1#ifndef SPACE_HPP
2#define SPACE_HPP
3
4/**
5 * \file space.hpp
6 * \brief N-body simulation space with gravitational interaction and collision
7 * response.
8 * \author Le Bars, Yoann
9 * \ingroup PhysicsCore
10 *
11 * This file defines the `Space` class, which manages the state and dynamics of
12 * all bodies in the simulation. The physics model includes:
13 *
14 * - **Gravitational Interaction**: All bodies exert a gravitational force on
15 * each other, calculated using Newton's law of universal gravitation.
16 *
17 * - **Collision Response**: Collisions between spherical bodies are handled
18 * using an impulse-based method that accounts for both linear and
19 * rotational motion, including friction.
20 *
21 * - **Integration**: The simulation state (position and velocity) is
22 * advanced using the **Velocity Verlet** algorithm, a time-reversible and
23 * energy-conserving semi-implicit Euler integration scheme.
24 *
25 * - **Broad-Phase Detection**: A k-d tree is used to efficiently find
26 * potentially colliding pairs of bodies, avoiding an O(N²) check.
27 *
28 * This file is part of the pure C++ benchmark.
29 *
30 * This program is free software: you can redistribute it and/or modify it
31 * under the terms of the GNU General Public License as published by the Free
32 * Software Foundation, either version 3 of the License, or (at your option)
33 * any later version.
34 *
35 * This program is distributed in the hope that it will be useful, but WITHOUT
36 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
37 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for
38 * more details.
39 *
40 * You should have received a copy of the GNU General Public License along with
41 * this program. If not, see <https://www.gnu.org/licenses/>.
42 */
43
44#include <omp.h>
45
46#include <chrono>
47#include <cstddef>
48#include <memory>
49#include <optional>
50#include <stdexcept>
51#include <vector>
52
53#include "body.hpp"
54#include "config.hpp"
55#include "constants.hpp"
56#include "kdtree.hpp"
57
58// Forward-declare the test fixture class from the global namespace.
59class SpaceTest;
60
61namespace Model {
62
63 /**
64 * \brief A data container to hold all relevant information for a single
65 * collision event.
66 */
68 /// \brief A proxy to the first body in the collision.
70 /// \brief A proxy to the second body in the collision.
72 /// \brief The normalised vector pointing from body 1 to body 2.
73 Vector3d nVec_;
74 /// \brief The vector from the centre of body 1 to the contact point.
75 Vector3d r1Vec_;
76 /// \brief The vector from the centre of body 2 to the contact point.
77 Vector3d r2Vec_;
78 /// \brief The relative velocity between the two bodies at the contact
79 /// point.
80 Vector3d vRel_;
81 /// \brief The penetration depth of the two bodies.
82 double overlap_;
83 };
84
85 /**
86 * \brief Class describing a space in which move several bodies.
87 * \ingroup PhysicsCore
88 */
89 class Space {
90 public:
91 /* Grant the test fixture from the global namespace access to protected
92 and private members. */
93 friend class ::SpaceTest;
94
95 /// \brief Default constructor.
96 Space() : dt_(0.), G_(G), epsilon_(EPSILON) {}
97
98 /**
99 * \brief Initialise the space.
100 *
101 * \param config The simulation configuration object containing all
102 * simulation parameters.
103 */
104 explicit Space(const Configuration::SimulationConfig& config);
105
106 /// \brief Destructor to clean up OpenMP locks.
107 ~Space();
108
109 /**
110 * \brief Number of bodies getter.
111 *
112 * \return The number of bodies in the space.
113 */
114 std::size_t n() const { return bodies_.size(); }
115
116 /**
117 * \brief Get a proxy to the i-th body.
118 *
119 * \param i Body index.
120 *
121 * \returns A mutable proxy object for the i-th body.
122 */
123 BodyProxy body(std::size_t i) { return bodies_[i]; }
124
125 /**
126 * \brief Get a const proxy to the i-th body.
127 *
128 * \param i Body index.
129 *
130 * \returns A const proxy object for the i-th body.
131 */
132 ConstBodyProxy body(std::size_t i) const { return bodies_[i]; }
133
134 /**
135 * \brief Get the time step used in the last completed simulation
136 * frame.
137 *
138 * \return The previous time step value.
139 */
140 [[nodiscard]] double getPreviousTimeStep() const {
141 return previous_dt_;
142 }
143
144 /**
145 * \brief Computes one full step of the simulation.
146 *
147 * This method orchestrates the main simulation loop, including:
148 * 1. Advancing positions and velocities (Velocity Verlet).
149 * 2. Detecting and resolving collisions.
150 * 3. Calculating new forces (gravity).
151 * \param iteration The current simulation iteration number.
152 */
153 void computeDynamics(std::size_t iteration);
154
155 /**
156 * \brief Calculates the kinetic and potential energy of the system.
157 * \return A tuple containing:
158 * 1. Translational Kinetic Energy
159 * 2. Rotational Kinetic Energy
160 * 3. Potential Energy
161 */
162 [[nodiscard]]
163 std::tuple<double, double, double> calculateSystemEnergy() const;
164
165 /**
166 * \brief Retrieves the current positions and orientations of all
167 * bodies.
168 * \return A tuple containing a vector of positions and a vector of
169 * quaternions.
170 */
171 [[nodiscard]] std::tuple<Vector3dVec, QuaterniondVec>
172 getAllBodyTransforms() const;
173
174 /**
175 * \brief Retrieves all data required for display for all bodies.
176 *
177 * \return A tuple containing vectors of positions, quaternions,
178 * torques, and angular accelerations.
179 */
180 [[nodiscard]] std::tuple<Vector3dVec, QuaterniondVec, Vector3dVec,
181 Vector3dVec>
182 getDataForDisplay() const;
183 void printProfilingReport() const;
184
185 // --- Configuration Getters/Setters ---
186 /**
187 * \brief Get the universal gravitational constant.
188 * \return The gravitational constant G.
189 */
190 [[nodiscard]] double getG() const { return G_; }
191
192 /**
193 * \brief Set the universal gravitational constant.
194 * \param g The new gravitational constant value.
195 */
196 void setG(double g) { G_ = g; }
197
198 /**
199 * \brief Get the numerical precision threshold.
200 * \return The epsilon value.
201 */
202 [[nodiscard]] double getEpsilon() const { return epsilon_; }
203
204 /**
205 * \brief Set the numerical precision threshold.
206 * \param eps The new epsilon value.
207 */
208 void setEpsilon(double eps) { epsilon_ = eps; }
209
210 /**
211 * \brief Get the coefficient of restitution.
212 * \return The coefficient of restitution.
213 */
214 [[nodiscard]] double getCoeffRestitution() const {
215 return coeffRestitution_;
216 }
217
218 /**
219 * \brief Set the coefficient of restitution.
220 * \param e The new coefficient of restitution (typically in [0, 1]).
221 */
223
224 /**
225 * \brief Get the coefficient of kinetic friction.
226 * \return The coefficient of kinetic friction.
227 */
228 [[nodiscard]] double getCoeffFriction() const { return coeffFriction_; }
229
230 /**
231 * \brief Set the coefficient of kinetic friction.
232 * \param mu The new coefficient of kinetic friction.
233 */
234 void setCoeffFriction(double mu) { coeffFriction_ = mu; }
235
236 /**
237 * \brief Get the coefficient of static friction.
238 * \return The coefficient of static friction.
239 */
240 [[nodiscard]] double getCoeffStaticFriction() const {
242 }
243
244 /**
245 * \brief Set the coefficient of static friction.
246 * \param mu_s The new coefficient of static friction.
247 */
248 void setCoeffStaticFriction(double mu_s) {
250 }
251
252 /**
253 * \brief Get the linear damping factor.
254 * \return The linear damping factor.
255 */
256 [[nodiscard]] double getLinearDamping() const { return linearDamping_; }
257
258 /**
259 * \brief Set the linear damping factor.
260 * \param damping The new linear damping factor (typically in [0, 1]).
261 */
262 void setLinearDamping(double damping) { linearDamping_ = damping; }
263
264 /**
265 * \brief Get the angular damping factor.
266 * \return The angular damping factor.
267 */
268 [[nodiscard]] double getAngularDamping() const {
269 return angularDamping_;
270 }
271
272 /**
273 * \brief Set the angular damping factor.
274 * \param damping The new angular damping factor (typically in [0, 1]).
275 */
276 void setAngularDamping(double damping) { angularDamping_ = damping; }
277
278 /**
279 * \brief Get the current time step.
280 * \return The current time step value.
281 */
282 [[nodiscard]] double getTimeStep() const { return dt_; }
283
284 /**
285 * \brief Set the current time step.
286 * \param dt The new time step value.
287 */
288 void setTimeStep(double dt) { dt_ = dt; }
289
290 /**
291 * \brief Get the positional correction factor.
292 * \return The positional correction factor.
293 */
294 [[nodiscard]] double getPositionalCorrectionFactor() const {
296 }
297
298 /**
299 * \brief Set the positional correction factor.
300 * \param factor The new positional correction factor.
301 */
302 void setPositionalCorrectionFactor(double factor) {
304 }
305
306 /**
307 * \brief Check if diagnostics are enabled.
308 * \return True if diagnostics are enabled, false otherwise.
309 */
310 [[nodiscard]] bool isDiagnosticsEnabled() const {
311 return diagnosticsEnabled_;
312 }
313
314 /**
315 * \brief Enable or disable diagnostics.
316 * \param enabled True to enable diagnostics, false to disable.
317 */
318 void setDiagnosticsEnabled(bool enabled) {
319 diagnosticsEnabled_ = enabled;
320 }
321
322 // --- Test API (protected but accessible via friend) ---
323 /**
324 * \brief Provides direct access to the Bodies container for testing.
325 * \return A reference to the internal Bodies object.
326 * \note This method is public but intended for testing only.
327 * Access is controlled via the friend declaration.
328 */
329 [[nodiscard]] Bodies& getBodiesForTest() { return bodies_; }
330
331 /**
332 * \brief Factory method to create a CollisionContext if two bodies are
333 * colliding.
334 *
335 * This method is part of the internal collision handling API and is
336 * exposed to tests for unit testing collision logic.
337 *
338 * \param b1 First body.
339 * \param b2 Second body.
340 * \return A CollisionContext object if a collision occurs, otherwise
341 * `std::nullopt`.
342 * \note This method is public but intended for testing only.
343 */
344 [[nodiscard]] std::optional<CollisionContext> createCollisionContext(
345 BodyProxy& b1, BodyProxy& b2);
346
347 /**
348 * \brief Calculates and applies the normal impulse (restitution) for a
349 * collision.
350 *
351 * This method is part of the internal collision handling API and is
352 * exposed to tests for unit testing collision logic.
353 *
354 * \param ctx The collision context.
355 * \return A tuple containing the total impulse vector, its magnitude,
356 * and the magnitude of the bias impulse component.
357 * \note This method is public but intended for testing only.
358 */
359 [[nodiscard]] std::tuple<Vector3d, double, double>
361
362 /**
363 * \brief Calculates the tangential (frictional) impulse for a
364 * collision.
365 *
366 * This method is part of the internal collision handling API and is
367 * exposed to tests for unit testing collision logic.
368 *
369 * \param ctx The collision context.
370 * \param jn The magnitude of the normal impulse, used for the friction
371 * limit.
372 * \param impulseN The normal impulse vector.
373 * \return The friction impulse vector.
374 * \note This method is public but intended for testing only.
375 */
376 [[nodiscard]] Vector3d applyFrictionImpulse(const CollisionContext& ctx,
377 double jn,
378 const Vector3d& impulseN);
379
380 /**
381 * \brief A greedy graph colouring algorithm to assign a colour
382 * (integer) to each body such that no two adjacent bodies (i.e.,
383 * colliding bodies) share the same colour.
384 *
385 * This method implements a parallel greedy graph colouring algorithm
386 * used to partition collision pairs into independent sets. Bodies with
387 * the same colour can have their collisions resolved in parallel
388 * without conflicts.
389 *
390 * \par Algorithm Overview
391 * The algorithm uses an iterative approach with two passes per
392 * iteration:
393 * - <b>Pass 1</b>: Each body (in parallel) assigns itself the
394 * smallest available colour not used by any of its neighbours.
395 * - <b>Pass 2</b>: Detect conflicts where adjacent bodies share the
396 * same colour (can occur due to race conditions in parallel
397 * execution).
398 *
399 * The algorithm iterates until no conflicts remain.
400 *
401 * \par Example
402 * Consider a collision graph with 4 bodies:
403 * \code
404 * graph[0] = {1, 2} // Body 0 collides with bodies 1 and 2
405 * graph[1] = {0, 3} // Body 1 collides with bodies 0 and 3
406 * graph[2] = {0} // Body 2 collides with body 0
407 * graph[3] = {1} // Body 3 collides with body 1
408 * \endcode
409 *
410 * After colouring, a possible result is:
411 * \code
412 * colors[0] = 0 // Bodies 0 and 3 can be processed in parallel
413 * colors[1] = 1 // Body 1 can be processed alone
414 * colors[2] = 1 // Bodies 1 and 2 can be processed in parallel
415 * colors[3] = 0 // (same as body 0)
416 * \endcode
417 *
418 * \par Thread Safety
419 * This method uses OpenMP locks to synchronize access to the shared
420 * `colors_` vector. Reads of neighbour colours are protected by locks
421 * to avoid data races, while writes are minimal and lock-protected.
422 *
423 * \param graph The collision graph where `graph[i]` contains
424 * indices of bodies colliding with `i`.
425 * The graph should be undirected (if body i collides with
426 * j, then j should appear in graph[i] and i in graph[j]).
427 *
428 * \note This method is public but intended for testing only.
429 * \note The algorithm is guaranteed to terminate as conflicts decrease
430 * with each iteration (worst-case: O(N) iterations for N bodies).
431 */
432 void colorGraph(const std::vector<std::vector<std::size_t>>& graph);
433
434 private:
435 /**
436 * \brief Handles collision detection and response between two
437 * bodies.
438 *
439 * \param b1 A proxy to the first body in the collision.
440 * \param b2 A proxy to the second body in the collision.
441 */
442 void handleCollision(BodyProxy& b1, BodyProxy& b2);
443
444 /**
445 * \brief Initialises the bodies in the simulation with random
446 * properties.
447 *
448 * \param n The number of bodies to create.
449 * \param dens The density of the bodies.
450 * \param seed The seed for random number generation.
451 */
452 void initializeBodies(std::size_t n, double dens, unsigned int seed);
453
454 /**
455 * \brief Resolves interpenetration by moving bodies apart along the
456 * collision normal.
457 */
459
460 /**
461 * \brief Calculates the effective mass of two colliding bodies in a
462 * given direction.
463 *
464 * \param ctx The collision context.
465 * \param direction_vec The direction vector (e.g., normal or
466 * tangent). \return The scalar effective mass.
467 */
468 [[nodiscard]] double getEffectiveMass(const CollisionContext& ctx,
469 const Vector3d& direction_vec);
470
471 // --- Member Variables ---
472 // --- Simulation Stages ---
473 /**
474 * \brief Logs system energy and checks for instability if diagnostics
475 * are enabled.
476 * \param iteration The current simulation iteration number.
477 */
478 void logSystemEnergy(std::size_t iteration);
479
480 /**
481 * \brief Builds the collision graph using broad-phase and
482 * narrow-phase detection.
483 *
484 * This method implements a two-phase collision detection algorithm
485 * to efficiently identify all colliding body pairs in the simulation.
486 * The result is stored in `collisionGraph_`, where `collisionGraph_[i]`
487 * contains the indices of all bodies colliding with body `i`.
488 *
489 * \par Algorithm Overview
490 * The algorithm uses a two-phase approach:
491 *
492 * <b>Phase 1: Broad-Phase Detection (Spatial Acceleration)</b>
493 * - Rebuilds the k-d tree with current body positions
494 * - For each body, performs a radius search to find all bodies within
495 * a safe distance (body radius + maximum body radius in system)
496 * - Uses thread-local storage to collect candidate pairs in parallel
497 * - Complexity: O(N log N) where N is the number of bodies
498 *
499 * <b>Phase 2: Narrow-Phase Detection (Exact Collision Test)</b>
500 * - For each candidate pair from Phase 1, performs an exact collision
501 * test using distance and radius comparison
502 * - Only pairs that actually overlap are added to the collision graph
503 * - Complexity: O(C) where C is the number of candidate pairs
504 *
505 * \par Example
506 * Consider a system with 3 bodies:
507 * \code
508 * Body 0: position (0, 0, 0), radius 1.0
509 * Body 1: position (1.5, 0, 0), radius 1.0 // Overlaps with body 0
510 * Body 2: position (5, 0, 0), radius 1.0 // No collision
511 * \endcode
512 *
513 * After execution:
514 * \code
515 * collisionGraph_[0] = {1} // Body 0 collides with body 1
516 * collisionGraph_[1] = {0} // Body 1 collides with body 0
517 * collisionGraph_[2] = {} // Body 2 has no collisions
518 * \endcode
519 *
520 * \par Thread Safety
521 * This method is fully parallelized using OpenMP:
522 * - Broad-phase searches are performed in parallel using thread-local
523 * storage to avoid race conditions
524 * - Candidate pairs are collected per-thread and merged sequentially
525 * after parallel execution
526 * - The final collision graph is built sequentially from thread-local
527 * results
528 *
529 * \par Performance Optimizations
530 * - Maximum radius is cached and only recalculated when bodies are
531 * added/removed (not every frame)
532 * - Results vectors are pre-allocated to avoid repeated allocations
533 * - k-d tree is rebuilt each frame to reflect position changes
534 * - Thread-local storage minimizes lock contention
535 *
536 * \note This method is called once per simulation step in
537 * `computeDynamics()`.
538 * \note The collision graph is undirected: if body i collides with j,
539 * then j appears in `collisionGraph_[i]` and i appears in
540 * `collisionGraph_[j]`.
541 */
542 void buildCollisionGraph();
543
544 /// \brief Resolves all detected collisions using graph colouring.
545 void resolveCollisions();
546
547 /// \brief Calculates and applies gravitational forces to all bodies.
548 void applyGravity();
549
550 /**
551 * \brief Determines the optimal time step for the next frame based on
552 * current velocities and accelerations.
553 */
555
556 /**
557 * \brief Applies linear and angular damping to all bodies.
558 *
559 * This method applies velocity damping to reduce energy over time,
560 * which helps stabilize the simulation and prevent numerical
561 * explosions.
562 */
563 void applyDamping();
564
565 /// \brief Time step (in s).
566 double dt_;
567
568 /// \brief Previous time step (in s).
570
571 /// \brief Universal gravitational constant (in m³ / kg /s²).
572 double G_;
573
574 /// \brief Numerical precision threshold.
575 double epsilon_;
576
577 /// \brief Coefficient of restitution for collisions.
579
580 /// \brief Coefficient of kinetic (sliding) friction.
582
583 /// \brief Coefficient of static friction.
585
586 /// \brief The percentage of overlap to correct in each frame.
588
589 /// \brief Damping factor for linear velocity.
591
592 /// \brief Damping factor for angular velocity.
594
595 /// \brief Flag indicating if diagnostic output is enabled.
597
598 /// \brief Frequency of energy log output (every N iterations).
599 std::size_t logFreq_;
600
601 /// \brief SoA container for all bodies in the simulation.
602 Bodies bodies_; // NOLINT
603
604 /// \brief k-d tree for broad-phase collision detection.
605 std::unique_ptr<KDTree> kdTree_; // NOLINT
606
607 /// \brief Adjacency list for the collision graph.
608 std::vector<std::vector<std::size_t>> collisionGraph_; // NOLINT
609
610 /**
611 * \brief Stores the colour of each body for parallel collision
612 * processing.
613 */
614 std::vector<int> colors_; // NOLINT
615
616 /// \brief OpenMP locks for thread-safe collision graph updates.
617 std::vector<omp_lock_t> graphLocks_;
618
619 /**
620 * \brief Thread-local storage for collision pairs found during
621 * broad-phase. This reduces lock contention on the global
622 * `collisionGraph_` during parallel processing.
623 */
624 std::vector<std::vector<std::pair<std::size_t, std::size_t>>>
626
627 /// \brief The total energy from the previous logged step.
628 double lastTotalEnergy_ = std::numeric_limits<double>::infinity();
629
630 /// \brief Cached maximum radius for collision detection optimization.
631 /// This is recalculated only when bodies are added/removed, not every
632 /// frame.
633 mutable double cached_max_radius_ = 0.0;
634
635 // --- Profiling Timers ---
636 /* These members accumulate the time spent in different simulation
637 stages. */
638 /// \brief Time spent in integration steps.
639 mutable double profIntegration_ = 0.0;
640
641 /// \brief Time spent in collision graph building.
642 mutable double profCollisionGraph_ = 0.0;
643
644 /// \brief Time spent in collision response.
645 mutable double profCollisionResponse_ = 0.0;
646
647 /// \brief Time spent in force calculation.
648 mutable double profForceCalculation_ = 0.0;
649 };
650}
651
652#endif // SPACE_HPP
SoA container for simulation bodies and proxies for AoS-like access.
Structure-of-Arrays (SoA) container for all bodies in the simulation.
Definition: body.hpp:274
std::size_t size() const
Get number of bodies in the simulation.
Definition: body.hpp:371
A proxy object that provides an AoS-like interface to a body stored in the Bodies SoA container.
Definition: body.hpp:141
A const proxy object that provides a read-only AoS-like interface to a body stored in the Bodies SoA ...
Definition: body.hpp:252
Class describing a space in which move several bodies.
Definition: space.hpp:89
Vector3d applyFrictionImpulse(const CollisionContext &ctx, double jn, const Vector3d &impulseN)
Calculates the tangential (frictional) impulse for a collision.
Definition: space.cpp:548
double G_
Universal gravitational constant (in m³ / kg /s²).
Definition: space.hpp:572
double profCollisionResponse_
Time spent in collision response.
Definition: space.hpp:645
~Space()
Destructor to clean up OpenMP locks.
Definition: space.cpp:202
double previous_dt_
Previous time step (in s).
Definition: space.hpp:569
double getTimeStep() const
Get the current time step.
Definition: space.hpp:282
double coeffStaticFriction_
Coefficient of static friction.
Definition: space.hpp:584
void initializeBodies(std::size_t n, double dens, unsigned int seed)
Initialises the bodies in the simulation with random properties.
Definition: space.cpp:209
void applyGravity()
Calculates and applies gravitational forces to all bodies.
Definition: space.cpp:1085
std::tuple< Vector3dVec, QuaterniondVec, Vector3dVec, Vector3dVec > getDataForDisplay() const
Retrieves all data required for display for all bodies.
Definition: space.cpp:400
double getPositionalCorrectionFactor() const
Get the positional correction factor.
Definition: space.hpp:294
std::vector< omp_lock_t > graphLocks_
OpenMP locks for thread-safe collision graph updates.
Definition: space.hpp:617
double getCoeffFriction() const
Get the coefficient of kinetic friction.
Definition: space.hpp:228
BodyProxy body(std::size_t i)
Get a proxy to the i-th body.
Definition: space.hpp:123
void resolveInterpenetration(const CollisionContext &ctx)
Resolves interpenetration by moving bodies apart along the collision normal.
Definition: space.cpp:485
void setCoeffRestitution(double e)
Set the coefficient of restitution.
Definition: space.hpp:222
double getCoeffStaticFriction() const
Get the coefficient of static friction.
Definition: space.hpp:240
void setLinearDamping(double damping)
Set the linear damping factor.
Definition: space.hpp:262
double getLinearDamping() const
Get the linear damping factor.
Definition: space.hpp:256
double getEffectiveMass(const CollisionContext &ctx, const Vector3d &direction_vec)
Calculates the effective mass of two colliding bodies in a given direction.
Definition: space.cpp:507
void logSystemEnergy(std::size_t iteration)
Logs system energy and checks for instability if diagnostics are enabled.
Definition: space.cpp:788
void resolveCollisions()
Resolves all detected collisions using graph colouring.
Definition: space.cpp:966
void setPositionalCorrectionFactor(double factor)
Set the positional correction factor.
Definition: space.hpp:302
void setTimeStep(double dt)
Set the current time step.
Definition: space.hpp:288
void colorGraph(const std::vector< std::vector< std::size_t > > &graph)
A greedy graph colouring algorithm to assign a colour (integer) to each body such that no two adjacen...
Definition: space.cpp:700
void computeDynamics(std::size_t iteration)
Computes one full step of the simulation.
Definition: space.cpp:908
double lastTotalEnergy_
The total energy from the previous logged step.
Definition: space.hpp:628
double linearDamping_
Damping factor for linear velocity.
Definition: space.hpp:590
double getPreviousTimeStep() const
Get the time step used in the last completed simulation frame.
Definition: space.hpp:140
double getG() const
Get the universal gravitational constant.
Definition: space.hpp:190
bool isDiagnosticsEnabled() const
Check if diagnostics are enabled.
Definition: space.hpp:310
std::vector< std::vector< std::size_t > > collisionGraph_
Adjacency list for the collision graph.
Definition: space.hpp:608
double profCollisionGraph_
Time spent in collision graph building.
Definition: space.hpp:642
void setG(double g)
Set the universal gravitational constant.
Definition: space.hpp:196
std::size_t logFreq_
Frequency of energy log output (every N iterations).
Definition: space.hpp:599
std::size_t n() const
Number of bodies getter.
Definition: space.hpp:114
void setCoeffStaticFriction(double mu_s)
Set the coefficient of static friction.
Definition: space.hpp:248
Space()
Default constructor.
Definition: space.hpp:96
Bodies bodies_
SoA container for all bodies in the simulation.
Definition: space.hpp:602
double profIntegration_
Time spent in integration steps.
Definition: space.hpp:639
bool diagnosticsEnabled_
Flag indicating if diagnostic output is enabled.
Definition: space.hpp:596
void setCoeffFriction(double mu)
Set the coefficient of kinetic friction.
Definition: space.hpp:234
double positionalCorrectionFactor_
The percentage of overlap to correct in each frame.
Definition: space.hpp:587
Bodies & getBodiesForTest()
Provides direct access to the Bodies container for testing.
Definition: space.hpp:329
void setEpsilon(double eps)
Set the numerical precision threshold.
Definition: space.hpp:208
ConstBodyProxy body(std::size_t i) const
Get a const proxy to the i-th body.
Definition: space.hpp:132
std::optional< CollisionContext > createCollisionContext(BodyProxy &b1, BodyProxy &b2)
Factory method to create a CollisionContext if two bodies are colliding.
Definition: space.cpp:642
double cached_max_radius_
Cached maximum radius for collision detection optimization. This is recalculated only when bodies are...
Definition: space.hpp:633
double angularDamping_
Damping factor for angular velocity.
Definition: space.hpp:593
double dt_
Time step (in s).
Definition: space.hpp:566
double profForceCalculation_
Time spent in force calculation.
Definition: space.hpp:648
void applyDamping()
Applies linear and angular damping to all bodies.
Definition: space.cpp:994
void updateAdaptiveTimeStep()
Determines the optimal time step for the next frame based on current velocities and accelerations.
Definition: space.cpp:1005
void setAngularDamping(double damping)
Set the angular damping factor.
Definition: space.hpp:276
double coeffRestitution_
Coefficient of restitution for collisions.
Definition: space.hpp:578
std::unique_ptr< KDTree > kdTree_
k-d tree for broad-phase collision detection.
Definition: space.hpp:605
std::tuple< Vector3d, double, double > applyRestitutionImpulse(CollisionContext &ctx)
Calculates and applies the normal impulse (restitution) for a collision.
Definition: space.cpp:523
std::vector< std::vector< std::pair< std::size_t, std::size_t > > > thread_local_collision_pairs_
Thread-local storage for collision pairs found during broad-phase. This reduces lock contention on th...
Definition: space.hpp:625
double getCoeffRestitution() const
Get the coefficient of restitution.
Definition: space.hpp:214
void setDiagnosticsEnabled(bool enabled)
Enable or disable diagnostics.
Definition: space.hpp:318
double coeffFriction_
Coefficient of kinetic (sliding) friction.
Definition: space.hpp:581
double epsilon_
Numerical precision threshold.
Definition: space.hpp:575
std::vector< int > colors_
Stores the colour of each body for parallel collision processing.
Definition: space.hpp:614
double getEpsilon() const
Get the numerical precision threshold.
Definition: space.hpp:202
void buildCollisionGraph()
Builds the collision graph using broad-phase and narrow-phase detection.
Definition: space.cpp:821
std::tuple< Vector3dVec, QuaterniondVec > getAllBodyTransforms() const
Retrieves the current positions and orientations of all bodies.
Definition: space.cpp:384
double getAngularDamping() const
Get the angular damping factor.
Definition: space.hpp:268
void handleCollision(BodyProxy &b1, BodyProxy &b2)
Handles collision detection and response between two bodies.
Definition: space.cpp:664
std::tuple< double, double, double > calculateSystemEnergy() const
Calculates the kinetic and potential energy of the system.
Definition: space.cpp:443
Constants for the model.
constexpr double EPSILON
Default computing precision value.
Definition: constants.hpp:32
constexpr double G
Default universal gravitational constant (in m³⋅kg⁻¹⋅s⁻²).
Definition: constants.hpp:35
A k-d tree wrapper for broad-phase collision detection using nanoflann.
A data container to hold all relevant information for a single collision event.
Definition: space.hpp:67
Vector3d vRel_
The relative velocity between the two bodies at the contact point.
Definition: space.hpp:80
BodyProxy & b1
A proxy to the first body in the collision.
Definition: space.hpp:69
BodyProxy & b2
A proxy to the second body in the collision.
Definition: space.hpp:71
Vector3d r1Vec_
The vector from the centre of body 1 to the contact point.
Definition: space.hpp:75
Vector3d r2Vec_
The vector from the centre of body 2 to the contact point.
Definition: space.hpp:77
Vector3d nVec_
The normalised vector pointing from body 1 to body 2.
Definition: space.hpp:73
double overlap_
The penetration depth of the two bodies.
Definition: space.hpp:82