AQUAgpusph 4.1.2
Loading...
Searching...
No Matches
Files | Macros | Functions
Lela

Files

file  DivSigma.cl
 Solid particles interaction with the boundary.
 
file  DivSigmaShepard.cl
 Renormalization of the deformation gradient.
 
file  EOS.cl
 Inverse EOS equation application.
 
file  GradU.cl
 Velocity gradient resulting from the interaction with the boundary.
 
file  GradUShepard.cl
 Renormalization of the deformation gradient.
 
file  InitP.cl
 Boundary pressure initialization.
 
file  Motion.cl
 Fixed boundary elements methods.
 
file  Stress.cl
 Fixed boundary elements methods.
 
file  ElasticBounce.cl
 The simplest boundary technique to assert the non-tresspasable boundary condition.
 
file  Corrector.cl
 Improved Euler time integration scheme corrector stage.
 
file  deltaSPH.cl
 delta-SPH for the linear elasticity module.
 
file  DivSigma.cl
 Solid particles interactions computation.
 
file  DivU.cl
 Velocity divergence computation.
 
file  GradU.cl
 Solid particles deformation gradient computation.
 
file  Rates.cl
 Velocity and density variation rates computation.
 
file  Shepard.cl
 Shepard renormalization factor for the Linear elasticity module.
 
file  Sigma.cl
 Stress tensor computation.
 
file  Sort.cl
 Sort all the particle variables by the cell indexes.
 
file  xSPH.cl
 delta-SPH methods, including the correction terms.
 

Macros

#define __DR_FACTOR__   1.5f
 The boundary elements effect is restricted to a quadrangular area of \( R \times R \), where \( R = DR_FACTOR \cdot \Delta r \).
 
#define __MIN_BOUND_DIST__   0.3f
 The elastic bounce is not tolerating that a particle becomes closer than this distance (multiplied by \( \Delta r \)).
 
#define __ELASTIC_FACTOR__   0.0f
 The amount of kinetic energy conserved in the interaction.
 
#define EXCLUDED_PARTICLE(index)   imove[index] != 2
 Restrict the aplication to the solid particles (imove=2)
 
#define EXCLUDED_PARTICLE(index)   imove[index] != 2
 Restrict the aplication to the solid particles (imove=2)
 

Functions

__kernel void entry (const __global int *imove, const __global vec *r, const __global vec *normal, const __global float *rho, const __global float *m, const __global matrix *sigma, __global vec *div_sigma, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells)
 Solid particles interaction with the boundary.
 
__kernel void entry (const __global int *imove, const __global float *shepard, __global vec *div_sigma, __global float *div_u, uint N)
 Renormalization of the deformation gradient.
 
__kernel void entry (const __global uint *iset, const __global int *imove, const __global float *p, __global float *rho, __constant float *refd, uint N, float cs, float p0)
 Inverse EOS to get the density from the pressure value.
 
__kernel void entry (const __global int *imove, const __global vec *r, const __global vec *u, const __global vec *normal, const __global float *m, __global matrix *grad_u, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells)
 Velocity gradient resulting from the interaction with the boundary.
 
__kernel void entry (const __global int *imove, const __global float *shepard, __global matrix *grad_u, uint N)
 Renormalization of the deformation gradient.
 
__kernel void entry (const __global int *imove, __global float *p, uint N, float p0)
 Pressure initialization.
 
__kernel void interpolation (const __global uint *iset, const __global int *imove, const __global vec *r, const __global float *m, const __global float *rho, __global vec *u, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells, uint BImotion_iset)
 Velocity interpolation at the boundary elements.
 
__kernel void shepard (const __global uint *iset, const __global int *imove, const __global float *shepard, __global vec *u, uint N, uint BImotion_iset)
 Velocity renormalization.
 
__kernel void euler (const __global uint *iset, const __global int *imove, __global vec *r, const __global vec *u, unsigned int N, float dt, uint BImotion_iset)
 Simple Euler time integration of the velocity.
 
