milupHPC documentation
  • src
  • sph
kernel.cu
Go to the documentation of this file.
1#include "../../include/sph/kernel.cuh"
2#include "../../include/cuda_utils/cuda_launcher.cuh"
3
4#define MIN_NUMBER_OF_INTERACTIONS_FOR_TENSORIAL_CORRECTION_TO_WORK 0
5
6__device__ SPH::SPH_kernel spiky_p = SPH::SmoothingKernel::spiky;
7__device__ SPH::SPH_kernel cubicSpline_p = SPH::SmoothingKernel::cubicSpline;
8__device__ SPH::SPH_kernel wendlandc2_p = SPH::SmoothingKernel::wendlandc2;
9__device__ SPH::SPH_kernel wendlandc4_p = SPH::SmoothingKernel::wendlandc4;
10__device__ SPH::SPH_kernel wendlandc6_p = SPH::SmoothingKernel::wendlandc6;
11
12namespace SPH {
13
14 __device__ void SmoothingKernel::spiky(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml) {
15
16 real r, q;
17 r = 0;
18 for (int d = 0; d < DIM; d++) {
19 r += dx[d] * dx[d];
20 dWdx[d] = 0;
21 }
22 r = cuda::math::sqrt(r);
23 *dWdr = 0;
24 *W = 0;
25 q = r/sml;
26
27#if DIM == 1
28 printf("Error, this kernel can only be used with DIM == 2,3\n");
29 //assert(0);
30#elif DIM == 2
31 if (q > 1) {
32 *W = 0;
33 } else if (q >= 0.0) {
34 *W = 10./(M_PI * sml * sml) * (1-q) * (1-q) * (1-q);
35 *dWdr = -30./(M_PI * sml * sml * sml) * (1-q) * (1-q);
36 }
37#else
38 if (q > 1) {
39 *W = 0;
40 } else if (q >= 0.0) {
41 *W = 15./(M_PI * sml * sml * sml) * (1-q) * (1-q) * (1-q);
42 *dWdr = -45/(M_PI * sml * sml * sml * sml) * (1-q) * (1-q);
43 }
44#endif
45 for (int d = 0; d < DIM; d++) {
46 dWdx[d] = *dWdr/r * dx[d];
47 }
48 }
49
50 __device__ void SmoothingKernel::cubicSpline(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml) {
51
52 real r, q, f;
53 r = 0;
54 for (int d = 0; d < DIM; d++) {
55 r += dx[d] * dx[d];
56 dWdx[d] = 0;
57 }
58 r = cuda::math::sqrt(r);
59 *dWdr = 0;
60 *W = 0;
61 q = r/sml;
62
63 f = 4./3. * 1./sml;
64#if DIM > 1
65 f = 40./(7 * M_PI) * 1./(sml * sml);
66#if DIM > 2
67 f = 8./M_PI * 1./(sml * sml * sml);
68#endif
69#endif
70 if (q > 1) {
71 *W = 0;
72 *dWdr = 0.0;
73 } else if (q > 0.5) {
74 *W = 2. * f * (1.-q) * (1.-q) * (1-q);
75 *dWdr = -6. * f * 1./sml * (1.-q) * (1.-q);
76 } else if (q <= 0.5) {
77 *W = f * (6. * q * q * q - 6. * q * q + 1.);
78 *dWdr = 6. * f/sml * (3 * q * q - 2 * q);
79 }
80 for (int d = 0; d < DIM; d++) {
81 dWdx[d] = *dWdr/r * dx[d];
82 }
83 }
84
85 // Wendland C2 from Dehnen & Aly 2012
86 __device__ void SmoothingKernel::wendlandc2(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml) {
87
88 real r, q;
89 r = 0;
90 for (int d = 0; d < DIM; d++) {
91 r += dx[d]*dx[d];
92 dWdx[d] = 0;
93 }
94 r = cuda::math::sqrt(r);
95 *dWdr = 0;
96 *W = 0;
97 if (r > sml) {
98 *W = 0;
99 } else {
100 q = r/sml;
101#if DIM == 1
102 *W = 5./(4. * sml) * (1-q) * (1-q) * (1-q) * (1+3*q) * (q < 1);
103 *dWdr = -15/(sml * sml) * q * (1-q) * (1-q) * (q < 1);
104#elif DIM == 2
105 *W = 7./(M_PI * sml * sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1+4 * q) * (q < 1);
106 *dWdr = -140./(M_PI * sml * sml * sml) * q * (1-q) * (1-q) * (1-q) * (q < 1);
107#else //DIM == 3
108 *W = 21./(2 * M_PI * sml * sml * sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1+4 * q) * (q < 1);
109 *dWdr = -210./(M_PI * sml * sml * sml * sml) * q * (1-q) * (1-q) * (1-q) * (q < 1);
110#endif
111 for (int d = 0; d < DIM; d++) {
112 dWdx[d] = *dWdr/r * dx[d];
113 }
114 }
115 }
116
117// Wendland C4 from Dehnen & Aly 2012
118 __device__ void SmoothingKernel::wendlandc4(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml) {
119
120 real r, q;
121 r = 0;
122 for (int d = 0; d < DIM; d++) {
123 r += dx[d]*dx[d];
124 dWdx[d] = 0;
125 }
126 r = cuda::math::sqrt(r);
127 *dWdr = 0;
128 *W = 0;
129
130 if (r > sml) {
131 *W = 0;
132 } else {
133 q = r/sml;
134#if DIM == 1
135 *W = 3./(2.*sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1+5*q+8*q*q) * (q < 1);
136 *dWdr = -21./(sml*sml) * q * (1-q) * (1-q) * (1-q) * (1-q) * (1+4*q) * (q < 1);
137#elif DIM == 2
138 *W = 9./(M_PI*sml*sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1.+6*q+35./3.*q*q) * (q < 1);
139 *dWdr = -54./(M_PI*sml*sml*sml) * (1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (1.-35.*q*q+105.*q*q*q) * (q< 1);
140#else //DIM == 3
141 *W = 495./(32.*M_PI*sml*sml*sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1.+6.*q+35./3.*q*q) * (q < 1);
142 *dWdr = -1485./(16.*M_PI*sml*sml*sml*sml) * (1-q) * (1-q) * (1-q) * (1-q) * (1-q) * (1.-35.*q*q+105.*q*q*q) * (q< 1);
143#endif
144 for (int d = 0; d < DIM; d++) {
145 dWdx[d] = *dWdr/r * dx[d];
146 }
147 }
148 }
149
150
151 // Wendland C6 from Dehnen & Aly 2012
152 __device__ void SmoothingKernel::wendlandc6(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml) {
153
154 real r, q;
155 r = 0;
156 for (int d = 0; d < DIM; d++) {
157 r += dx[d]*dx[d];
158 dWdx[d] = 0;
159 }
160 r = cuda::math::sqrt(r);
161 *dWdr = 0;
162 *W = 0;
163 if (r > sml) {
164 *W = 0;
165 } else {
166 q = r/sml;
167#if DIM == 1
168 *W = 55./(32.*sml) * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (1+7*q+19*q*q+21*q*q*q) * (q < 1);
169 *dWdr = -165./(16*sml*sml) * q * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (3+18*q+35*q*q) * (q < 1);
170#elif DIM == 2
171 *W = 78./(7.*M_PI*sml*sml) * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (1.+8.*q+25.*q*q+32*q*q*q) * (q < 1);
172
173 *dWdr = -1716./(7.*M_PI*sml*sml*sml) * q * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (1.+7*q+16*q*q) * (q < 1);
174#else // DIM == 3
175 *W = 1365./(64.*M_PI*sml*sml*sml) * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) * (1.+8.*q+25.*q*q+32*q*q*q) * (q < 1);
176 *dWdr = -15015./(32.*M_PI*sml*sml*sml*sml) * q * (1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q)*(1-q) *
177 (1.+7*q+16*q*q) * (q < 1);
178#endif
179 for (int d = 0; d < DIM; d++) {
180 dWdx[d] = *dWdr/r * dx[d];
181 }
182 }
183 }
184
185
186 CUDA_CALLABLE_MEMBER real fixTensileInstability(SPH_kernel kernel, Particles *particles, int p1, int p2) {
187
188 real hbar;
189 real dx[DIM];
190 real W;
191 real W2;
192 real dWdr;
193 real dWdx[DIM];
194
195 W = 0;
196 W2 = 0;
197 dWdr = 0;
198 for (int d = 0; d < DIM; d++) {
199 dx[d] = 0.0;
200 dWdx[d] = 0;
201 }
202 dx[0] = particles->x[p1] - particles->x[p2];
203#if DIM > 1
204 dx[1] = particles->y[p1] - particles->y[p2];
205#if DIM > 2
206 dx[2] = particles->z[p1] - particles->z[p2];
207#endif
208#endif
209
210 hbar = 0.5 * (particles->sml[p1] + particles->sml[p2]);
211 // calculate kernel for r and particle_distance
212 kernel(&W, dWdx, &dWdr, dx, hbar);
213 //TODO: matmean_particle_distance
214 //dx[0] = matmean_particle_distance[p_rhs.materialId[a]];
215 for (int d = 1; d < DIM; d++) {
216 dx[d] = 0;
217 }
218 kernel(&W2, dWdx, &dWdr, dx, hbar);
219
220 return W/W2;
221
222 }
223
224//#if (NAVIER_STOKES || BALSARA_SWITCH || INVISCID_SPH || INTEGRATE_ENERGY)
225#if BALSARA_SWITCH
226 namespace Kernel {
227 __global__ void CalcDivvandCurlv(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles) {
228
229 int i, inc, j, k, m, d, dd;
230 // absolute values of div v and curl v */
231 real divv;
232 real curlv[DIM];
233 real W, dWdr;
234 real Wj, dWdrj, dWdxj[DIM];
235 real dWdx[DIM], dx[DIM];
236 real sml;
237 real vi[DIM], vj[DIM];
238 real r;
239 inc = blockDim.x * gridDim.x;
240 for (i = threadIdx.x + blockIdx.x * blockDim.x; i < numParticles; i += inc) {
241 //if (EOS_TYPE_IGNORE == matEOS[p_rhs.materialId[i]] || p_rhs.materialId[i] == EOS_TYPE_IGNORE) {
242 // continue;
243 //}
244 k = particles->noi[i];
245 divv = 0;
246 for (m = 0; m < DIM; m++) {
247 curlv[m] = 0;
248 dWdx[m] = 0;
249 }
250 sml = particles->sml[i];
251 // interaction partner loop
252 for (m = 0; m < k; m++) {
253 j = interactions[i*MAX_NUM_INTERACTIONS + m];
254 // get the kernel values
255 #if VARIABLE_SML
256 sml = 0.5 *(particles->sml[i] + particles->sml[j]);
257 #endif
258 dx[0] = particles->x[i] - particles->x[j];
259 #if DIM > 1
260 dx[1] = particles->y[i] - particles->y[j];
261 #if DIM > 2
262 dx[2] = particles->z[i] - particles->z[j];
263 #endif
264 #endif
265
266
267 #if AVERAGE_KERNELS
268 kernel(&W, dWdx, &dWdr, dx, particles->sml[i]);
269 kernel(&Wj, dWdxj, &dWdrj, dx, particles->sml[j]);
270 # if SHEPARD_CORRECTION
271 //TODO: shephard correction
272 W /= particles->shepard_correction[i];
273 Wj /= particles->shepard_correction[j];
274 for (d = 0; d < DIM; d++) {
275 dWdx[d] /= p_rhs.shepard_correction[i];
276 dWdxj[d] /= p_rhs.shepard_correction[j];
277 }
278 # endif
279 W = 0.5 * (W + Wj);
280 for (d = 0; d < DIM; d++) {
281 dWdx[d] = 0.5 * (dWdx[d] + dWdxj[d]);
282 }
283 #else
284 kernel(&W, dWdx, &dWdr, dx, sml);
285 # if SHEPARD_CORRECTION
286 W /= p_rhs.shepard_correction[i];
287 for (d = 0; d < DIM; d++) {
288 dWdx[d] /= particles->shepard_correction[i];
289 }
290 # endif
291 #endif // AVERAGE_KERNELS
292
293 vi[0] = particles->vx[i];
294 vj[0] = particles->vx[j];
295 #if DIM > 1
296 vi[1] = particles->vy[i];
297 vj[1] = particles->vy[j];
298 #if DIM > 2
299 vi[2] = particles->vz[i];
300 vj[2] = particles->vz[j];
301 #endif
302 #endif
303 r = 0;
304 for (d = 0; d < DIM; d++) {
305 r += dx[d]*dx[d];
306 }
307 r = cuda::math::sqrt(r);
308 // divv
309 for (d = 0; d < DIM; d++) {
310 #if TENSORIAL_CORRECTION
311 for (dd = 0; dd < DIM; dd++) {
312 divv += particles->mass[j]/particles->rho[j] * (vj[d] - vi[d]) * particles->tensorialCorrectionMatrix[i*DIM*DIM+d*DIM+dd] * dWdx[dd];
313 }
314 #else
315 divv += particles->mass[j]/particles->rho[j] * (vj[d] - vi[d]) * dWdx[d];
316 #endif
317
318 }
319 /* curlv */
320 #if (DIM == 1 && BALSARA_SWITCH)
321 #error unset BALSARA SWITCH in 1D
322 #elif DIM == 2
323 // only one component in 2D
324 curlv[0] += particles->mass[j]/particles->rho[i] * ((vi[0] - vj[0]) * dWdx[1]
325 - (vi[1] - vj[1]) * dWdx[0]);
326 curlv[1] = 0;
327 #elif DIM == 3
328 curlv[0] += particles->mass[j]/particles->rho[i] * ((vi[1] - vj[1]) * dWdx[2]
329 - (vi[2] - vj[2]) * dWdx[1]);
330 curlv[1] += particles->mass[j]/particles->rho[i] * ((vi[2] - vj[2]) * dWdx[0]
331 - (vi[0] - vj[0]) * dWdx[2]);
332 curlv[2] += particles->mass[j]/particles->rho[i] * ((vi[0] - vj[0]) * dWdx[1]
333 - (vi[1] - vj[1]) * dWdx[0]);
334 #endif
335 }
336 for (d = 0; d < DIM; d++) {
337 //TODO: particles or particles_rhs: curlv and divv
338 particles->curlv[i*DIM+d] = curlv[d];
339 }
340 particles->divv[i] = divv;
341 }
342 }
343
344 namespace Launch {
345 real CalcDivvandCurlv(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles) {
346 ExecutionPolicy executionPolicy;
347 return cuda::launch(true, executionPolicy, ::SPH::Kernel::CalcDivvandCurlv, kernel, particles,
348 interactions, numParticles);
349 }
350 }
351
352 }
353#endif // (NAVIER_STOKES || BALSARA_SWITCH || INVISCID_SPH)
354
355#if ZERO_CONSISTENCY //SHEPARD_CORRECTION
356 // this adds zeroth order consistency but needs one more loop over all neighbours
357__global__ void shepardCorrection(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles) {
358
359 register int i, inc, j, m;
360 register real dr[DIM], h, dWdr;
361 inc = blockDim.x * gridDim.x;
362 real W, dWdx[DIM], Wj;
363 for (i = threadIdx.x + blockIdx.x * blockDim.x; i < numParticles; i += inc) {
364 real shepard_correction;
365 W = 0;
366 for (m = 0; m < DIM; m++) {
367 dr[m] = 0.0;
368 }
369 kernel(&W, dWdx, &dWdr, dr, particles->sml[i]);
370 shepard_correction = particles->mass[i]/particles->rho[i]*W;
371
372 for (m = 0; m < particles->noi[i]; m++) {
373 W = 0;
374 j = interactions[i*MAX_NUM_INTERACTIONS + m];
375 //if (EOS_TYPE_IGNORE == matEOS[p_rhs.materialId[j]] || p_rhs.materialId[j] == EOS_TYPE_IGNORE) {
376 // continue;
377 //}
378 dr[0] = particles->x[i] - particles->x[j];
379#if DIM > 1
380 dr[1] = particles->y[i] - particles->y[j];
381#if DIM > 2
382 dr[2] = particles->z[i] - particles->z[j];
383#endif
384#endif
385
386#if AVERAGE_KERNELS
387 kernel(&W, dWdx, &dWdr, dr, particles->sml[i]);
388 Wj = 0;
389 kernel(&Wj, dWdx, &dWdr, dr, particles->sml[j]);
390 W = 0.5*(W + Wj);
391#else
392 h = 0.5*(particles->sml[i] + particles->sml[j]);
393 kernel(&W, dWdx, &dWdr, dr, h);
394#endif
395
396 shepard_correction += particles->mass[j]/particles->rho[j]*W;
397 }
398 // TODO: particles or particles_rhs: shepard_correction
399 //particles->shepard_correction[i] = shepard_correction;
400 //printf("%g\n", shepard_correction);
401 }
402}
403#endif
404
405
406#if LINEAR_CONSISTENCY //TENSORIAL_CORRECTION
407 // this adds first order consistency but needs one more loop over all neighbours
408__global__ void tensorialCorrection(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles)
409{
410 register int i, inc, j, k, m;
411 register int d, dd;
412 int rv = 0;
413 inc = blockDim.x * gridDim.x;
414 register real r, dr[DIM], h, dWdr, tmp, f1, f2;
415 real W, dWdx[DIM];
416 real Wj, dWdxj[DIM];
417 real wend_f, wend_sml, q, distance;
418 for (i = threadIdx.x + blockIdx.x * blockDim.x; i < numParticles; i += inc) {
419 register real corrmatrix[DIM*DIM];
420 register real matrix[DIM*DIM];
421 for (d = 0; d < DIM*DIM; d++) {
422 corrmatrix[d] = 0;
423 matrix[d] = 0;
424 }
425 //if (EOS_TYPE_IGNORE == matEOS[p_rhs.materialId[i]] || p_rhs.materialId[i] == EOS_TYPE_IGNORE) {
426 // continue;
427 //}
428
429 k = particles->noi[i];
430
431 // loop over all interaction partner
432 for (m = 0; m < k; m++) {
433 j = interactions[i*MAX_NUM_INTERACTIONS+m];
434 //if (EOS_TYPE_IGNORE == matEOS[p_rhs.materialId[j]] || p_rhs.materialId[j] == EOS_TYPE_IGNORE) {
435 // continue;
436 //}
437 dr[0] = particles->x[i] - particles->x[j];
438#if DIM > 1
439 dr[1] = particles->y[i] - particles->y[j];
440#if DIM == 3
441 dr[2] = particles->z[i] - particles->z[j];
442
443 r = cuda::math::sqrt(dr[0]*dr[0]+dr[1]*dr[1]+dr[2]*dr[2]);
444#elif DIM == 2
445 r = cuda::math::sqrt(dr[0]*dr[0]+dr[1]*dr[1]);
446#endif
447#endif
448
449#if AVERAGE_KERNELS
450 kernel(&W, dWdx, &dWdr, dr, particles->sml[i]);
451 kernel(&Wj, dWdxj, &dWdr, dr, particles->sml[j]);
452# if SHEPARD_CORRECTION
453 W /= particles->shepard_correction[i];
454 Wj /= particles->shepard_correction[j];
455 for (d = 0; d < DIM; d++) {
456 dWdx[d] /= particles->shepard_correction[i];
457 dWdxj[d] /= particles->shepard_correction[j];
458 }
459 for (d = 0; d < DIM; d++) {
460 dWdx[d] = 0.5 * (dWdx[d] + dWdxj[d]);
461 }
462 W = 0.5 * (W + Wj);
463# endif
464
465
466#else
467 h = 0.5*(particles->sml[i] + particles->sml[j]);
468 kernel(&W, dWdx, &dWdr, dr, h);
469# if SHEPARD_CORRECTION
470 W /= particles->shepard_correction[i];
471 for (d = 0; d < DIM; d++) {
472 dWdx[d] /= particles->shepard_correction[i];
473 }
474# endif
475#endif // AVERAGE_KERNELS
476
477 for (d = 0; d < DIM; d++) {
478 for (dd = 0; dd < DIM; dd++) {
479 corrmatrix[d*DIM+dd] -= particles->mass[j]/particles->rho[j] * dr[d] * dWdx[dd];
480 }
481 }
482 } // end loop over interaction partners
483
484 rv = CudaUtils::invertMatrix(corrmatrix, matrix);
485 // if something went wrong during inversion, use identity matrix
486 if (rv < 0 || k < MIN_NUMBER_OF_INTERACTIONS_FOR_TENSORIAL_CORRECTION_TO_WORK) {
487 #if DEBUG_LINALG
488 if (threadIdx.x == 0) {
489 printf("could not invert matrix: rv: %d and k: %d\n", rv, k);
490 for (d = 0; d < DIM; d++) {
491 for (dd = 0; dd < DIM; dd++) {
492 printf("%e\t", corrmatrix[d*DIM+dd]);
493 }
494 printf("\n");
495 }
496 }
497 #endif
498 for (d = 0; d < DIM; d++) {
499 for (dd = 0; dd < DIM; dd++) {
500 matrix[d*DIM+dd] = 0.0;
501 if (d == dd)
502 matrix[d*DIM+dd] = 1.0;
503 }
504 }
505 }
506 for (d = 0; d < DIM*DIM; d++) {
507 // TODO: particles or particles_rhs: tensorialCorrectionMatrix
508 //particles->tensorialCorrectionMatrix[i*DIM*DIM+d] = matrix[d];
509
510 }
511 }
512}
513#endif
514}
515
ExecutionPolicy
Execution policy/instruction for CUDA kernel execution.
Definition: cuda_launcher.cuh:33
Particles
Particle(s) class based on SoA (Structur of Arrays).
Definition: particles.cuh:50
Particles::noi
integer * noi
(pointer to) number of interactions (array)
Definition: particles.cuh:118
Particles::x
real * x
(pointer to) x position (array)
Definition: particles.cuh:62
Particles::rho
real * rho
(pointer to) density (array)
Definition: particles.cuh:133
Particles::y
real * y
(pointer to) y position (array)
Definition: particles.cuh:70
Particles::sml
real * sml
(pointer to) smoothing length (array)
Definition: particles.cuh:113
Particles::mass
real * mass
(pointer to) mass (array)
Definition: particles.cuh:60
Particles::z
real * z
(pointer to) z position (array)
Definition: particles.cuh:78
Particles::vz
real * vz
(pointer to) z velocity (array)
Definition: particles.cuh:80
Particles::vx
real * vx
(pointer to) x velocity (array)
Definition: particles.cuh:64
Particles::vy
real * vy
(pointer to) y velocity (array)
Definition: particles.cuh:72
CUDA_CALLABLE_MEMBER
#define CUDA_CALLABLE_MEMBER
Definition: cuda_utilities.cuh:30
wendlandc6_p
__device__ SPH::SPH_kernel wendlandc6_p
Definition: kernel.cu:10
MIN_NUMBER_OF_INTERACTIONS_FOR_TENSORIAL_CORRECTION_TO_WORK
#define MIN_NUMBER_OF_INTERACTIONS_FOR_TENSORIAL_CORRECTION_TO_WORK
Definition: kernel.cu:4
spiky_p
__device__ SPH::SPH_kernel spiky_p
Definition: kernel.cu:6
cubicSpline_p
__device__ SPH::SPH_kernel cubicSpline_p
Definition: kernel.cu:7
wendlandc2_p
__device__ SPH::SPH_kernel wendlandc2_p
Definition: kernel.cu:8
wendlandc4_p
__device__ SPH::SPH_kernel wendlandc4_p
Definition: kernel.cu:9
CudaUtils::invertMatrix
__device__ int invertMatrix(real *m, real *inverted)
Invert matrix.
Definition: linalg.cu:279
Kernel
Definition: device_rhs.cuh:7
ProfilerIds::numParticles
const char *const numParticles
Definition: h5profiler.h:29
SPH::SmoothingKernel::wendlandc6
__device__ void wendlandc6(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml)
Definition: kernel.cu:152
SPH::SmoothingKernel::wendlandc2
__device__ void wendlandc2(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml)
Definition: kernel.cu:86
SPH::SmoothingKernel::wendlandc4
__device__ void wendlandc4(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml)
Definition: kernel.cu:118
SPH::SmoothingKernel::cubicSpline
__device__ void cubicSpline(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml)
Cubic spline kernel (Monaghan & Lattanzio 1985).
Definition: kernel.cu:50
SPH::SmoothingKernel::spiky
__device__ void spiky(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml)
Spiky kernel (Desbrun & Cani).
Definition: kernel.cu:14
SPH
SPH related functions and kernels.
Definition: density.cuh:23
SPH::fixTensileInstability
CUDA_CALLABLE_MEMBER real fixTensileInstability(SPH_kernel kernel, Particles *particles, int p1, int p2)
Calculates the kernel for the tensile instability fix (Monaghan 2000).
Definition: kernel.cu:186
SPH::SPH_kernel
void(* SPH_kernel)(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real h)
Function pointer to generic SPH kernel function.
Definition: kernel.cuh:26
cuda::math::sqrt
__device__ real sqrt(real a)
Square root of a floating point value.
Definition: cuda_utilities.cu:456
cuda::launch
real launch(bool timeKernel, const ExecutionPolicy &policy, void(*f)(Arguments...), Arguments... args)
CUDA execution wrapper function.
Definition: cuda_launcher.cuh:114
real
double real
Definition: parameter.h:15
MAX_NUM_INTERACTIONS
#define MAX_NUM_INTERACTIONS
Definition: parameter.h:106
DIM
#define DIM
Dimension of the problem.
Definition: parameter.h:38

milupHPC - src/sph/kernel.cu Source File
Generated on Wed Aug 31 2022 12:16:53 by Doxygen 1.9.3