PeriDyno 1.0.0
Loading...
Searching...
No Matches
SummationDensity.cu
Go to the documentation of this file.
1#include "SummationDensity.h"
2
3namespace dyno
4{
5 template<typename TDataType>
6 SummationDensity<TDataType>::SummationDensity()
7 : ParticleApproximation<TDataType>()
8 , m_factor(Real(1))
9 {
10 this->varRestDensity()->setValue(Real(1000));
11
12 auto callback = std::make_shared<FCallBackFunc>(
13 std::bind(&SummationDensity<TDataType>::calculateParticleMass, this));
14
15 this->varRestDensity()->attach(callback);
16 this->inSamplingDistance()->attach(callback);
17
18 this->inOther()->tagOptional(true);
19 //calculateParticleMass();
20 }
21
22 template<typename TDataType>
23 void SummationDensity<TDataType>::compute()
24 {
25 int p_num = this->inPosition()->getDataPtr()->size();
26 int n_num = this->inNeighborIds()->getDataPtr()->size();
27 if (p_num != n_num) {
28 Log::sendMessage(Log::Error, "The input array sizes of DensitySummation are not compatible!");
29 return;
30 }
31
32 if (this->outDensity()->size() != p_num) {
33 this->outDensity()->resize(p_num);
34 }
35
36 if (this->inOther()->isEmpty()) {
37 compute(
38 this->outDensity()->getData(),
39 this->inPosition()->getData(),
40 this->inNeighborIds()->getData(),
41 this->inSmoothingLength()->getData(),
42 m_particle_mass);
43 }
44 else {
45 compute(
46 this->outDensity()->getData(),
47 this->inPosition()->getData(),
48 this->inOther()->getData(),
49 this->inNeighborIds()->getData(),
50 this->inSmoothingLength()->getData(),
51 m_particle_mass);
52 }
53 }
54
55 template<typename Real, typename Coord, typename Kernel>
56 __global__ void SD_ComputeDensity(
57 DArray<Real> rhoArr,
58 DArray<Coord> posArr,
59 DArrayList<int> neighbors,
60 Real smoothingLength,
61 Real mass,
62 Kernel weight,
63 Real scale)
64 {
65 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
66 if (pId >= posArr.size()) return;
67
68 Real r;
69 Real rho_i = Real(0);
70 Coord pos_i = posArr[pId];
71 List<int>& list_i = neighbors[pId];
72 int nbSize = list_i.size();
73 for (int ne = 0; ne < nbSize; ne++)
74 {
75 int j = list_i[ne];
76 r = (pos_i - posArr[j]).norm();
77 rho_i += mass * weight(r, smoothingLength, scale);
78 }
79
80 rhoArr[pId] = rho_i;
81 }
82
83 template<typename Real, typename Coord, typename Kernel>
84 __global__ void SD_ComputeDensity(
85 DArray<Real> rhoArr,
86 DArray<Coord> posArr,
87 DArray<Coord> posQueried,
88 DArrayList<int> neighbors,
89 Real smoothingLength,
90 Real mass,
91 Kernel weight,
92 Real scale)
93 {
94 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
95 if (pId >= posArr.size()) return;
96
97 Real r;
98 Real rho_i = Real(0);
99 Coord pos_i = posArr[pId];
100 List<int>& list_i = neighbors[pId];
101 int nbSize = list_i.size();
102 for (int ne = 0; ne < nbSize; ne++)
103 {
104 int j = list_i[ne];
105 r = (pos_i - posQueried[j]).norm();
106 rho_i += mass * weight(r, smoothingLength, scale);
107 }
108
109 rhoArr[pId] = rho_i;
110 }
111
112 template<typename TDataType>
113 void SummationDensity<TDataType>::compute(
114 DArray<Real>& rho,
115 DArray<Coord>& pos,
116 DArrayList<int>& neighbors,
117 Real smoothingLength,
118 Real mass)
119 {
120 cuZerothOrder(rho.size(), this->varKernelType()->getDataPtr()->currentKey(), this->mScalingFactor,
121 SD_ComputeDensity,
122 rho,
123 pos,
124 neighbors,
125 smoothingLength,
126 mass);
127 }
128
129 template<typename TDataType>
130 void SummationDensity<TDataType>::compute(DArray<Real>& rho, DArray<Coord>& pos, DArray<Coord>& posQueried, DArrayList<int>& neighbors, Real smoothingLength, Real mass)
131 {
132 cuZerothOrder(rho.size(), this->varKernelType()->getDataPtr()->currentKey(), this->mScalingFactor,
133 SD_ComputeDensity,
134 rho,
135 pos,
136 posQueried,
137 neighbors,
138 smoothingLength,
139 mass);
140 }
141
142 template<typename TDataType>
143 void SummationDensity<TDataType>::calculateParticleMass()
144 {
145 Real rho_0 = this->varRestDensity()->getData();
146 Real d = this->inSamplingDistance()->getData();
147
148 m_particle_mass = d * d*d*rho_0;
149 }
150
151 DEFINE_CLASS(SummationDensity);
152}