__kernel void interpolation (const __global uint *iset, const __global int *imove, const __global vec *r, const __global float *m, const __global float *rho, __global float *p, __global matrix *S, const __global uint *icell, const __global uint *ihoc, __constant float *refd, uint N, uivec4 n_cells, vec g, uint BIstress_iset)
 Pressure and stress deviation interpolation at the boundary elements.
 
__kernel void shepard (const __global uint *iset, const __global int *imove, const __global float *shepard, __global float *p, __global matrix *S, uint N, uint BIstress_iset)
 Pressure and stress deviation renormalization.
 
__kernel void entry (const __global int *imove, const __global vec *r, const __global vec *normal, __global vec *u, __global vec *dudt, __global uint *icell, __global uint *ihoc, uint N, uivec4 n_cells, float dr, float dt)
 Performs the boundary effect on the fluid particles.
 
__kernel void entry (__global int *imove, __global matrix *S, __global matrix *dSdt, __global matrix *S_in, __global matrix *dSdt_in, unsigned int N, float dt)
 Improved Euler time integration scheme corrector stage.
 
__kernel void entry (const __global int *imove, const __global vec *r, const __global float *rho, const __global float *m, const __global matrix *sigma, __global vec *div_sigma, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells)
 Solid particles interactions computation.
 
__kernel void entry (const __global int *imove, const __global float *rho, const __global matrix *grad_u, __global float *div_u, unsigned int N)
 Velocity divergence computation.
 
__kernel void entry (const __global int *imove, const __global vec *r, const __global vec *u, const __global float *rho, const __global float *m, __global matrix *grad_u, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells)
 Solid particles deformation gradient computation.
 
__kernel void entry (const __global unsigned int *iset, const __global int *imove, const __global vec *div_sigma, const __global float *div_u, const __global matrix *grad_u, const __global matrix *S, __global vec *dudt, __global float *drhodt, __global matrix *dSdt, __constant float *shear_mod, unsigned int N, vec g)
 Velocity and density variation rates computation.
 
__kernel void entry (const __global float *p, const __global matrix *S, __global matrix *sigma, unsigned int N)
 Stress tensor computation.
 
__kernel void entry (const __global float *p_in, __global float *p, const __global matrix *S_in, __global matrix *S, const __global matrix *dSdt_in, __global matrix *dSdt, const __global unit *id_sorted, unsigned int N)
 Sort all the particle variables by the cell indexes.
 
__kernel void interpolation (const __global int *imove, const __global vec *r, const __global vec *u, const __global float *rho, const __global float *m, __global vec *xSPH_u, const __global uint *icell, const __global uint *ihoc, uint N, uivec4 n_cells)
 Velocity interpolation.
 
__kernel void xSPH (const __global int *imove, const __global vec *xSPH_u, const __global float *shepard, __global vec *u, unsigned int N, float xSPH_factor)
 xSPH velocity replacement.
 

Detailed Description

Macro Definition Documentation

◆ __DR_FACTOR__

#define __DR_FACTOR__   1.5f

The boundary elements effect is restricted to a quadrangular area of \( R \times R \), where \( R = DR_FACTOR \cdot \Delta r \).

◆ __ELASTIC_FACTOR__

#define __ELASTIC_FACTOR__   0.0f

The amount of kinetic energy conserved in the interaction.

A factor of 1 imply that the velocity of the particle will be preserved (except for the direction), while a factor of 0 imply that the particle will loss all its normal to the boundary velocity.

The tangential velocity is not affected.

◆ __MIN_BOUND_DIST__

#define __MIN_BOUND_DIST__   0.3f

The elastic bounce is not tolerating that a particle becomes closer than this distance (multiplied by \( \Delta r \)).

◆ EXCLUDED_PARTICLE [1/2]

#define EXCLUDED_PARTICLE (   index)    imove[index] != 2

Restrict the aplication to the solid particles (imove=2)

◆ EXCLUDED_PARTICLE [2/2]

#define EXCLUDED_PARTICLE (   index)    imove[index] != 2

Restrict the aplication to the solid particles (imove=2)

Function Documentation

◆ entry() [1/14]

