Private API

Peridynamics.failure_permit!Function
failure_permit!(body, set_name, fail_permit)

Set the failure permission for points of the set set_name of a body.

Arguments

  • body::AbstractBody: Body where the failure permission will be set.
  • set_name::Symbol: The name of a point set of this body.
  • fail_permit::Bool: If true, failure is allowed, and if false then no bonds of this point are allowed to break during the simulation.
Overwriting failure permission with `material!` and `failure_permit!`

The function material! calls failure_permit!, so if it is used afterwards, previously set failure permissions might be overwritten!

Throws

  • Error if the body does not contain a set with set_name.
source
Peridynamics.get_frac_paramsFunction
get_frac_params(::AbstractDamageModel, p::Dict{Symbol,Any}, δ::Float64, K::Float64)

Read or calculate the necessary fracture parameters for the specified damage model from the dictionary created with material!. This function has to be defined when creating a new damage model. Otherwise, a default method returns a empty named tuple (; ).

source
Peridynamics.set_failure_permissions!Function
set_failure_permissions!(body, set_name, params)

Grant or prohibit failure permission depending on the submitted fracture parameters by calling failure_permit!.

If fracture parameters are found, failure is allowed. If no fracture parameters are found, failure is not allowed.

source
Peridynamics.has_fractureFunction
has_fracture(mat, params)

Return true if at least one fracture parameter is set !=0 in params and the system therefore is supposed to have failure allowed or return false if not.

source
Peridynamics.pre_submission_checkFunction
pre_submission_check(body::Body; body_in_multibody_setup::Bool=false)
pre_submission_check(ms::AbstractMultibodySetup)

Check if necessary material parameters and conditions are defined when defining a Job.

source
Peridynamics.get_paramsetupFunction
get_paramsetup(body::AbstractBody, ::AbstractChunkHandler, ::SingleParamChunk)
get_paramsetup(body::AbstractBody, ch::AbstractChunkHandler, ::MultiParamChunk)

Return the parameters of a BodyChunk if only one parameter set is defined for the corresponding Body or the parameter handler if multiple parameter sets are defined.

source
Peridynamics.get_paramsFunction
get_params(paramhandler::ParameterHandler, point_id::Int)
get_params(params::AbstractPointParameters, ::Int)
get_params(chunk::BodyChunk, point_id::Int)

Return parameters of a specific point with index point_id of a Body with parameters params or parameter handler paramhandler or of the body chunk.

source
Peridynamics.find_pointsFunction
find_points(f, position)

Find all points whose positions meet function f.

The function f accepts only one positional argument and will be used in a findall call. Depending on the argument name, a different input will be processed. See point_set!.

source
Peridynamics.apply_precrack!Function
apply_precrack!(chunk, body, crack)

Apply the predefined crack of the body for the considered chunk by breaking the concerned bonds.

source
Peridynamics.point_sets_intersectFunction
point_sets_intersect(point_sets, key_a, key_b)

Return true if point sets key_a and key_b out of the dictionary point_sets have at least one common point, return false if they do not.

source
Peridynamics.velocity_databc!Function
velocity_databc!(body, data, set_name, dims)

Specifies velocity boundary conditions for points of the set set_name in body. The value of the boundary condition is assigned by reading the corresponding positions in the matrix data. Multiple dimensions can be handled at once.

Compatibility feature with other packages

This method / feature is used for compatibility with other packages developing with Peridynamics.jl. It is likely to change in the future, since the functionality of updating the values of the matrix during the simulation is not yet implemented. Consequently, at this stage, it is only available as a private API to facilitate future modifications and ensure easier implementation of changes.

Arguments

  • body::AbstractBody: Body the condition is specified on.
  • data::Matrix: A matrix of size length(dims) x n_points that contains the values of the boundary condition for each point in the body. But only the conditions of points contained in the set set_name are applied during the simulation! It should be noted, that the value of the data matrix is constant and currently cannot be updated during the simulation. The data matrix is not checked for NaN values, since this is handled in the apply_bc! function. If it contains NaN values, then these values are ignored.
  • set_name::Symbol: The name of a point set of this body. The condition applies only to the points in this set, even if the data matrix contains values for all points in the body.
  • dims::Vector{Union{Integer,Symbol}}: Vector containing the directions of the condition that should be applied, either specified as Symbol or integer.
    • x-direction: :x or 1
    • y-direction: :y or 2
    • z-direction: :z or 3
    It should not contain more than 3 elements, and the elements should be unique. The order of the elements does not matter, however it must match the values in the data matrix. So if the first column of the data matrix contains the values for the x-direction, then the first element of dims should be 1 or :x, and so on.

