cuIBM
A GPU-based Immersed Boundary Method code
Functions
kernels Namespace Reference

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...
 

Detailed Description

Contains all the custom-written CUDA kernels.

Contains all custom-written CUDA kernels.

namespace kernels

Function Documentation

__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 $ \hat{M} $.

Parameters
bc1array that contains the boundary conditions
Nnumber of u-velocity points on the boundary
nxnumber of cells in the x-direction
offsetindex in vector bc1 of the first element u-velocity point on the boundary
strideindex increment
dxcell-widths in the x-direction
dycell-widths in the y-direction
Ccoefficient from the Laplacian discretization at the boundary
bcboundary velocity
qflux vector
betalinear 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 $ \hat{M} $.

Parameters
bc1array that contains the boundary conditions
Nnumber of v-velocity points on the boundary
nxnumber of cells in the x-direction
numUnumber of u-velocity points in the domain
offsetindex in vector bc1 of the first element v-velocity point on the boundary
strideindex increment
dxcell-widths in the x-direction
dycell-widths in the y-direction
Ccoefficient from the Laplacian discretization at the boundary
bcboundary velocity
numUbcnumber of u-velocity points on the boundary
qflux vector
betalinear 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 $ \hat{M} $.

Parameters
bc1array that contains the boundary conditions
Nnumber of u-velocity points on the boundary
nxnumber of cells in the x-direction
offsetindex in vector bc1 of the first element u-velocity point on the boundary
strideindex increment
dxcell-widths in the x-direction
Ccoefficient from the Laplacian discretization at the boundary
bcboundary 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 $ \hat{M} $.

Parameters
bc1array that contains the boundary conditions
Nnumber of v-velocity points on the boundary
nxnumber of cells in the x-direction
numUnumber of u-velocity points in the domain
offsetindex in vector bc1 of the first element v-velocity point on the boundary
strideindex increment
dycell-widths in the y-direction
Ccoefficient from the Laplacian discretization at the boundary
bcboundary velocity
numUbcnumber 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 $ \hat{M} $.

Parameters
bc1array that contains the boundary conditions
Nnumber of v-velocity points on the boundary
nxnumber of cells in the x-direction
offsetindex in vector bc1 of the first element v-velocity point on the boundary
strideindex increment
dxcell-widths in the x-direction
Ccoefficient from the Laplacian discretization at the boundary
bcboundary velocity
timethe 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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity

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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity
bcBottombottom-boundary velocity
bcToptop-boundary velocity
bcLeftleft-boundary velocity
bcRightright-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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity
bcLeftleft-boundary velocity
bcRightright-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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity

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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity
bcBottombottom-boundary velocity
bcToptop-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.

Parameters
rnexplicit terms of the discretized momentum equation
Hexplicit convective terms of the discretized momentum equation
qvelcoity flux vector
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-increment
gammacoefficient of the convection term at the current time-step
zetacoefficient of the convection term at the previous time-step
alphacoefficient of the explicit diffusion term
nuviscosity
bcBottombottom-boundary velocity
bcToptop-boundary velocity
bcLeftleft-boundary velocity
bcRightright-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.

Parameters
FxYraw pointer to the vector storing the drag in the y-direction
qraw pointer to the vector storing all the fluxes
nuviscosity
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
Ix-index of the bottom-left corner cell of the control surface
Jy-index of the top-right corner cell of the control surface
ncxnumber of cells in the x-direction in the control volume
ncynumber 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.

Parameters
FxXraw pointer to the vector storing the drag in the x-direction
lambdaraw pointer to the vector storing all the pressure and Lagrangian forces
qraw pointer to the vector storing all the fluxes
nuviscosity
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
Ix-index of the bottom-left corner cell of the control surface
Jy-index of the top-right corner cell of the control surface
ncxnumber of cells in the x-direction in the control volume
ncynumber 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.

Parameters
FxUraw pointer to the vector storing the unsteady drag components
qraw pointer to the vector storing all the fluxes
qOldraw pointer to the vector sotring all the fluxes at the previous time-step
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
dttime increment
nxnumber of cells in the x-direction
nynumber of cells in the y-direcyion
Ix-index of the bottom-left cell of the control surface
Jy-index of the top-right cell of the control surface
ncxnumber of cells in the x-direction in the control volume
nycnumber of cells in the y-direction in the control volume

Definition at line 142 of file calculateForce.cu.

__global__ void kernels::fill_velB ( real velB,
real uB,
real vB,
int  totalPoints 
)

Stores an element of the u- and v- body-velocities into one single array.

Parameters
velBvector that contains both u- and v- velocities
uBu-velocity of body points (all bodies included)
vBv-velocity of body points (all bodies included)
totalPointsnumber 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.

Parameters
bc2array that contains boundary conditions
xminusleft-boundary velocities
xplusright-boundary velocities
dycell-widths in the x-direction
nxnumber of cells in the x-direction
nynumber 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.

Parameters
bc2array that contains boundary conditions
uBx-component of the body-velocity
vBy-component of the body-velcoity
totalPointsnumber of body-points (all bodies included)
nxnumber of cells in the x-direction
nynumber 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.

Parameters
bc2array that contains boundary conditions
yminusbottom-boundary velocities
yplustop-boundary velocities
dxcell-widths in the x-direction
nxnumber of cells in the x-direction
nynumber 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.

Parameters
MRowsrow index of elements of the mass matrix
MColscolumn index of elements of the mass matrix
MValsvalue of elements of the mass matrix
MinvRowsrow index of elements of the mass matrix inverse
MinvColscolumn index of elements of the mass matrix inverse
MinvValsvalue of elements of the mass matrix inverse
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-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.

