Automatic documentation
Table of contents on this page.
- Automatic documentation
Peridynamics functions
PeriDyn.__I
— Method__I(x::AbstractArray)
Returns x of the same type as of x (identity function). It is not a true identity function as it does not return a copy of x but returns x itself.
Arguments
x::AbstractArray
: input array
Returns
x
: x of the same type as of x (identity function)
PeriDyn.__O
— Method__O(x::AbstractArray)
Returns zero of the same type as of x.
Arguments
x::AbstractArray
: input array
Returns
zero(x)
: zero of the same type as of x
PeriDyn.__norm
— Method__norm(x::AbstractArray)
Returns magnitude of x ($\lVert x \rVert$).
Arguments
x::AbstractArray
: input array
Returns
y
: magnitude of x
PeriDyn.__unit
— Method__unit(x::AbstractArray)
Returns unit vector of x.
Arguments
x::AbstractArray
: input array
Returns
y
: unit vector of x
PeriDyn.cal_damage
— Methodcal_damage(env)
Calculates damage for each bond for a given environment.
Arguments
env
: environment
PeriDyn.cal_family!
— Methodcal_family!(family, x, horizon)
Calculates family for given positions and horizon.
Arguments
family
: familyx
: positionshorizon
: horizon
PeriDyn.cal_family
— Methodcal_family(x, horizon, max_neigh)
Calculates family for given positions and horizon.
Arguments
x
: positionshorizon
: horizonmax_neigh
: maximum number of neighbors
PeriDyn.cell_number
— Methodcell_number(i, j, k, N)
Calculates cell number for given i,j,k cell indices (when in a list).
Arguments
i
: cell index in x directionj
: cell index in y directionk
: cell index in z directionN
: number of cells in each direction
PeriDyn.dilatation!
— Methoddilatation!(theta, y, x, intact, family, volume, m, particle_size, horizon, device::Type{Val{:cuda}})
It gives dilatation as given ordinary state material model.
$\theta_i = \frac{3}{m_i} \sum_{j \in \mathcal{N}_i} \omega_{ij} \times \lVert X_{ij} \rVert \times e_{ij} \times vol_j \times s_{ij}$
where, $\omega_{ij}$ is influence function, $X_{ij}$ is $x_j - x_i$, $e_{ij}$ is bond extention, $vol_j$ is particle volume, $s_{ij}$is horizon correction factor.
Arguments
theta
: dilatationy
: deformed positionsx
: initial positionsintact
: intact bondsfamily
: family of particlesvolume
: particle volumem
: massparticle_size
: particle sizehorizon
: horizondevice::Type{Val{:cuda}}
: device type
Returns
nothing
:theta
is updated in-place
PeriDyn.dilatation!
— Methoddilatation!(theta, y, x, intact, family, volume, m, particle_size, horizon, device::Type{Val{:cpu}})
It gives dilatation as given ordinary state material model.
$\theta_i = \frac{3}{m_i} \sum_{j \in \mathcal{N}_i} \omega_{ij} \times \lVert X_{ij} \rVert \times e_{ij} \times vol_j \times s_{ij}$
where, $\omega_{ij}$ is influence function, $X_{ij}$ is $x_j - x_i$, $e_{ij}$ is bond extention, $vol_j$ is particle volume, $s_{ij}$is horizon correction factor.
Arguments
theta
: dilatationy
: deformed positionsx
: initial positionsintact
: intact bondsfamily
: family of particlesvolume
: particle volumem
: weighted volumeparticle_size
: particle sizehorizon
: horizondevice::Type{Val{:cpu}}
: device type
Returns
nothing
:theta
is updated in-place
PeriDyn.dilatation!
— Methoddilatation!(y, mat, device)
It gives dilatation as given by ordinary state material model. Calls dilatation!(theta, y, x, intact, family, volume, m, particle_size, horizon, device)
.
Arguments
y
: deformed positionsmat
: material modeldevice
: device type
Returns
nothing
: In-place operation
PeriDyn.dilatation
— Methoddilatation(y, x, intact, family, volume, m, particle_size, horizon, device::Type{Val{:cpu}})
It gives dilatation as given ordinary state material model.
$\theta_i = \frac{3}{m_i} \sum_{j \in \mathcal{N}_i} \omega_{ij} \times \lVert X_{ij} \rVert \times e_{ij} \times vol_j \times s_{ij}$
where, $\omega_{ij}$ is influence function, $X_{ij}$ is $x_j - x_i$, $e_{ij}$ is bond extention, $vol_j$ is particle volume, $s_{ij}$is horizon correction factor.
Arguments
y
: deformed positionsx
: initial positionsintact
: intact bondsfamily
: family of particlesvolume
: particle volumem
: massparticle_size
: particle sizehorizon
: horizondevice::Type{Val{:cpu}}
: device type
Returns
theta
: dilatation
PeriDyn.dilatation
— Methoddilatation(y, x, intact, family, volume, m, particle_size, horizon, device::Type{Val{:cuda}})
It gives dilatation as given ordinary state material model.
$\theta_i = \frac{3}{m_i} \sum_{j \in \mathcal{N}_i} \omega_{ij} \times \lVert X_{ij} \rVert \times e_{ij} \times vol_j \times s_{ij}$
where, $\omega_{ij}$ is influence function, $X_{ij}$ is $x_j - x_i$, $e_{ij}$ is bond extention, $vol_j$ is particle volume, $s_{ij}$is horizon correction factor.
Arguments
y
: deformed positionsx
: initial positionsintact
: intact bondsfamily
: family of particlesvolume
: particle volumem
: massparticle_size
: particle sizehorizon
: horizondevice::Type{Val{:cuda}}
: device type
Returns
theta
: dilatation
PeriDyn.get_cells
— Methodget_cells(x, horizon; max_part=30)
Calculates cells for given positions and horizon.
Arguments
x
: positionshorizon
: horizon
Keyword Arguments
max_part=30
: number of partitions in each direction
PeriDyn.horizon_correction
— Methodhorizon_correction(dr, ps, hr)
It gives horizon correction factor (It will give 1 as of now).
Arguments
dr
: $X_{ij}$ps
: particle sizehr
: horizon
Returns
s
: horizon correction factor
PeriDyn.influence_function
— Methodinfluence_function(dr)
It gives influence function as given ordinary state material model. Return $\frac{1}{ \lVert dr \rVert}$ as of now.
Arguments
dr
: $X_{ij}$
Returns
- $\omega_{ij}$: influence function
PeriDyn.loop_over_neighs!
— Methodloop_over_neighs!(y_, mat, fn_map, device; fn_reduce=(args...)->())
The function loops over neighbors and applies fn_map
to each neighbor and fn_reduce
to reduce the results.
The function fn_map
should have the following signature:
fn_map(mat, i, j, k, type1, type2,
xij, yij, extension, s,
wij, wji, items)
where,
mat
: material modeli
: particle indexj
: neighbor indexk
: family indextype1
: type of particlei
type2
: type of particlej
xij
:x_j - x_i
yij
:y_j - y_i
extension
:y_j - y_i
s
:s_ij
bond stretchwij
: influence functionw_ij
wji
: influence functionw_ji
items
: output ofget_items(mat)
function
The function fn_reduce
should have the following signature:
fn_reduce(mat, i, items)
where,
mat
: material modeli
: particle indexitems
: output ofget_items(mat)
function
Arguments
y_
: deformed positionsmat
: material modelfn_map
: function to be applied to each neighbordevice::Type{Val{:cuda}}
: device type
Keyword Arguments
fn_reduce=(args...)->()
: function to reduce the results
Returns
nothing
:y_
is updated in-place
PeriDyn.neigh_cells
— Methodneigh_cells(i, j, k, N)
Calculates neighbor cells for given i,j,k cell indices (when in a list).
Arguments
i
: cell index in x directionj
: cell index in y directionk
: cell index in z directionN
: number of cells in each direction
PeriDyn.weighted_volume
— Methodweighted_volume(x, volume, particle_size, family, horizon)
It gives weighted volume as given ordinary state material model.
$m_i = \sum_{j \in \mathcal{N}_i} \omega_{ij} \times \lVert X_{ij} \rVert^2 \times vol_j \times s_{ij}$
where,
- $\omega_{ij}$: influence function
- $X_{ij}$: $x_j - x_i$
- $vol_j$: particle volume
- $s_{ij}$: horizon correction
Arguments
x
: initial positionsvolume
: particle volumeparticle_size
: particle sizefamily
: family of particleshorizon
: horizon
Returns
m
: weighted volume
Simulation environment
PeriDyn.GeneralEnvironment
— TypeGeneralEnvironment <: SimulationEnvironment
Type for holding parameters for a simulation.
Fields
id::Int64
: ID of the environmenttype::AbstractArray{Int64,1}
: Type of each particlebid::AbstractArray{Int64,1}
: Block ID of each particleghost_particles::AbstractArray{Int64,1}
: Ghost particlesstate::Int64
: State of the environmenttime_step::Int64
: Time step of the environmentdt::T where T<:QF
: Time step sizey::AbstractArray{T,2} where T<:QF
: Position of each particlev::AbstractArray{T,2} where T<:QF
: Velocity of each particlef::AbstractArray{T,2} where T<:QF
: Force of each particlep::AbstractArray{T,2} where T<:QF
: Momentum of each particlevolume::AbstractArray{T,1} where T<:QF
: Volume of each particlemass::AbstractArray{T,1} where T<:QF
: Mass of each particledensity::AbstractArray{T,1} where T<:QF
: Density of each particleintact0::AbstractArray{Int64, 1}
: Intact particles informationneighs::AbstractArray{Int64,2}
: Neighbors of each particleboundary_conditions::AbstractArray{T, 1} where T
: Boundary conditionsshort_range_contact::AbstractArray{T, 1} where T
: Short range contactmaterial_blocks::AbstractArray{T, 1} where T
: Material blocksboundaries::Tuple
: BoundariesCollect!::Function
: Function for collecting dataParams::Dict{Symbol, Any}
: ParametersOut::Dict{Symbol, Any}
: Outputcprint::Function
: Function for printing
PeriDyn.SimulationEnvironment
— TypeSimulationEnvironment
Abstract type for holding parameters for a simulation.
PeriDyn.Env
— MethodEnv(id::Int64, materials,short_range_contact, boundary_conds, dt; state=2, bskin=0.5, units=false)
Constructor for GeneralEnvironment
type.
Arguments
id::Int64
: Environment IDmaterials::AbstractArray{Material, 1}
: Materialsshort_range_contact::AbstractArray{ContactModel, 1}
: Short range contact modelsboundary_conds::AbstractArray{BoundaryCondition, 1}
: Boundary conditionsdt::T where T<:QF
: Time stepstate::Int64
: State of the environmentbskin::T where T<:QF
: Boundary skinunits::Bool
: Units
Returns
env::GeneralEnvironment
: General environment
PeriDyn.printdata
— Methodprintdata(env)
Print (and save to log.txt) data of an environment during simulation.
If PRINT_DEFAULT_DATA
is set to true
, the following data will be printed:
envID
: Environment IDtime
: Simulation timepx
,py
,pz
: Total momentum in x, y, z directionsFx
,Fy
,Fz
: Total force in x, y, z directionsdamage
: Total damage- and the data printed by
cprint
function of the environment
Arguments
env
: Environment
PeriDyn.set_env_active!
— Methodset_env_active!(env)
Set environment state active.
PeriDyn.set_env_idel!
— Methodset_env_idel!(env)
Set environment state idel.
PeriDyn.set_env_inactive!
— Methodset_env_inactive!(env)
Set environment state inactive.
PeriDyn.set_ghost_particles!
— Methodset_ghost_particles!(env,ghost_particles)
Set ghost particles for an environment.
Material models
General material models and functions
PeriDyn.GeneralMaterial
— TypeGeneralMaterial
General peridynamics material type.
Fields
y::AbstractArray{Float64,2}
: Deformed position of the material.velocity::AbstractArray{Float64,2}
: Velocity of the material.x::AbstractArray{Float64,2}
: Position of the material.volume::AbstractArray{Float64,1}
: Volume of the material.type::AbstractArray{Int64,1}
: Type of the material.particle_size::Float64
: Particle size of the material.horizon::Float64
: Horizon of the material.family::VeryBigArray{Int64,2}
: Family of the material.intact::VeryBigArray{Bool, 2}
: Intact of the material.weighted_volume::AbstractArray{Float64,1}
: Weighted volume of the material.deformed::AbstractVector{Bool}
: Deformed of the material.skip_bb::Bool
: Skip bond based material (no idea why it is there).
PeriDyn.GeneralMaterial
— MethodGeneralMaterial(y0, v0, x, volume, type, horizon; max_neigh=100, particle_size=0, skip_bb=false)
Creates a GeneralMaterial
type.
Arguments
y0::AbstractArray{Float64,2}
: Initial Deformed position of the material.v0::AbstractArray{Float64,2}
: Initial velocity of the material.x::AbstractArray{Float64,2}
: Initial position of the material.volume::AbstractArray{Float64,1}
: Volume of the material.type::AbstractArray{Int64,1}
: Type of the material.horizon::Float64
: Horizon of the material.
Keyword Arguments
max_neigh::Int64
= 100: Maximum number of neighbors.particle_size::Float64
= 0: Particle size of the material.skip_bb::Bool
= false: Skip bond based material.
Returns
mat::GeneralMaterial
: General material type.
PeriDyn.GeneralMaterial
— MethodGeneralMaterial(items::Dict, args...; kwargs...)
Creates a GeneralMaterial
type.
Arguments
items::Dict
: Dictionary of the fields of theGeneralMaterial
type.args...
: Arguments of theGeneralMaterial
type.kwargs...
: Keyword arguments of theGeneralMaterial
type.
Returns
mat::GeneralMaterial
: General material type.
See also
PeriDyn.PeridynamicsMaterial
— TypePeridynamicsMaterial
Abstract Peridynamics material type.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(name, type, bid, gen, spc)
Creates a peridynamics material type with given name, type, block id, general material type and specific material type.
Arguments
name::String
: Name of the material.type::Union{UnitRange{Int64}, AbstractVector{Int64}}
: Type of the material.bid::Int64
: Block id of the material.gen::GeneralMaterial
: General material type.spc::SpecificMaterial
: Specific material type.
Returns
mat::PeridynamicsMaterial
: Peridynamics material type.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(bid, gen, spc; name="PeridynamicsMaterial")
Creates a peridynamics material type with given block id, general and specific material type. name
is optional.
Arguments
bid::Int64
: Block id of the material.gen::GeneralMaterial
: General material type.spc::SpecificMaterial
: Specific material type.name::String
: Name of the material.
Returns
mat::PeridynamicsMaterial
: Peridynamics material type.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(gen, spc; name="PeridynamicsMaterial")
Creates a peridynamics material type with given general and specific material type. name
is optional.
Arguments
gen::GeneralMaterial
: General material type.spc::SpecificMaterial
: Specific material type.name::String
: Name of the material.
Returns
mat::PeridynamicsMaterial
: Peridynamics material type.
PeriDyn.SpecificMaterial
— TypeSpecificMaterial
Abstract specific material type.
PeriDyn.force_density!
— Methodforce_density!(f, y, limits, mat::PeridynamicsMaterial)
Calculates the force density of the material. Calls force_density_T
with device
if all the particles are deformed.
Arguments
f::AbstractArray
: Force of the material.y::AbstractArray
: Deformed position of the material.limits::AbstractVector{Int64}
: Limits of the material.mat::PeridynamicsMaterial
: Peridynamics material type.
Returns
nothing
: In-place modification off
.
PeriDyn.force_density_T!
— Methodforce_density_T!(f, y, limits, mat::PeridynamicsMaterial, device::Symbol)
Calculates force density of material pair particles. Converts device
to Type{Val{:device}}
and calls force_density_T
.
Arguments
f::AbstractArray
: Force of the material.y::AbstractArray
: Deformed position of the material.limits::AbstractVector{Int64}
: Limits of the material.mat::PeridynamicsMaterial
: Peridynamics material type.device::Symbol
: Device
Returns
nothing
: In-place modification off
.
PeriDyn.force_density_T!
— Methodforce_density_T!(f, y, limits, mat::GeneralMaterial, device::Type{Val{:cpu}}; particles=nothing)
Calculates force density of material pair particles using CPU. Calls force_density_t_ij
for each pair of particles in the material. The force_density_t_ij
function is defined in the specific material type.
Arguments
f::AbstractArray
: Force of the material.y::AbstractArray
: Deformed position of the material.limits::AbstractVector{Int64}
: Limits of the material.mat::GeneralMaterial
: General material type.device::Type{Val{:cpu}}
: Deviceparticles::Union{Nothing, AbstractVector{Int64}}
: Particles to calculate the force.
Returns
nothing
: In-place modification off
.
PeriDyn.force_density_T!
— Methodforce_density_T!(f, y, limits, mat::GeneralMaterial, device::Type{Val{:cuda}}; particles=nothing)
Calculates force density of material pair particles using CUDA. Calls force_density_t_ij
for each pair of particles in the material. The force_density_t_ij
function is defined in the specific material type.
Arguments
f::AbstractArray
: Force of the material.y::AbstractArray
: Deformed position of the material.limits::AbstractVector{Int64}
: Limits of the material.mat::GeneralMaterial
: General material type.device::Type{Val{:cuda}}
: Deviceparticles::Union{Nothing, AbstractVector{Int64}}
: Particles to calculate the force.
Returns
nothing
: In-place modification off
.
PeriDyn.force_density_t_ij
— Methodforce_density_t_ij(mat, args...; kwargs...)
Force density function for peridynamics material. This is default implementation of force_density_t_ij
function
PeriDyn.init
— Methodinit(spc::T1, gen::T2) where {T1<:SpecificMaterial, T2<:GeneralMaterial}
Initialize specific material type. It uses information from GeneralMaterial
to create SpecificMaterial
type. This is default implementation of init
function which returns spc
without any change.
Arguments
spc::T1
: Specific material type.gen::T2
: General material type.
Returns
spc::T1
: Initialized specific material type.
PeriDyn.@newmaterial
— Macronewmaterial(name, fields...)
Creates a new material type. The macro creates two subtypes, one for the PeridynamicsMaterial
and one for the SpecificMaterial
with the given name. The name of the SpecificMaterial
is the given name with the suffix Specific
. The name of the PeridynamicsMaterial
is the given name with the suffix Material
.
For example, if the name is Elastic
, the macro creates two types ElasticMaterial
and ElasticSpecific
. The ElasticMaterial
is a subtype of PeridynamicsMaterial
and the ElasticSpecific
is a subtype of SpecificMaterial
. The ElasticSpecific
type is used to store the specific parameters of the material.
@newmaterial(Elastic,
"young_modulus::Matrix",
"poisson_ratio::Matrix",
"density::Vector")
The user is expected to define the force_density_t_ij
function for the ElasticSpecific
type.
Arguments
name::Symbol
: Name of the material.fields::Symbol
: Fields of the specific material.
Specific material models and functions
PeriDyn.SkipSpecific
— TypeSkipSpecific <: SpecificMaterial
Skip specific material type.
Fields
density::AbstractArray{<:QF, 1}
: Density.critical_stretch::AbstractArray{Float64, 2}
: Critical stretch.
Returns
SkipSpecific
: Skip specific material.
PeriDyn.BondBasedMaterial
— TypeBondBasedMaterial <: PeridynamicsMaterial
Bond based material type.
Fields
name::String
: Name of material.type::Array
: Type of material particles.bid::Int
: Material block id.general::GeneralMaterial
: General material type.specific::BondBasedSpecific
: Bond based specific material type.
PeriDyn.BondBasedSpecific
— TypeBondBasedSpecific
Bond based specific material type.
Fields
bond_stiffness::AbstractArray{T,2}
: Bond stiffness matrix.bulk_modulus::AbstractArray{T,2}
: Bulk modulus matrix.critical_stretch::AbstractArray{T,2}
: Critical stretch matrix.density::AbstractArray{T,1}
: Density vector.func::Function
: Bond force function.
Constructor
BondBasedSpecific(K::AbstractArray, critical_stretch::AbstractArray, density::AbstractArray; horizon=nothing, func=nothing)
BondBasedSpecific(K::Real, critical_stretch::Real, density::Real; kwargs...)
PeriDyn.BondBasedSpecific
— MethodBondBasedSpecific(K::Matrix, critical_stretch::Matrix, density::Vector; horizon=nothing, func=nothing)
Constructs BondBasedSpecific
material type.
Arguments
K::Vector
: Bulk modulus matrix.critical_stretch::Vector
: Critical stretch matrix.density::Vector
: Density vector.
Keyword Arguments
horizon::Real
: Horizon.func::Function
: Bond force function.
Returns
BondBasedSpecific
: Bond based specific material type.
PeriDyn.BondBasedSpecific
— MethodBondBasedSpecific(K::Real, critical_stretch::Real, density::Real; kwargs...)
Constructs BondBasedSpecific
material type.
Arguments
K::Real
: Bulk modulus matrix.critical_stretch::Real
: Critical stretch matrix.density::Real
: Density vector.
Keyword Arguments
horizon::Real
: Horizon.func::Function
: Bond force function.
Returns
BondBasedSpecific
: Bond based specific material type.
PeriDyn.BondBasedSpecific
— MethodBondBasedSpecific(K::Vector, critical_stretch::Vector, density::Vector; kwargs...)
Constructs BondBasedSpecific
material type.
Arguments
K::Vector
: Bulk modulus matrix.critical_stretch::Vector
: Critical stretch matrix.density::Vector
: Density vector.
Keyword Arguments
horizon::Real
: Horizon.func::Function
: Bond force function.
Returns
BondBasedSpecific
: Bond based specific material type.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(name, type, bid, gen, spc::BondBasedSpecific)
Constructs BondBasedMaterial
material type.
PeriDyn.bond_force
— Methodbond_force(s, bond_stiffness)
Bond force function.
Arguments
s::Real
: Bond stretch.bond_stiffness::Real
: Bond stiffness.
Returns
Real
: Bond force.
PeriDyn.force_density_t_ij
— Methodforce_density_t_ij(mat::BondBasedMaterial,
i, j, k, type1, type2,
xij, yij, extension, s,
wij, wji, items)
Calculates force density (actually acceleration) for bond based material type.
Arguments
mat::BondBasedMaterial
: Bond based material type.i::Int
: Index of particle $i$.j::Int
: Index of particle $j \in N_{i}$.k::Int
: Location index of particle $j$ infamily
andintact
.type1::Int
: Type of particle $i$.type2::Int
: Type of particle $j$.xij::Real
: $\lVert X_{ij} \rVert$yij::Real
: $\lVert Y_{ij} \rVert$extension::Real
: Bond extension.s::Real
: Bond stretch.wij::Real
: Influence function ($\omega_{ij}$) for pair $(i,j)$.wji::Real
: Influence function ($\omega_{ji}$) for pair $(j,i)$.items::Tuple
: Items for force calculation.
Returns
Real
: Force density.
PeriDyn.get_items
— Methodget_items(mat::BondBasedMaterial)
Returns items for force calculation.
Arguments
mat::BondBasedMaterial
: Bond based material type.
Returns
bond_stiffness::AbstractArray
: Bond stiffness.
PeriDyn.init
— Methodinit(spc::BondBasedSpecific, gen::GeneralMaterial)
Initializes BondBasedSpecific
material type.
Arguments
spc::BondBasedSpecific
: Bond based specific material type.gen::GeneralMaterial
: General material type.
Returns
BondBasedSpecific
: Bond based specific material type.
PeriDyn.out_of_loop!
— Methodout_of_loop!(y_, mat::BondBasedMaterial, device)
It is called once before the main loop for force calculation. It does nothing for BondBasedMaterial
.
PeriDyn.OrdinaryStateBasedMaterial
— TypeOrdinaryStateBasedMaterial <: PeridynamicsMaterial
Ordinary state based material type.
Fields
name::String
: Name of material.type::Array
: Type of material particles.bid::Int
: Material block id.general::GeneralMaterial
: General material type.specific::OrdinaryStateBasedSpecific
: Specific material type.
Returns
OrdinaryStateBasedMaterial
: Ordinary state based material.
PeriDyn.OrdinaryStateBasedSpecific
— TypeOrdinaryStateBasedSpecific <: SpecificMaterial
Ordinary state based specific material type.
Fields
bulk_modulus::AbstractArray{<:QF, 2}
: Bulk modulus.shear_modulus::AbstractArray{<:QF, 2}
: Shear modulus.critical_stretch::AbstractArray{Float64, 2}
: Critical stretch.density::AbstractArray{<:QF, 1}
: Density.theta::AbstractArray{<:QF, 1}
: Dilatation.
Returns
OrdinaryStateBasedSpecific
: Ordinary state based specific material.
PeriDyn.OrdinaryStateBasedSpecific
— MethodOrdinaryStateBasedSpecific(bulk_modulus, shear_modulus, critical_stretch, density)
Creates ordinary state based specific material type.
Arguments
bulk_modulus::AbstractVector{<:QF}
: Bulk modulus.shear_modulus::AbstractVector{<:QF}
: Shear modulus.critical_stretch::AbstractVector{Float64}
: Critical stretch.density::AbstractVector{<:QF}
: Density.
Returns
OrdinaryStateBasedSpecific
: Ordinary state based specific material.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(name, type, bid, gen, spc::OrdinaryStateBasedSpecific)
Creates ordinary state based material type.
Arguments
name::String
: Name of material.type::Int
: Type of material.bid::Int
: Bond type.gen::GeneralMaterial
: General material.spc::OrdinaryStateBasedSpecific
: Specific material.
Returns
OrdinaryStateBasedMaterial
: Ordinary state based material.
PeriDyn.force_density_t_ij
— Methodforce_density_t_ij(mat::OrdinaryStateBasedMaterial,
i, j, k, type1, type2,
xij, yij, extension, s,
wij, wji, items)
Calculates force density (actually acceleration) for ordinary state based material type.
Arguments
mat::BondBasedMaterial
: Ordinary state based material.i::Int
: Index of particle $i$.j::Int
: Index of particle $j \in N_{i}$.k::Int
: Location index of particle $j$ infamily
andintact
.type1::Int
: Type of particle $i$.type2::Int
: Type of particle $j$.xij::Real
: $\lVert X_{ij} \rVert$yij::Real
: $\lVert Y_{ij} \rVert$extension::Real
: Bond extension.s::Real
: Bond stretch.wij::Real
: Influence function ($\omega_{ij}$) for pair $(i,j)$.wji::Real
: Influence function ($\omega_{ji}$) for pair $(j,i)$.items::Tuple
: Items for force calculation.
Returns
Real
: Force density.
PeriDyn.get_items
— Methodget_items(mat::OrdinaryStateBasedMaterial)
Returns items for force calculation.
Arguments
mat::BondBasedMaterial
: Ordinary state based material.
Returns
Tuple
: Items for force calculation.
(m
, theta
, K
, G
) where m
is weighted volume, theta
is dilatation, K
is bulk modulus, G
is shear modulus.
PeriDyn.init
— Methodinit(spc::OrdinaryStateBasedSpecific, gen::GeneralMaterial)
Initializes ordinary state based specific material type. Set spc.theta
to zero.
Arguments
spc::OrdinaryStateBasedSpecific
: Specific material.gen::GeneralMaterial
: General material.
Returns
OrdinaryStateBasedSpecific
: Specific material.
PeriDyn.out_of_loop!
— Methodout_of_loop!(y_, mat::OrdinaryStateBasedMaterial, device)
Calculates dilatation for ordinary state based material type.
Arguments
y_::AbstractArray
: Deformed positions.mat::OrdinaryStateBasedMaterial
: Ordinary state based material.device::Type{Val{:device}}
: Device type.
Returns
nothing
PeriDyn.DruckerPrager
— TypeDruckerPrager <: PlasticFailureCriteria
Drucker-Prager failure criterion.
PeriDyn.ElastoPlasticSolidMaterial
— TypeElastoPlasticSolidMaterial <: PeridynamicsMaterial
Elastic-plastic solid material.
Fields
name::String
: Name of the materialtype::Array
: Type of the materialbid::Int
: Block idgeneral::GeneralMaterial
: General materialspecific::ElastoPlasticSolidSpecific
: Specific material
PeriDyn.ElastoPlasticSolidSpecific
— TypeElasticPlasticSolidMaterial <: PeridynamicsMaterial
Elastic-plastic solid material.
Fields
bulk_modulus::AbstractArray{<:QF,2}
: Bulk modulusshear_modulus::AbstractArray{<:QF,2}
: Shear moduluscritical_stretch::AbstractArray{Float64,2}
: Critical stretchdensity::AbstractArray{<:QF,1}
: Densityσγ::AbstractArray{<:QF,2}
: Yield stressψ::AbstractArray{<:QF,2}
: Yield stress valueedp::AbstractArray{<:QF,2}
: Plastic deviatoric straintd_norm2::AbstractArray{<:QF,1}
: td_norm2theta::AbstractArray{<:QF,1}
: Dilatationcriteria::PlasticFailureCriteria
: Plastic failure criteria
PeriDyn.ElastoPlasticSolidSpecific
— MethodElastoPlasticSolidSpecific(
bulk_modulus::AbstractArray{<:QF,1},
shear_modulus::AbstractArray{<:QF,1},
critical_stretch::Array{Float64,1},
density::AbstractArray{<:QF,1},
σγ::AbstractArray{<:QF,1};
criteria = VonMises())
Construct an ElastoPlasticSolidSpecific
material.
Arguments
bulk_modulus::AbstractArray{<:QF,1}
: Bulk modulusshear_modulus::AbstractArray{<:QF,1}
: Shear moduluscritical_stretch::Array{Float64,1}
: Critical stretchdensity::AbstractArray{<:QF,1}
: Densityσγ::AbstractArray{<:QF,1}
: Yield stresscriteria::PlasticFailureCriteria
: Plastic failure criteria
Returns
spc::ElastoPlasticSolidSpecific
: Specific material
PeriDyn.MohrCoulomb
— TypeMohrCoulomb <: PlasticFailureCriteria
Mohr-Coulomb failure criterion.
PeriDyn.PeridynamicsMaterial
— MethodPeridynamicsMaterial(name, type, bid, gen, spc::ElastoPlasticSolidSpecific)
Construct an ElastoPlasticSolidMaterial
.
Arguments
name::String
: Name of the materialtype::Array
: Type of the materialbid::Int
: Block idgen::GeneralMaterial
: General materialspc::ElastoPlasticSolidSpecific
: Specific material
Returns
mat::ElastoPlasticSolidMaterial
: Material
PeriDyn.PlasticFailureCriteria
— TypePlasticFailureCriteria
Abstract type for plastic failure criteria.
PeriDyn.Tresca
— TypeTresca <: PlasticFailureCriteria
Tresca failure criterion.
PeriDyn.VonMises
— TypeVonMises <: PlasticFailureCriteria
Von Mises failure criterion.
PeriDyn.effective_stress2
— Methodeffective_stress2(σ₁, σ₂, σ₃, criteria::VonMises)
Calculate the effective stress for Von Mises failure criterion.
\[σₑ^2 = ((σ₁ - σ₂)^2 + (σ₂ - σ₃)^2 + (σ₁ - σ₃)^2) / 2\]
Arguments
σ₁::Real
: First principal stress.σ₂::Real
: Second principal stress.σ₃::Real
: Third principal stress.criteria::VonMises
: Von Mises failure criterion.
Returns
Real
: Effective stress.
PeriDyn.force_density_t_ij
— Methodforce_density_t_ij(mat::ElastoPlasticSolidMaterial,
i, j, k, type1, type2,
xij, yij, extension, s,
wij, wji, items)
Calculate the force density.
Arguments
mat::BondBasedMaterial
: Bond based material type.i::Int
: Index of particle $i$.j::Int
: Index of particle $j \in N_{i}$.k::Int
: Location index of particle $j$ infamily
andintact
.type1::Int
: Type of particle $i$.type2::Int
: Type of particle $j$.xij::Real
: $\lVert X_{ij} \rVert$yij::Real
: $\lVert Y_{ij} \rVert$extension::Real
: Bond extension.s::Real
: Bond stretch.wij::Real
: Influence function ($\omega_{ij}$) for pair $(i,j)$.wji::Real
: Influence function ($\omega_{ji}$) for pair $(j,i)$.items::Tuple
: Items for force calculation.
Returns
Real
: Force density.
PeriDyn.get_items
— Methodget_items(mat::ElastoPlasticSolidMaterial)
Get the items of the material.
Arguments
mat::ElastoPlasticSolidMaterial
: Material
Returns
items
: Tuple of items
(weighted_volume
, theta
, td_norm2
, bulk_modulus
, shear_modulus
, σγ
, ψ
, edp
, horizon
, volume
, criteria
)
PeriDyn.init
— Methodinit(spc::ElastoPlasticSolidSpecific, gen::GeneralMaterial)
Initialize the specific material. It sets the yield stress value ψ
equal to $75/8π * \sigma_e^2 / \delta^5$ where $\sigma_e$ is the effective yield stress and $\delta$ is the horizon. It initializes the dilatation theta
, the plastic deviatoric strain edp
and the td_norm2
to zero.
Arguments
spc::ElastoPlasticSolidSpecific
: Specific materialgen::GeneralMaterial
: General material
Returns
spc::ElastoPlasticSolidSpecific
: Initialized specific material
PeriDyn.out_of_loop!
— Methodout_of_loop!(y_, mat::ElastoPlasticSolidMaterial, device)
Calculate the dilatation and the td_norm2 out of the loop.
Arguments
y_::AbstractArray{<:QF,2}
: Displacementmat::ElastoPlasticSolidMaterial
: Materialdevice
: Device
PeriDyn.td_norm2_ij!
— Methodtd_norm2_ij!(mat, i, j, k, type1, type2,
xij, yij, extension, s,
wij, wji, items)
Calculate the td_norm2.
Arguments
mat::BondBasedMaterial
: Bond based material type.i::Int
: Index of particle $i$.j::Int
: Index of particle $j \in N_{i}$.k::Int
: Location index of particle $j$ infamily
andintact
.type1::Int
: Type of particle $i$.type2::Int
: Type of particle $j$.xij::Real
: $\lVert X_{ij} \rVert$yij::Real
: $\lVert Y_{ij} \rVert$extension::Real
: Bond extension.s::Real
: Bond stretch.wij::Real
: Influence function ($\omega_{ij}$) for pair $(i,j)$.wji::Real
: Influence function ($\omega_{ji}$) for pair $(j,i)$.items::Tuple
: Items for force calculation.
Returns
- 'nothing': Update
td_norm2
in-place.
Contact models and functions
PeriDyn.ContactModel
— TypeContactModel
Abstract type for contact model.
PeriDyn.ContactModel11
— TypeContactModel11
Abstract type for contact model for a single material block.
PeriDyn.ContactModel12
— TypeContactModel12
Abstract type for contact model between two material blocks.
Base.show
— MethodPrint the details of a ContactModel object.
Arguments
io::IO
: The output IO stream.i::Union{ContactModel11, ContactModel12}
: The ContactModel object to print.
PeriDyn.ContactModel11_gcal
— MethodContactModel11_gcal(mat1, distanceX, max_neighs)
Calculate the parameters for ContactModel11
based on the input material, distance, and maximum neighbors.
Arguments:
mat1
: The material (typeContactModel11
).distanceX
: The distance factor.max_neighs
: The maximum number of neighbors.
Returns a tuple containing the calculated parameters for ContactModel11
.
PeriDyn.ContactModel12_gcal
— MethodContactModel12_gcal(mat1, mat2, distanceX, max_neighs)
Calculate the parameters for ContactModel12
based on the input materials, distance, and maximum neighbors.
Arguments:
mat1
: The first material (typeContactModel11
).mat2
: The second material (typeContactModel11
).distanceX
: The distance factor.max_neighs
: The maximum number of neighbors.
Returns a tuple containing the calculated parameters for ContactModel12
.
PeriDyn._cudaconvert
— Method_cudaconvert(x::T) where T <: Union{ContactModel11,ContactModel12}
Converts a single ContactModel11
or ContactModel12
object to a CUDA-compatible type.
PeriDyn._cudaconvert
— Method_cudaconvert(x::Vector{T}) where T <: Union{ContactModel11,ContactModel12}
Converts a vector of ContactModel11
or ContactModel12
objects to CUDA-compatible types.
PeriDyn.collision_box
— Methodcollision_box(x1::Array{Float64,2}, x2::Array{Float64,2}, skin::Float64)
Calculates collision box between two material blocks.
Arguments
x1
: Positions of material points (block 1)x2
: Positions of material points (block 2)skin
: Extra distance to consider (usually >= particle size)
Output
box_min
: Minimum position limits for overlapbox_max
: Maximum position limits for overlapifoverlap
: Boolean indicating if there is an overlap
PeriDyn.short_range_contact!
— Methodshort_range_contact!(y, f, type, bid, vol, RM::ContactModel11, device::Type{Val{:cpu}})
Updates (inplace) the contact acceleration of material points.
Arguments
y
: Positions of material pointsf
: Acceleration of material pointstype
: Type of material pointsbid
: BID valuesvol
: Volume valuesRM::ContactModel11
: Repulsion modeldevice::Type{Val{:cpu}}
: Device type for CPU acceleration
Output
- No return value. The function updates
f
in place.
PeriDyn.short_range_contact!
— Methodshort_range_contact!(y, f, type, bid, vol, RM::ContactModel11, device::Type{Val{:cuda}})
Updates (inplace) the contact acceleration of material points.
Arguments
y
: Positions of material pointsf
: Acceleration of material pointstype
: Type of material pointsbid
: BID valuesvol
: Volume valuesRM::ContactModel11
: Repulsion modeldevice::Type{Val{:cuda}}
: Device type for CUDA acceleration
Output
- No return value. The function updates
f
in place.
PeriDyn.short_range_contact!
— Methodshort_range_contact!(y, f, type, bid, vol, RM::ContactModel11)
Updates (inplace) the contact acceleration of material points.
Arguments
y
: Positions of material pointsf
: Acceleration of material pointstype
: Type of material pointsbid
: BID valuesvol
: Volume valuesRM::ContactModel11
: Repulsion model
Output
- No return value. The function updates
f
in place.
PeriDyn.short_range_contact!
— Methodshort_range_contact!(y, f, type, ContactModel)
Updates (inplace) the contact acceleration of material points.
Arguments
y
: Positions of material points.f
: Acceleration of material points.type
: Type of material points.ContactModel
: Repulsion model (see contacts.jl for more details).
Output
- None (Inplace update of f (acceleration)).
PeriDyn.short_range_contact!
— Methodshort_range_contact!(y, f, type, bid, vol, RM::ContactModel12)
Updates (inplace) the contact acceleration of material points.
Arguments
y
: Positions of material pointsf
: Acceleration of material pointstype
: Type of material pointsbid
: BID valuesvol
: Volume valuesRM::ContactModel12
: Repulsion model
Output
- No return value. The function updates
f
in place.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(neighbors, x, search_distance, equi_dist, family, intact, max_part)
Update the neighbor list for contact force calculation.
Arguments
neighbors
: Array storing the neighbor indicesx
: Positions of material pointssearch_distance
: Maximum search distance for neighborsequi_dist
: Equilibrium distance for contactfamily
: Array indicating the family relationship between material pointsintact
: Array indicating if the family relationship is intactmax_part
: Maximum number of particles in a cell
Output
- No return value. The function updates the
neighbors
array in place.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(y, type, RM::ContactModel11, device::Symbol; kwargs...)
Update neighbor list for contact force calculation (1-1 interaction) on a specific device.
Arguments
y
: Positions of material pointstype
: Type of material pointsRM::ContactModel11
: Repulsion modeldevice::Symbol
: Device type for accelerationkwargs...
: Additional keyword arguments
Output
- No return value. The function updates
RM.neighs
in place.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(y, type, RM::ContactModel11, device::Type{Val{:cpu}}; max_part=30)
Update the neighbor list for contact force calculation (1-1 interaction).
Arguments
y
: Positions of material pointstype
: Type of material pointsRM::ContactModel11
: Repulsion model for 1-1 interactiondevice::Type{Val{:cpu}}
: Device type (CPU)max_part=30
: Maximum number of particles in a cell
Output
- No return value. The function updates the neighbor list in the
RM
object.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(y, type, RM::ContactModel11, device::Type{Val{:cuda}}; max_part=30)
Update neighbor list for contact force calculation (1-1 interaction).
Arguments
y
: Positions of material pointstype
: Type of material pointsRM::ContactModel11
: Repulsion model for 1-1 interactiondevice::Type{Val{:cuda}}
: Device type (CUDA)max_part=30
: Maximum number of particles in a cell
Output
- No return value. The function updates the neighbor list in the
RM
object.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(y, type, RM::ContactModel11; kwargs...)
Update neighbor list for contact force calculation (1-1 interaction).
Arguments
y
: Positions of material pointstype
: Type of material pointsRM::ContactModel11
: Repulsion modelkwargs...
: Additional keyword arguments
Output
- No return value. The function updates
RM.neighs
in place.
PeriDyn.update_contact_neighs!
— Methodupdate_contact_neighs!(y, type, RM::ContactModel12; max_part=nothing)
Update neighbor list for contact force calculation (1-2 interaction).
Arguments
y
: Positions of material pointstype
: Type of material pointsRM::ContactModel12
: Repulsion modelmax_part=nothing
: Maximum number of particles (optional)
Output
- No return value. The function updates
RM.neighs
in place.
PeriDyn.NonLinearSpringContactModel11
— TypeNonLinearSpringContactModel11(exponent::Float64, stifness::Float64, type::AbstractVector{Int64}, bid::Int64, material::GeneralMaterial, name::String, equi_dist::Float64, distance::Float64, neighs::AbstractArray{Int64, 2}, max_neighs::Int64)
Nonlinear contact model for 1-1 material blocks.
PeriDyn.NonLinearSpringContactModel12
— TypeNonLinearSpringContactModel12(exponent::Float64, stifness::Float64, pair::Vector{AbstractVector{Int64}}, name::String, equi_dist::Float64, distance::Float64, neighs::AbstractArray{Int64, 2}, max_neighs::Int64)
Nonlinear contact model for 1-2 material blocks.
PeriDyn.NonLinearSpringContactModel
— MethodNonLinearSpringContactModel(exponent::Float64, stifness::Float64, mat1::PeridynamicsMaterial, mat2::PeridynamicsMaterial; distanceD=1.0, distanceX=3.0, max_neighs=50)
Constructs a nonlinear contact model for 1-2 material blocks.
Arguments:
exponent
: The exponent for the contact force calculation.stifness
: The stiffness coefficient for the contact force calculation.mat1
: The first PeridynamicsMaterial.mat2
: The second PeridynamicsMaterial.distanceD
: The distance factor for determining the search distance of neighbors within the same material.distanceX
: The distance factor for determining the search distance of neighbors between different materials.max_neighs
: The maximum number of neighbors to consider.
Returns: A NonLinearSpringContactModel12 object.
PeriDyn.NonLinearSpringContactModel
— MethodNonLinearSpringContactModel(exponent::Float64, stifness::Float64, mat1::PeridynamicsMaterial; distanceX=5, max_neighs=50)
Constructs a nonlinear contact model for 1-1 material blocks.
Arguments:
exponent
: The exponent for the contact force calculation.stifness
: The stiffness coefficient for the contact force calculation.mat1
: The PeridynamicsMaterial.distanceX
: The distance factor for determining the search distance of neighbors.max_neighs
: The maximum number of neighbors to consider.
Returns: A NonLinearSpringContactModel11 object.
PeriDyn.get_contact_force_fn
— Methodget_contact_force_fn(RepMod::NonLinearSpringContactModel11)
Returns a contact force function for a NonLinearSpringContactModel11 object.
Arguments:
RepMod
: The NonLinearSpringContactModel11 object.
Returns: A contact force function that takes a distance vector and returns the contact acceleration.
PeriDyn.LinearSpringContactModel
— MethodLinearSpringContactModel(stiffness::Float64, mat1::PeridynamicsMaterial, mat2::PeridynamicsMaterial; distanceX=5, max_neighs=50)
Constructs a linear contact model for 1-2 material blocks. The constructors internally call the NonLinearSpringContactModel
constructor with a nonlinearity parameter of 1 and the provided arguments.
Arguments:
stiffness
: The stiffness coefficient for the contact model.mat1
: PeridynamicsMaterial representing the first material block.mat2
: PeridynamicsMaterial representing the second material block.distanceX
: Optional keyword argument specifying the distance parameter (default: 5).max_neighs
: Optional keyword argument specifying the maximum number of neighbors (default: 50).
Returns:
- An instance of the LinearSpringContactModel.
PeriDyn.LinearSpringContactModel
— MethodLinearSpringContactModel(stiffness::Float64, mat1::PeridynamicsMaterial; distanceX=5, max_neighs=50)
Constructs a linear contact model for 1-1 material blocks. The constructors internally call the NonLinearSpringContactModel
constructor with a nonlinearity parameter of 1 and the provided arguments.
Arguments:
stiffness
: The stiffness coefficient for the contact model.mat1
: PeridynamicsMaterial representing the material block.distanceX
: Optional keyword argument specifying the distance parameter (default: 5).max_neighs
: Optional keyword argument specifying the maximum number of neighbors (default: 50).
Returns:
- An instance of the LinearSpringContactModel.
PeriDyn.ShortRangeContactModel
— MethodShortRangeContactModel(K, horizon, mat1::PeridynamicsMaterial, mat2::PeridynamicsMaterial; mfactor=1.0, kwargs...)
Create a short range contact model with the given stifness K
and horizon horizon
.
Arguments
K
: Stifness of the contact model.horizon
: Horizon of the contact model.mat1
: Material of the first body.mat2
: Material of the second body.
Keyword Arguments
mfactor
: Multiplication factor for the stifness. Default is1.0
.kwargs...
: Keyword arguments for theLinearSpringContactModel
constructor.
PeriDyn.ShortRangeContactModel
— MethodShortRangeContactModel(K, horizon, mat1::PeridynamicsMaterial; mfactor=1.0, kwargs...)
Create a short range contact model with the given stifness K
and horizon horizon
.
Arguments
K
: Stifness of the contact model.horizon
: Horizon of the contact model.mat1
: Material of the first body.
Keyword Arguments
mfactor
: Multiplication factor for the stifness. Default is1.0
.kwargs...
: Keyword arguments for theLinearSpringContactModel
constructor.
Boundary conditions
PeriDyn.BoundaryCondition
— TypeBoundaryCondition
Abstract type for boundary conditions.
PeriDyn.apply_bc!
— Methodapply_bc!(env, BC::BoundaryCondition, on::Symbol)
Apply the specified boundary condition BC
to the given environment env
on the specified aspect on
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The boundary condition to apply.on::Symbol
: The aspect on which the boundary condition is applied (:position
or:velocity
).
PeriDyn.apply_bc!
— Methodapply_bc!(env, BC::BoundaryCondition, ::Type{Val{:position}})
Apply the general boundary condition BC
to the position aspect of the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The general boundary condition to apply.
PeriDyn.apply_bc!
— Methodapply_bc!(env, BC::BoundaryCondition, ::Type{Val{:velocity}})
Apply the general boundary condition BC
to the velocity aspect of the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The general boundary condition to apply.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::BoundaryCondition)
Apply the boundary condition BC
to the given environment env
at the start of the simulation (t=0).
PeriDyn.check!
— Methodcheck!(env, BC::BoundaryCondition)
Check if the boundary condition BC
has changed and updates the BC
. Used for dynamic boundary conditions.
PeriDyn.@general_bc_p
— Macrogeneral_bc_p()
Macro for defining the common fields of a general boundary condition.
The following fields are defined:
bool
: Boolean array specifying the affected material points.dims
: Boolean array specifying the affected dimensions.last
: Array of the same type as the state vector specifying the last state of the affected material points.onlyatstart
: Flag indicating if the boundary condition is applied only at the start (default:false
).xF
: function for updating the position.vF
: function for updating the velocity.
PeriDyn.FixBC
— TypeFixBC
Struct representing the FixBC boundary condition.
Fields
bool
: Boolean array specifying the affected elements.dims
: Boolean array specifying the affected dimensions.last
: Last position of the affected elements.onlyatstart
: Flag indicating if the boundary condition is applied only at the start.xF
: function for updating the position.vF
: function for updating the velocity.
PeriDyn.FixBC
— MethodFixBC(bool; onlyatstart=false)
Construct a FixBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.
Keyword Arguments
dims
: Boolean array specifying the affected dimensions (default:[true, true, true]
).onlyatstart
: Flag indicating if the boundary condition is applied only at the start (default:false
).
Returns
A FixBC object representing the boundary condition.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::FixBC)
Apply the FixBC boundary condition at time 0 to the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The FixBC boundary condition to apply.
PeriDyn.ToFroBC
— TypeToFroBC
Struct representing the ToFroBC boundary condition.
Fields
bool
: Boolean array specifying the affected elements.dims
: Boolean array specifying the affected dimensions.last
: Last position of the affected elements.onlyatstart
: Flag indicating if the boundary condition is applied only at the start.xF
: function for updating the position.vF
: function for updating the velocity.rate
: Rate at which the elements move.direction
: Direction of movement.freq
: Frequency at which the direction of movement changes.applyafter
: Number of steps after which the frequency is applied.
PeriDyn.ToFroBC
— MethodToFroBC(bool, rate, freq; applyafter=0, onlyatstart=false)
Construct a ToFroBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.rate
: Rate at which the elements move.freq
: Frequency at which the direction of movement changes.
Keyword Arguments
dims
: Boolean array specifying the affected dimensions (default:[true, true, true]
).applyafter
: Number of steps after which the frequency is applied (default: 0).onlyatstart
: Flag indicating if the boundary condition is applied only at the start (default:false
).
Returns
A ToFroBC object representing the boundary condition.
PeriDyn.MoveBC
— MethodMoveBC(bool, rate; kwargs...)
Create a MoveBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.rate
: Rate at which the elements move.kwargs
: Additional keyword arguments passed toToFroBC
.
Returns
A MoveBC object representing the boundary condition. All MoveBC objects are ToFroBC objects with frequency Inf.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::ToFroBC)
Apply the ToFroBC boundary condition at time 0 to the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The ToFroBC boundary condition to apply.
PeriDyn.check!
— Methodcheck!(BC::ToFroBC, env)
Perform a check on the ToFroBC boundary condition.
Arguments
BC
: The ToFroBC boundary condition to check.env
: The environment associated with the boundary condition.
PeriDyn.DeltaScaleBC
— TypeDeltaScaleBC
Struct representing the DeltaScaleBC boundary condition.
Fields
bool
: Boolean array specifying the affected elements.dims
: Boolean array specifying the affected dimensions.last
: Last position of the affected elements.onlyatstart
: Flag indicating if the boundary condition is applied only at the start.xF
: function for updating the position.vF
: function for updating the velocity.
PeriDyn.DeltaScaleBC
— MethodDeltaScaleBC(bool, scale, fixpoint; onlyatstart=false)
Construct a DeltaScaleBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.scale
: Scale factor applied to the elements.
Keyword Arguments
dims
: Boolean array specifying the affected dimensions (default:[true, true, true]
).fixpoint
: Reference point used for scaling.onlyatstart
: Flag indicating if the boundary condition is applied only at the start (default:false
).
Returns
A DeltaScaleBC object representing the boundary condition.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::DeltaScaleBC)
Apply the DeltaScaleBC boundary condition at time 0 to the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The DeltaScaleBC boundary condition to apply.
PeriDyn.check!
— Methodcheck!(BC::DeltaScaleBC, env)
Perform a check on the DeltaScaleBC boundary condition.
Arguments
BC
: The DeltaScaleBC boundary condition to check.env
: The environment associated with the boundary condition.
PeriDyn.ScaleFixWaitBC
— TypeScaleFixWaitBC
Structure representing a ScaleFixWaitBC boundary condition.
Fields
bool
: Boolean array specifying the affected elements.dims
: Boolean array specifying the affected dimensions.last
: Last position of the affected elements.onlyatstart
: Flag indicating if the boundary condition is applied only at the start.xF
: function for updating the position.vF
: function for updating the velocity.checkF
: Function for checking if the boundary condition needs to be applied.
PeriDyn.ScaleFixWaitBC
— MethodScaleFixWaitBC(bool, scale, fixpoint, wait, scalebool; applyafter=0, onlyatstart=false)
Construct a ScaleFixWaitBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.scale
: Scale factor for the elements.fixpoint
: Fix point for the elements.wait
: Number of time steps to wait before applying the condition.scalebool
: Boolean array specifying the elements to be scaled.
Keyword Arguments
dims
: Boolean array specifying the affected dimensions (default:[true, true, true]
).applyafter
: Number of time steps after which the condition is applied (default: 0).onlyatstart
: Boolean indicating whether the condition is applied only at the start (default: false).
Returns
- An instance of ScaleFixWaitBC boundary condition.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::ScaleFixWaitBC)
Apply the ScaleFixWaitBC boundary condition at time 0.
Arguments
env
: Environment in which the condition is applied.BC
: Instance of ScaleFixWaitBC boundary condition.
PeriDyn.check!
— Methodcheck!(BC::ScaleFixWaitBC, env)
Check if the ScaleFixWaitBC boundary condition needs to be applied.
Arguments
BC
: Instance of ScaleFixWaitBC boundary condition.env
: Environment in which the condition is applied.
PeriDyn.ContainerBC
— TypeContainerBC
Struct representing the ContainerBC boundary condition.
Fields
bool
: Boolean array specifying the affected elements.dims
: Boolean array specifying the affected dimensions.last
: Last position of the affected elements.onlyatstart
: Flag indicating if the boundary condition is applied only at the start.xF
: function for updating the position.vF
: function for updating the velocity.
PeriDyn.ContainerBC
— MethodContainerBC(bool; limits=nothing, onlyatstart=false)
Construct a ContainerBC boundary condition.
Arguments
bool
: Boolean array specifying the affected elements.limits
: Limits of the container (default:nothing
).
Keyword Arguments
dims
: Boolean array specifying the affected dimensions (default:[true, true, true]
).onlyatstart
: Flag indicating if the boundary condition is applied only at the start (default:false
).
Returns
A ContainerBC object representing the boundary condition.
PeriDyn.apply_bc_at0!
— Methodapply_bc_at0!(env, BC::ContainerBC)
Apply the ContainerBC boundary condition at time 0 to the given environment env
.
Arguments
env
: The environment to which the boundary condition is applied.BC
: The ContainerBC boundary condition to apply.
Solvers and integrators
General solver functions
PeriDyn.DynamicSolver
— TypeDynamicSolver <: Solver
DynamicSolver is an abstract type for dynamic solvers.
PeriDyn.QuasiStaticSolver
— TypeQuasiStaticSolver <: Solver
QuasiStaticSolver is an abstract type for quasi-static solvers.
PeriDyn.Solver
— TypeSolver
Solver is an abstract type for solvers.
PeriDyn.apply_bc_at0
— Methodapply_bc_at0!(env, start_at)
Apply boundary conditions at t = t0.
Arguments
env
: the simulation environment.start_at
: Int64, the starting step.
PeriDyn.check_boundaries!
— Methodcheck_boundaries!(env)
Check if any material point is outside the defined boundaries. Do nothing.
Arguments
env
: the simulation environment.
PeriDyn.run!
— Methodrun!(envs, N::Int64, solver::Solver; filewrite_freq::Int64=10,
neigh_update_freq::Int64=1, average_prop_freq::Int64=1,
out_dir::String="datafile",
start_at::Int64=0, write_from::Int=0, ext::Symbol=:jld,
max_part=30)
Run simulation for a list of environments.
Arguments
envs
: Array{Environment}, the list of environments.N
: Int64, the number of steps.solver
: Solver, the solver.
Keyword Arguments
filewrite_freq
: Int64, the frequency of writing data files to disk. Default is 10.neigh_update_freq
: Int64, the frequency of updating neighbors. Default is 1.average_prop_freq
: Int64, the frequency of calculating average properties. Default is 1.out_dir
: String, the directory where the data files are saved. Default is "datafile".start_at
: Int64, the starting step. Default is 0.write_from
: Int, the starting index of the data files. Default is 0.ext
: Symbol, the extension of the data files. Default is :jld.max_part
: Int, the maximum number of particles in a neighborhood. Default is 30.
See also
PeriDyn.simulate!
— Methodsimulate!(args...; out_dir="datafile", append_date=true, kwargs...)
Simulate a list of environments. The last argument should be a solver. The args
and kwargs
are passed to run!
function. The out_dir
is the directory where the data files are saved. The append_date
is a boolean value. If true
, the date is appended to the out_dir
path.
See also
PeriDyn.update_acc!
— Methodupdate_acc!(env)
Updates acceleration of all the material points in a simulation environment.
Arguments
env
: the simulation environment.
See also
PeriDyn.update_contact_acc!
— Methodupdate_contact_acc!(env)
Update the forces (acceleration) due to contact.
Arguments
env
: the simulation environment.
See also
PeriDyn.update_mat_acc!
— Methodupdate_mat_acc!(env)
Update the forces (acceleration) due to material deformation.
\[\rho \ddot{u}(x_i, t) = \left[\sum_{j=1}^{N_j} \left\{T[x_i, t]\langle x_j-x_i \rangle - T[x_j, t]\langle x_i-x_j \rangle \right\}V_j + b(x_i, t)\right]\]
where $T$ is the force density, $V_j$ is the volume of the $j$th particle, $\rho$ is the density, $u$ is the displacement, and $b$ is the body force. The summation is over all the particles in the neighborhood given by $N_j$.
Arguments
env
: the simulation environment.
See also
PeriDyn.update_misc!
— Methodupdate_misc!(env)
Update the misc items such as momentum etc.
Arguments
env
: the simulation environment.
PeriDyn.update_neighs!
— Methodfunction update_neighs!(envs)
Updates neighbors of each material point for a list of simulation environments.
Arguments
envs
: Array{Environment}, the list of simulation environments.
See also
update_contact_neighs!
Quasi-static solver
PeriDyn.QSDrag
— TypeQSDrag
Quasi-static solver struct.
Fields
step_size
: Float64, the step size.drag
: Float64, the drag coefficient.max_iter
: Int64, the maximum number of iterations. Default is 1000.x_tol
: Float64, the tolerance for position. Default is 1.0e-3.f_tol
: Float64, the tolerance for force. Default is 1.0e-3.
Example
solver = QSDrag(1.0e-3, 1.0e-3)
See also
QSDrag
apply_solver!
minimize!
PeriDyn.QSDrag
— MethodQSDrag(step_size, drag; max_iter=1000, x_tol=1.0e-3, f_tol=1.0e-3)
Create a QSDrag solver.
Arguments
step_size
: Float64, the step size.drag
: Float64, the drag coefficient.
Keyword Arguments
max_iter
: Int64, the maximum number of iterations. Default is 1000.x_tol
: Float64, the tolerance for position. Default is 1.0e-3.f_tol
: Float64, the tolerance for force. Default is 1.0e-3.
Example
solver = QSDrag(1.0e-3, 1.0e-3)
PeriDyn.apply_solver!
— Methodapply_solver!(env, solver::QSDrag)
Apply the QSDrag solver to the environment.
Arguments
env
: GeneralEnvironment, the environment.solver
: QSDrag, the QSDrag solver.
Example
solver = QSDrag(1.0e-3, 1.0e-3)
apply_solver!(env, solver)
PeriDyn.minimize!
— Methodminimize!(env, solver::QSDrag)
Minimize the environment using the QSDrag solver. The QSDrag solver is a quasi-static solver that uses a drag force to minimize the energy of the system. The drag force is given by -λv(1 + k|v|)
, where λ
is the drag coefficient and k
is a constant. The constant k
is set to λ/10
by default. The drag force is applied to the particles in the system, and the particles are moved by Δy = 0.5FΔt^2 + vΔt
, where F
is the drag force, Δt
is the step size, and v
is the velocity of the particles. The particles are clipped to a maximum displacement of ps/10
, where ps
is the particle size, and the velocity is clipped to a maximum velocity of ps/10/Δt
. The solver is applied to the system until the maximum number of iterations is reached or the maximum displacement and force are below the given tolerances.
Arguments
env
: the simulation environment.solver
: QSDrag, the QSDrag solver.
Example
solver = QSDrag(1.0e-3, 1.0e-3)
minimize!(env, solver)
Dynamic solver
PeriDyn.DSVelocityVerlet
— TypeDSVelocityVerlet
The velocity verlet algorithm.
PeriDyn.apply_solver!
— Methodapply_solver!(env, solver::DSVelocityVerlet)
Apply the velocity verlet algorithm to the simulation environment.
Arguments
env
: the simulation environment.solver
: DSVelocityVerlet, the velocity verlet algorithm.
Example
solver = DSVelocityVerlet()
apply_solver!(env, solver)
PeriDyn.velocity_verlet_step!
— Methodvelocity_verlet_step!(env, solver::DSVelocityVerlet)
Apply one step of the velocity verlet algorithm to the simulation environment.
Steps:
- Update the velocity.
- Update the position.
- Apply the boundary conditions to the position.
- Update the acceleration.
- Update the velocity.
- Apply the boundary conditions to the velocity.
- Check the boundary conditions.
- Check the boundaries.
Arguments
env
: the simulation environment.solver
: DSVelocityVerlet, the velocity verlet algorithm.
Example
solver = DSVelocityVerlet()
velocity_verlet_step!(env, solver)
Input and output functions
PeriDyn.filepath_
— Methodfilepath_(file_prefix::String; append_date=false)
Returns the path to the output folder. If append_date
is true
, the date is appended.
Arguments
file_prefix
: String, the prefix of the output folder.
Keyword Arguments
append_date
: Bool, the boolean value to append date to the output folder. Default is
false
.
PeriDyn.load_from_file
— Methodload_from_file!(filename)
Loads data from a file into the specified environment.
Arguments
filename
: String, the name of the file to load data from.
Returns
env
: the simulation environment.
PeriDyn.print_data_file!
— Methodprint_data_file!(envs::Array{GeneralEnvironment}, file_prefix::String, i::Int64)
Prints data files to disk for a list of simulation environments.
Arguments
envs
: Array{GeneralEnvironment}, the list of simulation environments.file_prefix
: String, the prefix of the output folder.i
: Int64, the step number.
PeriDyn.save_state!
— Methodsave_state!(filename, env)
Save env
to disk.
Arguments
filename
: String, the filename.env
: GeneralEnvironment, the simulation environment.
Keyword Arguments
force
: Bool, the boolean value to force saving to data file format. Default isfalse
.
See also
save_state_ovito_bc!
PeriDyn.save_state_ovito_bc!
— Methodsave_state_ovito_bc!(filename, env)
Save env
to disk for ovito visualization. The boundary conditions are saved as type.
Arguments
filename
: String, the filename.env
: GeneralEnvironment, the simulation environment.
Keyword Arguments
force
: Bool, the boolean value to force saving to data file format. Default isfalse
.
See also
save_state!
PeriDyn.save_to_file
— Methodsave_to_file(env, filename)
Saves data from the specified environment to a file.
Arguments
env
: GeneralEnvironment, the environment to save.filename
: String, the name of the file to save data to.
PeriDyn.write_data
— Methodwrite_data(filename; kwargs...)
Writes the data file.
Arguments
filename
: String, the name of the file to be written.kwargs
: Keyword arguments, additional options for writing the file.
Note: This function supports writing files in the .data
and .jld2
formats.
PeriDyn.write_global_data
— Methodwrite_global_data(filename; kwargs...)
Writes the global data file.
Arguments
filename
: String, the name of the file to be written.kwargs
: Keyword arguments, additional options for writing the file.
Note: This function writes a data file in the .jld2
format.
PeriDyn.jld2array
— Methodjld2array(file, N; start=0, step=100)
Loads data from JLD files into an array.
Arguments
file
: String, the base name of the JLD files.N
: Int, the number of files to load.start
: Int, the index of the first file to load. Default is 0.step
: Int, the step size between files to load. Default is 100.
Returns
- Array, an array containing the data loaded from the JLD files.
Note: This function iterates over a range of files and loads data from each JLD file using the jldread
function.
PeriDyn.jld2ovito
— Methodjld2ovito(file, N; start=0, step=100)
Converts JLD files to Ovito data files.
Arguments
file
: String, the base name of the JLD files.N
: Int, the number of files to convert.start
: Int, the index of the first file to convert. Default is 0.step
: Int, the step size between files to convert. Default is 100.
Note: This function iterates over a range of files and converts each JLD file to an Ovito data file using the write_ovito
function.
PeriDyn.jldread
— Methodjldread(filename::String)
Reads data from a JLD file.
Arguments
filename
: String, the name of the JLD file to read.
Returns
- Dict, a dictionary containing the data read from the JLD file.
Note: This function uses the load
function from the JLD package to read the data from the file, and converts it to a dictionary format.
PeriDyn.write_ovito
— Methodwrite_ovito(filename::String; kwargs...)
Writes data to an Ovito data file.
Arguments
filename
: String, the name of the file to write.kwargs
: Keyword arguments, the data to write to the file.
Note: This function writes data in Ovito data file format, with each column specified by a keyword argument.
PeriDyn.write_ovito_cell_ids
— Methodwrite_ovito_cell_ids(filename::String, y::Matrix, horizon)
Writes cell IDs to an Ovito data file.
Arguments
filename
: String, the name of the file to write.y
: Matrix, the coordinates of the particles.horizon
: The horizon value.
Note: This function writes the cell IDs to an Ovito data file, based on the particle coordinates and horizon value.
Miscellaneous functions
Logging
PeriDyn.log_data
— Methodlog_data(; kwargs...)
Log the data to the log file (and print).
Arguments
kwargs
: Keyword arguments, the data to be logged.
PeriDyn.log_detail
— Methodlog_detail(x)
Log a detailed information message. Prints only if the log level is 4 or higher.
Arguments
x
: String, the information message.
PeriDyn.log_impinfo
— Methodlog_impinfo(x)
Log an important information message. Prints only if the log level is 2 or higher.
Arguments
x
: String, the information message.
PeriDyn.log_info
— Methodlog_info(x)
Log an information message.
Arguments
x
: String, the information message. Prints only if the log level is 3 or higher.
PeriDyn.log_warning
— Methodlog_warning(x)
Log a warning message. Prints only if the log level is 1 or higher.
Arguments
x
: String, the warning message.
PeriDyn.set_loglevel
— Methodset_loglevel(x)
Set the log level.
Arguments
x
: Int64, the log level.
PeriDyn.write_headers
— Methodwrite_headers(envs)
Write the column headers of the log file.
Arguments
envs
: Vector{Environment}, the simulation environments.
Macros
PeriDyn.get_ij
— Methodget_ij(j, i, x)
Get the difference of two vectors depending on the spatial dimensions.
$X_j - X_i$
Arguments
j
: Int64, the index of the first vector.i
: Int64, the index of the second vector.x
: Matrix, the matrix of the vectors.
PeriDyn.get_magnitude
— Methodget_magnitude(a)
Get the magnitude of a vector depending on the spatial dimensions.
Arguments
a
: Vector, the vector.
PeriDyn.map_reduce
— Methodmap_reduce(f, op, iter; init=0.0)
Map and reduce the operations to the vectors depending on the spatial dimensions.
Arguments
f
: Function, the function.op
: Function, the operation.iter
: Any, the iterator.
Keyword Arguments
init
: Any, the initial value.
PeriDyn.refresh
— Methodrefresh()
Refresh the functions which depends on macros.
PeriDyn.set_multi_threading
— Methodset_multi_threading(x)
Set the multi threading.
Arguments
x
: Bool, the multi threading.
PeriDyn.@_ij
— Macro_ij(j, i, x)
Get the difference of two vectors depending on the spatial dimensions.
$X_j - X_i$
Arguments
j
: Int64, the index of the first vector.i
: Int64, the index of the second vector.x
: Matrix, the matrix of the vectors.
PeriDyn.@_magnitude
— Macro_magnitude(a)
Get the magnitude of a vector depending on the spatial dimensions.
Arguments
a
: Vector, the vector.
PeriDyn.@applyops
— Macroapplyops(x)
Apply the operations to the vectors depending on the spatial dimensions.
Arguments
x
: String, the operations.
PeriDyn.@check_nan
— Macrocheck_nan(var, name)
Check if there is nan or inf in the variable.
Arguments
var
: Any, the variable.name
: String, the name of the variable.
PeriDyn.@def
— Macrodef(name, defination)
Define a macro.
Arguments
name
: Symbol, the name of the macro.defination
: Any, the defination of the macro.
PeriDyn.@gatherops
— Macrogatherops(x)
Gather the operations to the vectors depending on the spatial dimensions.
Arguments
x
: String, the operations.
PeriDyn.@map_reduce
— Macromap_reduce(f, op, iter, init)
Map and reduce the operations to the vectors depending on the spatial dimensions.
Arguments
f
: Function, the function.op
: Function, the operation.iter
: Any, the iterator.init
: Any, the initial value.
PeriDyn.@time_code
— Macrotime_code(ex, name)
Time the execution of the expression.
Arguments
ex
: Any, the expression.name
: String, the name of the expression.
Representation
PeriDyn.DPanel
— MethodDPanel(x; title_style="bold blue", title_justify=:center, kwargs...)
Display a variable in a panel.
Arguments
x
: Any, the text to be displayed.
Keyword Arguments
title_style
: String, the style of the title.title_justify
: Symbol, the justification of the title.kwargs...
: Any, keyword arguments to be passed toPanel
.
PeriDyn.SPanel
— MethodSPanel(x; title_style="bold green", title_justify=:center, kwargs...)
Display a variable in a panel.
Arguments
x
: Any, the text to be displayed.
Keyword Arguments
title_style
: String, the style of the title.title_justify
: Symbol, the justification of the title.kwargs...
: Any, keyword arguments to be passed toPanel
.
PeriDyn.array_repr
— Methodarray_repr(item)
Represent an array in nice format.
Arguments
item
: Any, the array.
PeriDyn.variable_color
— Methodvariable_color(x; kalar="#773399")
Color a variable.
Arguments
x
: Any, the variable to be colored.
Keyword Arguments
kalar
: String, the color to use.
PeriDyn.variable_txt
— Methodvariable_txt(item)
Represent a variable in nice format.
Arguments
item
: Any, the variable.
Utilities
PeriDyn.SymMat
— TypeSymMat
Create an symmetrical NxN matrix from a vector of length N(N+1)/2.
Arguments
v::Array
: Vector of length N(N+1)/2.N::Int
: Number of layers.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.SymMat
— MethodSymMat
Create an symmetrical NxN matrix from a vector of length N(N+1)/2.
Arguments
v::Array
: Vector of length N(N+1)/2.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.make_NN
— Methodmake_NN(layers::T; act = Flux.relu) where {T}
Create an neural network with given layers and activation function.
Arguments
layers::Tuple
: Tuple of layers.
Keyword Arguments
act::Function
: Activation function.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.make_NN
— Methodmake_NN(layers::T, N; act=Flux.relu) where {T}
Create an symmetrical NxN matrix from a vector of length N(N+1)/2.
Arguments
layers::Tuple
: Tuple of layers.N::Int
: Number of layers.
Keyword Arguments
act::Function
: Activation function.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.make_matrix
— Methodmake_matrix(S::Array{T,1})
Create an symmetrical NxN matrix from a vector of length N(N+1)/2.
Arguments
S::Array
: Vector of length N(N+1)/2.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.make_matrix_gm
— Methodmake_matrix(S::Array{T,1})
Create an symmetrical NxN matrix from a vector of length N(N+1)/2.
Arguments
S::Array
: Vector of length N(N+1)/2.
Returns
M::Matrix
: Symmetrical NxN matrix.
PeriDyn.make_vector
— Methodmake_vector(M::Array{T,2})
Create an array of upper triangle of an symmetrical NxN matrix.
Arguments
M::Matrix
: Symmetrical NxN matrix.
Returns
S::Array
: Vector of length N(N+1)/2.
PeriDyn.perturb
— Methodperturb(x::AbstractArray; e=1.0e-6)
Perturb an array with Gaussian noise.
Arguments
x::AbstractArray
: Array to be perturbed.
Keyword Arguments
e::Float64
: Standard deviation of the Gaussian noise.
Returns
x::AbstractArray
: Perturbed array.
PeriDyn.perturb
— Methodperturb(i::Int64, j::Int64; e=1.0e-6)
Create a Matrix of size i x j with Gaussian noise.
Arguments
i::Int64
: Number of rows.j::Int64
: Number of columns.
Keyword Arguments
e::Float64
: Standard deviation of the Gaussian noise.
Returns
x::AbstractArray
: Matrix of size i x j with Gaussian noise.
PeriDyn.repack!
— Methodrepack!(d::Dict, keys_, vals)
Repack a dictionary from its components inplace.
Arguments
d::Dict
: Dictionary to be repacked.keys_
: Keys of the dictionary.vals
: Values of the dictionary.
Returns
d::Dict
: Dictionary containing the components.
PeriDyn.repack
— Methodrepack(args...; keys_ = (:x, :v, :y, :volume, :type))
Repack a dictionary from its components.
Arguments
args...
: Components to be packed.
Keyword Arguments
keys_ = (:x, :v, :y, :volume, :type)
: Keys of the dictionary.
Returns
d::Dict
: Dictionary containing the components.
PeriDyn.unpack
— Methodunpack(d::Dict)
Unpack a dictionary into its components.
Arguments
d::Dict
: Dictionary to be unpacked.
Returns
x::Array{Float64, 1}
: x coordinates.v::Array{Float64, 1}
: v coordinates.y::Array{Float64, 1}
: y coordinates.volume::Array{Float64, 1}
: Volume of the mesh.type::Array{Int64, 1}
: Type of the mesh.