Poisson solver

SemiLagrangian.PoissonConstType
struct PoissonConst{T, N, Nsp, Nv, type, typeadd}
PoissonConst{T, Nsp, Nv}
PoissonConst(adv::Advection{T, Nsp, Nv, Nsum}; isfftbig=true)

Constant data for the computation of poisson coefficients

Arguments

  • adv::Advection{T, Nsp, Nv, Nsum, timeopt} : Advection constant data
  • isfftbig=true: if true compute the fttbig structure

Implementation

  • adv : Advection constant data
  • `v_k' : vector of vector of fourier coefficents for integration for each space dimension
  • fctv_k : Array of space dimensions of the inverse of the norm of fourier coefficients
  • pfftbig : Fourier data for space dimensions
source
SemiLagrangian.PoissonVarType
mutable struct PoissonVar{T, N, Nsp, Nv, type, typeadd} <: SemiLagrangian.AbstractExtDataAdv
PoissonVar{T, N, Nsp, Nv} <: AbstractExtDataAdv{T}
PoissonVar(pc::PoissonConst{T, N, Nsp, Nv})

mutable structure of variable data for the poisson computation

Arguments

  • pc::PoissonConst{T, N, Nsp, Nv} : poisson constant data

Implementation

  • pc::PoissonConst{T, Nsp, Nv} : poisson constant data
  • rho::Array{T, Nsp} : result of the compute_charge that is the sum along velocity dimensions
  • t_elfield::NTuple{Nsp,Array{Complex{T}, Nsp}} : electric fields initialized at each beginning of velocity advection subseries
source
SemiLagrangian.compute_charge!Method
compute_charge!(rho, t_mesh_v, f)
compute_charge!(rho, mesh_v, f)

Compute charge density from phase space distribution f.

ρ(x,t) = ∫ f(x,v,t) dv

Arguments

  • rho::Array{T,Nsp}: output result density array.
  • t_mesh_v::NTuple{Nv,UniformMesh{T}}: velocity mesh.
  • f::Array{T,Nsum}: distribution function array.
source
SemiLagrangian.compute_charge!Method
compute_charge!(self, advd)
compute_charge!( self::PoissonVar, advd::AdvectionData)

Compute charge density

ρ(x,t) = ∫ f(x,v,t) dv

Argument

  • self::AdvectionData : mutable structure of variables data.
source
SemiLagrangian.compute_eeMethod
compute_ee(self)
compute_ee(self::AdvectionData)

Compute electric energy || E(t,.) ||_L2

Argument

  • self::AdvectionData: advection data structure.
source
SemiLagrangian.compute_eeMethod
compute_ee(t_mesh_sp, t_elf)
compute_ee(t_mesh_sp, t_elf)

Compute electric energy || E(t,.) ||_L2

Arguments

  • t_mesh_sp::NTuple{N,UniformMesh{T}}: space mesh.
  • t_elf::NTuple{N,Array{T,N}}: electric field.
source
SemiLagrangian.compute_elfield!Method
compute_elfield!(self)
compute_elfield!( self:PoissonVar)

computation of electric field ∇.e = - ρ

Argument

  • self::PoissonVar : mutable structure of variables data.
source
SemiLagrangian.compute_elfield!Method
compute_elfield!(elf, mesh, rho)
compute_elfield!(elf::Array{T,1}, mesh::UniformMesh{T}, rho::Array{T,1}) where{T}

Computation of electric field of one dimension. ∇.e = - ρ

Argument

  • elf::Array{T,1}: output Vector.
  • mesh::UniformMesh{T} : mesh of the vector
  • rho::Array{T,1} : rho computed before
source
SemiLagrangian.compute_keMethod
compute_ke(self)
compute_ke(self::AdvectionData)

Compute kinetic Energy.

∫∫ v^2 f(x,v,t) dv dx

Arguments

  • self::AdvectionData: mutable structure of variables data.
source
SemiLagrangian.compute_keMethod
compute_ke(t_mesh_sp, t_mesh_v, f)
compute_ke(t_mesh_sp, t_mesh_v, f)

Compute kinetic Energy from phase space distribution f.

∫∫ v^2 f(x,v,t) dv dx

Arguments

  • t_mesh_sp::NTuple{Nsp, UniformMesh{T}}: space mesh.
  • t_mesh_v::NTuple{Nv, UniformMesh{T}}: velocity mesh.
  • f::Array{T,Nsum}: distribution function array.
source
SemiLagrangian.isvelocityMethod
isvelocity(pv, advd)
initcoef!(pv::PoissonVar{T, N,Nsp, Nv}, self::AdvectionData{T, N})

Implementation of the interface function that is called at the begining of each advection This is implementation for Vlasov-Poisson equation

source