Skip to content

Commit 8179611

Browse files
committed
introducing velocity profiles for circular emitters
1 parent 7bf5480 commit 8179611

File tree

8 files changed

+45
-16
lines changed

8 files changed

+45
-16
lines changed

SPlisHSPlasH/Emitter.cpp

Lines changed: 22 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -10,9 +10,10 @@
1010
using namespace SPH;
1111

1212
Emitter::Emitter(FluidModel* model, const unsigned int width, const unsigned int height, const Vector3r& pos,
13-
const Matrix3r& rotation, const Real velocity, const unsigned int type, const bool useBoundary)
13+
const Matrix3r& rotation, const Real velocity, const unsigned int type, const bool useBoundary,
14+
const unsigned int velocityProfile)
1415
: m_model(model), m_width(width), m_x(pos), m_rotation(rotation), m_velocity(velocity), m_type(type),
15-
m_useBoundary(useBoundary)
16+
m_useBoundary(useBoundary), m_velocityProfile(velocityProfile)
1617
{
1718
Simulation* sim = Simulation::getCurrent();
1819
m_depth = getDepth();
@@ -101,7 +102,7 @@ void Emitter::emitParticles(std::vector<unsigned int>& reusedParticles, unsigned
101102
const Vector3r axisDepth = m_rotation.col(0);
102103
const Vector3r axisHeight = m_rotation.col(1);
103104
const Vector3r axisWidth = m_rotation.col(2);
104-
Vector3r emitVel = m_velocity * axisDepth;
105+
Vector3r bulkEmitVel = m_velocity * axisDepth;
105106
Simulation* sim = Simulation::getCurrent();
106107
const Real particleRadius = sim->getParticleRadius();
107108
const Real particleDiameter = static_cast<Real>(2.0) * particleRadius;
@@ -128,7 +129,6 @@ void Emitter::emitParticles(std::vector<unsigned int>& reusedParticles, unsigned
128129
{
129130
m_model->setPosition(indexNextNewParticle, spawnPos);
130131
m_model->setParticleState(indexNextNewParticle, ParticleState::AnimatedByEmitter);
131-
m_model->setVelocity(indexNextNewParticle, emitVel);
132132
m_model->setObjectId(indexNextNewParticle, m_objectId);
133133

134134
indexNextNewParticle++;
@@ -149,9 +149,24 @@ void Emitter::emitParticles(std::vector<unsigned int>& reusedParticles, unsigned
149149

150150
if (insideEmitter) {
151151
// Advect ALL particles inside the emitter.
152-
// TODO: Doesn't get all particles within the boundary. Perhaps use getSizeExtraMargin()?
153-
m_model->setPosition(i, tempPos + timeStepSize * emitVel);
154-
m_model->setVelocity(i, emitVel);
152+
//
153+
// TODO: Doesn't get all particles within the rigid body. Perhaps use getSizeExtraMargin()?
154+
Vector3r localEmitVel;
155+
if (m_type == 1 && m_velocityProfile != 0) {
156+
// different velocity profiles for circular emitters
157+
const Vector3r localPos = tempPos - m_x;
158+
const Real r = sqrt(pow(axisWidth.dot(localPos), 2.0) + (pow(axisHeight.dot(localPos), 2.0)));
159+
const Vector3r maxEmitVel =
160+
(1.0 / (1.0 - (2.0 / (static_cast<Real>(m_velocityProfile) + 2.0)))) * bulkEmitVel;
161+
const Real relativeRadius = r / (m_size[1] / 2.0);
162+
localEmitVel = maxEmitVel * (1.0 - pow(relativeRadius, static_cast<Real>(m_velocityProfile)));
163+
}
164+
else {// box shaped emitter or uniform velocity profile
165+
localEmitVel = bulkEmitVel;
166+
}
167+
168+
m_model->setPosition(i, tempPos + timeStepSize * localEmitVel);
169+
m_model->setVelocity(i, localEmitVel);
155170
}
156171
if (!insideEmitter && m_model->getParticleState(i) == ParticleState::AnimatedByEmitter
157172
&& m_model->getObjectId(i) == m_objectId)
@@ -172,7 +187,6 @@ void Emitter::emitParticles(std::vector<unsigned int>& reusedParticles, unsigned
172187
// spawn a new particle
173188
m_model->setPosition(indexNextNewParticle, tempPos - axisDepth * particleDiameter * m_depth);
174189
m_model->setParticleState(indexNextNewParticle, ParticleState::AnimatedByEmitter);
175-
m_model->setVelocity(indexNextNewParticle, emitVel);
176190
m_model->setObjectId(indexNextNewParticle, m_objectId);
177191

178192
indexNextNewParticle++;

SPlisHSPlasH/Emitter.h

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,8 @@ namespace SPH
1212
{
1313
public:
1414
Emitter(FluidModel* model, const unsigned int width, const unsigned int height, const Vector3r& pos,
15-
const Matrix3r& rotation, const Real velocity, const unsigned int type = 0, const bool useBoundary = false);
15+
const Matrix3r& rotation, const Real velocity, const unsigned int type = 0, const bool useBoundary = false,
16+
const unsigned int velocityProfile = 0);
1617
virtual ~Emitter();
1718

1819
protected:
@@ -30,6 +31,7 @@ namespace SPH
3031
unsigned int m_emitCounter{0};
3132
unsigned int m_objectId;
3233
bool m_useBoundary;
34+
unsigned int m_velocityProfile;// 0: constant, 1: linear, 2: quadratic
3335

3436
FORCE_INLINE bool inBox(const Vector3r &x, const Vector3r &xBox, const Matrix3r &rotBox, const Vector3r &scaleBox)
3537
{

SPlisHSPlasH/EmitterSystem.cpp

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -84,9 +84,10 @@ void EmitterSystem::reset()
8484

8585
void EmitterSystem::addEmitter(const unsigned int width, const unsigned int height, const Vector3r& pos,
8686
const Matrix3r& rotation, const Real velocity, const unsigned int type,
87-
const bool useBoundary)
87+
const bool useBoundary, const unsigned int velocityProfile)
8888
{
89-
m_emitters.push_back(new Emitter(m_model, width, height, pos, rotation, velocity, type, useBoundary));
89+
m_emitters.push_back(
90+
new Emitter(m_model, width, height, pos, rotation, velocity, type, useBoundary, velocityProfile));
9091
}
9192

9293
void EmitterSystem::enableReuseParticles(const Vector3r &boxMin /*= Vector3r(-1, -1, -1)*/, const Vector3r &boxMax /*= Vector3r(1, 1, 1)*/)

SPlisHSPlasH/EmitterSystem.h

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -33,7 +33,8 @@ namespace SPH
3333
void enableReuseParticles(const Vector3r &boxMin = Vector3r(-1, -1, -1), const Vector3r &boxMax = Vector3r(1, 1, 1));
3434
void disableReuseParticles();
3535
void addEmitter(const unsigned int width, const unsigned int height, const Vector3r& pos, const Matrix3r& rotation,
36-
const Real velocity, const unsigned int type, const bool useBoundary = false);
36+
const Real velocity, const unsigned int type, const bool useBoundary = false,
37+
const unsigned int velocityProfile = 0);
3738
unsigned int numEmitters() const { return static_cast<unsigned int>(m_emitters.size()); }
3839
std::vector<Emitter*> &getEmitters() { return m_emitters; }
3940

SPlisHSPlasH/Utilities/SceneParameterObjects.cpp

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -139,6 +139,7 @@ int EmitterParameterObject::EMITTER_ROTANGLE = -1;
139139
int EmitterParameterObject::EMITTER_STARTTIME = -1;
140140
int EmitterParameterObject::EMITTER_ENDTIME = -1;
141141
int EmitterParameterObject::EMITTER_TYPE = -1;
142+
int EmitterParameterObject::EMITTER_VELOCITYPROFILE = -1;
142143

143144
void EmitterParameterObject::initParameters()
144145
{
@@ -185,6 +186,11 @@ void EmitterParameterObject::initParameters()
185186
EMITTER_USEBOUNDARY = createBoolParameter("useBoundary", "Use Boundary", &useBoundary);
186187
setGroup(EMITTER_USEBOUNDARY, "Emitter");
187188
setDescription(EMITTER_USEBOUNDARY, "Creates a boundary model around the emitter if defined");
189+
190+
EMITTER_VELOCITYPROFILE =
191+
createNumericParameter<unsigned int>("velocityProfile", "Velocity Profile", &velocityProfile);
192+
setGroup(EMITTER_VELOCITYPROFILE, "Emitter");
193+
setDescription(EMITTER_VELOCITYPROFILE, "Defines the velocity profile of the emitter.");
188194
}
189195

190196

SPlisHSPlasH/Utilities/SceneParameterObjects.h

Lines changed: 5 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -145,6 +145,7 @@ namespace Utilities
145145
Real emitEndTime;
146146
unsigned int type; // type: 0 = rectangular, 1 = circle
147147
bool useBoundary;
148+
unsigned int velocityProfile;
148149

149150
EmitterParameterObject()
150151
{
@@ -160,11 +161,12 @@ namespace Utilities
160161
emitEndTime = std::numeric_limits<Real>::max();
161162
type = 0;
162163
useBoundary = false;
164+
velocityProfile = 0;
163165
}
164166

165167
EmitterParameterObject(std::string id_, unsigned int width_, unsigned int height_, Vector3r x_, Real velocity_,
166168
Vector3r axis_, Real angle_, Real emitStartTime_, Real emitEndTime_, unsigned int type_,
167-
bool useBoundary_)
169+
bool useBoundary_, unsigned int velocityProfile_)
168170
{
169171
id = id_;
170172
width = width_;
@@ -177,6 +179,7 @@ namespace Utilities
177179
emitEndTime = emitEndTime_;
178180
type = type_;
179181
useBoundary = useBoundary_;
182+
velocityProfile = velocityProfile_;
180183
}
181184

182185
static int EMITTER_ID;
@@ -190,6 +193,7 @@ namespace Utilities
190193
static int EMITTER_ENDTIME;
191194
static int EMITTER_TYPE;
192195
static int EMITTER_USEBOUNDARY;
196+
static int EMITTER_VELOCITYPROFILE;
193197

194198
virtual void initParameters();
195199
};

Simulator/SimulatorBase.cpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1316,7 +1316,7 @@ void SimulatorBase::createEmitters()
13161316
}
13171317
Matrix3r rot = AngleAxisr(ed->angle, ed->axis.normalized()).toRotationMatrix();
13181318
model->getEmitterSystem()->addEmitter(ed->width, ed->height, ed->x, rot, ed->velocity, ed->type,
1319-
ed->useBoundary);
1319+
ed->useBoundary, ed->velocityProfile);
13201320

13211321
Emitter* emitter = model->getEmitterSystem()->getEmitters().back();
13221322
if (ed->useBoundary) {

pySPlisHSPlasH/UtilitiesModule.cpp

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -340,10 +340,11 @@ void UtilitiesModule(py::module m) {
340340
py::class_<Utilities::EmitterParameterObject, GenParam::ParameterObject>(m_sub_sub, "EmitterData")
341341
.def(py::init<>())
342342
.def(py::init<std::string, unsigned int, unsigned int, Vector3r, Real, Vector3r, Real, Real, Real, unsigned int,
343-
bool>(),
343+
bool, unsigned int>(),
344344
"id"_a = "Fluid", "width"_a = 5, "height"_a = 5, "x"_a = Vector3r::Zero(), "velocity"_a = 1,
345345
"axis"_a = Vector3r(1, 0, 0), "angle"_a = 0.0, "emitStartTime"_a = 0,
346-
"emitEndTime"_a = std::numeric_limits<Real>::max(), "type"_a = 0, "useBoundary"_a = false)
346+
"emitEndTime"_a = std::numeric_limits<Real>::max(), "type"_a = 0, "useBoundary"_a = false,
347+
"velocityProfile"_a = 0)
347348
.def_readwrite("id", &Utilities::EmitterParameterObject::id)
348349
.def_readwrite("width", &Utilities::EmitterParameterObject::width)
349350
.def_readwrite("height", &Utilities::EmitterParameterObject::height)

0 commit comments

Comments
 (0)