Reference
ParticleField
FLOWVPM.ParticleField — TypeParticleField(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 isFLOAT_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 withrVPMformulation, 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)
Formulations
FLOWVPM.rVPM — ConstantrVPMAlias 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)
FLOWVPM.cVPM — ConstantcVPMAlias for the classic VPM formulation.
FLOWVPM.formulation_tube_continuity — Constantformulation_tube_continuityAlias for mass conserving tube formulation. Enforces conservation of mass for a vortex tube (f=1/2, g=0)
FLOWVPM.formulation_tube_momentum — Constantformulation_tube_momentumAlias for momentum conserving tube formulation. Enforces conservation of momentum for a vortex tube (f=1/4, g=1/4)
Viscous Models
FLOWVPM.Inviscid — TypeInviscid()Creates an inviscid scheme with zero kinematic viscosity.
FLOWVPM.CoreSpreading — TypeCoreSpreading(nu, sgm0, zeta::Tzeta=zeta_fmm; <keyword arguments>)Creates a core spreading viscous scheme with the given parameters.
Arguments
nu::Real Kinematic viscosity.sgm0::Real Core size after reset.zeta::Function = zeta_fmmBasis function evaluation method.beta::Real = 1.5Maximum core size growth σ/σ_0.itmax::Int = 15Maximum number of RBF iterations.tol::Real = 1e-3RBF interpolation tolerance.iterror::Bool = trueThrow error if RBF didn't converge.verbose::Bool = falseVerbose on RBF interpolation.v_lvl::Int = 2Verbose printing tab level.debug::Bool = falsePrint verbose for debugging.rbf::Function = rbf_conjugategradientRBF function.rr0s::Array{R, 1} = zeros(R, 3)Initial field residuals.rrs::Array{R, 1} = zeros(R, 3)Current field residuals.prev_rrs::Array{R, 1} = zeros(R, 3)Previous field residuals.pAps::Array{R, 1} = zeros(R, 3)pAp product.alphas::Array{R, 1} = zeros(R, 3)Alpha coefficients.betas::Array{R, 1} = zeros(R, 3)flags::Array{Bool, 1} = zeros(Bool, 3)
FLOWVPM.ParticleStrengthExchange — TypeParticleStrengthExchange(nu; <keyword arguments>)Creates a particle strength exchange viscous scheme with the given parameters.
Arguments
nu::Real Kinematic viscosity.recalculate_vols::Bool = trueWhether to recalculate particle volumes.
FMM
FLOWVPM.FMM — Type`FMM(; p::Int=4, ncrit::Int=10, theta::Real=0.5, shrink_recenter::Bool=true,
relative_tolerance::Real=1e-3, absolute_tolerance::Real=1e-6,
autotune_p::Bool=true, autotune_ncrit::Bool=true,
autotune_reg_error::Bool=true, default_rho_over_sigma::Real=1.0)`Parameters for FMM solver.
Arguments
p: Order of multipole expansion.ncrit: Maximum number of particles per leaf.theta: Neighborhood criterion. This criterion defines the distance where the far field starts. The criterion is that if θ*r < R1+R2 the interaction between two cells is resolved through P2P, where r is the distance between cell centers, and R1 and R2 are each cell radius. This means that at θ=1, P2P is done only on cells that have overlap; at θ=0.5, P2P is done on cells that their distance is less than double R1+R2; at θ=0.25, P2P is done on cells that their distance is less than four times R1+R2; at θ=0, P2P is done on all cells.shrink_recenter: If true, shrink and recenter multipole expansions to account for nonzero particle radius.relative_tolerance: Relative error tolerance for FMM calls.absolute_tolerance: Absolute error tolerance fallback for FMM calls in case relative tolerance becomes too small.autotune_p: If true, automatically adjust p to optimize performance.autotune_ncrit: If true, automatically adjust ncrit to optimize performance.autotune_reg_error: If true, constrain regularization error in FMM calls.default_rho_over_sigma: Default value for ρ/σ in FMM calls (unused ifautotune_reg_erroris true).
Subfilter-Scale Models
FLOWVPM.noSFS — Constant`noSFS = NoSFS{FLOAT_TYPE}()`Alias for the no subfilter-scale model.
FLOWVPM.SFS_Cs_nobackscatter — Constant`SFS_Cs_nobackscatter = ConstantSFS(Estr_fmm; Cs=1.0, clippings=(clipping_backscatter,))`Alias for the Constant SFS model with no backscatter.
FLOWVPM.SFS_Cd_twolevel_nobackscatter — Constant`SFS_Cd_twolevel_nobackscatter = DynamicSFS(Estr_fmm, pseudo3level_beforeUJ, pseudo3level_positive_afterUJ; alpha=0.999, clippings=(clipping_backscatter,))`Alias for the Dynamic SFS model with two levels and no backscatter. This is the recommended SFS model for high fidelity modeling.
FLOWVPM.SFS_Cd_threelevel_nobackscatter — Constant`SFS_Cd_threelevel_nobackscatter = DynamicSFS(Estr_fmm, pseudo3level_beforeUJ, pseudo3level_positive_afterUJ; alpha=0.667, clippings=(clipping_backscatter,))`Alias for the Dynamic SFS model with three levels and no backscatter. This is similar to the two level version but uses a lower value of alpha (0.667).
Particle Kernels
FLOWVPM.singular — ConstantsingularAlias for the singular kernel.
FLOWVPM.gaussian — ConstantgaussianAlias for the Gaussian kernel.
FLOWVPM.gaussianerf — ConstantgaussianerfAlias for the Gaussian error function kernel.
FLOWVPM.winckelmans — ConstantwinckelmansAlias for the Winckelmans kernel.
Relaxation Schemes
FLOWVPM.norelaxation — ConstantnorelaxationAlias for the no relaxation scheme.
FLOWVPM.pedrizzetti — ConstantpedrizzettiAlias for the Pedrizzetti relaxation scheme.
FLOWVPM.correctedpedrizzetti — ConstantcorrectedpedrizzettiAlias for the corrected Pedrizzetti relaxation scheme. Is a modification to the pedrizzetti relaxation that preserves the vortex strength magnitude.
Time Integration
FLOWVPM.euler — Functioneuler(pfield::ParticleField, dt::Real; relax::Bool=false, custom_UJ=nothing)Convects the pfield by timestep dt using a forward Euler step.
Arguments
pfield::ParticleFieldThe particle field to integrate.dt::RealThe time step.relax::BoolWhether to apply relaxation (default: false).custom_UJOptional custom function for updating U and J.
FLOWVPM.rungekutta3 — FunctionSteps the field forward in time by dt in a third-order low-storage Runge-Kutta integration scheme. See Notebook entry 20180105.
rungekutta3(pfield::ParticleField, dt::Real; relax::Bool=false, custom_UJ=nothing)Steps the field forward in time by dt in a third-order low-storage Runge-Kutta integration scheme using the VPM reformulation. See Notebook entry 20180105 (RK integration) and notebook 20210104 (reformulation).
Arguments
pfield::ParticleFieldThe particle field to integrate.dt::R3The time step.relax::BoolWhether to apply relaxation (default: false).custom_UJOptional custom function for updating U and J.
UJ
FLOWVPM.UJ_fmm — FunctionUJ_fmm(pfield)
Calculates the velocity and Jacobian that the field exerts on itself through a fast-multipole approximation, saving U and J on the particles.
NOTE: This method accumulates the calculation on the properties U and J of every particle without previously emptying those properties.
FLOWVPM.UJ_direct — MethodUJ_direct(pfield)
Calculates the velocity and Jacobian that the field exerts on itself by direct particle-to-particle interaction, saving U and J on the particles.
NOTE: This method accumulates the calculation on the properties U and J of every particle without previously emptying those properties.