Parameters
MRowsrow index of elements of the mass matrix
MColscolumn index of elements of the mass matrix
MValsvalue of elements of the mass matrix
MinvRowsrow index of elements of the mass matrix inverse
MinvColscolumn index of elements of the mass matrix inverse
MinvValsvalue of elements of the mass matrix inverse
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dxcell-widths in the x-direction
dycell-widths in the y-direction
dttime-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.

Parameters
ARowsrows of the COO matrix A
AColscolumns of the COO matrix A
AValsvalues of the COO matrix A
MValsvalues of the COO matrix M
LRowsrows of the COO matrix L
LColscolumns of the COO matrix L
LValsvalues of the COO matrix A
ASizenumber of entries of the COO matrix A
alphaimplicit 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.

Parameters
ARowsrows of the COO matrix A
AColscolumns of the COO matrix A
AValsvalues of the COO matrix A
MValsvalues of the COO matrix M
LRowsrows of the COO matrix L
LColscolumns of the COO matrix L
LValsvalues of the COO matrix A
ASizenumber of entries of the COO matrix A
alphaimplicit coefficient of the diffusive scheme
tagsXtag to check if the node is next to an immersed boundary
tagsYtag 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.

Parameters
ERowsrow index of elements of the interpolation matrix
EColscolumn index of elements of the interpolation matrix
EValsvalue of elements of the interpolation matrix
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
xx-component of grid points
yy-component of grid points
dxcell-widths in the x-direction
totalPointsnumber of body points (all bodies included)
xBx-coordinate of body points (all bodies included)
yBy-coordinate of body points (all bodies included)
Ix-index of the cell in which the body point is located
Jy-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).

Parameters
QTRowsrow index of elements of the matrix
QTColscolumn index of elements of the matrix
QTValsvalue of elements of the matrix
nxnumber of cells in the x-direction
nynumber 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.

Parameters
FyYraw pointer to the vector storing the lift components in the y-direction
qraw pointer to the vector storing all the fluxes
lambdaraw pointer to the vector storing the pressure and Lagrangian forces
nuviscosity
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
nxnumber of cells in the x-direction
nynumber of cells in the y-direcyion
Ix-index of the bottom-left cell of the control surface
Jy-index of the top-right cell of the control surface
ncxnumber of cells in the x-direction in the control volume
nycnumber 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.

Parameters
FyXraw pointer to the vector storing the lift components in the x-direction
qraw pointer to the vector storing all the fluxes
nuviscosity
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
nxnumber of cells in the x-direction
nynumber of cells in the y-direcyion
Ix-index of the bottom-left cell of the control surface
Jy-index of the top-right cell of the control surface
ncxnumber of cells in the x-direction in the control volume
nycnumber 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.

Parameters
FyUraw pointer to the vector storing the unsteady lift components
qraw pointer to the vector storing all the fluxes
qOldraw pointer to the vector sotring all the fluxes at the previous time-step
dxraw pointer to the vector storing the cell widths in the x-direction
dyraw pointer to the vector storing the cell widths in the y-direction
dttime increment
nxnumber of cells in the x-direction
nynumber of cells in the y-direcyion
Ix-index of the bottom-left cell of the control surface
Jy-index of the top-right cell of the control surface
ncxnumber of cells in the x-direction in the control volume
nycnumber 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.

Parameters
QRowsrow index of elements of the matrix
QColscolumn index of elements of the matrix
QValsvalue of elements of the matrix
QSizesize of the matrix
tagstags 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.

Parameters
QTRowsrow index of elements of the matrix
QTColscolumn index of elements of the matrix
QTValsvalue of elements of the matrix
QTSizesize of the matrix
tagstags 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.

Parameters
QTRowsrow index of elements of the matrix
QTColscolumn index of elements of the matrix
QTValsvalue of elements of the matrix
ERowsrow index of elements of the interpolation matrix
EColscolumn index of elements of the interpolation matrix
EValsvalue of elements of the interpolation matrix
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
xx-component of grid points
yy-component of grid points
dxcell-widths in the x-direction
totalPointsnumber of body-points (all bodies included)
xBx-component of body-points (all bodies included)
yBy-component of body-points (all bodies included)
Ix-index of grid cells in which body points are located
Jy-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.

Parameters
rhs1vector to update
numUVnumber of velocity points in the domain
tagsvector 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.

Parameters
rhs1vector to update
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dttime-step size
dxgrid-spacings along a gridline in the x-direction
tagsvector used to differentiate regular points from forcing points
coeffscoefficients of interpolation
uvvelocity 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.

Parameters
rhs1vector to update
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dttime-step size
dxgrid-spacings along a gridline in the x-direction
tagsvector used to differentiate regular points from forcing points
coeffscoefficients of interpolation
coeffs2coefficients of interpolation
uvvelocity 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.

Parameters
rhs1vector to update
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dttime-step size
dygrid-spacings along a gridline in the y-direction
tagsvector used to differentiate regular points from forcing points
coeffscoefficients of interpolation
uvvelocity 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.

Parameters
rhs1vector to update
nxnumber of cells in the x-direction
nynumber of cells in the y-direction
dttime-step size
dygrid-spacings along a gridline in the y-direction
tagsvector used to differentiate regular points from forcing points
coeffscoefficients of interpolation
coeffs2coefficients of interpolation
uvvelocity vector

Definition at line 135 of file updateRHS1.cu.