cuIBM
A GPU-based Immersed Boundary Method code
|
Contains all the custom-written CUDA kernels. More...
Functions | |
__global__ void | dragLeftRight (real *FxX, real *q, real *lambda, real nu, real *dx, real *dy, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculates drag using a control-volume approach (left-right). More... | |
__global__ void | dragBottomTop (real *FxY, real *q, real nu, real *dx, real *dy, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculate drag using a control-volume approach (bottom-top). More... | |
__global__ void | dragUnsteady (real *FxU, real *q, real *qOld, real *dx, real *dy, real dt, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculate drag using a control-volume approach (unsteady). More... | |
__global__ void | liftLeftRight (real *FyX, real *q, real nu, real *dx, real *dy, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculate lift using a control-volume approach (left-right). More... | |
__global__ void | liftBottomTop (real *FyY, real *q, real *lambda, real nu, real *dx, real *dy, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculate lift using a control-volume approach (bottom-top). More... | |
__global__ void | liftUnsteady (real *FyU, real *q, real *qOld, real *dx, real *dy, real dt, int nx, int ny, int I, int J, int ncx, int ncy) |
Calculate lift using a control-volume approach (unsteady). More... | |
__global__ void | forceX (real *f, real *q, real *rn, int *tags, int nx, int ny, real *dx, real *dy, real dt, real alpha, real nu) |
Kernel not usable. More... | |
__global__ void | forceY () |
Kernel not usable. More... | |
__global__ void | generateA (int *ARows, int *ACols, real *AVals, real *MVals, int *LRows, int *LCols, real *LVals, int ASize, real alpha) |
Generates a block of the matrix resulting from implicit terms in the momentum equation. More... | |
__global__ void | generateADirectForcing (int *ARows, int *ACols, real *AVals, real *MVals, int *LRows, int *LCols, real *LVals, int ASize, real alpha, int *tags) |
Generates a block of the matrix resulting from implicit terms in the momentum equation for the direct forcing method. More... | |
__global__ void | bc1DirichletU (real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a Dirichlet-type condition. More... | |
__global__ void | bc1DirichletV (real *bc1, int N, int nx, int numU, int offset, int stride, real *dy, real C, real *bc, int numUbc) |
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given boundary with a Dirichlet-type condition. More... | |
__global__ void | bc1ConvectiveU (real *bc1, int N, int nx, int offset, int stride, real *dx, real *dy, real C, real *bc, real *q, real beta) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a convective-type condition. More... | |
__global__ void | bc1ConvectiveV (real *bc1, int N, int nx, int numU, int offset, int stride, real *dx, real *dy, real C, real *bc, int numUbc, real *q, real beta) |
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given boundary with a convective-type condition. More... | |
__global__ void | bc1SpecialU (real *bc1, int N, int nx, int offset, int stride, real *dx, real C, real *bc, real time) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a special-type condition. More... | |
__global__ void | fillBC2_v (real *bc2, real *yminus, real *yplus, real *dx, int nx, int ny) |
Computes inhomogeneous terms of the discrete divergence operator from the bottom and top boundaries at the v-velocity locations. More... | |
__global__ void | fillBC2_u (real *bc2, real *xminus, real *xplus, real *dy, int nx, int ny) |
Computes inhomogeneous terms of the discrete divergence operator from the left and right boundaries at the u-velocity locations. More... | |
__global__ void | fillBC2_uvB (real *bc2, real *uB, real *vB, int totalPoints, int nx, int ny) |
Computes inhomogeneous terms of the discrete divergence operator from the no-slip constraint at the body-point locations. More... | |
__global__ void | generateE (int *ERows, int *ECols, real *EVals, int nx, int ny, real *x, real *y, real *dx, int totalPoints, real *xB, real *yB, int *I, int *J) |
Computes elements of the interpolation matrix. More... | |
__global__ void | fillM_u (int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt) |
Computes an element of the mass matrix and its inverse for a x-velocity node. More... | |
__global__ void | fillM_v (int *MRows, int *MCols, real *MVals, int *MinvRows, int *MinvCols, real *MinvVals, int nx, int ny, real *dx, real *dy, real dt) |
Computes an element of the mass matrix and its inverse for a y-velocity node. More... | |
__global__ void | updateQ (int *QRows, int *QCols, real *QVals, int QSize, int *tags) |
Update the gradient operator at forcing nodes. More... | |
__global__ void | updateQT (int *QTRows, int *QTCols, real *QTVals, int QTSize, int *tags, real *coeffs) |
Update the divergence operator at forcing nodes. More... | |
void | generateQT (int *QTRows, int *QTCols, real *QTVals, int nx, int ny) |
Generates the divergence matrix (on the host). More... | |
__global__ void | updateQT (int *QTRows, int *QTCols, real *QTVals, int *ERows, int *ECols, real *EVals, int nx, int ny, real *x, real *y, real *dx, int totalPoints, real *xB, real *yB, int *I, int *J) |
Updates elements of the divergence matrix and the interpolation matrix. More... | |
__global__ void | convectionTermU (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms. More... | |
__global__ void | convectionTermV (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms. More... | |
__global__ void | convectionTermVBottomTop (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the bottom or top boundaries. More... | |
__global__ void | convectionTermUBottomTop (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop, real *bcLeft, real *bcRight) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the bottom or top boundaries. More... | |
__global__ void | convectionTermULeftRight (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcLeft, real *bcRight) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the left or right boundaries. More... | |
__global__ void | convectionTermVLeftRight (real *rn, real *H, real *q, int nx, int ny, real *dx, real *dy, real dt, real gamma, real zeta, real alpha, real nu, real *bcBottom, real *bcTop, real *bcLeft, real *bcRight) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the left or right boundaries. More... | |
__global__ void | fill_velB (real *velB, real *uB, real *vB, int totalPoints) |
Stores an element of the u- and v- body-velocities into one single array. More... | |
__global__ void | updateRHS1 (real *rhs1, int numUV, int *tags) |
Update the RHS vector of the velocity system at forcing nodes. More... | |
__global__ void | updateRHS1X (real *rhs1, int nx, int ny, real dt, real *dx, int *tags, real *coeffs, real *uv) |
Perform 1D linear interpolation at forcing points for the x-velocity. More... | |
__global__ void | updateRHS1Y (real *rhs1, int nx, int ny, real dt, real *dy, int *tags, real *coeffs, real *uv) |
Perform 1D linear interpolation at forcing points for the y-velocity. More... | |
__global__ void | updateRHS1X (real *rhs1, int nx, int ny, real dt, real *dx, int *tags, real *coeffs, real *coeffs2, real *uv) |
Perform 1D quadratic interpolation at forcing points for the x-velocity. More... | |
__global__ void | updateRHS1Y (real *rhs1, int nx, int ny, real dt, real *dy, int *tags, real *coeffs, real *coeffs2, real *uv) |
Perform 1D quadratic interpolation at forcing points for the y-velocity. More... | |
Contains all the custom-written CUDA kernels.
Contains all custom-written CUDA kernels.
namespace kernels
__global__ void kernels::bc1ConvectiveU | ( | real * | bc1, |
int | N, | ||
int | nx, | ||
int | offset, | ||
int | stride, | ||
real * | dx, | ||
real * | dy, | ||
real | C, | ||
real * | bc, | ||
real * | q, | ||
real | beta | ||
) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a convective-type condition.
Each element is also normalized with the corresponding diagonal element of the matrix .
bc1 | array that contains the boundary conditions |
N | number of u-velocity points on the boundary |
nx | number of cells in the x-direction |
offset | index in vector bc1 of the first element u-velocity point on the boundary |
stride | index increment |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
C | coefficient from the Laplacian discretization at the boundary |
bc | boundary velocity |
q | flux vector |
beta | linear interpolation coefficient (U*dt/dx) |
check if idx too high
Definition at line 98 of file generateBC1.cu.
__global__ void kernels::bc1ConvectiveV | ( | real * | bc1, |
int | N, | ||
int | nx, | ||
int | numU, | ||
int | offset, | ||
int | stride, | ||
real * | dx, | ||
real * | dy, | ||
real | C, | ||
real * | bc, | ||
int | numUbc, | ||
real * | q, | ||
real | beta | ||
) |
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given boundary with a convective-type condition.
Each element is also normalized with the corresponding diagonal element of the matrix .
bc1 | array that contains the boundary conditions |
N | number of v-velocity points on the boundary |
nx | number of cells in the x-direction |
numU | number of u-velocity points in the domain |
offset | index in vector bc1 of the first element v-velocity point on the boundary |
stride | index increment |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
C | coefficient from the Laplacian discretization at the boundary |
bc | boundary velocity |
numUbc | number of u-velocity points on the boundary |
q | flux vector |
beta | linear interpolation coefficient (U*dt/dx) |
check if idx too high
Definition at line 136 of file generateBC1.cu.
__global__ void kernels::bc1DirichletU | ( | real * | bc1, |
int | N, | ||
int | nx, | ||
int | offset, | ||
int | stride, | ||
real * | dx, | ||
real | C, | ||
real * | bc | ||
) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a Dirichlet-type condition.
Each element is also normalized with the corresponding diagonal element of the matrix .
bc1 | array that contains the boundary conditions |
N | number of u-velocity points on the boundary |
nx | number of cells in the x-direction |
offset | index in vector bc1 of the first element u-velocity point on the boundary |
stride | index increment |
dx | cell-widths in the x-direction |
C | coefficient from the Laplacian discretization at the boundary |
bc | boundary velocity |
check if idx too high
Definition at line 35 of file generateBC1.cu.
__global__ void kernels::bc1DirichletV | ( | real * | bc1, |
int | N, | ||
int | nx, | ||
int | numU, | ||
int | offset, | ||
int | stride, | ||
real * | dy, | ||
real | C, | ||
real * | bc, | ||
int | numUbc | ||
) |
Computes inhomogeneous term from the discrete Laplacian operator for the v-velocity at a given boundary with a Dirichlet-type condition.
Each element is also normalized with the corresponding diagonal element of the matrix .
bc1 | array that contains the boundary conditions |
N | number of v-velocity points on the boundary |
nx | number of cells in the x-direction |
numU | number of u-velocity points in the domain |
offset | index in vector bc1 of the first element v-velocity point on the boundary |
stride | index increment |
dy | cell-widths in the y-direction |
C | coefficient from the Laplacian discretization at the boundary |
bc | boundary velocity |
numUbc | number of u-velocity points on the boundary |
check if idx too high
Definition at line 66 of file generateBC1.cu.
__global__ void kernels::bc1SpecialU | ( | real * | bc1, |
int | N, | ||
int | nx, | ||
int | offset, | ||
int | stride, | ||
real * | dx, | ||
real | C, | ||
real * | bc, | ||
real | time | ||
) |
Computes inhomogeneous term from the discrete Laplacian operator for the u-velocity at a given boundary with a special-type condition.
It computes a sinusoidal displacement of the boundary in the x-direction. Each element is also normalized with the corresponding diagonal element of the matrix .
bc1 | array that contains the boundary conditions |
N | number of v-velocity points on the boundary |
nx | number of cells in the x-direction |
offset | index in vector bc1 of the first element v-velocity point on the boundary |
stride | index increment |
dx | cell-widths in the x-direction |
C | coefficient from the Laplacian discretization at the boundary |
bc | boundary velocity |
time | the current time |
check if idx too high
Definition at line 171 of file generateBC1.cu.
__global__ void kernels::convectionTermU | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu | ||
) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
transfer from global to shared memory
check bounds for convective term in the x-direction
< check if we compute globally
< check if element within block computes
X-component
Definition at line 38 of file generateRN.cu.
References BSZ.
__global__ void kernels::convectionTermUBottomTop | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu, | ||
real * | bcBottom, | ||
real * | bcTop, | ||
real * | bcLeft, | ||
real * | bcRight | ||
) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the bottom or top boundaries.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
bcBottom | bottom-boundary velocity |
bcTop | top-boundary velocity |
bcLeft | left-boundary velocity |
bcRight | right-boundary velocity |
Boundary Conditions for u at the Bottom ***********************************************************************************
Calculate the Diffusion Term and Velocity Components for the convection term
Convection Term
Calculate rn
Boundary conditions for u at the Top **************************************************************************************
Calculate the Diffusion Term and Velocity Components for the convection term
Convection Term
Calculate rn
Definition at line 291 of file generateRN.cu.
__global__ void kernels::convectionTermULeftRight | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu, | ||
real * | bcLeft, | ||
real * | bcRight | ||
) |
Computes u-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the left or right boundaries.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
bcLeft | left-boundary velocity |
bcRight | right-boundary velocity |
Boundary Conditions for u at the Left *************************************************************************************
Convection Term
Diffusion Term
Calculate rn
Boundary conditions for u at the Right ************************************************************************************
Convection Term
Diffusion Term
Calculate rn
Definition at line 400 of file generateRN.cu.
__global__ void kernels::convectionTermV | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu | ||
) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
transfer from global to shared memory
check bounds for convective term in the x-direction
< check if we compute globally
< check if element within block computes
Y-component
Definition at line 122 of file generateRN.cu.
References BSZ.
__global__ void kernels::convectionTermVBottomTop | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu, | ||
real * | bcBottom, | ||
real * | bcTop | ||
) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the bottom or top boundaries.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
bcBottom | bottom-boundary velocity |
bcTop | top-boundary velocity |
Boundary Conditions for v at the Bottom ***********************************************************************************
Convection Term
2nd order Adams-Bashforth
Diffusion Term
Calculate rn
Boundary conditions for v at the Top **************************************************************************************
Convection Term
2nd order Adams-Bashforth
Diffusion Term
Calculate rn
Definition at line 211 of file generateRN.cu.
__global__ void kernels::convectionTermVLeftRight | ( | real * | rn, |
real * | H, | ||
real * | q, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | gamma, | ||
real | zeta, | ||
real | alpha, | ||
real | nu, | ||
real * | bcBottom, | ||
real * | bcTop, | ||
real * | bcLeft, | ||
real * | bcRight | ||
) |
Computes v-component of the vector resulting from the explicit terms of the discretized momentum equation, and element of the explcit convective terms on the left or right boundaries.
rn | explicit terms of the discretized momentum equation |
H | explicit convective terms of the discretized momentum equation |
q | velcoity flux vector |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
gamma | coefficient of the convection term at the current time-step |
zeta | coefficient of the convection term at the previous time-step |
alpha | coefficient of the explicit diffusion term |
nu | viscosity |
bcBottom | bottom-boundary velocity |
bcTop | top-boundary velocity |
bcLeft | left-boundary velocity |
bcRight | right-boundary velocity |
Boundary Conditions for v at the Left *************************************************************************************
Convection Term
Diffusion Term
Calculate rn
Boundary Conditions for v at the Right ************************************************************************************
Convection Term
Diffusion Term
Calculate rn
Definition at line 480 of file generateRN.cu.
__global__ void kernels::dragBottomTop | ( | real * | FxY, |
real * | q, | ||
real | nu, | ||
real * | dx, | ||
real * | dy, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculate drag using a control-volume approach (bottom-top).
Evaluate the contribution from the bottom and top parts of the control surface.
FxY | raw pointer to the vector storing the drag in the y-direction |
q | raw pointer to the vector storing all the fluxes |
nu | viscosity |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
I | x-index of the bottom-left corner cell of the control surface |
J | y-index of the top-right corner cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
ncy | number of cells in the y-direction in the control volume |
Definition at line 87 of file calculateForce.cu.
__global__ void kernels::dragLeftRight | ( | real * | FxX, |
real * | q, | ||
real * | lambda, | ||
real | nu, | ||
real * | dx, | ||
real * | dy, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculates drag using a control-volume approach (left-right).
Evaluate the contribution from the left and right parts of the control surface.
FxX | raw pointer to the vector storing the drag in the x-direction |
lambda | raw pointer to the vector storing all the pressure and Lagrangian forces |
q | raw pointer to the vector storing all the fluxes |
nu | viscosity |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
I | x-index of the bottom-left corner cell of the control surface |
J | y-index of the top-right corner cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
ncy | number of cells in the y-direction in the control volume |
Definition at line 41 of file calculateForce.cu.
__global__ void kernels::dragUnsteady | ( | real * | FxU, |
real * | q, | ||
real * | qOld, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculate drag using a control-volume approach (unsteady).
Evaluate the unsteady contribution of the control volume.
FxU | raw pointer to the vector storing the unsteady drag components |
q | raw pointer to the vector storing all the fluxes |
qOld | raw pointer to the vector sotring all the fluxes at the previous time-step |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
dt | time increment |
nx | number of cells in the x-direction |
ny | number of cells in the y-direcyion |
I | x-index of the bottom-left cell of the control surface |
J | y-index of the top-right cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
nyc | number of cells in the y-direction in the control volume |
Definition at line 142 of file calculateForce.cu.
Stores an element of the u- and v- body-velocities into one single array.
velB | vector that contains both u- and v- velocities |
uB | u-velocity of body points (all bodies included) |
vB | v-velocity of body points (all bodies included) |
totalPoints | number of body points (all bodies included) |
Definition at line 26 of file generateVelB.cu.
__global__ void kernels::fillBC2_u | ( | real * | bc2, |
real * | xminus, | ||
real * | xplus, | ||
real * | dy, | ||
int | nx, | ||
int | ny | ||
) |
Computes inhomogeneous terms of the discrete divergence operator from the left and right boundaries at the u-velocity locations.
bc2 | array that contains boundary conditions |
xminus | left-boundary velocities |
xplus | right-boundary velocities |
dy | cell-widths in the x-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
Definition at line 52 of file generateBC2.cu.
__global__ void kernels::fillBC2_uvB | ( | real * | bc2, |
real * | uB, | ||
real * | vB, | ||
int | totalPoints, | ||
int | nx, | ||
int | ny | ||
) |
Computes inhomogeneous terms of the discrete divergence operator from the no-slip constraint at the body-point locations.
bc2 | array that contains boundary conditions |
uB | x-component of the body-velocity |
vB | y-component of the body-velcoity |
totalPoints | number of body-points (all bodies included) |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
Definition at line 74 of file generateBC2.cu.
__global__ void kernels::fillBC2_v | ( | real * | bc2, |
real * | yminus, | ||
real * | yplus, | ||
real * | dx, | ||
int | nx, | ||
int | ny | ||
) |
Computes inhomogeneous terms of the discrete divergence operator from the bottom and top boundaries at the v-velocity locations.
bc2 | array that contains boundary conditions |
yminus | bottom-boundary velocities |
yplus | top-boundary velocities |
dx | cell-widths in the x-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
Definition at line 30 of file generateBC2.cu.
__global__ void kernels::fillM_u | ( | int * | MRows, |
int * | MCols, | ||
real * | MVals, | ||
int * | MinvRows, | ||
int * | MinvCols, | ||
real * | MinvVals, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt | ||
) |
Computes an element of the mass matrix and its inverse for a x-velocity node.
MRows | row index of elements of the mass matrix |
MCols | column index of elements of the mass matrix |
MVals | value of elements of the mass matrix |
MinvRows | row index of elements of the mass matrix inverse |
MinvCols | column index of elements of the mass matrix inverse |
MinvVals | value of elements of the mass matrix inverse |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
Definition at line 33 of file generateM.cu.
__global__ void kernels::fillM_v | ( | int * | MRows, |
int * | MCols, | ||
real * | MVals, | ||
int * | MinvRows, | ||
int * | MinvCols, | ||
real * | MinvVals, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt | ||
) |
Computes an element of the mass matrix and its inverse for a y-velocity node.
MRows | row index of elements of the mass matrix |
MCols | column index of elements of the mass matrix |
MVals | value of elements of the mass matrix |
MinvRows | row index of elements of the mass matrix inverse |
MinvCols | column index of elements of the mass matrix inverse |
MinvVals | value of elements of the mass matrix inverse |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dx | cell-widths in the x-direction |
dy | cell-widths in the y-direction |
dt | time-increment |
Definition at line 69 of file generateM.cu.
__global__ void kernels::forceX | ( | real * | f, |
real * | q, | ||
real * | rn, | ||
int * | tags, | ||
int | nx, | ||
int | ny, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
real | alpha, | ||
real | nu | ||
) |
Kernel not usable.
transfer from global to shared memory
check bounds for convective term in the x-direction
< check if we compute globally
< check if element within block computes
X-component
Definition at line 299 of file calculateForce.cu.
References BSZ.
Referenced by NSWithBody< memoryType >::calculateForce(), bodies< memoryType >::initialise(), and DirectForcingSolver< memoryType >::writeData().
__global__ void kernels::forceY | ( | ) |
Kernel not usable.
Definition at line 354 of file calculateForce.cu.
Referenced by NSWithBody< memoryType >::calculateForce(), and bodies< memoryType >::initialise().
__global__ void kernels::generateA | ( | int * | ARows, |
int * | ACols, | ||
real * | AVals, | ||
real * | MVals, | ||
int * | LRows, | ||
int * | LCols, | ||
real * | LVals, | ||
int | ASize, | ||
real | alpha | ||
) |
Generates a block of the matrix resulting from implicit terms in the momentum equation.
It assembles the matrix A
as a combination of the Laplacian matrix L
and the mass matrix M
. A = M-alpha*L The parameter alpha is the coefficient of the implicit part of the diffusion term. It is 1 for a backward Euler scheme, 0.5 for a Crank-Nicolson scheme, and 0 for a fully explicit scheme.
ARows | rows of the COO matrix A |
ACols | columns of the COO matrix A |
AVals | values of the COO matrix A |
MVals | values of the COO matrix M |
LRows | rows of the COO matrix L |
LCols | columns of the COO matrix L |
LVals | values of the COO matrix A |
ASize | number of entries of the COO matrix A |
alpha | implicit coefficient of the diffusive scheme |
Definition at line 38 of file generateA.cu.
Referenced by NavierStokesSolver< memoryType >::assembleMatrices().
__global__ void kernels::generateADirectForcing | ( | int * | ARows, |
int * | ACols, | ||
real * | AVals, | ||
real * | MVals, | ||
int * | LRows, | ||
int * | LCols, | ||
real * | LVals, | ||
int | ASize, | ||
real | alpha, | ||
int * | tags | ||
) |
Generates a block of the matrix resulting from implicit terms in the momentum equation for the direct forcing method.
It assembles the matrix A
as a combination of the Laplacian matrix L
and the mass matrix M
. The parameter alpha is the coefficient of the implicit part of the diffusion term. It is 1 for a backward Euler scheme, 0.5 for a Crank-Nicolson scheme, and 0 for a fully explicit scheme. The left-hand side matrix A is set up as M-alpha*L, where M is the mass matrix, and L the Laplacian matrix. But in the case of the direct forcing method, some rows are determined by interpolation relations, and the rows of L are modified appropriately. For these rows alone, the rows of A are given by M-L.
ARows | rows of the COO matrix A |
ACols | columns of the COO matrix A |
AVals | values of the COO matrix A |
MVals | values of the COO matrix M |
LRows | rows of the COO matrix L |
LCols | columns of the COO matrix L |
LVals | values of the COO matrix A |
ASize | number of entries of the COO matrix A |
alpha | implicit coefficient of the diffusive scheme |
tagsX | tag to check if the node is next to an immersed boundary |
tagsY | tag to check if the node is next to an immersed boundary |
Definition at line 77 of file generateA.cu.
__global__ void kernels::generateE | ( | int * | ERows, |
int * | ECols, | ||
real * | EVals, | ||
int | nx, | ||
int | ny, | ||
real * | x, | ||
real * | y, | ||
real * | dx, | ||
int | totalPoints, | ||
real * | xB, | ||
real * | yB, | ||
int * | I, | ||
int * | J | ||
) |
Computes elements of the interpolation matrix.
ERows | row index of elements of the interpolation matrix |
ECols | column index of elements of the interpolation matrix |
EVals | value of elements of the interpolation matrix |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
x | x-component of grid points |
y | y-component of grid points |
dx | cell-widths in the x-direction |
totalPoints | number of body points (all bodies included) |
xB | x-coordinate of body points (all bodies included) |
yB | y-coordinate of body points (all bodies included) |
I | x-index of the cell in which the body point is located |
J | y-index of the cell in which the body point is located |
Definition at line 73 of file generateE.cu.
References deltaDeviceE().
void kernels::generateQT | ( | int * | QTRows, |
int * | QTCols, | ||
real * | QTVals, | ||
int | nx, | ||
int | ny | ||
) |
Generates the divergence matrix (on the host).
QTRows | row index of elements of the matrix |
QTCols | column index of elements of the matrix |
QTVals | value of elements of the matrix |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
QT is an (np + 2*nb) x nuv matrix
Generate the GT part
Definition at line 110 of file generateQT.cu.
Referenced by NavierStokesSolver< memoryType >::assembleMatrices(), NavierStokesSolver< memoryType >::generateQT(), and TairaColoniusSolver< memoryType >::generateQT().
__global__ void kernels::liftBottomTop | ( | real * | FyY, |
real * | q, | ||
real * | lambda, | ||
real | nu, | ||
real * | dx, | ||
real * | dy, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculate lift using a control-volume approach (bottom-top).
Evaluate the contribution from the bottom and top parts of the control surface.
FyY | raw pointer to the vector storing the lift components in the y-direction |
q | raw pointer to the vector storing all the fluxes |
lambda | raw pointer to the vector storing the pressure and Lagrangian forces |
nu | viscosity |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direcyion |
I | x-index of the bottom-left cell of the control surface |
J | y-index of the top-right cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
nyc | number of cells in the y-direction in the control volume |
Definition at line 231 of file calculateForce.cu.
__global__ void kernels::liftLeftRight | ( | real * | FyX, |
real * | q, | ||
real | nu, | ||
real * | dx, | ||
real * | dy, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculate lift using a control-volume approach (left-right).
Evaluate the contribution from the left and right parts of the control surface.
FyX | raw pointer to the vector storing the lift components in the x-direction |
q | raw pointer to the vector storing all the fluxes |
nu | viscosity |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
nx | number of cells in the x-direction |
ny | number of cells in the y-direcyion |
I | x-index of the bottom-left cell of the control surface |
J | y-index of the top-right cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
nyc | number of cells in the y-direction in the control volume |
Definition at line 177 of file calculateForce.cu.
__global__ void kernels::liftUnsteady | ( | real * | FyU, |
real * | q, | ||
real * | qOld, | ||
real * | dx, | ||
real * | dy, | ||
real | dt, | ||
int | nx, | ||
int | ny, | ||
int | I, | ||
int | J, | ||
int | ncx, | ||
int | ncy | ||
) |
Calculate lift using a control-volume approach (unsteady).
Evaluate the unsteady contribution of the control volume.
FyU | raw pointer to the vector storing the unsteady lift components |
q | raw pointer to the vector storing all the fluxes |
qOld | raw pointer to the vector sotring all the fluxes at the previous time-step |
dx | raw pointer to the vector storing the cell widths in the x-direction |
dy | raw pointer to the vector storing the cell widths in the y-direction |
dt | time increment |
nx | number of cells in the x-direction |
ny | number of cells in the y-direcyion |
I | x-index of the bottom-left cell of the control surface |
J | y-index of the top-right cell of the control surface |
ncx | number of cells in the x-direction in the control volume |
nyc | number of cells in the y-direction in the control volume |
Definition at line 278 of file calculateForce.cu.
__global__ void kernels::updateQ | ( | int * | QRows, |
int * | QCols, | ||
real * | QVals, | ||
int | QSize, | ||
int * | tags | ||
) |
Update the gradient operator at forcing nodes.
QRows | row index of elements of the matrix |
QCols | column index of elements of the matrix |
QVals | value of elements of the matrix |
QSize | size of the matrix |
tags | tags of regular and forcing nodes |
Definition at line 65 of file generateQT.cu.
Referenced by DirectForcingSolver< memoryType >::generateQT().
__global__ void kernels::updateQT | ( | int * | QTRows, |
int * | QTCols, | ||
real * | QTVals, | ||
int | QTSize, | ||
int * | tags, | ||
real * | coeffs | ||
) |
Update the divergence operator at forcing nodes.
QTRows | row index of elements of the matrix |
QTCols | column index of elements of the matrix |
QTVals | value of elements of the matrix |
QTSize | size of the matrix |
tags | tags of regular and forcing nodes |
coeffs |
Definition at line 87 of file generateQT.cu.
Referenced by TairaColoniusSolver< memoryType >::generateQT(), and TairaColoniusSolver< memoryType >::updateSolverState().
__global__ void kernels::updateQT | ( | int * | QTRows, |
int * | QTCols, | ||
real * | QTVals, | ||
int * | ERows, | ||
int * | ECols, | ||
real * | EVals, | ||
int | nx, | ||
int | ny, | ||
real * | x, | ||
real * | y, | ||
real * | dx, | ||
int | totalPoints, | ||
real * | xB, | ||
real * | yB, | ||
int * | I, | ||
int * | J | ||
) |
Updates elements of the divergence matrix and the interpolation matrix.
QTRows | row index of elements of the matrix |
QTCols | column index of elements of the matrix |
QTVals | value of elements of the matrix |
ERows | row index of elements of the interpolation matrix |
ECols | column index of elements of the interpolation matrix |
EVals | value of elements of the interpolation matrix |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
x | x-component of grid points |
y | y-component of grid points |
dx | cell-widths in the x-direction |
totalPoints | number of body-points (all bodies included) |
xB | x-component of body-points (all bodies included) |
yB | y-component of body-points (all bodies included) |
I | x-index of grid cells in which body points are located |
J | y-index of grid cells in which body points are located |
Definition at line 184 of file generateQT.cu.
References deltaDeviceQT().
__global__ void kernels::updateRHS1 | ( | real * | rhs1, |
int | numUV, | ||
int * | tags | ||
) |
Update the RHS vector of the velocity system at forcing nodes.
Elements associated with a forcing node are set to 0.
rhs1 | vector to update |
numUV | number of velocity points in the domain |
tags | vector used to differentiate regular points from forcing points |
Definition at line 30 of file updateRHS1.cu.
Referenced by DirectForcingSolver< memoryType >::assembleRHS1().
__global__ void kernels::updateRHS1X | ( | real * | rhs1, |
int | nx, | ||
int | ny, | ||
real | dt, | ||
real * | dx, | ||
int * | tags, | ||
real * | coeffs, | ||
real * | uv | ||
) |
Perform 1D linear interpolation at forcing points for the x-velocity.
rhs1 | vector to update |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dt | time-step size |
dx | grid-spacings along a gridline in the x-direction |
tags | vector used to differentiate regular points from forcing points |
coeffs | coefficients of interpolation |
uv | velocity vector |
Definition at line 54 of file updateRHS1.cu.
__global__ void kernels::updateRHS1X | ( | real * | rhs1, |
int | nx, | ||
int | ny, | ||
real | dt, | ||
real * | dx, | ||
int * | tags, | ||
real * | coeffs, | ||
real * | coeffs2, | ||
real * | uv | ||
) |
Perform 1D quadratic interpolation at forcing points for the x-velocity.
rhs1 | vector to update |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dt | time-step size |
dx | grid-spacings along a gridline in the x-direction |
tags | vector used to differentiate regular points from forcing points |
coeffs | coefficients of interpolation |
coeffs2 | coefficients of interpolation |
uv | velocity vector |
Definition at line 108 of file updateRHS1.cu.
__global__ void kernels::updateRHS1Y | ( | real * | rhs1, |
int | nx, | ||
int | ny, | ||
real | dt, | ||
real * | dy, | ||
int * | tags, | ||
real * | coeffs, | ||
real * | uv | ||
) |
Perform 1D linear interpolation at forcing points for the y-velocity.
rhs1 | vector to update |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dt | time-step size |
dy | grid-spacings along a gridline in the y-direction |
tags | vector used to differentiate regular points from forcing points |
coeffs | coefficients of interpolation |
uv | velocity vector |
Definition at line 80 of file updateRHS1.cu.
__global__ void kernels::updateRHS1Y | ( | real * | rhs1, |
int | nx, | ||
int | ny, | ||
real | dt, | ||
real * | dy, | ||
int * | tags, | ||
real * | coeffs, | ||
real * | coeffs2, | ||
real * | uv | ||
) |
Perform 1D quadratic interpolation at forcing points for the y-velocity.
rhs1 | vector to update |
nx | number of cells in the x-direction |
ny | number of cells in the y-direction |
dt | time-step size |
dy | grid-spacings along a gridline in the y-direction |
tags | vector used to differentiate regular points from forcing points |
coeffs | coefficients of interpolation |
coeffs2 | coefficients of interpolation |
uv | velocity vector |
Definition at line 135 of file updateRHS1.cu.