MPS-Basic
Loading...
Searching...
No Matches
mps.cpp
Go to the documentation of this file.
1#include "mps.hpp"
2
3#include "particle.hpp"
4#include "weight.hpp"
5
6#include <queue>
7
8// This include is for checking if the pressure calculator is explicit.
9// This is needed because we need to update the pressure again when using EMPS.
11
12using std::cerr;
13using std::endl;
14
16 const Input& input,
17 std::unique_ptr<PressureCalculator::Interface>&& pressureCalculator,
18 std::unique_ptr<SurfaceDetector::Interface>&& surfaceDetector
19) {
20 this->settings = input.settings;
21 this->domain = input.settings.domain;
22 this->particles = input.particles;
23 this->pressureCalculator = std::move(pressureCalculator);
24 this->surfaceDetector = std::move(surfaceDetector);
26
30}
31
34 calGravity();
37
39 collision();
40
43 auto pressures = pressureCalculator->calc(particles);
44 for (auto& particle : particles) {
45 particle.pressure = pressures[particle.id];
46 }
47
51
52 // Update pressure again when using EMPS
53 if (auto explicitPressureCalculator = dynamic_cast<PressureCalculator::Explicit*>(pressureCalculator.get())) {
54 auto pressures = explicitPressureCalculator->calc(particles);
55 for (auto& particle : particles) {
56 particle.pressure = pressures[particle.id];
57 }
58 }
59
60 calCourant();
61}
62
64#pragma omp parallel for
65 for (auto& p : particles) {
66 if (p.type == ParticleType::Fluid) {
67 p.acceleration += settings.gravity;
68
69 } else {
70 p.acceleration.setZero();
71 }
72 }
73}
74
75void MPS::calViscosity(const double& re) {
76 double n0 = refValuesForLaplacian.n0;
77 double lambda = refValuesForLaplacian.lambda;
78 double a = (settings.kinematicViscosity) * (2.0 * settings.dim) / (n0 * lambda);
79
80#pragma omp parallel for
81 for (auto& pi : particles) {
82 if (pi.type != ParticleType::Fluid)
83 continue;
84
85 Eigen::Vector3d viscosityTerm = Eigen::Vector3d::Zero();
86
87 for (auto& neighbor : pi.neighbors) {
88 auto& pj = particles[neighbor.id];
89
90 if (neighbor.distance < settings.re_forLaplacian) {
91 double w = weight(neighbor.distance, re);
92 viscosityTerm += (pj.velocity - pi.velocity) * w;
93 }
94 }
95
96 viscosityTerm *= a;
97 pi.acceleration += viscosityTerm;
98 }
99}
100
102#pragma omp parallel for
103 for (auto& p : particles) {
104 if (p.type == ParticleType::Fluid) {
105 p.velocity += p.acceleration * settings.dt;
106 p.position += p.velocity * settings.dt;
107 }
108 p.acceleration.setZero();
109 }
110}
111
113 for (auto& pi : particles) {
114 if (pi.type != ParticleType::Fluid)
115 continue;
116
117 for (auto& neighbor : pi.neighbors) {
118 Particle& pj = particles[neighbor.id];
119 if (pj.type == ParticleType::Fluid && pj.id >= pi.id)
120 continue;
121
122 if (neighbor.distance < settings.collisionDistance) {
123
124 double invMi = pi.inverseDensity();
125 double invMj = pj.inverseDensity();
126 double mass = 1.0 / (invMi + invMj);
127
128 Eigen::Vector3d normal = (pj.position - pi.position).normalized();
129 double relVel = (pj.velocity - pi.velocity).dot(normal);
130 double impulse = 0.0;
131 if (relVel < 0.0)
132 impulse = -(1 + settings.coefficientOfRestitution) * relVel * mass;
133 pi.velocity -= impulse * invMi * normal;
134 pj.velocity += impulse * invMj * normal;
135
136 double depth = settings.collisionDistance - neighbor.distance;
137 double positionImpulse = depth * mass;
138 pi.position -= positionImpulse * invMi * normal;
139 pj.position += positionImpulse * invMj * normal;
140
141 // cerr << "WARNING: Collision between particles " << pi.id << " and " << pj.id << " occurred."
142 // << endl;
143 }
144 }
145 }
146}
147
148void MPS::calNumberDensity(const double& re) {
149#pragma omp parallel for
150 for (auto& pi : particles) {
151 pi.numberDensity = 0.0;
152
153 if (pi.type == ParticleType::Ghost)
154 continue;
155
156 for (auto& neighbor : pi.neighbors)
157 pi.numberDensity += weight(neighbor.distance, re);
158 }
159}
160
162#pragma omp parallel for
163 for (auto& pi : particles) {
164 if (pi.type == ParticleType::Ghost || pi.type == ParticleType::DummyWall) {
165 pi.boundaryCondition = FluidState::Ignored;
166
167 } else { // Fluid particles
168 if (surfaceDetector->isFreeSurface(particles, pi)) {
169 pi.boundaryCondition = FluidState::FreeSurface;
170 } else {
171 pi.boundaryCondition = FluidState::Inner;
172 }
173 }
174 }
175}
176
177void MPS::setMinimumPressure(const double& re) {
178#pragma omp parallel for
179 for (auto& p : particles) {
180 p.minimumPressure = p.pressure;
181 }
182
183 for (auto& pi : particles) {
184 if (pi.type == ParticleType::Ghost || pi.type == ParticleType::DummyWall)
185 continue;
186
187 for (auto& neighbor : pi.neighbors) {
188 Particle& pj = particles[neighbor.id];
190 continue;
191 if (pj.id > pi.id)
192 continue;
193
194 if (neighbor.distance < re) {
195 pi.minimumPressure = std::min(pi.minimumPressure, pj.pressure);
196 pj.minimumPressure = std::min(pj.minimumPressure, pi.pressure);
197 }
198 }
199 }
200}
201
202void MPS::calPressureGradient(const double& re) {
203 double a = settings.dim / refValuesForGradient.n0;
204
205#pragma omp parallel for
206 for (auto& pi : particles) {
207 if (pi.type != ParticleType::Fluid)
208 continue;
209
210 Eigen::Vector3d grad = Eigen::Vector3d::Zero();
211 for (auto& neighbor : pi.neighbors) {
212 Particle& pj = particles[neighbor.id];
214 continue;
215
216 if (neighbor.distance < re) {
217 double w = weight(neighbor.distance, re);
218 // double dist2 = pow(neighbor.distance, 2);
219 double dist2 = (pj.position - pi.position).squaredNorm();
220 double pij = (pj.pressure - pi.minimumPressure) / dist2;
221 grad += (pj.position - pi.position) * pij * w;
222 }
223 }
224 grad *= a;
225 pi.acceleration -= grad / pi.density;
226 }
227}
228
230#pragma omp parallel for
231 for (auto&& p : particles) {
232 if (p.type == ParticleType::Fluid) {
233 p.velocity += p.acceleration * settings.dt;
234 p.position += p.acceleration * settings.dt * settings.dt;
235 }
236
237 p.acceleration.setZero();
238 }
239}
240
242 courant = 0.0;
243
244 for (auto& pi : particles) {
245 if (pi.type != ParticleType::Fluid)
246 continue;
247
248 double iCourant = (pi.velocity.norm() * settings.dt) / settings.particleDistance;
249 courant = std::max(courant, iCourant);
250 }
251
253 cerr << "ERROR: Courant number is larger than CFL condition. Courant = " << courant << endl;
254 }
255}
void setMinimumPressure(const double &re)
set minimum pressure for pressure gradient calculation
Definition mps.cpp:177
void stepForward()
Definition mps.cpp:32
void calViscosity(const double &re)
calculate viscosity term of Navier-Stokes equation
Definition mps.cpp:75
std::unique_ptr< PressureCalculator::Interface > pressureCalculator
Interface for pressure calculation.
Definition mps.hpp:34
void moveParticleUsingPressureGradient()
move particles in correction step
Definition mps.cpp:229
Settings settings
Settings for the simulation.
Definition mps.hpp:27
double courant
Maximum courant number among all particles.
Definition mps.hpp:36
void calPressureGradient(const double &re)
calculate pressure gradient term
Definition mps.cpp:202
MPS()=default
RefValues refValuesForLaplacian
Reference values for the simulation ( , )
Definition mps.hpp:29
void moveParticle()
move particles in prediction step
Definition mps.cpp:101
void setBoundaryCondition()
set boundary condition of pressure Poisson equation
Definition mps.cpp:161
RefValues refValuesForNumberDensity
Reference values for the simulation ( , )
Definition mps.hpp:28
void calCourant()
calculate Courant number
Definition mps.cpp:241
void collision()
calculate collision between particles when they are too close
Definition mps.cpp:112
NeighborSearcher neighborSearcher
Neighbor searcher for neighbor search.
Definition mps.hpp:47
std::unique_ptr< SurfaceDetector::Interface > surfaceDetector
Interface for free surface detection.
Definition mps.hpp:48
void calNumberDensity(const double &re)
calculate number density of each particle
Definition mps.cpp:148
RefValues refValuesForGradient
Reference values for the simulation ( , )
Definition mps.hpp:30
void calGravity()
calculate gravity term
Definition mps.cpp:63
Particles particles
Particles in the simulation.
Definition mps.hpp:31
Domain domain
Domain of the simulation.
Definition mps.hpp:32
void setNeighbors(Particles &particles)
Class for particle in MPS method.
Definition particle.hpp:47
double inverseDensity() const
calculate inverse of density for collision
Definition particle.cpp:12
Eigen::Vector3d position
position of the particle
Definition particle.hpp:55
double pressure
pressure of the particle
Definition particle.hpp:58
int id
index of the particle
Definition particle.hpp:50
double minimumPressure
minimum pressure of the particle
Definition particle.hpp:64
Eigen::Vector3d velocity
velocity of the particle
Definition particle.hpp:56
ParticleType type
type of the particle
Definition particle.hpp:51
int size() const
Get the number of particles.
Definition particles.cpp:21
Class for explicit pressure calculation.
Definition explicit.hpp:14
Struct for reference values of MPS method.
Definition refvalues.hpp:8
double n0
reference value of number density for source term of pressure Poisson equation
Definition refvalues.hpp:10
double lambda
coefficient for laplacian
Definition refvalues.hpp:11
@ Inner
inner fluid particle
@ FreeSurface
free surface particle
@ Ignored
Ghost or dummy.
@ Ghost
Ghost particle (outside of the domain, not used for calculation)
@ DummyWall
Dummy wall particle (pressure is not calculated)
@ Fluid
Fluid particle.
Represents the input data for MPS simulation.
Definition input.hpp:12
Particles particles
Initial particles arrangement in the simulation.
Definition input.hpp:14
Settings settings
Settings for the simulation.
Definition input.hpp:13
Eigen::Vector3d gravity
Gravity.
Definition settings.hpp:32
double reMax
Maximum of effective radius.
Definition settings.hpp:57
double collisionDistance
Distance for collision detection.
Definition settings.hpp:50
double particleDistance
Initial distance between particles.
Definition settings.hpp:18
int dim
Dimension of the simulation.
Definition settings.hpp:17
double re_forNumberDensity
Effective radius for number density.
Definition settings.hpp:54
double kinematicViscosity
Kinematic viscosity.
Definition settings.hpp:28
double cflCondition
CFL condition.
Definition settings.hpp:22
double coefficientOfRestitution
Coefficient of restitution.
Definition settings.hpp:51
double re_forGradient
Effective radius for gradient.
Definition settings.hpp:55
double dt
Time step.
Definition settings.hpp:19
Domain domain
domain of the simulation
Definition settings.hpp:25
double re_forLaplacian
Effective radius for Laplacian.
Definition settings.hpp:56
double weight(double dis, double re)
Wight function for MPS method presented by Koshizuka and Oka, 1996
Definition weight.cpp:5