Throws

  • Errors if the body does not contain a set with set_name.
  • Errors if the directions are not correctly specified.
  • Errors if the dimensions of the data matrix are incorrect.
source
Peridynamics.forcedensity_databc!Function
forcedensity_databc!(body, data, set_name, dims)

Specifies forcedensity boundary conditions for points of the set set_name in body. The value of the boundary condition is assigned by reading the corresponding positions in the matrix data. Multiple dimensions can be handled at once.

Compatibility feature with other packages

This method / feature is used for compatibility with other packages developing with Peridynamics.jl. It is likely to change in the future, since the functionality of updating the values of the matrix during the simulation is not yet implemented. Consequently, at this stage, it is only available as a private API to facilitate future modifications and ensure easier implementation of changes.

Arguments

  • body::AbstractBody: Body the condition is specified on.
  • data::Matrix: A matrix of size length(dims) x n_points that contains the values of the boundary condition for each point in the body. But only the conditions of points contained in the set set_name are applied during the simulation! It should be noted, that the value of the data matrix is constant and currently cannot be updated during the simulation. The data matrix is not checked for NaN values, since this is handled in the apply_bc! function. If it contains NaN values, then these values are ignored.
  • set_name::Symbol: The name of a point set of this body. The condition applies only to the points in this set, even if the data matrix contains values for all points in the body.
  • dims::Vector{Union{Integer,Symbol}}: Vector containing the directions of the condition that should be applied, either specified as Symbol or integer.
    • x-direction: :x or 1
    • y-direction: :y or 2
    • z-direction: :z or 3
    It should not contain more than 3 elements, and the elements should be unique. The order of the elements does not matter, however it must match the values in the data matrix. So if the first column of the data matrix contains the values for the x-direction, then the first element of dims should be 1 or :x, and so on.

Throws

  • Errors if the body does not contain a set with set_name.
  • Errors if the directions are not correctly specified.
  • Errors if the dimensions of the data matrix are incorrect.
source
Peridynamics.InterfaceErrorType
InterfaceError

A type for a customized error that is thrown when a material model is not implemented correctly.

Fields

  • type::DataType: Type that is used.
  • func::String: Function that is used.
source
Peridynamics.HaloExchangeType
HaloExchange

A type used for communication between body chunks. This type is used with both MPI and multithreaded simulations. Note that the tag is only used with MPI and has the value 0 in multithreaded simulations.

Fields

  • tag::Int: Tag used for the MPI sending and receiving commands.
  • src_chunk_id::Int: Index of the chunk that sends information.
  • dest_chunk_id::Int: Index of the chunk that receives information.
  • src_idxs::Vector{Int}: Indices of the points in the source chunk that send information.
  • dest_idxs::Vector{Int}: Indices of the points in the destination chunk that receive information.
source
Peridynamics.JobOptionsType
JobOptions{F,V}

A type that contains the options of a job.

Type Parameters

  • F: Type for fields of simulation.
  • V: Type for basename of vtk-files.

Fields

  • export_allowed::Bool: Specify if data is exported for the job.
  • root::String: Path of the folder where all data is saved.
  • vtk::String: Path of the folder where vtk-files are saved.
  • logfile::String: Complete path of the logfile.
  • freq::Int: Frequency of exported time steps.
  • fields::F: Exported fields of the job.
  • vtk_filebase::V: Basename of exported vtk-files.
source
Peridynamics.MPIHaloInfoType
MPIHaloInfo

A type providing information about the halo points of body chunks.

