Particle Field

FLOWVPM.ParticleFieldType
ParticleField(maxparticles::Int, R=FLOAT_TYPE; <keyword arguments>)

Create a new particle field with maxparticles particles. The particle field is created with the default values for the other parameters.

Arguments

  • maxparticles::Int: Maximum number of particles in the field.
  • R=FLOAT_TYPE: Type of the particle field. Default is FLOAT_TYPE.

Keyword Arguments

  • formulation::Formulation=ReformulatedVPM{FLOAT_TYPE}(0, 1/5): VPM governing equations. Default is reformulated to conserve mass for a tube and angular momentum for a sphere.
  • viscous::ViscousScheme=Inviscid(): Viscous scheme. Note that with rVPM formulation, artificial viscosity is not needed for numerical stability, as is common in VPM.
  • np::Int=0: Number of particles currently in the field.
  • nt::Int=0: Current time step number.
  • t::Real=0: Current simulation time.
  • transposed::Bool=true: If true, the transposed scheme of the stretching term is used (recommended for stability).
  • fmm::FMM: Fast multipole method tuning and auto-tuning settings.
  • M::Array{R, 1}=zeros(R, 4): Auxilliary memory for computations. Should not be modified for most purposes.
  • SFS::SubFilterScale=NoSFS{FLOAT_TYPE}(): Subfilter-scale turbulence model.
  • kernel::Kernel=Kernel(zeta_gauserf, g_gauserf, dgdr_gauserf, g_dgdr_gauserf): Regularization scheme. Default is Gaussian smoothing of the vorticity field.
  • UJ=UJ_fmm: Method used to compute the $N$-body problem. Default uses the fast multipole method to achieve $O(N)$ complexity.
  • Uinf::Function=(t) -> [0.0,0.0,0.0]: Uniform freestream velocity function Uinf(t).
  • relaxation::Relaxation=Relaxation(relax_pedrizzetti, 1, 0.3): Relaxation scheme to re-align the vorticity field to be divergence-free.
  • integration::Tintegration=rungekutta3: Time integration scheme. Default is a Runge-Kutta 3rd order, low-memory scheme.
  • useGPU::Int: Run on GPU if >0, CPU if 0. Default is 0. (Experimental and does not accelerate SFS calculations)
source
FLOWVPM.rVPMConstant
rVPM

Alias for the reformulated VPM formulation. Enforces conservation of mass and momentum for a spherical fluid element and is the default formulation (f=0, g=1/5)

source
FLOWVPM.add_particleFunction

add_particle(pfield::ParticleField, X, Gamma, sigma; <keyword arguments>)

Add a particle to the field.

Arguments

  • pfield::ParticleField : Particle field to add the particle to.
  • X : Position of the particle.
  • Gamma : Strength of the particle.
  • sigma : Smoothing radius of the particle.
  • vol : Volume of the particle. Default is 0.
  • circulation : Circulation of the particle. Default is 1.
  • C : SFS parameter of the particle. Default is 0.
  • static : If true, the particle is static. Default is false.
source

add_particle(pfield::ParticleField, P)

Add a copy of Particle P to the field.

source
FLOWVPM.remove_particleFunction

remove_particle(pfield::ParticleField, i)

Remove the i-th particle in the field. This is done by moving the last particle that entered the field into the memory slot of the target particle. To remove particles sequentally, you will need to go from the last particle back to the first one (see documentation of get_particleiterator for an example).

source
FLOWVPM.get_npFunction
`get_np(pfield::ParticleField)`

Returns current number of particles in the field.
source
FLOWVPM.get_particleiteratorFunction
`get_particleiterator(pfield::ParticleField; start_i=1, end_i=np)`

Return an iterator over particles that can be used as follows

julia> # Initiate particle field
       pfield = FLOWVPM.ParticleField(10);

julia> # Add particles
       for i in 1:4
           FLOWVPM.add_particle(pfield, (i*10^0, i*10^1, i*10^2), zeros(3), 1.0)
       end

julia> # Iterate over particles
       for P in FLOWVPM.get_particleiterator(pfield)
           println(P.var[1:3])
       end
[1.0, 10.0, 100.0]
[2.0, 20.0, 200.0]
[3.0, 30.0, 300.0]
[4.0, 40.0, 400.0]
source