-
Notifications
You must be signed in to change notification settings - Fork 9
/
kernel.h
161 lines (157 loc) · 3.97 KB
/
kernel.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
/**
* @author Christoph Schaefer [email protected]
*
* @section LICENSE
* Copyright (c) 2019 Christoph Schaefer
*
* This file is part of miluphcuda.
*
* miluphcuda is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* miluphcuda is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with miluphcuda. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef _KERNEL_H
#define _KERNEL_H
#include "parameter.h"
/**
* @brief Calculates b-spline kernel and derivatives for an interaction
*
*/
__device__ void fastKernelvalueAndDerivative(
double &W, double &dWdx
#if DIM > 1
, double &dWdy
#endif
#if DIM == 3
, double &dWdz
#endif
, double &dWdr, int particle1, int particle2,
double dx
#if DIM > 1
, double dy
#endif
#if DIM == 3
, double dz
#endif
);
/// function pointer to generic SPH kernel function
/// cuda allows function pointers since Fermi architecture.
/// however, we need a typedef to set the function pointers
typedef void (*SPH_kernel) (double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief *The* standard cubic b-spline.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void cubic_spline(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Spiky kernel.
* @details only implemented for 2D and 3D.
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void spiky(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Quartic spline kernel.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void quartic_spline(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Quintic spline kernel.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void quintic_spline(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Wendland C2 kernel.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void wendlandc2(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Wendland C4 kernel.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void wendlandc4(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Wendland C6 kernel.
*
* @param W
* @param dWdx
* @param dWdr
* @param dx
* @param h
* @return __device__
*/
__device__ void wendlandc6(double *W, double dWdx[DIM], double *dWdr, double dx[DIM], double h);
/**
* @brief Calculates the kernel for the tensile instability fix following Monaghan 2000.
*
* @param particle1
* @param particle2
* @return __device__
*/
__device__ double fixTensileInstability( int particle1, int particle2);
/**
* @brief Calculates the tensorial correction factors for linear consistency.
*
* @param interactions
* @return __global__
*/
__global__ void tensorialCorrection(int *interactions);
/**
* @brief Calculates the zeroth order corrections for the kernel sum.
*
* @param interactions
* @return __global__
*/
__global__ void shepardCorrection(int *interactions);
/**
* @brief Calculates \f$ \nabla \cdot \vec{v} \f$ and \f$ \nabla \times \vec{v} \f$
*
* @param interactions
* @return __global__
*/
__global__ void CalcDivvandCurlv(int *interactions);
#endif