Libary
Peridynamics.BodySetup
Peridynamics.BondBasedMaterial
Peridynamics.Contact
Peridynamics.ExportSettings
Peridynamics.ForceDensityBC
Peridynamics.PDContactAnalysis
Peridynamics.PDSingleBodyAnalysis
Peridynamics.PointCloud
Peridynamics.PosDepVelBC
Peridynamics.PreCrack
Peridynamics.TimeDiscretization
Peridynamics.VelocityBC
Peridynamics.VelocityIC
Peridynamics.calc_stable_user_timestep
Peridynamics.read_inp
Peridynamics.sphere_radius
Peridynamics.submit
Types
Peridynamics.PointCloud
— TypePointCloud
Peridynamic spatial discretization with material points defining a point cloud.
Fields
n_points::Int
: number of material pointsposition::Matrix{Float64}
: coordinates of points in reference configurationvolume::Vector{Float64}
: material point volumesfailure_flag::BitVector
: if failure of point is possible: element=true
radius::Vector{Float64}
: radius of the material point sphere
PointCloud(position, volume[, point_sets])
Create a PointCloud
by specifying position and volume of all material points. The radius of the point sphere is derived from the volume with the function sphere_radius
and the failure_flag
set to true
for all points.
Arguments
position::Matrix{Float64}
: the position of the material points ($3 \times N$-matrix for $N$ material points)volume::Vector{Float64}
: the volume of the material pointspoint_sets::Dict{String,Vector{Int}}=Dict{String,Vector{Int}}()
: optional point sets
Returns
PointCloud
: point cloud with specified material point position and volume
Examples
Creating a PointCloud
with 4 manually defined points:
julia> position = [
0.0 1.0 0.0 0.0
0.0 0.0 1.0 0.0
0.0 0.0 0.0 1.0
]
3×4 Matrix{Float64}:
0.0 1.0 0.0 0.0
0.0 0.0 1.0 0.0
0.0 0.0 0.0 1.0
julia> volume = [1.0, 1.0, 1.0, 1.0]
4-element Vector{Float64}:
1.0
1.0
1.0
1.0
julia> pc = PointCloud(position, volume)
4-points PointCloud
julia> pc.failure_flag
4-element BitVector:
1
1
1
1
julia> pc.radius
4-element Vector{Float64}:
0.6203504908994001
0.6203504908994001
0.6203504908994001
0.6203504908994001
PointCloud(W, L, H, Δx; center_x=0, center_y=0, center_z=0)
Generate a uniformly distributed PointCloud
with width W
, length L
, height H
and point spacing Δx
. Optional keyword arguments provide the possibility to set the center.
Arguments
W::Real
: width of thePointCloud
-block in x-directionL::Real
: length of thePointCloud
-block in y-directionH::Real
: height of thePointCloud
-block in z-directionΔx::Real
: point spacing in x-, y- and z- directioncenter_x::Real=0
: x-coordinate of thePointCloud
-block centercenter_y::Real=0
: y-coordinate of thePointCloud
-block centercenter_z::Real=0
: z-coordinate of thePointCloud
-block center
Returns
PointCloud
: point cloud with withW
, lengthL
, heightH
and point spacingΔx
Examples
Cube with side length 1 and point spacing $Δx = 0.1$:
julia> PointCloud(1, 1, 1, 0.1)
1000-points PointCloud
Peridynamics.BondBasedMaterial
— TypeBondBasedMaterial <: AbstractPDMaterial
Bond based peridynamic material model.
Fields
δ::Float64
: horizonrho::Float64
: densityE::Float64
: young's modulusnu::Float64
: poisson ratioG::Float64
: shear modulusK::Float64
: bulk modulusbc::Float64
: bond constantGc::Float64
: critical energy release rateεc::Float64
: critical bond strain
BondBasedMaterial(; horizon::Real, rho::Real, E::Real, Gc::Real)
Specify a material only with horizon
, density rho
, Young's modulus E
and critical energy release rate Gc
.
Keywords
horizon::Real
: horizonrho::Real
:densityE::Real
: young's modulusGc::Real
: critical energy release rate
Peridynamics.PreCrack
— TypePreCrack(point_id_set_a::Vector{Int}, point_id_set_b::Vector{Int})
Definition of an preexisting crack in the model. Points in point_id_set_a
cannot have interactions with points in point_id_set_b
.
Fields
point_id_set_a::Vector{Int}
: first point-id setpoint_id_set_b::Vector{Int}
: second point-id set
Peridynamics.TimeDiscretization
— TypeTimeDiscretization
Time discretization type for setting the number of timesteps and the timestep Δt
.
Fields
n_timesteps::Int
: number of time stepsΔt::Float64
: constant time stepalg::Symbol
: algorithm used for time integration. Possible values::verlet
: Velocity verlet algorithm for explicit time integration:dynrelax
: Adaptive dynamic relaxation for quasistatic time integration
TimeDiscretization(n_timesteps::Int[, Δt::Real]; alg::Symbol=:verlet)
Arguments
n_timesteps::Int
: number of time stepsΔt::Real
: optional specified time step
Keywords
alg::Symbol
: optional specification of algorithm used for time integration. Possible values::verlet
(default): Velocity verlet algorithm for explicit time integration:dynrelax
: Adaptive dynamic relaxation for quasistatic time integration
Peridynamics.ExportSettings
— TypeExportSettings
Export settings.
Fields
path::String
: path where results will be savedexportfreq::Int
: export frequency, will export everyexportfreq
-th timestepresultfile_prefix::String
: prefix of the result-filenamelogfile::String
: name of logfileexportflag::Bool
: disable export for a simulation where saved results are not needed
ExportSettings([path::String, freq::Int])
Create ExportSettings
only by path
and freq
. If no arguments are specified, the exportflag
will be set to false
and export disabled.
Arguments
path::String
: path where results will be savedfreq::Int
: export frequency
Peridynamics.BodySetup
— TypeBodySetup
Setup of multiple bodies for PDContactAnalysis
.
Fields:
pc::PointCloud
: point cloudmat::AbstractPDMaterial
: material modelprecracks::Vector{PreCrack}
: predefined cracksbcs::Vector{<:AbstractBC}
: boundary conditionsics::Vector{<:AbstractIC}
: initial conditionscalc_timestep::Bool
: use body for time step calculation
Peridynamics.Contact
— TypeContact
Contact definition.
Fields
body_id_set::Tuple{Int,Int}
: bodies for which contact is definedsearch_radius::Float64
: search radius for contact definition to activatespring_constant::Float64
: spring constant used for contact force calculation, default value:spring_constant = 1.0e12
Peridynamics.PDSingleBodyAnalysis
— TypePDSingleBodyAnalysis{T<:AbstractPDMaterial} <: AbstractPDAnalysis
Peridynamic single body analysis.
Fields
name::String
: simulation namepc::PointCloud
: point cloudmat::T
: material model for the bodyprecracks::Vector{PreCrack}
: predefined cracksbcs::Vector{<:AbstractBC}
: boundary conditionsics::Vector{<:AbstractIC}
: initial conditionstd::TimeDiscretization
: time discretizationes::ExportSettings
: export settings
Peridynamics.PDContactAnalysis
— TypePDContactAnalysis <: AbstractPDAnalysis
Peridynamic contact analysis.
Fields
name::String
: simulation namebody_setup::Vector{BodySetup}
: bodies used in this simulationn_bodies::Int
: number of bodiescontact::Vector{Contact}
: contact definitionstd::TimeDiscretization
: time discretizationes::ExportSettings
: export settings
Peridynamics.VelocityBC
— TypeVelocityBC <: AbstractBC
Velocity boundary condition. The value of the velocity is calculated every time step with the function fun
and applied to the dimension dim
of the points specified by point_id_set
.
Fields
fun::Function
: functionf(t)
with current timet
as argument, that calculates the value of the velocity for each timesteppoint_id_set::Vector{Int}
: point-id set with all points for which the boundary condition gets applied every timestepdim::Int
: dimension on which the boundary condition gets applied to. Possible values:- x-direction:
dim=1
- y-direction:
dim=2
- z-direction:
dim=3
- x-direction:
Examples
The constant velocity $v = 0.1$ in $y$-direction gets applied to the first $10$ points:
julia> VelocityBC(t -> 0.1, 1:10, 2)
VelocityBC(var"#7#8"(), [1, 2, 3, 4, 5, 6, 7, 8, 9, 10], 2)
Peridynamics.VelocityIC
— TypeVelocityIC <: AbstractIC
Velocity initial condition. The value val
of the velocity is applied as initial condition to the dimension dim
of the points specified by point_id_set
.
Fields
val::Float64
: value of the velocitypoint_id_set::Vector{Int}
: point-id set with all points for which the initial condition gets applied todim::Int
: dimension on which the initial condition gets applied to. Possible values:- x-direction:
dim=1
- y-direction:
dim=2
- z-direction:
dim=3
- x-direction:
Peridynamics.PosDepVelBC
— TypePosDepVelBC <: AbstractBC
Position dependent velocity boundary condition. The value of the force density is calculated every time step with the function fun
and applied to the dimension dim
of the points specified by point_id_set
.
Fields
fun::Function
: functionf(x, y, z, t)
with x-, y-, z-position of the point and current timet
as arguments, calculates the value of the force density for each timestep dependent of the point positionpoint_id_set::Vector{Int}
: point-id set with all points for which the boundary condition gets applied todim::Int
: dimension on which the initial condition gets applied to. Possible values:- x-direction:
dim=1
- y-direction:
dim=2
- z-direction:
dim=3
- x-direction:
Peridynamics.ForceDensityBC
— TypeForceDensityBC <: AbstractBC
Force density boundary condition. The value of the force density is calculated every time step with the function fun
and applied to the dimension dim
of the points specified by point_id_set
.
Fields
fun::Function
: functionf(t)
with current timet
as argument, that calculates the value of the force density for each timesteppoint_id_set::Vector{Int}
: point-id set with all points for which the boundary condition gets applied every timestepdim::Int
: dimension on which the boundary condition gets applied to. Possible values:- x-direction:
dim=1
- y-direction:
dim=2
- z-direction:
dim=3
- x-direction:
Functions
Peridynamics.calc_stable_user_timestep
— Functioncalc_stable_user_timestep(pc::PointCloud, mat::AbstractPDMaterial, Sf::Float64=0.7)
Function to determine the stable timestep for the specified point cloud.
Arguments
pc::PointCloud
: point cloudmat::AbstractPDMaterial
: material modelSf::Float64
: safety factor for time step, default valueSf = 0.7
Returns
Float64
: stable user timestepΔt
Peridynamics.read_inp
— Functionread_inp(file::String)
Read Abaqus .inp-file and convert meshes to PointCloud
objects with the help of the AbaqusReader.jl
package. Every element is converted to a point. The center of the element becomes the position of the point and the element volume becomes the point volume. Element sets defined in Abaqus are converted to corresponding point sets.
Currently supported mesh elements: [:Tet4, :Hex8]
Arguments
file::String
: path to Abaqus .inp-file
Returns
PointCloud
: generated point cloud with element volume as point volume
Examples
The specimen of the TensileTest.jl
example script:
julia> pointcloud = read_inp("examples/models/TensileTestMesh.inp")
[ Info: 21420 nodes found
[ Info: Parsing elements. Type: C3D8. Topology: Hex8
[ Info: Creating elset bottom
[ Info: Creating elset top
16900-points PointCloud
julia> pointcloud.point_sets
Dict{String, Vector{Int64}} with 2 entries:
"bottom" => [11701, 11702, 11703, 11704, 11705, 11706, 11707, 11708, 117…
"top" => [7501, 7502, 7503, 7504, 7505, 7506, 7507, 7508, 7509, 7510 …
Peridynamics.sphere_radius
— Functionsphere_radius(vol::T) where {T<:Real}
Calculate the radius $r$ of the sphere by equation
\[r = \sqrt[3]{\frac{3 \; V}{4 \; \pi}}\]
with specified sphere volume $V$.
Arguments
vol::T where {T<:Real}
: volume $V$ of the sphere
Returns
Float64
: radius $r$ of sphere with volumevol
Peridynamics.submit
— Functionsubmit(sim::T) where {T<:AbstractPDAnalysis}
Submit a simulation job to determine the results of the specified analysis.
Arguments
sim::T where {T<:AbstractPDAnalysis}
: simulation job