__kernel void entry ( __global int *  imove,
__global matrix S,
__global matrix dSdt,
__global matrix S_in,
__global matrix dSdt_in,
unsigned int  N,
float  dt 
)

Improved Euler time integration scheme corrector stage.

Time integration is based in the following quasi-second order Predictor-Corrector integration scheme:

\( \S_{n+1} = \S_{n} + \Delta t \left. \frac{\mathrm{d}S}{\mathrm{d}t} \right\vert_{n+1/2} + \frac{\Delta t}{2} \left( \left. \frac{\mathrm{d}S}{\mathrm{d}t} \right\vert_{n + 1/2} - \left. \frac{\mathrm{d}S}{\mathrm{d}t} \right\vert_{n - 1/2} \right) \)

Parameters
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
isetSet of particles index.
SDeviatory stress \( S_{n+1} \).
dSdtDeviatory stress rate of change \( \left. \frac{d S}{d t} \right\vert_{n+1} \).
S_inDeviatory stress \( S_{n+1/2} \).
dSdt_inDeviatory stress rate of change \( \left. \frac{d S}{d t} \right\vert_{n+1/2} \).
NNumber of particles.
dtTime step \( \Delta t \).
See also
lela/Predictor.cl

◆ entry() [2/14]

__kernel void entry ( const __global float *  p,
const __global matrix S,
__global matrix sigma,
unsigned int  N 
)

Stress tensor computation.

\[ \sigma = p \mathcal{I} + S \]

See also
https://en.wikipedia.org/wiki/Linear_elasticity
https://en.wikipedia.org/wiki/Hooke's_law
Parameters
pPressure \( p \).
SDeviatory stress \( S \).
sigmaStress tensor \( \sigma \).
NNumber of particles.

◆ entry() [3/14]

__kernel void entry ( const __global float *  p_in,
__global float *  p,
const __global matrix S_in,
__global matrix S,
const __global matrix dSdt_in,
__global matrix dSdt,
const __global unit id_sorted,
unsigned int  N 
)

Sort all the particle variables by the cell indexes.

Parameters
p_inUnsorted pressure \( p \).
pSorted pressure \( p \).
S_inUnsorted Deviatory stress \( S \).
SSorted Deviatory stress \( S \).
dSdt_inUnsorted Deviatory stress rate of change \( \frac{d S}{d t} \).
dSdtSorted Deviatory stress rate of change \( \frac{d S}{d t} \). \( \left. \frac{d S}{d t} \right\vert_{n+1} \).
id_sortedPermutations list from the unsorted space to the sorted one.
NNumber of particles.

◆ entry() [4/14]

__kernel void entry ( const __global int *  imove,
__global float *  p,
uint  N,
float  p0 
)

Pressure initialization.

Since the pressure field is resulting from the density field (applying the EOS), letting free the pressure at the boundary may be dangerous due to the value must be unassigned. Therefore the pressure field will be initialized as the background pressure, letting the user to don't sdpecifically assign a pressure value without crashing the simulation.

Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
pPressure \( p \).
p0Background pressure \( p_0 \).
NTotal number of particles and boundary elements.
BImotion_isetSet of particles affected

◆ entry() [5/14]

__kernel void entry ( const __global int *  imove,
const __global float *  rho,
const __global matrix grad_u,
__global float *  div_u,
unsigned int  N 
)

Velocity divergence computation.

Since we already computed the gradient of the velocity, we can get the divergence as the sum of the diagonal components:

\[ \nabla \cdot \mathbf{u} = \sum_i^d \left(\nabla \mathbf{u}\right)_ii \]

Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors (ignored by this preset).
  • imove < 0 for boundary elements/particles.
rhoDensity \( \rho \).
grad_uGradient of the deformation \( \nabla \mathbf{r}^{*} \).
div_uVelocity divergence \( \rho \nabla \cdot \mathbf{u} \).
NNumber of particles.

◆ entry() [6/14]

__kernel void entry ( const __global int *  imove,
const __global float *  shepard,
__global matrix grad_u,
uint  N 
)

