milupHPC documentation
  • include
  • sph
kernel.cuh
Go to the documentation of this file.
1
8#ifndef MILUPHPC_KERNEL_CUH
9#define MILUPHPC_KERNEL_CUH
10
11#include "../particles.cuh"
12#include "../parameter.h"
13#include "../helper.cuh"
14
15//#include <boost/mpi.hpp>
16#include <assert.h>
17#include "../parameter.h"
18#include "../cuda_utils/linalg.cuh"
19
20namespace SPH {
21
22
26 typedef void (*SPH_kernel)(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real h);
27
28 namespace SmoothingKernel {
29
56 __device__ void spiky(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml);
57
86 __device__ void cubicSpline(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml);
87
116 __device__ void wendlandc2(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml);
117
146 __device__ void wendlandc4(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml);
147
176 __device__ void wendlandc6(real *W, real dWdx[DIM], real *dWdr, real dx[DIM], real sml);
177 }
178
179
180 //TODO: implement:
192 CUDA_CALLABLE_MEMBER real fixTensileInstability(SPH_kernel kernel, Particles *particles, int p1, int p2);
193
194#if (NAVIER_STOKES || BALSARA_SWITCH || INVISCID_SPH || INTEGRATE_ENERGY)
195
199 namespace Kernel {
200 __global__ void CalcDivvandCurlv(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles);
201 namespace Launch {
202 real CalcDivvandCurlv(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles);
203 }
204 }
205#endif
206
207#if ZERO_CONSISTENCY //SHEPARD_CORRECTION
208
212 __global__ void shepardCorrection(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles);
213#endif
214
215#if LINEAR_CONSISTENCY //TENSORIAL_CORRECTION
216
220 __global__ void tensorialCorrection(SPH_kernel kernel, Particles *particles, int *interactions, int numParticles);
221#endif
222
223}
224
225
226#endif //MILUPHPC_KERNEL_CUH
Particles
Particle(s) class based on SoA (Structur of Arrays).
Definition: particles.cuh:50
CUDA_CALLABLE_MEMBER
#define CUDA_CALLABLE_MEMBER
Definition: cuda_utilities.cuh:30
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
real
double real
Definition: parameter.h:15
DIM
#define DIM
Dimension of the problem.
Definition: parameter.h:38

milupHPC - include/sph/kernel.cuh Source File
Generated on Wed Aug 31 2022 12:16:52 by Doxygen 1.9.3