Fields

  • point_ids::Dict{Int,Vector{Int}}: Indices of local and halo points of each chunk.
  • halos_points::Dict{Int,Vector{Int}}: Indices of halo points of each chunk.
  • localizers::Dict{Int,Dict{Int,Int}}: Localizes global indices to local indices in the corresponding chunks.
  • hidxs_by_src::Dict{Int,Dict{Int,Vector{Int}}}: Dict specifying for each chunk the indices of halo points but depending on the chunk they belong to.
source
Peridynamics.MPIBodyDataHandlerType
MPIBodyDataHandler

A type for handling information and communication for a body chunk in MPI.

Type Parameters

  • Sys: Type of system of the body.
  • M: Type of the material model.
  • P: Material parameter type.
  • S: Storage type of the system.
  • Bufs: Type for buffers.

Fields

  • chunk::BodyChunk{Sys,M,P,S}: Body chunk of the body.
  • n_halo_fields::Int: Number of chunks in communication with this chunk due to existing halo points.
  • lth_exs_send::Vector{HaloExchange}: Local-to-halo-exchanges of the body chunk that send information.
  • lth_exs_recv::Vector{HaloExchange}: Local-to-halo-exchanges of the body chunk that receive information.
  • htl_exs_send::Vector{HaloExchange}: Halo-to-local-exchanges of the body chunk that send information.
  • htl_exs_recv::Vector{HaloExchange}: Halo-to-local-exchanges of the body chunk that receive information.
  • lth_send_bufs::Vector{Bufs}: Buffers for local-to-halo-exchanges that send information.
  • lth_recv_bufs::Vector{Bufs}: Buffers for local-to-halo-exchanges that receive information.
  • htl_send_bufs::Vector{Bufs}: Buffers for halo-to-local-exchanges that send information.
  • htl_recv_bufs::Vector{Bufs}: Buffers for halo-to-local-exchanges that receive information.
  • field_to_buf::Dict{Symbol,Int}: Dict specifying the index in the buffer of each field of the simulation.
  • lth_reqs::Vector{Vector{MPI.Request}}: Pre-allocated buffer for requests for local-to-halo-exchanges.
  • htl_reqs::Vector{Vector{MPI.Request}}: Pre-allocated buffer for requests for halo-to-local-exchanges.
source
Peridynamics.ParameterHandlerType
ParameterHandler

A type used to manage multiple point parameters defined for the same body. It is used to assign different point parameters to the points of a body.

Type Parameters

  • P<:AbstractPointParameters: Point parameter type.

Fields

  • parameters::Vector{P}: All parameter sets defined in the simulation.
  • point_mapping::Vector{Int}: Vector assigning the related parameter set to each material point.
source
Peridynamics.ThreadsBodyDataHandlerType
ThreadsBodyDataHandler

A type for handling all data of a body in multithreading simulations.

Type Parameters

  • Sys: Type of system of the body.
  • M: Type of the material model.
  • P: Material parameter type.
  • S: Storage type of the system.

Fields

  • n_chunks::Int: Number of chunks of the body.
  • chunks::Vector{BodyChunk{Sys,M,P,S}}: All body chunks of the body.
  • lth_exs::Vector{Vector{HaloExchange}}: All local-to-halo-exchanges of each body chunk of the body.
  • htl_exs::Vector{Vector{HaloExchange}}: All halo-to-local-exchanges of each body chunk of the body.
source
Peridynamics.ThreadsMultibodyDataHandlerType
ThreadsMultibodyDataHandler

A type for handling all data of multiple bodies in multithreading simulations.

Type Parameters

  • BDH: Body data handler type.
  • PC: Position cache type.
  • VC: Volume cache type.

Fields

  • n_bodies::Int: Number of bodies in the simulation.
  • body_dhs::BDH: Tuple containing all body data handlers.
  • body_names::Vector{Symbol}: Names of the bodies.
  • body_idxs::Dict{Symbol,Int}: Names of bodies assigned to their indices.
  • srf_contacts::Vector{ShortRangeForceContact}: All short range force contacts of this simulation.
  • position_caches::PC: Positions of all points of all bodies (should be of type Vector{Matrix{Float64}}).
  • volume_caches::VC: Volumes of all points of all bodies (should be of type Vector{Vector{Float64}}).
