PeriDyno 1.0.0
Loading...
Searching...
No Matches
DualParticleFluid.cu
Go to the documentation of this file.
1#include "DualParticleFluid.h"
2//DataType
3#include "Auxiliary/DataSource.h"
4
5//Collision
6#include "Collision/NeighborPointQuery.h"
7
8//ParticleSystem
9#include "ParticleSystem/Module/ImplicitViscosity.h"
10#include "ParticleSystem/Module/ParticleIntegrator.h"
11
12//DualParticleSystem
13#include "Module/DualParticleIsphModule.h"
14
15
16
17namespace dyno
18{
19 __global__ void DPS_AttributeReset(
20 DArray<Attribute> att
21 )
22 {
23 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
24 if (pId >= att.size()) return;
25
26 att[pId].setFluid();
27 att[pId].setDynamic();
28 }
29
30 template<typename TDataType>
31 DualParticleFluid<TDataType>::DualParticleFluid()
32 : DualParticleFluid<TDataType>::DualParticleFluid(2)
33 {
34 }
35
36
37 template<typename TDataType>
38 DualParticleFluid<TDataType>::DualParticleFluid(int key)
39 : ParticleFluid<TDataType>()
40 {
41 this->varVirtualParticleSamplingStrategy()->getDataPtr()->setCurrentKey(key);
42 this->varReshuffleParticles()->setValue(false);
43
44 this->animationPipeline()->clear();
45
46 auto smoothingLength = std::make_shared<FloatingNumber<TDataType>>();
47 this->animationPipeline()->pushModule(smoothingLength);
48 smoothingLength->varValue()->setValue(Real(0.012));
49
50 auto samplingDistance = std::make_shared<FloatingNumber<TDataType>>();
51 this->animationPipeline()->pushModule(smoothingLength);
52 samplingDistance->varValue()->setValue(Real(0.005));
53
54 if (key == EVirtualParticleSamplingStrategy::SpatiallyAdaptiveStrategy)
55 {
56 auto m_adaptive_virtual_position = std::make_shared<VirtualSpatiallyAdaptiveStrategy<TDataType>>();
57 this->statePosition()->connect(m_adaptive_virtual_position->inRPosition());
58 m_adaptive_virtual_position->varSamplingDistance()->setValue(Real(0.005)); /**Virtual particle radius*/
59 m_adaptive_virtual_position->varCandidatePointCount()->getDataPtr()->setCurrentKey(VirtualSpatiallyAdaptiveStrategy<TDataType>::neighbors_33);
60 vpGen = m_adaptive_virtual_position;
61 }
62 else if (key == EVirtualParticleSamplingStrategy::ParticleShiftingStrategy)
63 {
64 auto m_virtual_particle_shifting = std::make_shared<VirtualParticleShiftingStrategy<TDataType >>();
65 this->stateFrameNumber()->connect(m_virtual_particle_shifting->inFrameNumber());
66 this->stateFrameNumber()->connect(m_virtual_particle_shifting->inFrameNumber());
67 this->statePosition()->connect(m_virtual_particle_shifting->inRPosition());
68 this->stateTimeStep()->connect(m_virtual_particle_shifting->inTimeStep());
69 this->animationPipeline()->pushModule(m_virtual_particle_shifting);
70 vpGen = m_virtual_particle_shifting;
71 }
72 else if (key == EVirtualParticleSamplingStrategy::ColocationStrategy)
73 {
74 auto m_virtual_equal_to_Real = std::make_shared<VirtualColocationStrategy<TDataType >>();
75 this->statePosition()->connect(m_virtual_equal_to_Real->inRPosition());
76 this->animationPipeline()->pushModule(m_virtual_equal_to_Real);
77 vpGen = m_virtual_equal_to_Real;
78 }
79
80 this->animationPipeline()->pushModule(vpGen);
81 vpGen->outVirtualParticles()->connect(this->stateVirtualPosition());
82
83 auto m_nbrQuery = std::make_shared<NeighborPointQuery<TDataType>>();
84 smoothingLength->outFloating()->connect(m_nbrQuery->inRadius());
85 this->statePosition()->connect(m_nbrQuery->inPosition());
86 this->animationPipeline()->pushModule(m_nbrQuery);
87
88 auto m_rv_nbrQuery = std::make_shared<NeighborPointQuery<TDataType>>();
89 smoothingLength->outFloating()->connect(m_rv_nbrQuery->inRadius());
90 this->statePosition()->connect(m_rv_nbrQuery->inOther());
91 vpGen->outVirtualParticles()->connect(m_rv_nbrQuery->inPosition());
92 this->animationPipeline()->pushModule(m_rv_nbrQuery);
93
94 auto m_vr_nbrQuery = std::make_shared<NeighborPointQuery<TDataType>>();
95 smoothingLength->outFloating()->connect(m_vr_nbrQuery->inRadius());
96 this->statePosition()->connect(m_vr_nbrQuery->inPosition());
97 vpGen->outVirtualParticles()->connect(m_vr_nbrQuery->inOther());
98 this->animationPipeline()->pushModule(m_vr_nbrQuery);
99
100 auto m_vv_nbrQuery = std::make_shared<NeighborPointQuery<TDataType>>();
101 smoothingLength->outFloating()->connect(m_vv_nbrQuery->inRadius());
102 vpGen->outVirtualParticles()->connect(m_vv_nbrQuery->inPosition());
103 this->animationPipeline()->pushModule(m_vv_nbrQuery);
104
105 auto m_dualIsph = std::make_shared<DualParticleIsphModule<TDataType>>();
106 smoothingLength->outFloating()->connect(m_dualIsph->varSmoothingLength());
107 this->stateTimeStep()->connect(m_dualIsph->inTimeStep());
108 this->statePosition()->connect(m_dualIsph->inRPosition());
109 vpGen->outVirtualParticles()->connect(m_dualIsph->inVPosition());
110 this->stateVelocity()->connect(m_dualIsph->inVelocity());
111 //this->stateParticleAttribute()->connect(m_dualIsph->inParticleAttribute());
112 //this->stateBoundaryNorm()->connect(m_dualIsph->inBoundaryNorm());
113 m_nbrQuery->outNeighborIds()->connect(m_dualIsph->inNeighborIds());
114 m_rv_nbrQuery->outNeighborIds()->connect(m_dualIsph->inRVNeighborIds());
115 m_vr_nbrQuery->outNeighborIds()->connect(m_dualIsph->inVRNeighborIds());
116 m_vv_nbrQuery->outNeighborIds()->connect(m_dualIsph->inVVNeighborIds());
117 this->stateTimeStep()->connect(m_dualIsph->inTimeStep());
118 this->animationPipeline()->pushModule(m_dualIsph);
119
120 auto m_integrator = std::make_shared<ParticleIntegrator<TDataType>>();
121 this->stateTimeStep()->connect(m_integrator->inTimeStep());
122 this->statePosition()->connect(m_integrator->inPosition());
123 this->stateVelocity()->connect(m_integrator->inVelocity());
124 this->stateParticleAttribute()->connect(m_integrator->inAttribute());
125 this->animationPipeline()->pushModule(m_integrator);
126
127 auto m_visModule = std::make_shared<ImplicitViscosity<TDataType>>();
128 m_visModule->varViscosity()->setValue(Real(0.3));
129 this->stateTimeStep()->connect(m_visModule->inTimeStep());
130 smoothingLength->outFloating()->connect(m_visModule->inSmoothingLength());
131 this->stateTimeStep()->connect(m_visModule->inTimeStep());
132 this->statePosition()->connect(m_visModule->inPosition());
133 this->stateVelocity()->connect(m_visModule->inVelocity());
134 m_nbrQuery->outNeighborIds()->connect(m_visModule->inNeighborIds());
135 this->animationPipeline()->pushModule(m_visModule);
136 }
137
138
139 template<typename TDataType>
140 DualParticleFluid<TDataType>::~DualParticleFluid()
141 {
142
143 }
144
145 template<typename TDataType>
146 void DualParticleFluid<TDataType>::resetStates()
147 {
148 this->ParticleFluid<TDataType>::resetStates();
149
150 auto ptSet = this->statePointSet()->getDataPtr();
151 if(ptSet != nullptr)
152 {
153 auto pts = ptSet->getPoints();
154 this->stateBoundaryNorm()->resize(pts.size());
155 this->stateParticleAttribute()->resize(pts.size());
156
157 cuExecute(pts.size(), DPS_AttributeReset,
158 this->stateParticleAttribute()->getData());
159
160 this->stateBoundaryNorm()->getDataPtr()->reset();
161 }
162
163 if (this->stateVirtualPointSet()->isEmpty())
164 {
165 this->stateVirtualPointSet()->allocate();
166 }
167
168 if (!this->stateVirtualPosition()->isEmpty())
169 {
170 auto virtualPoints = this->stateVirtualPointSet()->getDataPtr();
171 virtualPoints->setPoints(this->stateVirtualPosition()->getData());
172 }
173 else
174 {
175 auto virtualPoints = this->stateVirtualPointSet()->getDataPtr();
176 virtualPoints->clear();
177 }
178 }
179
180 template<typename TDataType>
181 void DualParticleFluid<TDataType>::preUpdateStates()
182 {
183 this->varReshuffleParticles()->setValue(false);
184 this->ParticleFluid<TDataType>::preUpdateStates();
185
186 this->stateBoundaryNorm()->resize(this->statePosition()->size());
187 this->stateBoundaryNorm()->reset();
188 this->stateParticleAttribute()->resize(this->statePosition()->size());
189
190 cuExecute(this->statePosition()->size(), DPS_AttributeReset,
191 this->stateParticleAttribute()->getData());
192 }
193
194
195 template<typename TDataType>
196 void DualParticleFluid<TDataType>::postUpdateStates()
197 {
198 this->ParticleSystem<TDataType>::postUpdateStates();
199
200 if (!this->stateVirtualPosition()->isEmpty())
201 {
202 auto virtualPoints = this->stateVirtualPointSet()->getDataPtr();
203 virtualPoints->setPoints(this->stateVirtualPosition()->getData());
204 }
205 else
206 {
207 auto virtualPoints = this->stateVirtualPointSet()->getDataPtr();
208 virtualPoints->clear();
209 }
210 }
211
212 DEFINE_CLASS(DualParticleFluid);
213}
214
215