ConservationLawsParticles.jl
Documentation for ConservationLawsParticles.jl.
Particle method for 1D conservation laws.
This package implements the deterministic particle schemes described in the article
E. Radici, F. Stra: Entropy solutions of mildly singular non-local scalar conservation laws with congestion via deterministic particle method. SIAM Journal on Mathematical Analysis 55.3 (2023), pp. 2001-2041. DOI: https://doi.org/10.1137/21M1462994. arXiv: https://arxiv.org/abs/2107.10760.
You can cite the article as
@article{Radici-Stra-2023,
author = {Radici, Emanuela and Stra, Federico},
title = {Entropy solutions of mildly singular nonlocal scalar conservation
laws with congestion via deterministic particle methods},
journal = {SIAM Journal on Mathematical Analysis},
volume = {55},
number = {3},
pages = {2001-2041},
year = {2023},
issn = {0036-1410, 1095-7154},
doi = {10.1137/21M1462994},
url = {https://epubs.siam.org/doi/10.1137/21M1462994},
eprint = {2107.10760},
eprinttype = {arxiv},
eprintclass = {math.AP},
keywords = {entropy solutions, nonlinear mobility, nonlocal interaction,
deterministic particle method, Lagrangian discretization,
numerical scheme},
addendum = {MSC: 35L65, 65M12, 65M75},
}The convergence rate of the scheme is studied (among other results) in the follow-up article
@article{Marconi-Radici-Stra-2024,
author = {Marconi, Elio and Radici, Emanuela and Stra, Federico},
title = {Stability of quasi-entropy solutions of non-local scalar
conservation laws},
journal = {Calculus of Variations and Partial Differential Equations},
volume = {64},
number = {3},
pages = {},
year = {2024},
date = {2024-11-16},
issn = {1432-0835},
doi = {10.1007/s00526-024-02848-9},
url = {https://link.springer.com/article/10.1007/s00526-024-02848-9},
eprint = {2211.02450},
eprinttype = {arxiv},
eprintclass = {math.AP},
keywords = {entropy solutions, conservation laws},
addendum = {MSC: 35L65, 65M12},
}At these links you can find a BibLaTeX file and a PDF file with the compiled bibliography.
Library
ConservationLawsParticles.AbstractModelConservationLawsParticles.DiffusiveIntegratedModelConservationLawsParticles.DiffusiveSampledModelConservationLawsParticles.HyperbolicModelConservationLawsParticles.IntegratedInteractionConservationLawsParticles.IntegratedModelConservationLawsParticles.ParabolicModelConservationLawsParticles.SampledInteractionConservationLawsParticles.SampledModelConservationLawsParticles.compute_interactionConservationLawsParticles.diffusionConservationLawsParticles.eachspeciesConservationLawsParticles.external_velocityConservationLawsParticles.gaussian_particlesConservationLawsParticles.integrated_interactionConservationLawsParticles.interactionConservationLawsParticles.make_velocitiesConservationLawsParticles.make_velocityConservationLawsParticles.mobilityConservationLawsParticles.num_speciesConservationLawsParticles.pwc_densitiesConservationLawsParticles.pwc_densities!ConservationLawsParticles.pwc_densityConservationLawsParticles.pwc_densityConservationLawsParticles.sampled_interactionConservationLawsParticles.speciesConservationLawsParticles.@time_independent
Public
ConservationLawsParticles.AbstractModel — Typeabstract type AbstractModel{S}Abstract type representing a particle model.
ConservationLawsParticles.DiffusiveIntegratedModel — TypeDiffusiveIntegratedModel((V₁, ...), ((W₁₁, ...), ...), (mob₁, ...), (D₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - integrated interactions
Wᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ, - diffusions
Dᵢ.
See also DiffusiveSampledModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = DiffusiveIntegratedModel(
(V, V),
((W_attr, W_rep),
(W_rep, W_attr)),
(mobρ, mobσ),
(Diffusion(1), nothing));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> velocities_diff(x, model, 0.) .|> v -> round(v; digits=2)
([6.03, -0.76, 0.63, -6.33], [23.26, 0.9, 0.71, -8.77, -55.23])ConservationLawsParticles.DiffusiveSampledModel — TypeDiffusiveSampledModel((V₁, ...), ((W′₁₁, ...), ...), (mob₁, ...), (D₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - sampled interactions
W′ᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ, - diffusions
Dᵢ.
See also DiffusiveIntegratedModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = DiffusiveSampledModel(
(V, V),
((Wprime_attr, Wprime_rep),
(Wprime_rep, Wprime_attr)),
(mobρ, mobσ),
(Diffusion(1), nothing));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> velocities_diff(x, model, 0.) .|> v -> round(v; digits=2)
([5.7, -0.81, 0.34, -6.64], [22.41, 1.12, 1.52, -7.87, -54.54])ConservationLawsParticles.HyperbolicModel — TypeHyperbolicModel((V₁, ...), ((I₁₁, ...), ...), (mob₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - integrated interactions
Iᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ.
See also ParabolicModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = HyperbolicModel(
(V, V),
((SampledInteraction(Wprime_attr), IntegratedInteraction(W_rep)),
(IntegratedInteraction(W_rep), SampledInteraction(Wprime_attr))),
(mobρ, mobσ));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> abstract_velocities(x, model, 0.) .|> v -> round(v; digits=2)
([6.27, 0.26, -0.38, -6.56], [22.92, 0.82, 0.71, -8.68, -54.89])ConservationLawsParticles.IntegratedInteraction — TypeIntegratedInteraction(W)Represents the integrated interaction induced by the function W.
The interaction can be evaluated with compute_interaction.
See also SampledInteraction.
ConservationLawsParticles.IntegratedModel — TypeIntegratedModel((V₁, ...), ((W₁₁, ...), ...), (mob₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - integrated interactions
Wᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ.
See also SampledModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = IntegratedModel(
(V, V),
((W_attr, W_rep),
(W_rep, W_attr)),
(mobρ, mobσ));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> velocities(x, model, 0.) .|> v -> round(v; digits=2)
([6.62, 0.31, -0.43, -6.91], [23.26, 0.9, 0.71, -8.77, -55.23])ConservationLawsParticles.ParabolicModel — TypeParabolicModel((V₁, ...), ((I₁₁, ...), ...), (mob₁, ...), (D₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - interactions
Iᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ, - diffusions
Dᵢ.
See also HyperbolicModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = ParabolicModel(
(V, V),
((SampledInteraction(Wprime_attr), IntegratedInteraction(W_rep)),
(IntegratedInteraction(W_rep), SampledInteraction(Wprime_attr))),
(mobρ, mobσ),
(Diffusion(1), nothing));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> abstract_velocities(x, model, 0.) .|> v -> round(v; digits=2)
([5.68, -0.8, 0.68, -5.97], [22.92, 0.82, 0.71, -8.68, -54.89])ConservationLawsParticles.SampledInteraction — TypeSampledInteraction(W′)Represents the sampled interaction induced by the function W′.
The interaction can be evaluated with compute_interaction.
See also IntegratedInteraction.
ConservationLawsParticles.SampledModel — TypeSampledModel((V₁, ...), ((W′₁₁, ...), ...), (mob₁, ...))Represents a particles system with:
- external velocities
Vᵢ, - sampled interactions
W′ᵢⱼ(this is the effect of the speciesjon the speciesi), - mobilities
mobᵢ.
See also IntegratedModel.
Examples
julia> using ConservationLawsParticles.Examples, RecursiveArrayTools
julia> model = SampledModel(
(V, V),
((Wprime_attr, Wprime_rep),
(Wprime_rep, Wprime_attr)),
(mobρ, mobσ));
julia> x = ArrayPartition(gaussian_particles(2, 4), collect(range(-3, 4, length=5)));
julia> round.(x; digits=2)
([-2.0, -0.3, 0.3, 2.0], [-3.0, -1.25, 0.5, 2.25, 4.0])
julia> velocities(x, model, 0.) .|> v -> round(v; digits=2)
([6.29, 0.25, -0.72, -7.23], [22.41, 1.12, 1.52, -7.87, -54.54])ConservationLawsParticles.compute_interaction — Functioncompute_interaction(t, x, interaction, ys[, dens_diff])Evaluates the interaction generated by the particles ys at (t, x).
interaction can be a SampledInteraction, an IntegratedInteraction, or any other form of interaction that provides this interface.
See also SampledInteraction, IntegratedInteraction, sampled_interaction, integrated_interaction.
ConservationLawsParticles.diffusion — Methoddiffusion(mod::AbstractModel, i::Integer) -> Any
Returns the diffusion associated to the species i.
ConservationLawsParticles.eachspecies — Methodeachspecies(mod::AbstractModel) -> Any
Returns an iterator over the indices of the species of the model.
ConservationLawsParticles.external_velocity — Methodexternal_velocity(mod::AbstractModel, i::Integer) -> Any
Returns the external velocity field associated to the species i.
ConservationLawsParticles.gaussian_particles — Methodgaussian_particles(width::Real, number::Integer)Generates number particles defining a Gaussian distribution with variance 1/2 truncated to the interval [-width, width].
ConservationLawsParticles.integrated_interaction — Functionintegrated_interaction([t,] x, W, ys[, dens_diff])Computes $-(W' * \rho)(t, x) = -(W * \rho')(t,x)$, where $\rho$ is the piecewise-constant density associated to the particles ys.
It t is omitted, then W(x) is assumed independent of time.
To ensure the correctness of the computation, dens_diff must coincide with diff(pwc_density(ys)). It can be pre-computed and passed explicitly to allow reuse (as an optimization).
See also sampled_interaction, compute_interaction.
ConservationLawsParticles.interaction — Methodinteraction(
mod::AbstractModel,
i::Integer,
j::Integer
) -> Any
Returns the interaction exerted on the species i by the species j.
ConservationLawsParticles.mobility — Methodmobility(mod::AbstractModel, i::Integer) -> Any
Returns the mobility associated to the species i.
ConservationLawsParticles.num_species — Methodnum_species(mod::AbstractModel) -> Any
Returns the number of species of the model.
ConservationLawsParticles.pwc_densities — Functionpwc_densities(xs::AbstractVector{<:Real}...)
pwc_densities(xs::Tuple{Vararg{AbstractVector{<:Real}}})Given a list of n families of particles of lengths l₁, …, lₙ, returns an n-tuple whose s-th entry is an n×2×lₛ array densₛ.
This array contains the information about the densities of all the species around the particles of the s-th species. This array is indexed as [o,side,i], where o∈{1, …, n} is the index of the other species whose density we want to know, side is either 1 if we want to know the left density or 2 for the right density, and i∈{1, …, lₛ} is the index of the particle.
In other words, densₛ[o,side,i] is the density of the o-th species on the left (side=1) or right (side=2) of the i-th particle of the s-th species.
See also pwc_densities! for an inplace version.
Examples
julia> pwc_densities([0,1,2], [0,2,3])
([0.0 0.5; 0.0 0.25;;; 0.5 0.5; 0.25 0.25;;; 0.5 0.0; 0.25 0.5], [0.0 0.5; 0.0 0.25;;; 0.5 0.0; 0.25 0.5;;; 0.0 0.0; 0.5 0.0])ConservationLawsParticles.pwc_densities! — Functionpwc_densities!(dens, xs::AbstractVector{<:Real}...)
pwc_densities!(dens, xs::Tuple{Vararg{AbstractVector{<:Real}}})This is an inplace version of pwc_densities that writes the result in dens.
See the documentation of pwc_densities for an explanation.
ConservationLawsParticles.pwc_density — Methodpwc_density(x::AbstractVector, y::AbstractVector)The returned x_dens is indexed as x_dens[x_or_y, left_or_right, i].
ConservationLawsParticles.pwc_density — Methodpwc_density(x::AbstractVector)Returns the piecewise constant probability density from the quantile particle positions.
Let the quantile particles be at positions x₀, x₁, …, xₙ (remember that in Julia they are actually x[1], x[2], ..., x[n], x[n+1]). The returned array R is such that R[1] = R[n+2] = 0 and R[i] = 1 / (n * (x[i] - x[i-1])) for the intermediate indices (this formula holds also for the first and last entry if we assume x[0] = -∞ and x[n+2] = ∞).
Examples
julia> pwc_density([0, 1, 3])
4-element Vector{Float64}:
0.0
0.5
0.25
0.0ConservationLawsParticles.sampled_interaction — Functionsampled_interaction([t,] x, W′, ys)Computes $-(W' * \dot\rho)(t, x)$, which is the sampled approximation of $-(W' * \rho)(t, x)$, where $\rho$ is the piecewise-constant density associated to the particles ys.
It t is omitted, then W′(x) is assumed independent of time.
See also integrated_interaction, compute_interaction.
ConservationLawsParticles.species — Methodspecies(state, i::Integer) -> Any
Returns the particles associated to the species i.
ConservationLawsParticles.@time_independent — MacroAutomatically define a method which takes time as first argument and discards it.
Examples
The definition
@time_independent V(x) = -x^3is equivalent to
V(x) = -x^3
V(t, x) = V(x)This also works with more than one argument, for instance
@time_independent V(x₁, x₂) = x₁ * x₂is equivalent to
V(x₁, x₂) = x₁ * x₂
V(t, x₁, x₂) = V(x₁, x₂)Private
ConservationLawsParticles.make_velocities — Methodmake_velocities((V₁, ...), ((W′₁₁, ...), ...), (mob₁, ...))Creates a function velocities(dx, x, p, t) that computes the velocity of the particles under the influence of the external velocities Vᵢ, mutual interactions W′ᵢⱼ and congestions given by mobᵢ.
See also make_velocity.
ConservationLawsParticles.make_velocity — Methodmake_velocity(V::Function, Wprime::Function, mobility::Function)Creates a function velocity(dx, x, p, t) that computes the velocity of the particles under the influence of an external velocity V, a mutual interaction Wprime and the congestion given by mobility.
See also make_velocities.