source
Peridynamics.BodyChunkType
BodyChunk{System, Material, Params, Storage}

A type that contains all data of a body chunk. For parallel simulations, the body is divided into multiple chunks. Each BodyChunk instance contains all necessary information for the simulation on this specific chunk. This type is used for multithreading and MPI.

Type Parameters

  • System<:AbstractSystem: Type of the system.
  • Material<:AbstractMaterial: Type of the material model of the system.
  • Params<:AbstractParameterSetup: Material parameters of the points in the body chunk.
  • Storage<:AbstractStorage: Storage of all information that changes during the simulation.

Fields

  • body_name::Symbol: Name of the body in multibody simulations.
  • system::System: System with all information that is known before the simulation.
  • mat::Material: Material model of the system.
  • paramsetup::Params: Material parameters of the points in the body chunk.
  • storage::Storage: Storage of all information that changes during the simulation.
  • psets::Dict{Symbol,Vector{Int}}: Point sets of the chunk with local indices.
  • sdbcs::Vector{SingleDimBC}: Single dimension boundary conditions.
  • pdsdbcs::Vector{PosDepSingleDimBC}: Position dependent single dimension boundary conditions.
  • cells::Vector{MeshCell{VTKCellType,Tuple{Int64}}}: Cells for vtk export.
source
Peridynamics.BondType
Bond

Type that describes a bond of two points in a peridynamics body.

Fields

  • neighbor::Int: The index of the neighbor point with which the bond is formed.
  • length::Float64: The length of the bond.
  • fail_permit::Bool: Describes whether failure is allowed or not for this bond.
source
Peridynamics.BondSystemType
BondSystem{Correction}

A type for a system for all peridynamic formulations that work with just bonds of two points.

Type Parameters

  • Correction<:AbstractCorrection: Applied surface correction.

Fields

  • position::Matrix{Float64}: Positions of all points of the system.
  • volume::Vector{Float64}: Volumes of the points of the system.
  • bonds::Vector{Bond}: Vector containing all bonds of the bond system.
  • n_neighbors::Vector{Int}: Number of neighbors for each point of the system.
  • bond_ids::Vector{UnitRange{Int}}: Range of the bonds vector containing bonds of considered point.
  • correction::Correction: Applied surface correction.
  • chunk_handler::ChunkHandler: Type to handle the chunks for the simulation. See ChunkHandler.
source
Peridynamics.ChunkHandlerType
ChunkHandler

A type to handle a body chunk and its communication to other chunks.

Fields

  • n_loc_points::Int: Number of local points that belong to the body chunk.
  • point_ids::Vector{Int}: Indices of all local and halo points of the chunk.
  • loc_points::UnitRange{Int}: Indices of local points of the chunk.
  • halo_points::Vector{Int}: Indices of halo points of the chunk.
  • hidxs_by_src::Dict{Int,UnitRange{Int}}: Dict specifying the indices of halo Points depending on the body chunk they belong to.
  • localizer::Dict{Int,Int}: Localizes global indices to local indices in this chunk.
source
Peridynamics.PointDecompositionType
PointDecomposition

A type that describes how a body is divided into multiple body chunks.

Fields

  • n_chunks::Int: Number of body chunks.
  • decomp::Vector{UnitRange{Int}}: Indices of the points belonging to each chunk.
  • point_src::Dict{Int,Int}: Dict that assigns all point indices to the chunk they belong to.
source
Peridynamics.TwoNeighborInteractionType
TwoNeighborInteraction

Type for two-neighbor interactions.

Fields

  • oni_j::Int: One-neighbor interaction of considered point with point j.
  • oni_k::Int: One-neighbor interaction of considered point with point k.
  • surface::Float64: Surface spread by this two-neighbor interaction.
source
Peridynamics.ThreeNeighborInteractionType
ThreeNeighborInteraction

Type for three-neighbor interactions.

Fields

  • oni_j::Int: One-neighbor interaction of considered point with point j.
  • oni_k::Int: One-neighbor interaction of considered point with point k.
  • oni_l::Int: One-neighbor interaction of considered point with point l.
  • volume::Float64: Volume spread by this three-neighbor interaction.
