18#ifndef MILUPHPC_PARTICLES_CUH
19#define MILUPHPC_PARTICLES_CUH
151#if VARIABLE_SML || INTEGRATE_SML
175#if SOLID || NAVIER_STOKES
191 real *alpha_jutzi_old;
212 real *dalpha_epspordt;
221 real *shepardCorrection;
224#if LINEAR_CONSISTENCY
226 real *tensorialCorrectionMatrix;
246 real *damage_porjutzi;
248 real *ddamage_porjutzidt;
416#if VARIABLE_SML || INTEGRATE_SML
440#if SOLID || NAVIER_STOKES
479 real *epsilon_v,
real *depsilon_vdt);
489#if LINEAR_CONSISTENCY
725#if VARIABLE_SML || INTEGRATE_SML
732 __global__
void setVariableSML(
Particles *particles,
real *dsmldt);
744 __global__
void setSMLCorrection(
Particles *particles,
real *sml_omega);
747 void setSMLCorrection(
Particles *particles,
real *sml_omega);
778#if SOLID || NAVIER_STOKES
785 __global__
void setSolidNavierStokes(
Particles *particles,
real *sigma);
793 void setSolidNavierStokes(
Particles *particles,
real *sigma);
803 __global__
void setArtificialStress(
Particles *particles,
real *R);
838 real *epsilon_v,
real *depsilon_vdt);
863 real *epsilon_v,
real *depsilon_vdt);
873 __global__
void setZeroConsistency(
Particles *particles,
real *shepardCorrection);
881 void setZeroConsistency(
Particles *particles,
real *shepardCorrection);
884#if LINEAR_CONSISTENCY
891 __global__
void setLinearConsistency(
Particles *particles,
real *tensorialCorrectionMatrix);
899 void setLinearConsistency(
Particles *particles,
real *tensorialCorrectionMatrix);
941 __global__
void setPalphaPorosity(
Particles *particles,
real *damage_porjutzi,
real *ddamage_porjutzidt);
950 void setPalphaPorosity(
Particles *particles,
real *damage_porjutzi,
real *ddamage_porjutzidt);
1019#if VARIABLE_SML || INTEGRATE_SML
1065#if VARIABLE_SML || INTEGRATE_SML
1151#if VARIABLE_SML || INTEGRATE_SML
CUDA_CALLABLE_MEMBER void setSML(real *sml)
CUDA_CALLABLE_MEMBER IntegratedParticles()
CUDA_CALLABLE_MEMBER ~IntegratedParticles()
CUDA_CALLABLE_MEMBER void reset(integer index)
idInteger * uid
unique identifier
CUDA_CALLABLE_MEMBER void setIntegrateDensity(real *drhodt)
CUDA_CALLABLE_MEMBER void set(idInteger *uid, real *rho, real *e, real *dedt, real *p, real *cs, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az)
Particle(s) class based on SoA (Structur of Arrays).
integer * noi
(pointer to) number of interactions (array)
real * e
(pointer to) internal energy (array)
CUDA_CALLABLE_MEMBER real weightedEntry(integer index, Entry::Name entry)
integer * materialId
(pointer to) material identifier (array)
real * x
(pointer to) x position (array)
real * rho
(pointer to) density (array)
integer * level
(pointer to) level of the (pseudo-)particles
idInteger * uid
(pointer to) unique identifier (array)
real * drhodt
(pointer to) time derivative of density (array)
CUDA_CALLABLE_MEMBER void set(integer *numParticles, integer *numNodes, real *mass, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az, integer *level, idInteger *uid, integer *materialId, real *sml, integer *nnl, integer *noi, real *e, real *dedt, real *cs, real *rho, real *p)
Setter (DIM = 3) assigning variables to pointer members.
real * g_ay
(pointer to) y gravitational acceleration (array)
real * ay
(pointer to) y acceleration (array)
integer * nodeType
(pointer to) node type
integer * nnl
(pointer to) near(est) neighbor list (array)
real * az
(pointer to) z acceleration (array)
real * p
(pointer to) pressure (array)
real * cs
(pointer to) sound of speed (array)
real * y
(pointer to) y position (array)
real * ax
(pointer to) x acceleration (array)
CUDA_CALLABLE_MEMBER void reset(integer index)
Reset (specific) entries.
real * u
energy (kinetic + gravitational for now)
real * g_ax
(pointer to) x gravitational acceleration (array)
CUDA_CALLABLE_MEMBER Particles()
CUDA_CALLABLE_MEMBER ~Particles()
real * sml
(pointer to) smoothing length (array)
real * dedt
(pointer to) time derivative of internal energy (array)
real * mass
(pointer to) mass (array)
real * z
(pointer to) z position (array)
integer * numParticles
number of particles
integer * numNodes
number of nodes
real * vz
(pointer to) z velocity (array)
CUDA_CALLABLE_MEMBER real distance(integer index_1, integer index_2)
Distance of two particles.
CUDA_CALLABLE_MEMBER void setLeapfrog(real *ax_old, real *ay_old, real *az_old, real *g_ax_old, real *g_ay_old, real *g_az_old)
real * muijmax
(pointer) to max(mu_ij) (array) needed for artificial viscosity and determining timestp
CUDA_CALLABLE_MEMBER void setGravity(real *g_ax, real *g_ay, real *g_az)
Constructor (DIM = 3) assigning gravitational acceleration to member variables.
real * g_az
(pointer to) z gravitational acceleration (array)
real * vx
(pointer to) x velocity (array)
CUDA_CALLABLE_MEMBER void setNodeType(integer *nodeType)
Setter for node type.
real * vy
(pointer to) y velocity (array)
CUDA_CALLABLE_MEMBER void setArtificialViscosity(real *muijmax)
Setter for artificial viscosity (entry).
CUDA_CALLABLE_MEMBER void setIntegrateDensity(real *drhodt)
Setter in dependence of INTEGRATE_DENSITY.
CUDA_CALLABLE_MEMBER void setU(real *u)
Setter for energy.
#define CUDA_CALLABLE_MEMBER
void setSML(IntegratedParticles *integratedParticles, real *sml)
void setIntegrateDensity(IntegratedParticles *integratedParticles, real *drhodt)
void set(IntegratedParticles *integratedParticles, idInteger *uid, real *rho, real *e, real *dedt, real *p, real *cs, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az)
__global__ void set(IntegratedParticles *integratedParticles, idInteger *uid, real *rho, real *e, real *dedt, real *p, real *cs, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az)
__global__ void setSML(IntegratedParticles *integratedParticles, real *sml)
__global__ void setIntegrateDensity(IntegratedParticles *integratedParticles, real *drhodt)
void setU(Particles *particles, real *u)
void setIntegrateDensity(Particles *particles, real *drhodt)
void setGravity(Particles *particles, real *g_ax, real *g_ay, real *g_az)
void setLeapfrog(Particles *particles, real *ax_old, real *ay_old, real *az_old, real *g_ax_old, real *g_ay_old, real *g_az_old)
void set(Particles *particles, integer *numParticles, integer *numNodes, real *mass, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az, integer *level, idInteger *uid, integer *materialId, real *sml, integer *nnl, integer *noi, real *e, real *dedt, real *cs, real *rho, real *p)
real test(Particles *particles, bool time=false)
void setNodeType(Particles *particles, integer *nodeType)
real info(Particles *particles, integer n, integer m, integer k)
real check4nans(Particles *particles, integer n)
void setArtificialViscosity(Particles *particles, real *muijmax)
__global__ void set(Particles *particles, integer *numParticles, integer *numNodes, real *mass, real *x, real *y, real *z, real *vx, real *vy, real *vz, real *ax, real *ay, real *az, integer *level, idInteger *uid, integer *materialId, real *sml, integer *nnl, integer *noi, real *e, real *dedt, real *cs, real *rho, real *p)
__global__ void setArtificialViscosity(Particles *particles, real *muijmax)
__global__ void setNodeType(Particles *particles, integer *nodeType)
__global__ void setIntegrateDensity(Particles *particles, real *drhodt)
__global__ void setGravity(Particles *particles, real *g_ax, real *g_ay, real *g_az)
__global__ void info(Particles *particles, integer n, integer m, integer k)
__global__ void setLeapfrog(Particles *particles, real *ax_old, real *ay_old, real *az_old, real *g_ax_old, real *g_ay_old, real *g_az_old)
__global__ void check4nans(Particles *particles, integer n)
__global__ void setU(Particles *particles, real *u)
__global__ void test(Particles *particles)
Particle class related functions and kernels.
const char *const numParticles
Settings via preprocessor directives, typedefs, constants, structs.