PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
Loading...
Searching...
No Matches
Integrator.cuh
Go to the documentation of this file.
1
10
11#ifndef PBAT_GPU_IMPL_VBD_INTEGRATOR_H
12#define PBAT_GPU_IMPL_VBD_INTEGRATOR_H
13
14#include "pbat/gpu/Aliases.h"
15#include "pbat/gpu/impl/common/Buffer.cuh"
16#include "pbat/gpu/impl/contact/VertexTriangleMixedCcdDcd.cuh"
17#include "pbat/sim/vbd/Data.h"
18#include "pbat/sim/vbd/Enums.h"
19// clang-format off
25#include "Kernels.cuh"
26// clang-format on
27
28#include <cuda/api/stream.hpp>
29#include <string_view>
30#include <vector>
31
32namespace pbat::gpu::impl::vbd {
33
39{
40 public:
44
50 Integrator(Data const& data);
57 void Step(GpuScalar dt, GpuIndex iterations, GpuIndex substeps = GpuIndex{1});
70 void TracedStep(
71 GpuScalar dt,
72 GpuIndex iterations,
73 GpuIndex substeps,
74 GpuIndex t,
75 std::string_view dir = ".");
82 Eigen::Vector<GpuScalar, 3> const& min,
83 Eigen::Vector<GpuScalar, 3> const& max);
88 void SetBlockSize(GpuIndex blockSize);
89
113 void UpdateActiveSet();
122 virtual void Solve(kernels::BackwardEulerMinimization& bdf, GpuIndex iterations);
134 void TracedSolve(
136 GpuIndex iterations,
137 GpuIndex t,
138 GpuIndex s,
139 std::string_view dir);
152 void UpdateBdfState(GpuScalar sdt);
153
154 public:
157
158 Eigen::Vector<GpuScalar, 3> mWorldMin;
159 Eigen::Vector<GpuScalar, 3> mWorldMax;
169
173
178
185
192
198
202
206 cuda::stream_t mStream;
207
215};
216
217} // namespace pbat::gpu::impl::vbd
218
219#endif // PBAT_GPU_IMPL_VBD_INTEGRATOR_H
VBD kernels.
Definition Buffer.cuh:21
Vertex-triangle contact detection using a mixed CCD/DCD approach.
Definition VertexTriangleMixedCcdDcd.cuh:18
void UpdateActiveSet()
Computes active contact pairs in the current configuration.
Definition Integrator.cu:250
virtual void Solve(kernels::BackwardEulerMinimization &bdf, GpuIndex iterations)
VBD to iterate on the BDF minimization problem.
Definition Integrator.cu:275
EInitializationStrategy mInitializationStrategy
Strategy to use to determine the initial BCD iterate.
Definition Integrator.cuh:204
common::Buffer< GpuScalar, 3 > mInertialTargetPositions
Inertial target for vertex positions.
Definition Integrator.cuh:171
common::Buffer< GpuScalar, 3 > mVelocities
Current vertex velocities.
Definition Integrator.cuh:175
GpuScalar mSmoothFrictionRelativeVelocityThreshold
Definition Integrator.cuh:196
GpuScalar mCollisionPenalty
Collision penalty coefficient.
Definition Integrator.cuh:194
void SetSceneBoundingBox(Eigen::Vector< GpuScalar, 3 > const &min, Eigen::Vector< GpuScalar, 3 > const &max)
Set the bounding box over the scene used for spatial partitioning.
Definition Integrator.cu:150
common::Buffer< GpuScalar > mShapeFunctionGradients
4x3x|# elements| shape function gradients
Definition Integrator.cuh:182
contact::VertexTriangleMixedCcdDcd cd
Collision detector.
Definition Integrator.cuh:161
common::Buffer< GpuIndex > mVertexTetrahedronPrefix
Vertex-tetrahedron adjacency list's prefix sum.
Definition Integrator.cuh:187
common::Buffer< GpuIndex > mPadj
Partition vertices.
Definition Integrator.cuh:201
void ComputeInertialTargets(GpuScalar sdt, GpuScalar sdt2)
Compute the inertial target positions of the BDF1 objective.
Definition Integrator.cu:190
common::Buffer< GpuScalar, 3 > x
Vertex positions.
Definition Integrator.cuh:155
void InitializeActiveSet(GpuScalar dt)
Detect candidate contact pairs for the current time step.
Definition Integrator.cu:163
GpuIndex mActiveSetUpdateFrequency
Active set update frequency.
Definition Integrator.cuh:160
cuda::stream_t mStream
Cuda stream on which this VBD instance will run.
Definition Integrator.cuh:206
void UpdateBdfState(GpuScalar sdt)
Update the BDF state (i.e. positions and velocities) after solving the BDF minimization.
Definition Integrator.cu:329
Eigen::Vector< GpuScalar, 3 > mWorldMin
World AABB min.
Definition Integrator.cuh:158
void RunVbdIteration(kernels::BackwardEulerMinimization &bdf)
Run a single iteration of the VBD's BDF minimization.
Definition Integrator.cu:303
common::Buffer< GpuScalar, 3 > mVelocitiesAtT
Previous vertex velocities.
Definition Integrator.cuh:174
common::Buffer< GpuScalar > XVA
|x.Size()| array of vertex areas for contact response
Definition Integrator.cuh:167
common::Buffer< GpuScalar, 3 > xb
Write buffer for positions which handles data races.
Definition Integrator.cuh:172
common::Buffer< GpuScalar > mQuadratureWeights
Definition Integrator.cuh:179
common::Buffer< GpuScalar > FA
|F.Size()| array of triangle areas for contact response
Definition Integrator.cuh:168
Integrator(Data const &data)
Construct Integrator from data.
Definition Integrator.cu:22
common::Buffer< GpuIndex, 4 > T
Tetrahedral mesh elements.
Definition Integrator.cuh:156
GpuScalar mDetHZero
Numerical zero for hessian determinant check.
Definition Integrator.cuh:184
void TracedStep(GpuScalar dt, GpuIndex iterations, GpuIndex substeps, GpuIndex t, std::string_view dir=".")
Execute one simulation step with tracing to disk.
Definition Integrator.cu:105
common::Buffer< GpuIndex > fc
Definition Integrator.cuh:163
GpuIndexVectorX mPptr
Definition Integrator.cuh:199
void Step(GpuScalar dt, GpuIndex iterations, GpuIndex substeps=GpuIndex{1})
Execute one simulation step.
Definition Integrator.cu:82
common::Buffer< GpuScalar, 3 > mExternalAcceleration
External acceleration.
Definition Integrator.cuh:176
common::Buffer< GpuIndex > mVertexTetrahedronLocalVertexIndices
Definition Integrator.cuh:190
common::Buffer< GpuScalar > mMass
Lumped mass matrix diagonals.
Definition Integrator.cuh:177
common::Buffer< GpuIndex > mVertexTetrahedronNeighbours
Vertex-tetrahedron adjacency list's neighbour list.
Definition Integrator.cuh:189
pbat::sim::vbd::Data Data
Data type for VBD.
Definition Integrator.cuh:43
kernels::BackwardEulerMinimization BdfDeviceParameters(GpuScalar dt, GpuScalar dt2)
Obtain Backward Euler time stepping minimization problem for device code.
Definition Integrator.cu:349
void TracedSolve(kernels::BackwardEulerMinimization &bdf, GpuIndex iterations, GpuIndex t, GpuIndex s, std::string_view dir)
VBD to iterate on the BDF minimization problem.
Definition Integrator.cu:284
GpuIndex mGpuThreadBlockSize
Number of threads per CUDA thread block.
Definition Integrator.cuh:205
common::Buffer< GpuScalar > mLameCoefficients
2x|# elements| 1st and 2nd Lame parameters
Definition Integrator.cuh:183
Eigen::Vector< GpuScalar, 3 > mWorldMax
World AABB max.
Definition Integrator.cuh:159
pbat::sim::vbd::EInitializationStrategy EInitializationStrategy
Enum type for initialization strategy.
Definition Integrator.cuh:41
void SetBlockSize(GpuIndex blockSize)
Set the number of threads per GPU block for per-vertex elastic energy computation.
Definition Integrator.cu:158
common::Buffer< GpuScalar, 3 > mPositionsAtT
Previous vertex positions.
Definition Integrator.cuh:170
GpuScalar mRayleighDamping
Rayleigh damping coefficient.
Definition Integrator.cuh:193
GpuScalar mFrictionCoefficient
Coefficient of friction.
Definition Integrator.cuh:195
void InitializeBcdSolution(GpuScalar sdt, GpuScalar sdt2)
Compute starting point of BCD minimization.
Definition Integrator.cu:218
Type aliases for GPU code.
Vertex Block Descent (VBD) algorithms.
Definition AndersonIntegrator.cu:13
EInitializationStrategy
Initialization strategies for the VBD time step minimization.
Definition Enums.h:9
float GpuScalar
Scalar type for GPU code.
Definition Aliases.h:19
Eigen::Vector< GpuIndex, Eigen::Dynamic > GpuIndexVectorX
Index vector type for GPU code.
Definition Aliases.h:28
std::int32_t GpuIndex
Index type for GPU code.
Definition Aliases.h:20
Device-side BFD1 minimization problem.
Definition Kernels.cuh:40
VBD simulation configuration.
Definition Data.h:15