source
Peridynamics.InteractionSystemType
InteractionSystem

A peridynamic system type that is mainly designed for continuum-kinematics-inspired peridynamics [JMS19].

Fields

  • position::Matrix{Float64}: Positions of all points of the system.
  • one_nis::Vector{Bond}: Vector containing all one-neighbor interactions (bonds) of the system.
  • two_nis::Vector{TwoNeighborInteraction}: Vector containing all two-neighbor interactions of the system.
  • three_nis::Vector{ThreeNeighborInteraction}: Vector containing all three-neighbor interactions of the system.
  • volume::Vector{Float64}: Volumes of the points of the system.
  • volume_one_nis::Vector{Float64}: Effective volumes of one-neighbor interactions.
  • volume_two_nis::Vector{Float64}: Effective volumes of two-neighbor interactions.
  • volume_three_nis::Vector{Float64}: Effective volumes of three-neighbor interactions.
  • n_one_nis::Vector{Int}: Number of one-neighbor interactions for each point of the system.
  • n_two_nis::Vector{Int}: Number of two-neighbor interactions for each point of the system.
  • n_three_nis::Vector{Int}: Number of three-neighbor interactions for each point of the system.
  • one_ni_idxs::Vector{UnitRange{Int}}: Range of the one-neighbor interactions vector containing interactions of considered point.
  • two_ni_idxs::Vector{UnitRange{Int}}: Range of the two-neighbor interactions vector containing interactions of considered point.
  • three_ni_idxs::Vector{UnitRange{Int}}: Range of the three-neighbor interactions vector containing interactions of considered point.
  • chunk_handler::ChunkHandler: Type to handle the chunks for the simulation. See ChunkHandler.
source
Peridynamics.PointSetsPreCrackType
PointSetsPreCrack

Type describing a predefined crack in a peridynamic body.

Fields

  • set_a::Symbol: Point set containing points on one side of the crack.
  • set_b::Symbol: Point set with points on other side of the crack.
  • filter_bonds::Bool: If true, the involved bonds are filtered out so no damage is present at the beginning of the simulation. Else, all involved bonds are marked broken from the beginning.
source
Peridynamics.StandardPointParametersType
StandardPointParameters

Type containing the material parameters for a standard peridynamics model using the bond-based, ordinary state-based or non-ordinary state-based correspondence formulation of peridynamics.

Fields

  • δ::Float64: Horizon.
  • rho::Float64: Density.
  • E::Float64: Young's modulus.
  • nu::Float64: Poisson's ratio.
  • G::Float64: Shear modulus.
  • K::Float64: Bulk modulus.
  • λ::Float64: 1st Lamé parameter.
  • μ::Float64: 2nd Lamé parameter.
  • Gc::Float64: Critical energy release rate.
  • εc::Float64: Critical strain.
  • bc::Float64: Bond constant.
source
Peridynamics.SingleDimBCType
SingleDimBC{F}

Type for a boundary condition in a single dimension for a peridynamic simulation.

Type Parameters

  • F<:Function: Time dependent function which describes the boundary condition.

Fields

  • fun::F: Time dependent function which describes the boundary condition.
  • field::Symbol: Field of the condition (e.g. velocity, force density).
  • point_set::Symbol: Point set on which the condition is applied.
  • dim::UInt8: Dimension in which the condition is applied.
source
Peridynamics.PosDepSingleDimBCType
PosDepSingleDimBC{F}

Type for a position dependent boundary condition in a single dimension for a peridynamic simulation.

Type Parameters

  • F<:Function: Position and time dependent function which describes the boundary condition.

Fields

  • fun::F: Position and time dependent function which describes the boundary condition.
  • field::Symbol: Field of the condition (e.g. velocity, force density).
  • point_set::Symbol: Point set on which the condition is applied.
  • dim::UInt8: Dimension in which the condition is applied.
source
Peridynamics.CKIPointParametersType
CKIPointParameters

Type containing the material parameters for a continuum-kinematics-inspired peridynamics model.