Renormalization of the deformation gradient.

The main drawback of the boundary integrals formulation is the requirement of the renormalization of the computed differentiqal operators, which is destroying several conservation properties.

Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
shepardShepard term \( \gamma(\mathbf{x}) = \int_{\Omega} W(\mathbf{y} - \mathbf{x}) \mathrm{d}\mathbf{x} \).
grad_uGradient of the velocity \( \nabla \mathbf{u} \).
NTotal number of particles and boundary elements.
Here is the call graph for this function:

◆ entry() [7/14]

__kernel void entry ( const __global int *  imove,
const __global float *  shepard,
__global vec div_sigma,
__global float *  div_u,
uint  N 
)

Renormalization of the deformation gradient.

The main drawback of the boundary integrals formulation is the requirement of the renormalization of the computed differentiqal operators, which is destroying several conservation properties.

Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
shepardShepard term \( \gamma(\mathbf{x}) = \int_{\Omega} W(\mathbf{y} - \mathbf{x}) \mathrm{d}\mathbf{x} \).
div_sigmaDivergence of the stress tensor \( \frac{\nabla \cdot \sigma}{rho} \).
div_uVelocity divergence \( \rho \nabla \cdot \mathbf{u} \).
NTotal number of particles and boundary elements.
Here is the call graph for this function:

◆ entry() [8/14]

__kernel void entry ( const __global int *  imove,
const __global vec r,
const __global float *  rho,
const __global float *  m,
const __global matrix sigma,
__global vec div_sigma,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells 
)

Solid particles interactions computation.

