Automatic documentation

Table of contents on this page.

Peridynamics functions

PeriDyn.__IMethod
__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)
source
PeriDyn.__OMethod
__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
source
PeriDyn.__normMethod
__norm(x::AbstractArray)

Returns magnitude of x ($\lVert x \rVert$).

Arguments

  • x::AbstractArray: input array

Returns

  • y: magnitude of x
source
PeriDyn.__unitMethod
__unit(x::AbstractArray)

Returns unit vector of x.

Arguments

  • x::AbstractArray: input array

Returns

  • y: unit vector of x
source
PeriDyn.cal_damageMethod
cal_damage(env)

Calculates damage for each bond for a given environment.

Arguments

  • env: environment
source
PeriDyn.cal_family!Method
cal_family!(family, x, horizon)

Calculates family for given positions and horizon.

Arguments

  • family: family
  • x: positions
  • horizon: horizon
source
PeriDyn.cal_familyMethod
cal_family(x, horizon, max_neigh)

Calculates family for given positions and horizon.

Arguments

  • x: positions
  • horizon: horizon
  • max_neigh: maximum number of neighbors
source
PeriDyn.cell_numberMethod
cell_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 direction
  • j: cell index in y direction
  • k: cell index in z direction
  • N: number of cells in each direction
source
PeriDyn.dilatation!Method
dilatation!(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: dilatation
  • y: deformed positions
  • x: initial positions
  • intact: intact bonds
  • family: family of particles
  • volume: particle volume
  • m: mass
  • particle_size: particle size
  • horizon: horizon
  • device::Type{Val{:cuda}}: device type

Returns

  • nothing: theta is updated in-place
source
PeriDyn.dilatation!Method
dilatation!(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: dilatation
  • y: deformed positions
  • x: initial positions
  • intact: intact bonds
  • family: family of particles
  • volume: particle volume
  • m: weighted volume
  • particle_size: particle size
  • horizon: horizon
  • device::Type{Val{:cpu}}: device type

Returns

  • nothing: theta is updated in-place
source
PeriDyn.dilatation!Method
dilatation!(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 positions
  • mat: material model
  • device: device type

Returns

  • nothing: In-place operation
source
PeriDyn.dilatationMethod
dilatation(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 positions
  • x: initial positions
  • intact: intact bonds
  • family: family of particles
  • volume: particle volume
  • m: mass
  • particle_size: particle size
  • horizon: horizon
  • device::Type{Val{:cpu}}: device type

Returns

  • theta: dilatation
source
PeriDyn.dilatationMethod
dilatation(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 positions
  • x: initial positions
  • intact: intact bonds
  • family: family of particles
  • volume: particle volume
  • m: mass
  • particle_size: particle size
  • horizon: horizon
  • device::Type{Val{:cuda}}: device type

Returns

  • theta: dilatation
source
PeriDyn.get_cellsMethod
get_cells(x, horizon; max_part=30)

Calculates cells for given positions and horizon.

Arguments

  • x: positions
  • horizon: horizon

Keyword Arguments

  • max_part=30: number of partitions in each direction
source
PeriDyn.horizon_correctionMethod
horizon_correction(dr, ps, hr)

It gives horizon correction factor (It will give 1 as of now).

Arguments

  • dr: $X_{ij}$
  • ps: particle size
  • hr: horizon

Returns

  • s: horizon correction factor
source
PeriDyn.influence_functionMethod
influence_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
source
PeriDyn.loop_over_neighs!Method
loop_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 model
  • i: particle index
  • j: neighbor index
  • k: family index
  • type1: type of particle i
  • type2: type of particle j
  • xij: x_j - x_i
  • yij: y_j - y_i
  • extension: y_j - y_i
  • s: s_ij bond stretch
  • wij: influence function w_ij
  • wji: influence function w_ji
  • items: output of get_items(mat) function

The function fn_reduce should have the following signature:

fn_reduce(mat, i, items)

where,

  • mat: material model
  • i: particle index
  • items: output of get_items(mat) function

Arguments

  • y_: deformed positions
  • mat: material model
  • fn_map: function to be applied to each neighbor
  • device::Type{Val{:cuda}}: device type

Keyword Arguments

  • fn_reduce=(args...)->(): function to reduce the results

Returns

  • nothing: y_ is updated in-place
source
PeriDyn.neigh_cellsMethod
neigh_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 direction
  • j: cell index in y direction
  • k: cell index in z direction
  • N: number of cells in each direction
source
PeriDyn.weighted_volumeMethod
weighted_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 positions
  • volume: particle volume
  • particle_size: particle size
  • family: family of particles
  • horizon: horizon

Returns

  • m: weighted volume
source

Simulation environment

PeriDyn.GeneralEnvironmentType
GeneralEnvironment <: SimulationEnvironment

Type for holding parameters for a simulation.

Fields

  • id::Int64: ID of the environment
  • type::AbstractArray{Int64,1}: Type of each particle
  • bid::AbstractArray{Int64,1}: Block ID of each particle
  • ghost_particles::AbstractArray{Int64,1}: Ghost particles
  • state::Int64: State of the environment
  • time_step::Int64: Time step of the environment
  • dt::T where T<:QF: Time step size
  • y::AbstractArray{T,2} where T<:QF: Position of each particle
  • v::AbstractArray{T,2} where T<:QF: Velocity of each particle
  • f::AbstractArray{T,2} where T<:QF: Force of each particle
  • p::AbstractArray{T,2} where T<:QF: Momentum of each particle
  • volume::AbstractArray{T,1} where T<:QF: Volume of each particle
  • mass::AbstractArray{T,1} where T<:QF: Mass of each particle
  • density::AbstractArray{T,1} where T<:QF: Density of each particle
  • intact0::AbstractArray{Int64, 1}: Intact particles information
  • neighs::AbstractArray{Int64,2}: Neighbors of each particle
  • boundary_conditions::AbstractArray{T, 1} where T: Boundary conditions
  • short_range_contact::AbstractArray{T, 1} where T: Short range contact
  • material_blocks::AbstractArray{T, 1} where T: Material blocks
  • boundaries::Tuple: Boundaries
  • Collect!::Function: Function for collecting data
  • Params::Dict{Symbol, Any}: Parameters
  • Out::Dict{Symbol, Any}: Output
  • cprint::Function: Function for printing
source
PeriDyn.EnvMethod
Env(id::Int64, materials,short_range_contact, boundary_conds, dt; state=2, bskin=0.5, units=false)

Constructor for GeneralEnvironment type.

Arguments

  • id::Int64: Environment ID
  • materials::AbstractArray{Material, 1}: Materials
  • short_range_contact::AbstractArray{ContactModel, 1}: Short range contact models
  • boundary_conds::AbstractArray{BoundaryCondition, 1}: Boundary conditions
  • dt::T where T<:QF: Time step
  • state::Int64: State of the environment
  • bskin::T where T<:QF: Boundary skin
  • units::Bool: Units

Returns

  • env::GeneralEnvironment: General environment
source
PeriDyn.printdataMethod
printdata(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 ID
  • time: Simulation time
  • px, py, pz: Total momentum in x, y, z directions
  • Fx, Fy, Fz: Total force in x, y, z directions
  • damage: Total damage
  • and the data printed by cprint function of the environment

Arguments

  • env: Environment
source

Material models

General material models and functions

PeriDyn.GeneralMaterialType
GeneralMaterial

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).
source
PeriDyn.GeneralMaterialMethod
GeneralMaterial(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.
source
PeriDyn.PeridynamicsMaterialMethod
PeridynamicsMaterial(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.
source
PeriDyn.PeridynamicsMaterialMethod
PeridynamicsMaterial(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.
source
PeriDyn.PeridynamicsMaterialMethod
PeridynamicsMaterial(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.
source
PeriDyn.force_density!Method
force_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 of f.
source
PeriDyn.force_density_T!Method
force_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 of f.
source
PeriDyn.force_density_T!Method
force_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}}: Device
  • particles::Union{Nothing, AbstractVector{Int64}}: Particles to calculate the force.

Returns

  • nothing: In-place modification of f.
source
PeriDyn.force_density_T!Method
force_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}}: Device
  • particles::Union{Nothing, AbstractVector{Int64}}: Particles to calculate the force.

Returns

  • nothing: In-place modification of f.
source
PeriDyn.force_density_t_ijMethod
force_density_t_ij(mat, args...; kwargs...)

Force density function for peridynamics material. This is default implementation of force_density_t_ij function

source
PeriDyn.initMethod
init(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.
source
PeriDyn.@newmaterialMacro
newmaterial(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.
source

Specific material models and functions

PeriDyn.SkipSpecificType
SkipSpecific <: SpecificMaterial

Skip specific material type.

Fields

  • density::AbstractArray{<:QF, 1}: Density.
  • critical_stretch::AbstractArray{Float64, 2}: Critical stretch.

Returns

  • SkipSpecific: Skip specific material.
source
PeriDyn.BondBasedMaterialType
BondBasedMaterial <: 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.
source
PeriDyn.BondBasedSpecificType
BondBasedSpecific

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...)
source
PeriDyn.BondBasedSpecificMethod
BondBasedSpecific(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.
source
PeriDyn.BondBasedSpecificMethod
BondBasedSpecific(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.
source
PeriDyn.BondBasedSpecificMethod
BondBasedSpecific(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.
source
PeriDyn.bond_forceMethod
bond_force(s, bond_stiffness)

Bond force function.

Arguments

  • s::Real: Bond stretch.
  • bond_stiffness::Real: Bond stiffness.

Returns

  • Real: Bond force.
source
PeriDyn.force_density_t_ijMethod
force_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$ in family and intact.
  • 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.
source
PeriDyn.get_itemsMethod
get_items(mat::BondBasedMaterial)

Returns items for force calculation.

Arguments

  • mat::BondBasedMaterial: Bond based material type.

Returns

  • bond_stiffness::AbstractArray: Bond stiffness.
source
PeriDyn.initMethod
init(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.
source
PeriDyn.out_of_loop!Method
out_of_loop!(y_, mat::BondBasedMaterial, device)

It is called once before the main loop for force calculation. It does nothing for BondBasedMaterial.

source
PeriDyn.OrdinaryStateBasedMaterialType
OrdinaryStateBasedMaterial <: 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.
source
PeriDyn.OrdinaryStateBasedSpecificType
OrdinaryStateBasedSpecific <: 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.
source
PeriDyn.OrdinaryStateBasedSpecificMethod
OrdinaryStateBasedSpecific(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.
source
PeriDyn.PeridynamicsMaterialMethod
PeridynamicsMaterial(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.
source
PeriDyn.force_density_t_ijMethod
force_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$ in family and intact.
  • 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.
source
PeriDyn.get_itemsMethod
get_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.

source
PeriDyn.initMethod
init(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.
source
PeriDyn.out_of_loop!Method
out_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
source
PeriDyn.ElastoPlasticSolidMaterialType
ElastoPlasticSolidMaterial <: PeridynamicsMaterial

Elastic-plastic solid material.

Fields

  • name::String: Name of the material
  • type::Array: Type of the material
  • bid::Int: Block id
  • general::GeneralMaterial: General material
  • specific::ElastoPlasticSolidSpecific: Specific material
source
PeriDyn.ElastoPlasticSolidSpecificType
ElasticPlasticSolidMaterial <: PeridynamicsMaterial

Elastic-plastic solid material.

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
  • σγ::AbstractArray{<:QF,2}: Yield stress
  • ψ::AbstractArray{<:QF,2}: Yield stress value
  • edp::AbstractArray{<:QF,2}: Plastic deviatoric strain
  • td_norm2::AbstractArray{<:QF,1}: td_norm2
  • theta::AbstractArray{<:QF,1}: Dilatation
  • criteria::PlasticFailureCriteria: Plastic failure criteria
source
PeriDyn.ElastoPlasticSolidSpecificMethod
ElastoPlasticSolidSpecific(
    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 modulus
  • shear_modulus::AbstractArray{<:QF,1}: Shear modulus
  • critical_stretch::Array{Float64,1}: Critical stretch
  • density::AbstractArray{<:QF,1}: Density
  • σγ::AbstractArray{<:QF,1}: Yield stress
  • criteria::PlasticFailureCriteria: Plastic failure criteria

Returns

  • spc::ElastoPlasticSolidSpecific: Specific material
source
PeriDyn.PeridynamicsMaterialMethod
PeridynamicsMaterial(name, type, bid, gen, spc::ElastoPlasticSolidSpecific)

Construct an ElastoPlasticSolidMaterial.

Arguments

  • name::String: Name of the material
  • type::Array: Type of the material
  • bid::Int: Block id
  • gen::GeneralMaterial: General material
  • spc::ElastoPlasticSolidSpecific: Specific material

Returns

  • mat::ElastoPlasticSolidMaterial: Material
source
PeriDyn.effective_stress2Method
effective_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.
source
PeriDyn.force_density_t_ijMethod
force_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$ in family and intact.
  • 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.
source
PeriDyn.get_itemsMethod
get_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)

source
PeriDyn.initMethod
init(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 material
  • gen::GeneralMaterial: General material

Returns

  • spc::ElastoPlasticSolidSpecific: Initialized specific material
source
PeriDyn.out_of_loop!Method
out_of_loop!(y_, mat::ElastoPlasticSolidMaterial, device)

Calculate the dilatation and the td_norm2 out of the loop.

Arguments

  • y_::AbstractArray{<:QF,2}: Displacement
  • mat::ElastoPlasticSolidMaterial: Material
  • device: Device
source
PeriDyn.td_norm2_ij!Method
td_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$ in family and intact.
  • 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.
source

Contact models and functions

Base.showMethod

Print the details of a ContactModel object.

Arguments

  • io::IO: The output IO stream.
  • i::Union{ContactModel11, ContactModel12}: The ContactModel object to print.
source
PeriDyn.ContactModel11_gcalMethod
ContactModel11_gcal(mat1, distanceX, max_neighs)

Calculate the parameters for ContactModel11 based on the input material, distance, and maximum neighbors.

Arguments:

  • mat1: The material (type ContactModel11).
  • distanceX: The distance factor.
  • max_neighs: The maximum number of neighbors.

Returns a tuple containing the calculated parameters for ContactModel11.

source
PeriDyn.ContactModel12_gcalMethod
ContactModel12_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 (type ContactModel11).
  • mat2: The second material (type ContactModel11).
  • distanceX: The distance factor.
  • max_neighs: The maximum number of neighbors.

Returns a tuple containing the calculated parameters for ContactModel12.

source
PeriDyn._cudaconvertMethod
_cudaconvert(x::T) where T <: Union{ContactModel11,ContactModel12}

Converts a single ContactModel11 or ContactModel12 object to a CUDA-compatible type.

source
PeriDyn._cudaconvertMethod
_cudaconvert(x::Vector{T}) where T <: Union{ContactModel11,ContactModel12}

Converts a vector of ContactModel11 or ContactModel12 objects to CUDA-compatible types.

source
PeriDyn.collision_boxMethod
collision_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 overlap
  • box_max: Maximum position limits for overlap
  • ifoverlap: Boolean indicating if there is an overlap
source
PeriDyn.short_range_contact!Method
short_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 points
  • f: Acceleration of material points
  • type: Type of material points
  • bid: BID values
  • vol: Volume values
  • RM::ContactModel11: Repulsion model
  • device::Type{Val{:cpu}}: Device type for CPU acceleration

Output

  • No return value. The function updates f in place.
source
PeriDyn.short_range_contact!Method
short_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 points
  • f: Acceleration of material points
  • type: Type of material points
  • bid: BID values
  • vol: Volume values
  • RM::ContactModel11: Repulsion model
  • device::Type{Val{:cuda}}: Device type for CUDA acceleration

Output

  • No return value. The function updates f in place.
source
PeriDyn.short_range_contact!Method
short_range_contact!(y, f, type, bid, vol, RM::ContactModel11)

Updates (inplace) the contact acceleration of material points.

Arguments

  • y: Positions of material points
  • f: Acceleration of material points
  • type: Type of material points
  • bid: BID values
  • vol: Volume values
  • RM::ContactModel11: Repulsion model

Output

  • No return value. The function updates f in place.
source
PeriDyn.short_range_contact!Method
short_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)).
source
PeriDyn.short_range_contact!Method
short_range_contact!(y, f, type, bid, vol, RM::ContactModel12)

Updates (inplace) the contact acceleration of material points.

Arguments

  • y: Positions of material points
  • f: Acceleration of material points
  • type: Type of material points
  • bid: BID values
  • vol: Volume values
  • RM::ContactModel12: Repulsion model

Output

  • No return value. The function updates f in place.
source
PeriDyn.update_contact_neighs!Method
update_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 indices
  • x: Positions of material points
  • search_distance: Maximum search distance for neighbors
  • equi_dist: Equilibrium distance for contact
  • family: Array indicating the family relationship between material points
  • intact: Array indicating if the family relationship is intact
  • max_part: Maximum number of particles in a cell

Output

  • No return value. The function updates the neighbors array in place.
source
PeriDyn.update_contact_neighs!Method
update_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 points
  • type: Type of material points
  • RM::ContactModel11: Repulsion model
  • device::Symbol: Device type for acceleration
  • kwargs...: Additional keyword arguments

Output

  • No return value. The function updates RM.neighs in place.
source
PeriDyn.update_contact_neighs!Method
update_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 points
  • type: Type of material points
  • RM::ContactModel11: Repulsion model for 1-1 interaction
  • device::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.
source
PeriDyn.update_contact_neighs!Method
update_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 points
  • type: Type of material points
  • RM::ContactModel11: Repulsion model for 1-1 interaction
  • device::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.
source
PeriDyn.update_contact_neighs!Method
update_contact_neighs!(y, type, RM::ContactModel11; kwargs...)

Update neighbor list for contact force calculation (1-1 interaction).

Arguments

  • y: Positions of material points
  • type: Type of material points
  • RM::ContactModel11: Repulsion model
  • kwargs...: Additional keyword arguments

Output

  • No return value. The function updates RM.neighs in place.
source
PeriDyn.update_contact_neighs!Method
update_contact_neighs!(y, type, RM::ContactModel12; max_part=nothing)

Update neighbor list for contact force calculation (1-2 interaction).

Arguments

  • y: Positions of material points
  • type: Type of material points
  • RM::ContactModel12: Repulsion model
  • max_part=nothing: Maximum number of particles (optional)

Output

  • No return value. The function updates RM.neighs in place.
source
PeriDyn.NonLinearSpringContactModel11Type
NonLinearSpringContactModel11(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.

source
PeriDyn.NonLinearSpringContactModel12Type
NonLinearSpringContactModel12(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.

source
PeriDyn.NonLinearSpringContactModelMethod
NonLinearSpringContactModel(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.

source
PeriDyn.NonLinearSpringContactModelMethod
NonLinearSpringContactModel(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.

source
PeriDyn.get_contact_force_fnMethod
get_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.

source
PeriDyn.LinearSpringContactModelMethod
LinearSpringContactModel(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.
source
PeriDyn.LinearSpringContactModelMethod
LinearSpringContactModel(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.
source
PeriDyn.ShortRangeContactModelMethod
ShortRangeContactModel(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 is 1.0.
  • kwargs...: Keyword arguments for the LinearSpringContactModel constructor.
source
PeriDyn.ShortRangeContactModelMethod
ShortRangeContactModel(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 is 1.0.
  • kwargs...: Keyword arguments for the LinearSpringContactModel constructor.
source

Boundary conditions

PeriDyn.apply_bc!Method
apply_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).
source
PeriDyn.apply_bc!Method
apply_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.
source
PeriDyn.apply_bc!Method
apply_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.
source
PeriDyn.apply_bc_at0!Method
apply_bc_at0!(env, BC::BoundaryCondition)

Apply the boundary condition BC to the given environment env at the start of the simulation (t=0).

source
PeriDyn.check!Method
check!(env, BC::BoundaryCondition)

Check if the boundary condition BC has changed and updates the BC. Used for dynamic boundary conditions.

source
PeriDyn.@general_bc_pMacro
general_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.
source
PeriDyn.FixBCType
FixBC

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.
source
PeriDyn.FixBCMethod
FixBC(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.

source
PeriDyn.apply_bc_at0!Method
apply_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.
source
PeriDyn.ToFroBCType
ToFroBC

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.
source
PeriDyn.ToFroBCMethod
ToFroBC(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.

source
PeriDyn.MoveBCMethod
MoveBC(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 to ToFroBC.

Returns

A MoveBC object representing the boundary condition. All MoveBC objects are ToFroBC objects with frequency Inf.

source
PeriDyn.apply_bc_at0!Method
apply_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.
source
PeriDyn.check!Method
check!(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.
source
PeriDyn.DeltaScaleBCType
DeltaScaleBC

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.
source
PeriDyn.DeltaScaleBCMethod
DeltaScaleBC(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.

source
PeriDyn.apply_bc_at0!Method
apply_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.
source
PeriDyn.check!Method
check!(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.
source
PeriDyn.ScaleFixWaitBCType
ScaleFixWaitBC

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.
source
PeriDyn.ScaleFixWaitBCMethod
ScaleFixWaitBC(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.
source
PeriDyn.apply_bc_at0!Method
apply_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.
source
PeriDyn.check!Method
check!(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.
source
PeriDyn.ContainerBCType
ContainerBC

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.
source
PeriDyn.ContainerBCMethod
ContainerBC(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.

source
PeriDyn.apply_bc_at0!Method
apply_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.
source

Solvers and integrators

General solver functions

PeriDyn.apply_bc_at0Method
apply_bc_at0!(env, start_at)

Apply boundary conditions at t = t0.

Arguments

  • env: the simulation environment.
  • start_at: Int64, the starting step.
source
PeriDyn.check_boundaries!Method
check_boundaries!(env)

Check if any material point is outside the defined boundaries. Do nothing.

Arguments

  • env: the simulation environment.
source
PeriDyn.run!Method
run!(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

simulate!

source
PeriDyn.simulate!Method
simulate!(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

run!

source
PeriDyn.update_mat_acc!Method
update_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

force_density!

source
PeriDyn.update_misc!Method
update_misc!(env)

Update the misc items such as momentum etc.

Arguments

  • env: the simulation environment.
source
PeriDyn.update_neighs!Method
function 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!
source

Quasi-static solver

PeriDyn.QSDragType
QSDrag

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!
source
PeriDyn.QSDragMethod
QSDrag(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)
source
PeriDyn.apply_solver!Method
apply_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)
source
PeriDyn.minimize!Method
minimize!(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)
source

Dynamic solver

PeriDyn.apply_solver!Method
apply_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)
source
PeriDyn.velocity_verlet_step!Method
velocity_verlet_step!(env, solver::DSVelocityVerlet)

Apply one step of the velocity verlet algorithm to the simulation environment.

Steps:

  1. Update the velocity.
  2. Update the position.
  3. Apply the boundary conditions to the position.
  4. Update the acceleration.
  5. Update the velocity.
  6. Apply the boundary conditions to the velocity.
  7. Check the boundary conditions.
  8. Check the boundaries.

Arguments

  • env: the simulation environment.
  • solver: DSVelocityVerlet, the velocity verlet algorithm.

Example

solver = DSVelocityVerlet()
velocity_verlet_step!(env, solver)
source

Input and output functions

PeriDyn.filepath_Method
filepath_(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.

source
PeriDyn.load_from_fileMethod
load_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.
source
PeriDyn.print_data_file!Method
print_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.
source
PeriDyn.save_state!Method
save_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 is false.

See also

  • save_state_ovito_bc!
source
PeriDyn.save_state_ovito_bc!Method
save_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 is false.

See also

  • save_state!
source
PeriDyn.save_to_fileMethod
save_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.
source
PeriDyn.write_dataMethod
write_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.

source
PeriDyn.write_global_dataMethod
write_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.

source
PeriDyn.jld2arrayMethod
jld2array(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.

source
PeriDyn.jld2ovitoMethod
jld2ovito(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.

source
PeriDyn.jldreadMethod
jldread(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.

source
PeriDyn.write_ovitoMethod
write_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.

source
PeriDyn.write_ovito_cell_idsMethod
write_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.

source

Miscellaneous functions

Logging

PeriDyn.log_dataMethod
log_data(; kwargs...)

Log the data to the log file (and print).

Arguments

  • kwargs: Keyword arguments, the data to be logged.
source
PeriDyn.log_detailMethod
log_detail(x)

Log a detailed information message. Prints only if the log level is 4 or higher.

Arguments

  • x: String, the information message.
source
PeriDyn.log_impinfoMethod
log_impinfo(x)

Log an important information message. Prints only if the log level is 2 or higher.

Arguments

  • x: String, the information message.
source
PeriDyn.log_infoMethod
log_info(x)

Log an information message.

Arguments

  • x: String, the information message. Prints only if the log level is 3 or higher.
source
PeriDyn.log_warningMethod
log_warning(x)

Log a warning message. Prints only if the log level is 1 or higher.

Arguments

  • x: String, the warning message.
source
PeriDyn.write_headersMethod
write_headers(envs)

Write the column headers of the log file.

Arguments

  • envs: Vector{Environment}, the simulation environments.
source

Macros

PeriDyn.get_ijMethod
get_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.
source
PeriDyn.get_magnitudeMethod
get_magnitude(a)

Get the magnitude of a vector depending on the spatial dimensions.

Arguments

  • a: Vector, the vector.
source
PeriDyn.map_reduceMethod
map_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.
source
PeriDyn.@_ijMacro
_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.
source
PeriDyn.@_magnitudeMacro
_magnitude(a)

Get the magnitude of a vector depending on the spatial dimensions.

Arguments

  • a: Vector, the vector.
source
PeriDyn.@applyopsMacro
applyops(x)

Apply the operations to the vectors depending on the spatial dimensions.

Arguments

  • x: String, the operations.
source
PeriDyn.@check_nanMacro
check_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.
source
PeriDyn.@defMacro
def(name, defination)

Define a macro.

Arguments

  • name: Symbol, the name of the macro.
  • defination: Any, the defination of the macro.
source
PeriDyn.@gatheropsMacro
gatherops(x)

Gather the operations to the vectors depending on the spatial dimensions.

Arguments

  • x: String, the operations.
source
PeriDyn.@map_reduceMacro
map_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.
source
PeriDyn.@time_codeMacro
time_code(ex, name)

Time the execution of the expression.

Arguments

  • ex: Any, the expression.
  • name: String, the name of the expression.
source

Representation

PeriDyn.DPanelMethod
DPanel(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 to Panel.
source
PeriDyn.SPanelMethod
SPanel(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 to Panel.
source
PeriDyn.variable_colorMethod
variable_color(x; kalar="#773399")

Color a variable.

Arguments

  • x: Any, the variable to be colored.

Keyword Arguments

  • kalar: String, the color to use.
source

Utilities

PeriDyn.SymMatType
SymMat

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.
source
PeriDyn.SymMatMethod
SymMat

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.
source
PeriDyn.make_NNMethod
make_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.
source
PeriDyn.make_NNMethod
make_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.
source
PeriDyn.make_matrixMethod
make_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.
source
PeriDyn.make_matrix_gmMethod
make_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.
source
PeriDyn.make_vectorMethod
make_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.
source
PeriDyn.perturbMethod
perturb(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.
source
PeriDyn.perturbMethod
perturb(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.
source
PeriDyn.repack!Method
repack!(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.
source
PeriDyn.repackMethod
repack(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.
source
PeriDyn.unpackMethod
unpack(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.
source