Fields

  • δ::Float64: Horizon.
  • rho::Float64: Density.
  • E::Float64: Young's modulus.
  • nu::Float64: Poisson's ratio.
  • G::Float64: Shear modulus.
  • K::Float64: Bulk modulus.
  • λ::Float64: 1st Lamé parameter.
  • μ::Float64: 2nd Lamé parameter.
  • Gc::Float64: Critical energy release rate.
  • εc::Float64: Critical strain.
  • C1::Float64: Material constant for one-neighbor interactions.
  • C2::Float64: Material constant for two-neighbor interactions.
  • C3::Float64: Material constant for three-neighbor interactions.
source
Peridynamics.BACPointParametersType
BACPointParameters

Type containing the material parameters for a peridynamics model using the bond-associated correspondence formulation of Chen and Spencer.

Fields

  • δ::Float64: Horizon.
  • δb::Float64: Bond-associated horizon.
  • rho::Float64: Density.
  • E::Float64: Young's modulus.
  • nu::Float64: Poisson's ratio.
  • G::Float64: Shear modulus.
  • K::Float64: Bulk modulus.
  • λ::Float64: 1st Lamé parameter.
  • μ::Float64: 2nd Lamé parameter.
  • Gc::Float64: Critical energy release rate.
  • εc::Float64: Critical strain.
  • bc::Float64: Bond constant.
source
Peridynamics.SingleDimICType
SingleDimIC

Type for an initial condition in a single dimension for a peridynamic simulation.

Fields

  • value::Float64: Value of the condition.
  • field::Symbol: Field of the condition (e.g. velocity, force density).
  • point_set::Symbol: Point set on which the condition is applied.
  • dim::UInt8: Dimension in which the condition is applied.
source
Peridynamics.PosDepSingleDimICType
PosDepSingleDimIC{F}

Type for a position dependent initial condition in a single dimension for a peridynamic simulation.

Type Parameters

  • F<:Function: Position dependent function which describes the initial condition.

Fields

  • fun::F: Position dependent function which describes the initial condition.
  • field::Symbol: Field of the condition (e.g. velocity, force density).
  • point_set::Symbol: Point set on which the condition is applied.
  • dim::UInt8: Dimension in which the condition is applied.
source
Peridynamics.ShortRangeForceContactType
ShortRangeForceContact

A type for contact simulations with the short range forces algorithm.

Type Parameters

  • N: Neighborhood search object used by PointNeighbors.jl.

Fields

  • body_id_a::Symbol: Index of a body of the contact.
  • body_id_b::Symbol: Index of a body of the contact.
  • radius::Float64: Search radius for contact.
  • penalty_factor::Float64: Penalty factor for the contact simulation.
  • nhs::N: Neighborhood search object used by PointNeighbors.jl.
source
Peridynamics.@storageMacro
@storage material storage
@storage material solver storage

A macro for automatic creation of a storage type. The macro then maps a specified material type with the storage structure. When no solver is specified, this storage is always used independently of the solver type. Specifying a solver allows the usage of customized storages for a solver.

The following macros can be used before field definitions inside the storage struct:

  • @pointfield: Specifies a point field (data that each point has). Most of the time a Matrix (column i is a vectorial quantity of the local point i) or a Vector (entry i is a scalar value of local point i).
  • @lthfield: Specifies a point field with local-to-halo exchange during the simulation (creates a buffer and automatically updates this field when running the local-to-halo update functions).
  • @htlfield: Specifies a point field with halo-to-local exchange during the simulation, similar to @lthfield.

Example

Example definition of the storage for the bond-based material: ````julia @storage BBMaterial struct BBStorage <: AbstractStorage @lthfield position::Matrix{Float64} @pointfield displacement::Matrix{Float64} @pointfield velocity::Matrix{Float64} @pointfield velocityhalf::Matrix{Float64} @pointfield velocityhalfold::Matrix{Float64} @pointfield acceleration::Matrix{Float64} @pointfield bint::Matrix{Float64} @pointfield bintold::Matrix{Float64} @pointfield bext::Matrix{Float64} @pointfield densitymatrix::Matrix{Float64} @pointfield damage::Vector{Float64} bondactive::Vector{Bool} @pointfield nactive_bonds::Vector{Int} end

source