Compute the differential operators involved in the numerical scheme, taking into account just the solid-solid interactions.

Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors (ignored by this module).
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
rhoDensity \( \rho \).
mMass \( m \).
sigmaStress tensor \( \sigma \).
div_sigmaDivergence of the stress tensor \( \frac{\nabla \cdot \sigma}{rho} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction.
Here is the call graph for this function:

◆ entry() [9/14]

__kernel void entry ( const __global int *  imove,
const __global vec r,
const __global vec normal,
__global vec u,
__global vec dudt,
__global uint *  icell,
__global uint *  ihoc,
uint  N,
uivec4  n_cells,
float  dr,
float  dt 
)

Performs the boundary effect on the fluid particles.

Parameters
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
normalNormal \( \mathbf{n} \).
uVelocity \( \mathbf{u} \).
dudtVelocity rate of change \( \left. \frac{d \mathbf{u}}{d t} \right\vert_{n+1} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction
drDistance between particle \( \Delta r \).
dtTime step \( \Delta t \).

◆ entry() [10/14]

__kernel void entry ( const __global int *  imove,
const __global vec r,
const __global vec normal,
const __global float *  rho,
const __global float *  m,
const __global matrix sigma,
__global vec div_sigma,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells 
)

Solid particles interaction with the boundary.

Compute the differential operators involved in the numerical scheme, taking into account just the solid-solid interactions.

Parameters
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
normalNormal \( \mathbf{n} \).
rhoDensity \( \rho \).
mArea of the boundary element \( s \).
sigmaStress tensor \( \sigma \).
div_sigmaDivergence of the stress tensor \( \frac{\nabla \cdot \sigma}{rho} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction.
Here is the call graph for this function:

◆ entry() [11/14]

__kernel void entry ( const __global int *  imove,
const __global vec r,
const __global vec u,
const __global float *  rho,
const __global float *  m,
__global matrix grad_u,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells 
)

Solid particles deformation gradient computation.

Compute the gradient of the deformation vector:

\[ \nabla \mathbf{r}^{*} = \mathbf{r}^{*} \otimes \nabla \]

See also
https://en.wikipedia.org/wiki/Matrix_calculus
https://en.wikipedia.org/wiki/Outer_product
Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors (ignored by this preset).
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
uVelocity \( \mathbf{u} \).
rhoDensity \( \rho \).
mMass \( m \).
grad_uGradient of the velocity \( \nabla \mathbf{u} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction
Here is the call graph for this function:

◆ entry() [12/14]

__kernel void entry ( const __global int *  imove,
const __global vec r,
const __global vec u,
const __global vec normal,
const __global float *  m,
__global matrix grad_u,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells 
)

Velocity gradient resulting from the interaction with the boundary.

Compute the gradient of the velocity:

\[ \nabla \mathbf{u} = \mathbf{u} \otimes \nabla \]

See also
https://en.wikipedia.org/wiki/Matrix_calculus
https://en.wikipedia.org/wiki/Outer_product
Parameters
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors (ignored by this preset).
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
uVelocity \( \mathbf{u} \).
normalNormal \( \mathbf{n} \).
mArea of the boundary element \( s \).
grad_uGradient of the velocity \( \nabla \mathbf{u} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction
Here is the call graph for this function:

◆ entry() [13/14]

__kernel void entry ( const __global uint *  iset,
const __global int *  imove,
const __global float *  p,
__global float *  rho,
__constant float *  refd,
uint  N,
float  cs,
float  p0 
)

Inverse EOS to get the density from the pressure value.

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
pPressure \( p \).
rhoDensity \( \rho \).
refdDensity of reference of the fluid \( \rho_0 \).
NTotal number of particles and boundary elements.
csSpeed of sound \( c_s \).
p0Background pressure \( p_0 \).

◆ entry() [14/14]

__kernel void entry ( const __global unsigned int *  iset,
const __global int *  imove,
const __global vec div_sigma,
const __global float *  div_u,
const __global matrix grad_u,
const __global matrix S,
__global vec dudt,
__global float *  drhodt,
__global matrix dSdt,
__constant float *  shear_mod,
unsigned int  N,
vec  g 
)

Velocity and density variation rates computation.

The mass conservation and momentum equations are applied from the already computed differential operators:

  • \( \frac{\mathrm{d} \mathbf{u}}{\mathrm{d} t} = - \frac{\nabla \cdot \sigma}{rho} + \mathbf{g}\)
  • \( \frac{\mathrm{d} \rho}{\mathrm{d} t} = - \rho \nabla \cdot \mathbf{u}\)
Parameters
isetSet of particles index.
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors (ignored by this preset).
  • imove < 0 for boundary elements/particles.
div_sigmaDivergence of the stress tensor \( \frac{\nabla \cdot \sigma}{rho} \).
div_uVelocity divergence \( \rho \nabla \cdot \mathbf{u} \).
grad_uVelocity gradient \( \nabla \mathbf{u} \).
SDeviatory stress \( S \).
dudtVelocity rate of change \( \frac{d \mathbf{u}}{d t} \).
drhodtDensity rate of change \( \frac{d \rho}{d t} \).
dSdtDeviatory stress rate of change \( \frac{d S}{d t} \).
shear_modShear modulus \( \mu \).
NNumber of particles.
gGravity acceleration \( \mathbf{g} \).

◆ euler()

__kernel void euler ( const __global uint *  iset,
const __global int *  imove,
__global vec r,
const __global vec u,
unsigned int  N,
float  dt,
uint  BImotion_iset 
)

Simple Euler time integration of the velocity.

Since the velocity is resulting from the interpolation of the solid particle, it does not make any sense to integrate it with an improved Euler scheme.

\( \mathbf{u}_{n+1} = \mathbf{u}_{n} + \Delta t \left. \frac{\mathrm{d} \mathbf{u}}{\mathrm{d}t} \right\vert_{n+1/2} \right) \)

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
uVelocity \( \mathbf{u} \).
NNumber of particles.
dtTime step \( \Delta t \).
BImotion_isetSet of particles affected

◆ interpolation() [1/3]

__kernel void interpolation ( const __global int *  imove,
const __global vec r,
const __global vec u,
const __global float *  rho,
const __global float *  m,
__global vec xSPH_u,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells 
)

Velocity interpolation.

Parameters
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
uVelocity \( \mathbf{u} \).
rhoDensity \( \rho \).
mMass \( m \).
xSPH_uInterpolated velocity (Shepard renormalization is not applied yet).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction
Here is the call graph for this function:

◆ interpolation() [2/3]

__kernel void interpolation ( const __global uint *  iset,
const __global int *  imove,
const __global vec r,
const __global float *  m,
const __global float *  rho,
__global float *  p,
__global matrix S,
const __global uint *  icell,
const __global uint *  ihoc,
__constant float *  refd,
uint  N,
uivec4  n_cells,
vec  g,
uint  BIstress_iset 
)

Pressure and stress deviation interpolation at the boundary elements.

The values are computed using just the fluid information. The resulting interpolated values are not renormalized yet.

Just the elements with the flag imove = -3 are considered boundary elements.

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
mMass \( m \).
rhoDensity \( \rho \).
pPressure \( p \).
SDeviatory stress \( S \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
refdDensity of reference of the fluid \( \rho_0 \).
NNumber of particles.
n_cellsNumber of cells in each direction
gGravity acceleration \( \mathbf{g} \).
BIstress_isetSet of particles affected
Here is the call graph for this function:

◆ interpolation() [3/3]

__kernel void interpolation ( const __global uint *  iset,
const __global int *  imove,
const __global vec r,
const __global float *  m,
const __global float *  rho,
__global vec u,
const __global uint *  icell,
const __global uint *  ihoc,
uint  N,
uivec4  n_cells,
uint  BImotion_iset 
)

Velocity interpolation at the boundary elements.

The values are computed using just the solid information. The resulting interpolated values are not renormalized yet.

Just the elements with the flag imove = -3 are considered boundary elements.

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
rPosition \( \mathbf{r} \).
mMass \( m \).
rhoDensity \( \rho \).
uVelocity \( \mathbf{u} \).
icellCell where each particle is located.
ihocHead of chain for each cell (first particle found).
NNumber of particles.
n_cellsNumber of cells in each direction
BImotion_isetSet of particles affected
Here is the call graph for this function:

◆ shepard() [1/2]

__kernel void shepard ( const __global uint *  iset,
const __global int *  imove,
const __global float *  shepard,
__global float *  p,
__global matrix S,
uint  N,
uint  BIstress_iset 
)

Pressure and stress deviation renormalization.

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
shepardShepard term \( \gamma(\mathbf{x}) = \int_{\Omega} W(\mathbf{y} - \mathbf{x}) \mathrm{d}\mathbf{x} \).
pPressure \( p \).
SDeviatory stress \( S \).
NTotal number of particles and boundary elements.
BIstress_isetSet of particles affected
Here is the call graph for this function:

◆ shepard() [2/2]

__kernel void shepard ( const __global uint *  iset,
const __global int *  imove,
const __global float *  shepard,
__global vec u,
uint  N,
uint  BImotion_iset 
)

Velocity renormalization.

Parameters
isetSet of particles index.
imoveMoving flags.
  • imove = 2 for regular solid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
shepardShepard term \( \gamma(\mathbf{x}) = \int_{\Omega} W(\mathbf{y} - \mathbf{x}) \mathrm{d}\mathbf{x} \).
uVelocity \( \mathbf{u} \).
NTotal number of particles and boundary elements.
BImotion_isetSet of particles affected
Here is the call graph for this function:

◆ xSPH()

__kernel void xSPH ( const __global int *  imove,
const __global vec xSPH_u,
const __global float *  shepard,
__global vec u,
unsigned int  N,
float  xSPH_factor 
)

xSPH velocity replacement.

Parameters
imoveMoving flags.
  • imove > 0 for regular fluid particles.
  • imove = 0 for sensors.
  • imove < 0 for boundary elements/particles.
xSPH_uInterpolated velocity (Shepard renormalization is not applied yet).
shepardShepard term \( \gamma(\mathbf{x}) = \int_{\Omega} W(\mathbf{y} - \mathbf{x}) \mathrm{d}\mathbf{x} \).
uVelocity \( \mathbf{u} \).
NNumber of particles.
xSPH_factorRatio between the instantaneous velocity, and the interpolated one.
Here is the call graph for this function: