Peridynamics.jl public API
Material models
Peridynamics.BBMaterial
— TypeBBMaterial()
BBMaterial{Correction}()
A material type used to assign the material of a Body
with the standard bond-based formulation of peridynamics.
Possible correction methods are:
NoCorrection
: No correction is applied. (default)EnergySurfaceCorrection
: The energy based surface correction method of Le and Bobaru (2018) is applied.
Examples
julia> mat = BBMaterial()
BBMaterial{NoCorrection}()
julia> mat = BBMaterial{EnergySurfaceCorrection}()
BBMaterial{EnergySurfaceCorrection}()
BBMaterial{Correction}
Material type for the bond-based peridynamics formulation.
Type Parameters
Correction
: A correction algorithm type. See the constructor docs for more informations.
Allowed material parameters
When using material!
on a Body
with BBMaterial
, then the following parameters are allowed:
horizon::Float64
: Radius of point interactionsrho::Float64
: DensityE::Float64
: Young's modulusGc::Float64
: Critical energy release rateepsilon_c::Float64
: Critical strain
In bond-based peridynamics, the Poisson's ratio is limited to 1/4 for 3D simulations. Therefore the specification of this keyword is not allowed when using material!
, as it is hardcoded to nu = 1/4
.
Allowed export fields
When specifying the fields
keyword of Job
for a Body
with BBMaterial
, the following fields are allowed:
position::Matrix{Float64}
: Position of each pointdisplacement::Matrix{Float64}
: Displacement of each pointvelocity::Matrix{Float64}
: Velocity of each pointvelocity_half::Matrix{Float64}
: Velocity parameter for Verlet time solveracceleration::Matrix{Float64}
: Acceleration of each pointb_int::Matrix{Float64}
: Internal force density of each pointb_ext::Matrix{Float64}
: External force density of each pointdamage::Vector{Float64}
: Damage of each pointn_active_bonds::Vector{Int}
: Number of intact bonds of each point
Peridynamics.OSBMaterial
— TypeOSBMaterial()
OSBMaterial{Correction}()
A material type used to assign the material of a Body
with the ordinary state-based formulation of peridynamics.
Possible correction methods are:
NoCorrection
: No correction is applied (default)EnergySurfaceCorrection
: The energy based surface correction method of Le and Bobaru (2018) is applied
Examples
julia> mat = OSBMaterial()
OSBMaterial{NoCorrection}()
julia> mat = OSBMaterial{EnergySurfaceCorrection}()
OSBMaterial{EnergySurfaceCorrection}()
OSBMaterial{Correction}
Material type for the ordinary state-based peridynamics formulation.
Type Parameters
Correction
: A correction algorithm type. See the constructor docs for more informations.
Allowed material parameters
When using material!
on a Body
with OSBMaterial
, then the following parameters are allowed:
horizon::Float64
: Radius of point interactionsrho::Float64
: DensityE::Float64
: Young's modulusnu::Float64
: Poisson's ratioGc::Float64
: Critical energy release rateepsilon_c::Float64
: Critical strain
Allowed export fields
When specifying the fields
keyword of Job
for a Body
with OSBMaterial
, the following fields are allowed:
position::Matrix{Float64}
: Position of each pointdisplacement::Matrix{Float64}
: Displacement of each pointvelocity::Matrix{Float64}
: Velocity of each pointvelocity_half::Matrix{Float64}
: Velocity parameter for Verlet time solveracceleration::Matrix{Float64}
: Acceleration of each pointb_int::Matrix{Float64}
: Internal force density of each pointb_ext::Matrix{Float64}
: External force density of each pointdamage::Vector{Float64}
: Damage of each pointn_active_bonds::Vector{Int}
: Number of intact bonds of each point
Peridynamics.NOSBMaterial
— TypeNOSBMaterial(; maxdmg, maxjacobi, corr)
A material type used to assign the material of a Body
with the local continuum consistent (correspondence) formulation of non-ordinary state-based peridynamics.
Keywords
maxdmg::Float64
: Maximum value of damage a point is allowed to obtain. If this value is exceeded, all bonds of that point are broken because the deformation gradient would then possibly containNaN
values. (default:0.95
)maxjacobi::Float64
: Maximum value of the Jacobi determinant. If this value is exceeded, all bonds of that point are broken. (default:1.03
)corr::Float64
: Correction factor used for zero-energy mode stabilization. The stabilization algorithm of Silling (2017) is used. (default:100.0
)
This formulation is known to be not suitable for fracture simultations without stabilization of the zero-energy modes. Therefore be careful when doing fracture simulations and try out different paremeters for maxdmg
, maxjacobi
, and corr
.
Examples
julia> mat = NOSBMaterial()
NOSBMaterial(maxdmg=0.95, maxjacobi=1.03, corr=100.0)
NOSBMaterial
Material type for the local continuum consistent (correspondence) formulation of non-ordinary state-based peridynamics.
Fields
maxdmg::Float64
: Maximum value of damage a point is allowed to obtain. See the constructor docs for more informations.maxjacobi::Float64
: Maximum value of the Jacobi determinant. See the constructor docs for more informations.corr::Float64
: Correction factor used for zero-energy mode stabilization. See the constructor docs for more informations.
Allowed material parameters
When using material!
on a Body
with NOSBMaterial
, then the following parameters are allowed:
horizon::Float64
: Radius of point interactionsrho::Float64
: DensityE::Float64
: Young's modulusnu::Float64
: Poisson's ratioGc::Float64
: Critical energy release rateepsilon_c::Float64
: Critical strain
Allowed export fields
When specifying the fields
keyword of Job
for a Body
with NOSBMaterial
, the following fields are allowed:
position::Matrix{Float64}
: Position of each pointdisplacement::Matrix{Float64}
: Displacement of each pointvelocity::Matrix{Float64}
: Velocity of each pointvelocity_half::Matrix{Float64}
: Velocity parameter for Verlet time solveracceleration::Matrix{Float64}
: Acceleration of each pointb_int::Matrix{Float64}
: Internal force density of each pointb_ext::Matrix{Float64}
: External force density of each pointdamage::Vector{Float64}
: Damage of each pointn_active_bonds::Vector{Int}
: Number of intact bonds of each point
Peridynamics.CKIMaterial
— TypeCKIMaterial()
A material type used to assign the material of a Body
with the continuum-kinematics-inspired peridynamics fomulation.
Examples
julia> mat = CKIMaterial()
CKIMaterial()
CKIMaterial
Material type for the continuum-kinematics-inspired peridynamics framework.
Allowed material parameters
When using material!
on a Body
with CKIMaterial
, then the following parameters are allowed:
horizon::Float64
: Radius of point interactionsrho::Float64
: DensityE::Float64
: Young's modulusnu::Float64
: Poisson's ratioGc::Float64
: Critical energy release rateepsilon_c::Float64
: Critical strainC1::Float64
: One-neighbor interaction parameter (default:0.0
)C2::Float64
: Two-neighbor interaction parameter (default:0.0
)C3::Float64
: Two-neighbor interaction parameter (default:0.0
)
If any of the interaction parameters is used with material!
, the the Young's modulus and Poisson's ratio are ignored and only the specified interaction parameters will influence the force density calculated from that interaction.
If no interaction parameter is specified, then the Young's modulus and Poisson's ratio are used to calculate these parameters accordingly to Ekiz, Steinmann, and Javili (2022).
Allowed export fields
When specifying the fields
keyword of Job
for a Body
with CKIMaterial
, the following fields are allowed:
position::Matrix{Float64}
: Position of each pointdisplacement::Matrix{Float64}
: Displacement of each pointvelocity::Matrix{Float64}
: Velocity of each pointvelocity_half::Matrix{Float64}
: Velocity parameter for Verlet time solveracceleration::Matrix{Float64}
: Acceleration of each pointb_int::Matrix{Float64}
: Internal force density of each pointb_ext::Matrix{Float64}
: External force density of each pointdamage::Vector{Float64}
: Damage of each pointn_active_one_nis::Vector{Int}
: Number of intact one-neighbor interactions of each point
System related types
Peridynamics.NoCorrection
— TypeNoCorrection
A correction handler for materials that use the bond system. If NoCorrection
is used, then no correction will be applied.
See also BBMaterial
, OSBMaterial
for further information on how to use the correction type.
Peridynamics.EnergySurfaceCorrection
— TypeEnergySurfaceCorrection
A correction handler for materials that use the bond system. If EnergySurfaceCorrection
is used, then the energy based surface correction method of Le and Bobaru (2018) is used.
See also BBMaterial
, OSBMaterial
for further information on how to use the correction type.
Discretization
Peridynamics.Body
— TypeBody(material, position, volume)
Body(material, inp_file)
Constructs a Body
for a peridynamics simulation.
Arguments
material::AbstractMaterial
: The material which is defined for the whole body.position::AbstractMatrix
: A3×n
matrix with the point position of then
points.volume::AbstractVector
: A vector with the volume of each point.inp_file::AbstractString
: An Abaqus input file containing meshes, imported withread_inp
.
Throws
- Errors if the number of points is not larger than zero
- Errors if
position
is not a3×n
matrix and has the same length asvolume
- Errors if
position
orvolume
containNaN
values
Example
julia> Body(BBMaterial(), rand(3, 10), rand(10))
10-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
10-point set `all_points`
Please note that the fields are intended for internal use only. They are not part of the public API of Peridynamics.jl, and thus can be altered (or removed) at any time without it being considered a breaking change.
Body{Material,PointParameters}
Type Parameters
Material <: AbstractMaterial
: Type of the specified material modelPointParameters <: AbstractPointParameters
: Type of the point parameters
Fields
mat::Material
: The material formulation.n_points::Int
: The number of points that in the body.position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.fail_permit::Vector{Bool}
: A vector that describes if failure is allowed for each point.point_sets::Dict{Symbol,Vector{Int}}
: A dictionary containing point sets.point_params::Vector{PointParameters}
: A vector containing all point parameters.params_map::Vector{Int}
: A vector that maps parameters inpoint_params
to each point.single_dim_bcs::Vector{SingleDimBC}
: A vector with boundary conditions on a single dimension.posdep_single_dim_bcs::Vector{PosDepSingleDimBC}
: A vector with position dependent boundary conditions on a single dimension.single_dim_ics::Vector{SingleDimIC}
: A vector with initial conditions on a single dimension.posdep_single_dim_ics::Vector{PosDepSingleDimIC}
: A vector with position dependent initial conditions on a single dimension.point_sets_precracks::Vector{PointSetsPreCrack}
: A vector with predefined point set cracks.
Peridynamics.MultibodySetup
— TypeMultibodySetup(body_pairs...)
Setup for a peridynamic simulation with multiple bodies.
Arguments
body_pairs::Pair{Symbol,<:AbstractBody}
: Pairs of:body_name => body_object
. The name of the body has to be specified as a Symbol.
Throws
- Errors if less than 2 bodies are defined
Examples
julia> sphere = Body(BBMaterial(), pos_sphere, vol_sphere)
280-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
280-point set `all_points`
julia> plate = Body(BBMaterial(), pos_plate, vol_plate)
25600-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
25600-point set `all_points`
julia> ms = MultibodySetup(:sphere => sphere, :plate => plate)
25880-point MultibodySetup:
280-point Body{BBMaterial{NoCorrection}} with name `sphere`
25600-point Body{BBMaterial{NoCorrection}} with name `plate`
Please note that the fields are intended for internal use only. They are not part of the public API of Peridynamics.jl, and thus can be altered (or removed) at any time without it being considered a breaking change.
MultibodySetup{Bodies}
Type Parameters
Bodies <: Tuple
: All types of the different bodies in the multibody setup.
Fields
bodies::Bodies
: A Tuple containing all the bodies.body_names::Vector{Symbol}
: All body names.body_idxs::Dict{Symbol,Int}
: A Dict to get the body index with the body name.srf_contacts::Vector{ShortRangeForceContact}
: All short range force contacts.
Peridynamics.point_set!
— Functionpoint_set!(body, set_name, points)
point_set!(fun, body, set_name)
Add a point set to a Body
. The points of the set can be either specified directly with the points::AbstractVector
argument, or as the result of the filter function fun
. By default, a body already contains a point set with the name :all_points
, containg a set with all points.
Arguments
body::AbstractBody
:Body
where the set will be added.set_name::Symbol
: Name of the point set.points::AbstractVector
: Some vector containing the point indices of the set. The indices have to be in bounds with theposition
andvolume
ofbody
.fun::Function
: Function for filtering points. This function accepts only one positional argument and will be used in afindall
call. Depending on the argument name, a different input will be processed:x
: The function will receive the x-coordinate of each point inposition
ofbody
:points = findall(fun, @view(position[1, :]))
y
: The function will receive the y-coordinate of each point inposition
ofbody
:points = findall(fun, @view(position[2, :]))
z
: The function will receive the z-coordinate of each point inposition
ofbody
:points = findall(fun, @view(position[3, :]))
p
: The function will receive the a vector containing each dimension of each point inposition
ofbody
:points = findall(fun, eachcol(position))
Throws
- Errors if a point set with the same
set_name
already exists. - Errors if
points
are not in bounds withposition
andvolume
of thebody
.
Examples
Add a point set to body
with all points that have a x-corrdinate larger than zero:
julia> point_set!(x -> x > 0, body, :larger_than_zero)
julia> point_sets(body)
Dict{Symbol, Vector{Int64}} with 2 entries:
:larger_than_zero => [6, 7, 8, 9, 10, 16, 17, 18, 19, 20 … 9…
:all_points => [1, 2, 3, 4, 5, 6, 7, 8, 9, 10 … 991, 9…
Add a point set to body
with all points that are positioned inside a sphere with radius r
around the center. Note that the do
-syntax can be used, as fun
is the first argument of point_set!
:
julia> point_set!(body, :inside_sphere) do p
sqrt(p[1]^2 + p[2]^2 + p[3]^2) ≤ r
end
julia> point_sets(body)
Dict{Symbol, Vector{Int64}} with 2 entries:
:larger_than_zero => [6, 7, 8, 9, 10, 16, 17, 18, 19, 20 … 9…
:inside_sphere => [1, 2, 3, 4, 5, 6, 7, 8, 9, 10 … 991, 9…
Peridynamics.point_sets
— Functionpoint_sets(body)
Returns all point sets of body
.
Arguments
body::AbstractBody
:Body
.
Example
julia> body = Body(BBMaterial(), rand(3,100), rand(100))
100-point Body{BBMaterial{NoCorrection}}:
100-point set `all_points`
julia> point_set!(body, :set_a, 1:10) # first ten points
julia> Peridynamics.point_sets(body)
Dict{Symbol, Vector{Int64}} with 2 entries:
:all_points => [1, 2, 3, 4, 5, 6, 7, 8, 9, 10 … 91, 92, 93, 94, 95, 96, 97, 98, 9…
:set_a => [1, 2, 3, 4, 5, 6, 7, 8, 9, 10]
Peridynamics.failure_permit!
— Functionfailure_permit!(body, fail_permit)
failure_permit!(body, set_name, fail_permit)
Set the failure permission for points of a body
. By default, failure is allowed for all points. If no set_name
is specified, then the permission fail_permit
is set for all points of the 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
: Iftrue
, failure is allowed, and iffalse
then no bonds of this point are allowed to break during the simulation.
Throws
- Errors if the body does not contain a set with
set_name
.
Examples
julia> failure_permit!(body, false)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
1000 points with no failure permission
Peridynamics.material!
— Functionmaterial!(body, set_name; kwargs...)
material!(body; kwargs...)
Assign material point parameters to points of body
. If no set_name
is specified, then the parameters will be set for all points of the body.
Arguments
body::AbstractBody
:Body
.set_name::Symbol
: The name of a point set of this body.
Keywords
Allowed keywords depend on the selected material model. Please look at the documentation of the material you specified when creating the body. The default material keywords are:
horizon::Float64
: Radius of point interactionsrho::Float64
: DensityE::Float64
: Young's modulusnu::Float64
: Poisson's ratioGc::Float64
: Critical energy release rateepsilon_c::Float64
: Critical strain
Throws
- Errors if a kwarg is not eligible for specification with the body material.
Example
julia> material!(body; horizon=3.0, E=2.1e5, rho=8e-6, Gc=2.7)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
1 point parameter(s):
Parameters BBMaterial: δ=3.0, E=210000.0, nu=0.25, rho=8.0e-6, Gc=2.7
Peridynamics.velocity_bc!
— Functionvelocity_bc!(fun, body, set_name, dim)
Specifies velocity boundary conditions for points of the set set_name
in body
. The value of the boundary condition is calculated with the function fun
at every time step.
Arguments
fun::Function
: Condition function for the calculation of a value, should return aFloat64
. If the condition function returns aNaN
, then this value is ignored, which can be used to turn conditions off after a specified period of time. This function accepts one ore two positional arguments and is aware of the argument names. Possible arguments and names:fun(t)
: The function will receive the current timet
at every time step. This makes it possible to specify conditions that change over time.fun(p, t)
: This function will be processed for every point ofset_name
and receives the reference position of a point asSVector{3}
and the current timet
at every time step. This makes it possible to specify conditions that also depend on the position of a point.
body::AbstractBody
:Body
the condition is specified on.set_name::Symbol
: The name of a point set of this body.dim::Union{Integer,Symbol}
: Direction of the condition, either specified as Symbol or integer.- x-direction:
:x
or1
- y-direction:
:y
or2
- z-direction:
:z
or3
- x-direction:
Throws
- Errors if the body does not contain a set with
set_name
. - Errors if the direction is not correctly specified.
- Errors if function is not suitable as condition function and has the wrong arguments.
Example
julia> velocity_bc!(t -> 2.0, body, :all_points, 1)
julia> velocity_bc!((p,t) -> p[1] * t, body, :all_points, :y)
julia> velocity_bc!(t -> t > 0.00001 ? 1.0 : NaN, body, :all_points, :z)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
3 boundary condition(s):
BC on velocity: point_set=all_points, dim=1
BC on velocity: point_set=all_points, dim=3
Pos.-dep. BC on velocity: point_set=all_points, dim=2
Peridynamics.velocity_ic!
— Functionvelocity_ic!(body, set_name, dim, value)
velocity_ic!(fun, body, set_name, dim)
Specifies velocity initial conditions for points of the set set_name
in body
. The value
of the initial condition is specified before time integration. If a function fun
is specified, then the value is with that function.
Arguments
body::AbstractBody
:Body
the condition is specified on.set_name::Symbol
: The name of a point set of this body.dim::Union{Integer,Symbol}
: Direction of the condition, either specified as Symbol or integer.- x-direction:
:x
or1
- y-direction:
:y
or2
- z-direction:
:z
or3
- x-direction:
value::Real
: Value that is specified before time integration.fun::Function
: Condition function for the calculation of a value, should return aFloat64
. If the condition function returns aNaN
, then this value is ignored, which can be used to turn off the condition for a specified position. This function accepts one ore two positional arguments and is aware of the argument names. Possible arguments and names:fun(p)
: The function will receive the reference positionp
of a point asSVector{3}
.
Throws
- Errors if the body does not contain a set with
set_name
. - Errors if the direction is not correctly specified.
- Errors if function is not suitable as condition function and has the wrong arguments.
Example
julia> velocity_ic!(body, :all_points, :x, -100.0)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
1 initial condition(s):
IC on velocity: point_set=all_points, dim=1
Peridynamics.forcedensity_bc!
— Functionforcedensity_bc!(fun, body, set, dim)
Specifies force density boundary conditions for points of the set set_name
in body
. The value of the boundary condition is calculated with the function fun
at every time step.
Arguments
fun::Function
: Condition function for the calculation of a value, should return aFloat64
. If the condition function returns aNaN
, then this value is ignored, which can be used to turn conditions off after a specified period of time. This function accepts one ore two positional arguments and is aware of the argument names. Possible arguments and names:fun(t)
: The function will receive the current timet
at every time step. This makes it possible to specify conditions that change over time.fun(p, t)
: This function will be processed for every point ofset_name
and receives the reference position of a point asSVector{3}
and the current timet
at every time step. This makes it possible to specify conditions that also depend on the position of a point.
body::AbstractBody
:Body
the condition is specified on.set_name::Symbol
: The name of a point set of this body.dim::Union{Integer,Symbol}
: Direction of the condition, either specified as Symbol or integer.- x-direction:
:x
or1
- y-direction:
:y
or2
- z-direction:
:z
or3
- x-direction:
Throws
- Errors if the body does not contain a set with
set_name
. - Errors if the direction is not correctly specified.
- Errors if function is not suitable as condition function and has the wrong arguments.
Example
julia> forcedensity_bc!(t -> 8000.0, body, :all_points, :x)
julia> forcedensity_bc!((p,t) -> p[1] * t, body, :all_points, :y)
julia> forcedensity_bc!(t -> t > 0.00001 ? 8000.0 : NaN, body, :all_points, :z)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
3 boundary condition(s):
BC on force density: point_set=all_points, dim=1
BC on force density: point_set=all_points, dim=3
Pos.-dep. BC on force density: point_set=all_points, dim=2
Peridynamics.precrack!
— Functionprecrack!(body, set_a, set_b; update_dmg=true)
Creates a crack between two point sets by prohibiting interaction between points of different point sets. The points in set_a
are not allowed to interact with points in set_b
.
Arguments
body::AbstractBody
:Body
.set_a::Symbol
: The name of a point set of this body.set_b::Symbol
: The name of a point set of this body.
Keywords
update_dmg::Bool
: Iftrue
, the material points involved in the predefined crack are initially damaged. Iffalse
, the bonds involved are deleted and the material points involved with the predefined crack are not damaged in the reference results. (default:true
)
Throws
- Errors if the body does not contain sets with name
set_a
andset_b
. - Errors if the point sets intersect and a point is included in both sets.
Example
julia> point_set!(body, :a, 1:2)
julia> point_set!(body, :b, 3:4)
julia> precrack!(body, :a, :b)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
3 point set(s):
1000-point set `all_points`
2-point set `a`
2-point set `b`
1 predefined crack(s)
Peridynamics.contact!
— Functioncontact!(multibody_setup, name_body_a, name_body_b; kwargs...)
Defines a short range force contact between body name_body_a
and name_body_b
in the MultibodySetup
multibody_setup
.
Arguments
multibody_setup::MultibodySetup
:MultibodySetup
.name_body_a::Symbol
: The name of a body in this multibody setup.name_body_b::Symbol
: The name of a body in this multibody setup.
Keywords
radius::Float64
: Contact search radius. If a the distance of a point in bodyname_body_a
and a point in bodyname_body_b
is lower than this radius, a contact force is calculated. This radius should be in the order of the point spacing of a point cloud.penalty_factor::Float64
: Penalty factor for the short range force contact algorithm. (default:1e12
)
Throws
- Errors if
multibody_setup
does not contain bodies with namename_body_a
andname_body_b
. - Errors if the keyword
radius
is not specified orradius ≤ 0
. - Errors if
penalty_factor ≤ 0
.
Examples
julia> ms = MultibodySetup(:a => body_a, :b => body_b)
2000-point MultibodySetup:
1000-point Body{BBMaterial{NoCorrection}} with name `b`
1000-point Body{BBMaterial{NoCorrection}} with name `b`
julia> contact!(ms, :a, :b; radius=0.001)
julia> ms
2000-point MultibodySetup:
1000-point Body{BBMaterial{NoCorrection}} with name `b`
1000-point Body{BBMaterial{NoCorrection}} with name `b`
2 short range force contact(s)
Peridynamics.uniform_box
— Functionuniform_box(lx, ly, lz, ΔX0; kwargs...)
Creates a grid of uniformly distributed points in a cuboid with lengths lx
, ly
and lz
and point spacing ΔX0
.
Arguments
lx::Real
: Length in x-dimension.ly::Real
: Length in y-dimension.lz::Real
: Length in z-dimension.ΔX0::Real
: Spacing of the points.
Keywords
center
: The coordinates of the center of the cuboid. Default:(0, 0, 0)
Returns
position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.
Examples
julia> position, volume = uniform_box(10, 10, 10, 2);
julia> position
3×125 Matrix{Float64}:
-4.0 -2.0 0.0 2.0 4.0 -4.0 -2.0 … 0.0 2.0 4.0 -4.0 -2.0 0.0 2.0 4.0
-4.0 -4.0 -4.0 -4.0 -4.0 -2.0 -2.0 2.0 2.0 2.0 4.0 4.0 4.0 4.0 4.0
-4.0 -4.0 -4.0 -4.0 -4.0 -4.0 -4.0 4.0 4.0 4.0 4.0 4.0 4.0 4.0 4.0
julia> volume
125-element Vector{Int64}:
8
8
8
8
⋮
8
8
8
8
Peridynamics.uniform_sphere
— Functionuniform_sphere(diameter, ΔX0; kwargs...)
Creates a grid of uniformly distributed points in a sphere with a specific diameter
and the point spacing ΔX0
. Due to the uniform point spacings, edges on the surface of the sphere can occur.
Arguments
diameter::Real
: Diameter of the sphere.ΔX0::Real
: Spacing of the points.
Keywords
center
: The coordinates of the center of the sphere. Default:(0, 0, 0)
Returns
position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.
Examples
julia> position, volume = uniform_sphere(10, 2);
julia> position
3×81 Matrix{Float64}:
-2.0 0.0 2.0 -2.0 0.0 2.0 -2.0 … 0.0 2.0 -2.0 0.0 2.0 -2.0 0.0 2.0
-2.0 -2.0 -2.0 0.0 0.0 0.0 2.0 -2.0 -2.0 0.0 0.0 0.0 2.0 2.0 2.0
-4.0 -4.0 -4.0 -4.0 -4.0 -4.0 -4.0 4.0 4.0 4.0 4.0 4.0 4.0 4.0 4.0
julia> volume
81-element Vector{Int64}:
8
8
8
8
8
⋮
8
8
8
8
8
Peridynamics.uniform_cylinder
— Functionuniform_cylinder(diameter::Real, height::Real, ΔX0::Real; kwargs...)
Creates a grid of uniformly distributed points in a cylindrical shape with a specific height
in z-dimension, a diameter
and the point spacing ΔX0
. Due to the uniform point spacings, edges on the surface of the cylinder can occur.
Arguments
diameter::Real
: Diameter of the cylinder.height::Real
: Height of the cylinder.ΔX0::Real
: Spacing of the points.
Keywords
center
: The coordinates of the center of the cylinder. Default:(0, 0, 0)
Returns
position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.
Examples
julia> position, volume = uniform_cylinder(5, 10, 2);
julia> position
3×20 Matrix{Float64}:
-1.0 1.0 -1.0 1.0 -1.0 1.0 -1.0 … 1.0 -1.0 1.0 -1.0 1.0 -1.0 1.0
-1.0 -1.0 1.0 1.0 -1.0 -1.0 1.0 -1.0 1.0 1.0 -1.0 -1.0 1.0 1.0
-4.0 -4.0 -4.0 -4.0 -2.0 -2.0 -2.0 2.0 2.0 2.0 4.0 4.0 4.0 4.0
julia> volume
20-element Vector{Int64}:
8
8
8
8
8
⋮
8
8
8
8
8
Peridynamics.round_sphere
— Functionround_sphere(diameter, ΔX0; kwargs...)
Creates a grid of points distributed in a smooth sphere without edges on the surface with a specific diameter
and the point spacing ΔX0
. Internally, some parts of TrixiParticles.jl
were copied and adapted for this function.
Arguments
diameter::Real
: Diameter of the sphere.ΔX0::Real
: Spacing of the points.
Keywords
center
: The coordinates of the center of the sphere. Default:(0, 0, 0)
Returns
position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.
Examples
julia> position, volume = round_sphere(10, 2);
julia> position
3×63 Matrix{Float64}:
0.0 1.74586 0.539501 -1.41243 -1.41243 … -1.95518 -0.48289 1.35303 0.0
0.0 0.0 1.66041 1.02619 -1.02619 -0.941566 -2.11568 -1.69664 0.0
2.0 0.97569 0.97569 0.97569 0.97569 -3.36017 -3.36017 -3.36017 -4.0
julia> volume
63-element Vector{Float64}:
8.311091676163475
8.311091676163475
8.311091676163475
8.311091676163475
8.311091676163475
⋮
8.311091676163475
8.311091676163475
8.311091676163475
8.311091676163475
Peridynamics.round_cylinder
— Functionround_cylinder(diameter::Real, height::Real, ΔX0::Real; kwargs...)
Creates a grid of points distributed in a cylindrical shape with a specific height
in z-dimension, a diameter
and the point spacing ΔX0
. Due to a concentric point distribution in the x-y-plane, there are no sharp edges that appear on the surface of the cylinder.
Arguments
diameter::Real
: Diameter of the cylinder.height::Real
: Height of the cylinder.ΔX0::Real
: Spacing of the points.
Keywords
center
: The coordinates of the center of the cylinder. Default:(0, 0, 0)
Returns
position::Matrix{Float64}
: A3×n_points
matrix with the position of the points.volume::Vector{Float64}
: A vector with the volume of each point.
Examples
julia> position, volume = round_cylinder(5, 10, 2);
julia> position
3×30 Matrix{Float64}:
1.5 0.463525 -1.21353 -1.21353 … -1.21353 -1.21353 0.463525
0.0 1.42658 0.881678 -0.881678 0.881678 -0.881678 -1.42658
-5.0 -5.0 -5.0 -5.0 5.0 5.0 5.0
julia> volume
30-element Vector{Int64}:
13.089969389957473
13.089969389957473
13.089969389957473
13.089969389957473
13.089969389957473
⋮
13.089969389957473
13.089969389957473
13.089969389957473
13.089969389957473
Peridynamics.n_points
— Functionn_points(body)
Returns the total number of points in a body.
Arguments
body::Body
:Body
.
Returns
n_points::Int
: The number of points in the body.
Examples
julia> body = Body(BBMaterial(), pos, vol)
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
julia> n_points(body)
1000
n_points(multibody_setup)
Returns the total number of points in a multibody setup.
Arguments
multibody_setup::MultibodySetup
:MultibodySetup
.
Returns
n_points::Int
: The sum of all points from all bodies in the multibody setup.
Examples
julia> ms = MultibodySetup(:a => body_a, :b => body_b)
2000-point MultibodySetup:
1000-point Body{BBMaterial{NoCorrection}} with name `a`
1000-point Body{BBMaterial{NoCorrection}} with name `b`
julia> n_points(ms)
2000
Preprocessing & simulation setup
Peridynamics.AbaqusMeshConverter.read_inp
— Functionread_inp(file::String)
Read Abaqus .inp-file and convert meshes to a point cloud 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
position::Matrix{Float64}
: Point position (midpoint of every element)volume::Vector{Float64}
: Point volume (volume of every element)point_sets
: Element sets defined in the .inp-file
Peridynamics.mpi_isroot
— Functionmpi_isroot()
Helper function that returns a bool indicating if a process is the MPI root process. It can be safely used even for multithreading simulations, as it is always true
if the package is started in a normal Julia environment which is not started by MPI.
Peridynamics.force_mpi_run!
— Functionforce_mpi_run!()
Helper function to force the usage of the MPI backend. After this function is called, all following simulations will use MPI.
Peridynamics.force_threads_run!
— Functionforce_threads_run!()
Helper function to force the usage of the multithreading backend. After this function is called, all following simulations will use multithreading.
Peridynamics.enable_mpi_timers!
— Functionenable_mpi_timers!()
Helper function to enable timers defined with the TimerOutputs
for simulations with the MPI backend. The results of the timers then will be exported into the specified path of a Job
. By default, not timers will be used with MPI simulations. It can be safely used with multithreading.
See also disable_mpi_timers!
.
Peridynamics.disable_mpi_timers!
— Functiondisable_mpi_timers!()
Helper function to disable timers defined with the TimerOutputs
for simulations with the MPI backend. It is mainly used to reset the behaviour after a call of enable_mpi_timers!
. It can be safely used with multithreading.
Peridynamics.enable_mpi_progress_bars!
— Functionenable_mpi_progress_bars!()
Helper function to enable progress bars with MPI simulations on a personal computer. After this function is called, progress bars are beeing shown with MPI simulations like with multithreading simulations. This behavior can be reset to default with reset_mpi_progress_bars!
.
Progress bars are by default disabled with MPI simulations, because they can really mess up with the output files produced by a HPC system. Therefore, a warning is shown as a reminder to reset this behaviour before submitting a job to a cluster!
Peridynamics.reset_mpi_progress_bars!
— Functionreset_mpi_progress_bars!()
After this function is called, progress bars are again disabled on MPI simulations (standard setting). This will reset the behavior after a call of enable_mpi_progress_bars!
.
Peridynamics.@mpitime
— Macro@mpitime expression
Time the expression
if the mpi rank is zero. Lowers to:
if mpi_isroot()
@time expression
else
expression
end
See also: mpi_isroot
.
Peridynamics.@mpiroot
— Macro@mpiroot [option] expression
Run the code if the mpi rank is zero. Lowers to something similar as:
if mpi_isroot()
expression
end
Options
:wait
: All MPI ranks will wait until the root rank finishes evaluatingexpression
.
See also: mpi_isroot
.
Solving
Peridynamics.VelocityVerlet
— TypeVelocityVerlet(; kwargs...)
Time integration solver for the Velocity Verlet algorithm. Specify either the number of steps or the time the simulation should cover.
Keywords
time::Real
: The total time the simulation will cover. If this keyword is specified, the keywordsteps
is no longer allowed. (optional)steps::Int
: Number of calculated time steps. If this keyword is specified, the keywordtime
is no longer allowed. (optional)stepsize::Real
: Manually specify the size of the time step. (optional)safety_factor::Real
: Safety factor for step size to ensure stability. (default:0.7
)
Keep in mind that manually specifying the critical time step is dangerous! If the specified time step is too high and the CFL condition no longer holds, the simulation will give wrong results and maybe crash!
Throws
- Errors if both
time
andsteps
are specified as keywords. - Errors if neither
time
norsteps
are specified as keywords. - Errors if
safety_factor < 0
orsafety_factor > 1
.
Example
julia> VelocityVerlet(steps=2000)
VelocityVerlet:
n_steps 2000
safety_factor 0.7
julia> VelocityVerlet(time=0.001)
VelocityVerlet:
end_time 0.001
safety_factor 0.7
julia> VelocityVerlet(steps=2000, stepsize=0.0001)
┌ Warning: stepsize specified! Please be sure that the CFD-condition holds!
└ @ Peridynamics ~/Code/Peridynamics.jl/src/time_solvers/velocity_verlet.jl:66
VelocityVerlet:
n_steps 2000
Δt 0.0001
safety_factor 0.7
Peridynamics.DynamicRelaxation
— TypeDynamicRelaxation(; kwargs...)
Time integration solver for the adaptive dynamic relaxation algorithm used for quasi-static simulations.
Keywords
steps::Int
: Number of calculated time steps. If this keyword is specified, the keywordtime
is no longer allowed.stepsize::Real
: Manually specify the size of the time step. (default:1.0
)damping_factor::Real
: Damping factor to increase the value in the mass matrix. (default:1.0
)
Throws
- Errors if
steps < 0
. - Errors if
stepsize < 0
. - Errors if
damping_factor < 0
.
Example
julia> DynamicRelaxation(steps=1000)
DynamicRelaxation:
n_steps 1000
Δt 1
Λ 1
julia> DynamicRelaxation(steps=1000, damping_factor=2)
DynamicRelaxation:
n_steps 1000
Δt 1
Λ 2
Peridynamics.Job
— TypeJob(spatial_setup, time_solver; kwargs...)
A type that contains all the information necessary for a peridynamic simulation. You can submit
a Job
to start the simulation.
Arguments
spatial_setup
: ABody
orMultibodySetup
.time_solver
:VelocityVerlet
orDynamicRelaxation
.
Keywords
path::String
: Path to store results. If it does not exist yet it will be created during the simulation. (optional)freq::Int
: Output frequency of result files. A output file will be written everyfreq
-th time step. (default:10
)fields
: Fields that should be exported to output files. Allowed keywords depend on the selected material model. Please look at the documentation of the material you specified when creating the body. (default:(:displacement, :damage)
)If
spatial_setup
is aBody
, thefields
keyword can be of the form:fields::Symbol
: A symbol specifying a single output field.fields::NTuple{N,Symbol} where N
: A Tuple specifying multiple output fields.fields::Vector{Symbol}
: A Vector specifying multiple output fields.
If
spatial_setup
is aMultibodySetup
, thefields
keyword can also be specified for every body separately:fields::Dict{Symbol,T}
: A Dictionary containing the fields separately for every body.T
is here every possible type of thefields
keyword that can be used for a single body.
If no keyword is specified when creating a Job
, then no files will be exported.
Example
julia> job = Job(multibody_setup, verlet_solver; path="my_results/sim1")
Job:
spatial_setup 25880-point MultibodySetup
time_solver VelocityVerlet(n_steps=2000, safety_factor=0.7)
options export_allowed=true, freq=10
Peridynamics.submit
— Functionsubmit(job::Job; quiet=false)
Run the simulation by submitting the job.
Arguments
job::Job
: Job that contains all defined parameters and conditions.
Keywords
quiet::Bool
: Iftrue
, no outputs are printed in the terminal. (default:false
)
Postprocessing
Peridynamics.VtkReader.read_vtk
— Functionread_vtk(file::AbstractString)
Read vtu or pvtu file containing simulation results of a time step.
Arguments
file::String
: Path to VTK file in vtu or pvtu format
Returns
Dict{String, VecOrMat{Float64}}
: Simulation results as a dictionary
Examples
julia> read_vtk("results/fragmenting_cylinder/vtk/timestep_000520.pvtu")
Dict{Symbol, VecOrMat{Float64}} with 4 entries:
:position => [0.0263309 0.027315 … 0.0293543 0.030339; 0.000292969 0.000294475…
:displacement => [0.00583334 0.00581883 … 0.00585909 0.00584271; -0.000162852 -0.0…
:damage => [0.616071, 0.569343, 0.528571, 0.463415, 0.438776, 0.553571, 0.56…
:time => [9.69363e-5]
Peridynamics.process_each_export
— Functionprocess_each_export(f, vtk_path; kwargs...)
process_each_export(f, job; kwargs...)
A function for postprocessing every exported file. This function works with multithreading and MPI and determines the backend exactly like the submit
function.
Arguments
f::Function
: The processing function with signaturef(r0, r, id)
.vtk_path::AbstractString
: A path that should contain the export results of a simulation.job::Job
: A job object. The path of the VTK files will then be processed from the job options.
Keywords
serial::Bool
: Iftrue
, all results will be processed in the correct order of the time steps and on a single thread, cf. the MPI root rank.