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.
Keywords
dmgmodel::AbstractDamageModel: Damage model defining the fracture behavior. (default:CriticalStretch())
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.DM: A damage model type.
Allowed material parameters
When using material! on a Body with BBMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters
E::Float64: Young's modulus.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
In bond-based peridynamics, the Poisson's ratio is limited to 1/4 for 3D simulations. Therefore, only one additional elastic parameter is required. Optionally, the specification of a second keyword is allowed, if the parameter combination results in 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 point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_bonds::Vector{Int}: Number of intact bonds of each point.strain_energy_density::Vector{Float64}: Strain energy density of each point.
Peridynamics.DHBBMaterial — TypeDHBBMaterial()
DHBBMaterial{Correction}()A material type used to assign the material of a Body with the dual-horizon bond-based formulation of peridynamics.
Keywords
dmgmodel::AbstractDamageModel: Damage model defining the fracture behavior. (default:CriticalStretch())
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 = DHBBMaterial()
DHBBMaterial{NoCorrection}()
julia> mat = DHBBMaterial{EnergySurfaceCorrection}()
DHBBMaterial{EnergySurfaceCorrection}()DHBBMaterial{Correction}Material type for the dual-horizon bond-based peridynamics formulation.
Type Parameters
Correction: A correction algorithm type. See the constructor docs for more informations.DM: A damage model type.
Allowed material parameters
When using material! on a Body with DHBBMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters:
E::Float64: Young's modulus.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
In bond-based peridynamics, the Poisson's ratio is limited to 1/4 for 3D simulations. Therefore, only one additional elastic parameter is required. Optionally, the specification of a second keyword is allowed, if the parameter combination results in nu = 1/4.
Allowed export fields
When specifying the fields keyword of Job for a Body with DHBBMaterial, the following fields are allowed:
position::Matrix{Float64}: Position of each point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_bonds::Vector{Int}: Number of intact bonds of each point.strain_energy_density::Vector{Float64}: Strain energy density of each point.
Peridynamics.OSBMaterial — TypeOSBMaterial(; kernel, dmgmodel)
OSBMaterial{Correction}(; kernel, dmgmodel)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.
Keywords
kernel::Function: Kernel function used for weighting the interactions between points.
(default:linear_kernel)dmgmodel::AbstractDamageModel: Damage model defining the damage behavior.
(default:CriticalStretch())
Examples
julia> mat = OSBMaterial()
OSBMaterial{NoCorrection}(dmgmodel=CriticalStretch())
julia> mat = OSBMaterial{EnergySurfaceCorrection}()
OSBMaterial{EnergySurfaceCorrection}(dmgmodel=CriticalStretch())OSBMaterial{Correction,K,DM}Material type for the ordinary state-based peridynamics formulation.
Type Parameters
Correction: A correction algorithm type. See the constructor docs for more informations.K: A kernel function type. See the constructor docs for more informations.DM: A damage model type. See the constructor docs for more informations.
Fields
kernel::Function: Kernel function used for weighting the interactions between points.dmgmodel::AbstractDamageModel: Damage model defining the damage behavior. See the constructor docs for more informations.
Allowed material parameters
When using material! on a Body with OSBMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters:
E::Float64: Young's modulus.nu::Float64: Poisson's ratio.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
Note that exactly two elastic parameters are required to specify a material. Please choose two out of the six allowed elastic parameters.
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 point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_bonds::Vector{Int}: Number of intact bonds of each point.strain_energy_density::Vector{Float64}: Strain energy density of each point.
Peridynamics.CMaterial — TypeCMaterial(; kernel, model, zem, dmgmodel, maxdmg)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
kernel::Function: Kernel function used for weighting the interactions between points.
(default:linear_kernel)
The following kernels can be used:model::AbstractConstitutiveModel: Constitutive model defining the material behavior.
(default:SaintVenantKirchhoff())
The following models can be used:zem::AbstractZEMStabilization: Algorithm of zero-energy mode stabilization.
(default:ZEMSilling)
The following algorithms can be used:dmgmodel::AbstractDamageModel: Damage model defining the damage behavior.
(default:CriticalStretch)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 containNaNvalues.
(default:0.85)
This formulation is known to be not suitable for fracture simulations without stabilization of the zero-energy modes. Therefore be careful when doing fracture simulations and try out different parameters for maxdmg and zem.
Examples
julia> mat = CMaterial()
CMaterial{SaintVenantKirchhoff, ZEMSilling, typeof(linear_kernel), CriticalStretch}(maxdmg=0.85)CMaterial{CM,ZEM,K,DM}Material type for the local continuum consistent (correspondence) formulation of non-ordinary state-based peridynamics.
Type Parameters
CM: A constitutive model type. See the constructor docs for more informations.ZEM: A zero-energy mode stabilization type. See the constructor docs for more informations.K: A kernel function type. See the constructor docs for more informations.DM: A damage model type. See the constructor docs for more informations.
Fields
kernel::Function: Kernel function used for weighting the interactions between points. See the constructor docs for more informations.model::AbstractConstitutiveModel: Constitutive model defining the material behavior. See the constructor docs for more informations.zem::AbstractZEMStabilization: Zero-energy mode stabilization. See the constructor docs for more informations.dmgmodel::AbstractDamageModel: Damage model defining the damage behavior. See the constructor docs for more informations.maxdmg::Float64: Maximum value of damage a point is allowed to obtain. See the constructor docs for more informations.
Allowed material parameters
When using material! on a Body with CMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters:
E::Float64: Young's modulus.nu::Float64: Poisson's ratio.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
Note that exactly two elastic parameters are required to specify a material. Please choose two out of the six allowed elastic parameters.
Allowed export fields
When specifying the fields keyword of Job for a Body with CMaterial, the following fields are allowed:
position::Matrix{Float64}: Position of each point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_bonds::Vector{Int}: Number of intact bonds of each point.cauchy_stress::Matrix{Float64}: Cauchy stress tensor of each point.von_mises_stress::Vector{Float64}: Von Mises stress of each point.hydrostatic_stress::Vector{Float64}: Hydrostatic stress of each point.strain_energy_density::Vector{Float64}: Strain energy density of each point.
Peridynamics.CRMaterial — TypeCRMaterial(; kernel, model, zem, dmgmodel, maxdmg)The same as the CMaterial but with rotation of the stress tensor for large deformation simulations, therefore not all models are supported.
Supported models:
SaintVenantKirchhoffLinearElastic
Please take a look at the CMaterial docs for more information about the material!
Peridynamics.BACMaterial — TypeBACMaterial(; kernel, model, dmgmodel, maxdmg)A material type used to assign the material of a Body with the bond-associated correspondence formulation of Chen and Spencer (2019).
Keywords
kernel::Function: Kernel function used for weighting the interactions between points. (default:linear_kernel)model::AbstractConstitutiveModel: Constitutive model defining the material behavior. (default:SaintVenantKirchhoff())dmgmodel::AbstractDamageModel: Damage model defining the fracture behavior. (default:CriticalStretch())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 containNaNvalues. (default:0.85)
Examples
julia> mat = BACMaterial()
BACMaterial{SaintVenantKirchhoff, typeof(linear_kernel), CriticalStretch}()BACMaterial{CM,K,DM}Material type for the bond-associated correspondence formulation of Chen and Spencer (2019).
Type Parameters
CM: A constitutive model type. See the constructor docs for more informations.K: A kernel function type. See the constructor docs for more informations.DM: A damage model type.
Fields
kernel::Function: Kernel function used for weighting the interactions between points.model::AbstractConstitutiveModel: Constitutive model defining the material behavior.dmgmodel::AbstractDamageModel: Damage model defining the fracture behavior.maxdmg::Float64: Maximum value of damage a point is allowed to obtain. See the constructor docs for more informations.
Allowed material parameters
When using material! on a Body with BACMaterial, then the following parameters are allowed:
horizon::Float64: Radius of point interactions.rho::Float64: Density.E::Float64: Young's modulus.nu::Float64: Poisson's ratio.Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
Allowed export fields
When specifying the fields keyword of Job for a Body with BACMaterial, the following fields are allowed:
position::Matrix{Float64}: Position of each point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_bonds::Vector{Int}: Number of intact bonds of each point.
Peridynamics.CKIMaterial — TypeCKIMaterial(; dmgmodel)A material type used to assign the material of a Body with the continuum-kinematics-inspired peridynamics formulation.
Keywords
dmgmodel::AbstractDamageModel: Damage model defining the damage behavior.
(default:CriticalStretch())
Examples
julia> mat = CKIMaterial()
CKIMaterial(dmgmodel=CriticalStretch())CKIMaterial{DM}Material type for the continuum-kinematics-inspired peridynamics framework.
Type Parameters
DM: A damage model type. See the constructor docs for more informations.
Fields
dmgmodel::AbstractDamageModel: Damage model defining the damage behavior. See the constructor docs for more informations.
Allowed material parameters
When using material! on a Body with CKIMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters:
E::Float64: Young's modulus.nu::Float64: Poisson's ratio.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
Interaction parameters:
C1::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 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).
Note that exactly two elastic parameters are required to specify a material. Please choose two out of the six allowed elastic parameters.
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 point.displacement::Matrix{Float64}: Displacement of each point.velocity::Matrix{Float64}: Velocity of each point.velocity_half::Matrix{Float64}: Velocity parameter for Verlet time solver.acceleration::Matrix{Float64}: Acceleration of each point.b_int::Matrix{Float64}: Internal force density of each point.b_ext::Matrix{Float64}: External force density of each point.damage::Vector{Float64}: Damage of each point.n_active_one_nis::Vector{Int}: Number of intact one-neighbor interactions of each point.
Peridynamics.RKCMaterial — TypeRKCMaterial(; kernel, model, dmgmodel, maxdmg, reprkernel, regfactor)A material type used to assign the material of a Body with a reproducing kernel peridynamics (correspondence) formulation with bond-associated quadrature integration at the center of the bonds.
Keywords
kernel::Function: Kernel function used for weighting the interactions between points.
(default:cubic_b_spline_kernel)
The following kernels can be used:model::AbstractConstitutiveModel: Constitutive model defining the material behavior.
(default:SaintVenantKirchhoff())
The following models can be used:dmgmodel::AbstractDamageModel: Damage model defining the damage behavior.
(default:CriticalStretch)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 containNaNvalues.
(default:1.0)reprkernel::Symbol: A kernel function used for the reproducing kernel approximation of the moment matrix. This kernel is used to calculate the moment matrix and the gradient weights, which are used to approximate the deformation gradient.
(default::C1)
The following kernels can be used::C1: Linear reproducing kernel basis [x, y, z] with first-order accuracy for correspondence formulation. This is the default kernel.:RK1: First-order reproducing kernel basis [1, x, y, z] with constant term for enhanced stability in uniform deformation fields.:RK2: Second-order reproducing kernel basis [1, x, y, z, x², y², z²] with diagonal quadratic terms for improved accuracy in curved deformation fields.:PD2: Second-order peridynamic basis [x, y, z, x², xy, xz, y², yz, z²] with full quadratic terms but without constant term.
regfactor::Float64: A regularization factor used to stabilize the inversion of the moment matrix. The product ofregfactor * δ^3is used for the regularization withinvreg. The regularization helps to ensure numerical stability by mitigating issues such as ill-conditioning or singularity in the matrix operations.
(default:1e-13)
Examples
julia> mat = RKCMaterial()
RKCMaterial{SaintVenantKirchhoff, typeof(linear_kernel), CriticalStretch}(maxdmg=1.0)RKCMaterial{CM,K,DM}Material type for a reproducing kernel peridynamics (correspondence) formulation with bond-associated quadrature integration at the center of the bonds.
Type Parameters
CM: A constitutive model type. See the constructor docs for more informations.K: A kernel function type. See the constructor docs for more informations.DM: A damage model type. See the constructor docs for more informations.
Fields
kernel::Function: Kernel function used for weighting the interactions between points. See the constructor docs for more informations.model::AbstractConstitutiveModel: Constitutive model defining the material behavior. See the constructor docs for more informations.dmgmodel::AbstractDamageModel: Damage model defining the damage behavior. See the constructor docs for more informations.maxdmg::Float64: Maximum value of damage a point is allowed to obtain. See the constructor docs for more informations.reprkernel::Symbol: A kernel function used for the reproducing kernel approximation. See the constructor docs for more informations.regfactor::Float64: Regularization factor. See the constructor docs for more informations.
Allowed material parameters
When using material! on a Body with RKCMaterial, then the following parameters are allowed: Material parameters:
horizon::Float64: Radius of point interactionsrho::Float64: Density
Elastic parameters:
E::Float64: Young's modulusnu::Float64: Poisson's ratioG::Float64: Shear modulusK::Float64: Bulk moduluslambda::Float64: 1st Lamé parametermu::Float64: 2nd Lamé parameter
Fracture parameters:
Gc::Float64: Critical energy release rateepsilon_c::Float64: Critical strain
Note that exactly two elastic parameters are required to specify a material. Please choose two out of the six allowed elastic parameters.
Allowed export fields
When specifying the fields keyword of Job for a Body with RKCMaterial, 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 pointcauchy_stress::Matrix{Float64}: Cauchy stress tensor of each pointvon_mises_stress::Vector{Float64}: Von Mises stress of each pointhydrostatic_stress::Vector{Float64}: Hydrostatic stress of each point.strain_energy_density::Vector{Float64}: Strain energy density of each point.
Peridynamics.RKCRMaterial — TypeRKCRMaterial(; kernel, model, dmgmodel, maxdmg, reprkernel, regfactor)The same as the RKCMaterial but with rotation of the stress tensor for large deformation simulations, therefore not all models are supported.
Supported models:
SaintVenantKirchhoffLinearElastic
Please take a look at the RKCMaterial docs for more information about the material!
System or material related types
Peridynamics.CriticalStretch — TypeCriticalStretchA damage model based on the stretch of the bond. The bond is considered to be broken if the stretch exceeds a critical value. The critical value can be defined via the fracture energy Gc or the critical stretch εc using the material! function. The damage model is defined globally for the whole body as part of the material.
Peridynamics.NoCorrection — TypeNoCorrectionA 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 — TypeEnergySurfaceCorrectionA 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.
Peridynamics.ZEMSilling — TypeZEMSilling(; Cs)Zero-energy mode stabilization algorithm of Silling (2017). This is necessary for the correspondence formulation to stabilize the zero-energy modes. See also CMaterial on how to use this stabilization algorithm.
Keywords
Cs::Real: Stabilization factor. (default:100.0)
Peridynamics.ZEMWan — TypeZEMWan()Zero-energy mode stabilization algorithm of Wan et al. (2019), which is an improvement to Silling's algorithm that does not require a stabilization parameter. See also CMaterial on how to use this stabilization algorithm.
Peridynamics.LinearElastic — TypeLinearElasticLinear elastic constitutive model that can be specified when using a CMaterial and BACMaterial.
The strain energy density $\Psi$ is given by
\[\Psi = \frac{1}{2} \lambda \, \mathrm{tr}(\boldsymbol{E})^2 + \mu \, \mathrm{tr}(\boldsymbol{E} \cdot \boldsymbol{E}) \; ,\]
with the first and second Lamé parameters $\lambda$ and $\mu$, and the Green-Lagrange strain tensor $\boldsymbol{E}$
\[\boldsymbol{E} = \frac{1}{2} \left( \boldsymbol{F}^{\top} \boldsymbol{F} - \boldsymbol{I} \right) \; .\]
The first Piola-Kirchhoff stress $\boldsymbol{P}$ is given by
\[\begin{aligned} \boldsymbol{S} &= \mathbb{C} : \boldsymbol{E} \; , \\ \boldsymbol{P} &= \boldsymbol{F} \, \boldsymbol{S} \; , \end{aligned}\]
with the deformation gradient $\boldsymbol{F}$, the elastic stiffness tensor $\mathbb{C}$, and the second Piola-Kirchhoff stress $\boldsymbol{S}$.
Peridynamics.NeoHooke — TypeNeoHookeCompressible Neo-Hookean hyperelastic constitutive model that can be specified when using a CMaterial and BACMaterial.
The strain energy density $\Psi$ is given by
\[\Psi = \frac{1}{2} \mu \left( I_1 - 3 \right) - \mu \log(J) + \frac{1}{2} \lambda \log(J)^2 \; ,\]
with the first invariant $I_1 = \mathrm{tr}(\boldsymbol{C})$ of the right Cauchy-Green deformation tensor $\boldsymbol{C} = \boldsymbol{F}^{\top} \boldsymbol{F}$, the Jacobian $J = \mathrm{det}(\boldsymbol{F})$, and the first and second Lamé parameters $\lambda$ and $\mu$.
The first Piola-Kirchhoff stress $\boldsymbol{P}$ is given by
\[\begin{aligned} \boldsymbol{C} &= \boldsymbol{F}^{\top} \boldsymbol{F} \; , \\ \boldsymbol{S} &= \mu \left( \boldsymbol{I} - \boldsymbol{C}^{-1} \right) + \lambda \log(J) \boldsymbol{C}^{-1} \; , \\ \boldsymbol{P} &= \boldsymbol{F} \, \boldsymbol{S} \; , \end{aligned}\]
with the deformation gradient $\boldsymbol{F}$ and the second Piola-Kirchhoff stress $\boldsymbol{S}$.
Reference
Treloar, L. R. G. (1943). "The elasticity of a network of long-chain molecules—II." Transactions of the Faraday Society, 39, 241–246. DOI: 10.1039/TF9433900241
Peridynamics.NeoHookePenalty — TypeNeoHookePenaltyCompressible Neo-Hookean hyperelastic model with a penalty-type volumetric formulation, suitable for modeling rubber-like and biological materials. Can be specified when using a CMaterial and BACMaterial.
The strain energy density $\Psi$ is given by
\[\Psi = \frac{1}{2} G \left( \bar{I}_1 - 3 \right) + \frac{K}{8} \left( J^2 + J^{-2} - 2 \right) \; ,\]
with the modified first invariant $\bar{I}_1 = I_1 J^{-2/3}$ where $I_1 = \mathrm{tr}(\boldsymbol{C})$ is the first invariant of the right Cauchy-Green deformation tensor $\boldsymbol{C} = \boldsymbol{F}^{\top} \boldsymbol{F}$, the Jacobian $J = \mathrm{det}(\boldsymbol{F})$, the shear modulus $G$, and the bulk modulus $K$.
The first Piola-Kirchhoff stress $\boldsymbol{P}$ is given by
\[\begin{aligned} \boldsymbol{C} &= \boldsymbol{F}^{\top} \boldsymbol{F} \; , \\ \boldsymbol{S} &= G \left( \boldsymbol{I} - \frac{1}{3} \mathrm{tr}(\boldsymbol{C}) \boldsymbol{C}^{-1} \right) J^{-\frac{2}{3}} + \frac{K}{4} \left( J^2 - J^{-2} \right) \boldsymbol{C}^{-1} \; , \\ \boldsymbol{P} &= \boldsymbol{F} \, \boldsymbol{S} \; , \end{aligned}\]
with the deformation gradient $\boldsymbol{F}$ and the second Piola-Kirchhoff stress $\boldsymbol{S}$.
This model uses a penalty-type volumetric term $\Psi_{\mathrm{vol}} = \frac{K}{8}(J^2 + J^{-2} - 2)$, which is computationally efficient and widely used in commercial finite element codes for nearly-incompressible materials. The term penalizes volume changes from the reference configuration ($J = 1$).
This differs from other volumetric formulations such as:
- Standard Neo-Hookean (logarithmic): $\Psi_{\mathrm{vol}} = -\mu \ln J + \frac{\lambda}{2}\ln^2 J$
- Simo-Miehe (polyconvex): $\Psi_{\mathrm{vol}} = \frac{K}{4}(J^2 - 1 - 2\ln J)$
Error handling
If the Jacobian $J$ is smaller than the machine precision eps() or a NaN, the strain energy density and first Piola-Kirchhoff stress tensor are defined as zero: $\Psi = 0$ and $\boldsymbol{P} = \boldsymbol{0}$.
Peridynamics.SaintVenantKirchhoff — TypeSaintVenantKirchhoffSaint-Venant-Kirchhoff constitutive model that can be specified when using a CMaterial and BACMaterial.
The strain energy density $\Psi$ is given by
\[\Psi = \frac{1}{2} \lambda \, \mathrm{tr}(\boldsymbol{E})^2 + \mu \, \mathrm{tr}(\boldsymbol{E} \cdot \boldsymbol{E}) \; ,\]
with the first and second Lamé parameters $\lambda$ and $\mu$, and the Green-Lagrange strain tensor
\[\boldsymbol{E} = \frac{1}{2} \left( \boldsymbol{F}^{\top} \boldsymbol{F} - \boldsymbol{I} \right) \; .\]
The first Piola-Kirchhoff stress $\boldsymbol{P}$ is given by
\[\begin{aligned} \boldsymbol{S} &= \lambda \, \mathrm{tr}(\boldsymbol{E}) \, \boldsymbol{I} + 2 \mu \boldsymbol{E} \; , \\ \boldsymbol{P} &= \boldsymbol{F} \, \boldsymbol{S} \; , \end{aligned}\]
with the deformation gradient $\boldsymbol{F}$ and the second Piola-Kirchhoff stress $\boldsymbol{S}$.
This model is equivalent to the LinearElastic model, both using the same strain energy density function.
Peridynamics.const_one_kernel — Functionconst_one_kernel(δ, L)A kernel function $\omega$ (also called influence function) used for weighting the bonds in a family. The kernel function is simply defined as a constant value 1:
\[\omega = 1\]
Peridynamics.linear_kernel — Functionlinear_kernel(δ, L)A linear kernel function $\omega$ (also called influence function) used for weighting the bonds in a family. The kernel function is defined as
\[\omega = \frac{\delta}{L} \; ,\]
with the horizon $\delta$ and the initial bond length $L$.
Peridynamics.cubic_b_spline_kernel — Functioncubic_b_spline_kernel(δ, L)A cubic B-spline kernel function $\omega$ used for weighting the bonds in a family. The kernel function is defined as
\[\begin{aligned} \xi &= \frac{L}{\delta} \; , \\ \omega &= \left\{ \begin{array}{ll} \frac{2}{3} - 4 \xi^2 + 4 \xi^3 & \quad \text{if} \; 0 < \xi \leq 0.5 \; , \\[3pt] \frac{4}{3} - 4 \xi + 4 \xi^2 - \frac{4}{3} \xi^3 & \quad \text{if} \; 0.5 < \xi \leq 1 \; , \\[3pt] 0 & \quad \text{else} \; , \end{array} \right. \end{aligned}\]
with the horizon $\delta$ and the initial bond length $L$.
Peridynamics.cubic_b_spline_kernel_norm — Functioncubic_b_spline_kernel_norm(δ, L)A cubic B-spline kernel function $\omega$ used for weighting the bonds in a family. The kernel function is defined as
\[\begin{aligned} \xi &= \frac{L}{\delta} \; , \\ \omega &= \frac{8}{\pi \, \delta^3} \cdot \left\{ \begin{array}{ll} 1 - 6 \xi^2 + 6 \xi^3 & \quad \text{if} \; 0 < \xi \leq 0.5 \; , \\[3pt] 2 (1 - \xi)^3 & \quad \text{if} \; 0.5 < \xi \leq 1 \; , \\[3pt] 0 & \quad \text{else} \; , \end{array} \right. \end{aligned}\]
with the horizon $\delta$ and the initial bond length $L$. This kernel is properly normalized to satisfy the condition $\int_{\mathcal{H}(X)} \omega(\Delta X) dV' = 1$.
Discretization
Peridynamics.Body — TypeBody(material, position, volume)
Body(material, inp_file)Construct a Body for a peridynamics simulation.
Arguments
material::AbstractMaterial: The material which is defined for the whole body. Available material models:BBMaterial: Bond-based peridynamics.DHBBMaterial: Dual-horizon bond-based peridynamics.OSBMaterial: Ordinary state-based peridynamics, also called linear peridynamic solid (LPS).CMaterial: Correspondence formulation.CRMaterial: Correspondence formulation with stress rotation for objectivity enforcement.RKCMaterial: Reproducing kernel peridynamics with bond-associated higher-order integration.RKCRMaterial: Reproducing kernel peridynamics with bond-associated higher-order integration with stress rotation for objectivity enforcement.BACMaterial: Bond-associated correspondence formulation of Chen and Spencer.CKIMaterial: Continuum-kinematics-inspired peridynamics.
position::AbstractMatrix: A3×nmatrix with the point position of thenpoints.volume::AbstractVector: A vector with the volume of each point.inp_file::AbstractString: An Abaqus input file containing meshes, imported withread_inp.
Throws
- Error if the number of points is not larger than zero.
- Error if
positionis not a3×nmatrix and has the same length asvolume. - Error if
positionorvolumecontainNaNvalues.
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 model.PointParameters <: 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_pointsmatrix 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 different point parameter instances of the body. Each point can have its ownPointParametersinstance.params_map::Vector{Int}: A vector that maps each point index to a parameter instance inpoint_params.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.data_bcs::Vector{DataBC}: A vector with data boundary conditions.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
- Error 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:Bodywhere 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 thepositionandvolumeofbody.fun::Function: Function for filtering points. This function accepts only one positional argument and will be used in afindallcall. Depending on the argument name, a different input will be processed:x: The function will receive the x-coordinate of each point inpositionofbody:points = findall(fun, @view(position[1, :]))y: The function will receive the y-coordinate of each point inpositionofbody:points = findall(fun, @view(position[2, :]))z: The function will receive the z-coordinate of each point inpositionofbody:points = findall(fun, @view(position[3, :]))p: The function will receive the a vector containing each dimension of each point inpositionofbody:points = findall(fun, eachcol(position))
Throws
- Error if a point set with the same
set_namealready exists. - Error if
pointsare not in bounds withpositionandvolumeof 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)Return 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.no_failure! — Functionno_failure!(body::AbstractBody, set_name::Symbol)
no_failure!(body::AbstractBody)Disallow failure for all points of the point set set_name of the body. If no set_name is specified, failure is prohibited for the whole body.
Arguments
body::AbstractBody:Bodyfor which failure is prohibited.set_name::Symbol: The name of a point set of this body.
The function material! sets failure permissions due to the provided input parameters, so if it is used afterwards, previously set failure prohibitions might be overwritten!
Throws
- Error if the body does not contain a set with
set_name.
Examples
julia> no_failure!(body)
julia> body
1000-point Body{BBMaterial{NoCorrection}}:
1 point set(s):
1000-point set `all_points`
1000 points with failure prohibitedPeridynamics.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:
Material parameters:
horizon::Float64: Radius of point interactions.rho::Float64: Density.
Elastic parameters:
E::Float64: Young's modulus.nu::Float64: Poisson's ratio.G::Float64: Shear modulus.K::Float64: Bulk modulus.lambda::Float64: 1st Lamé parameter.mu::Float64: 2nd Lamé parameter.
Fracture parameters:
Gc::Float64: Critical energy release rate.epsilon_c::Float64: Critical strain.
To enable fracture in a simulation, define one of the allowed fracture parameters. If none are defined, fracture is disabled.
Throws
- Error 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.7Peridynamics.velocity_bc! — Functionvelocity_bc!(fun, body, set_name, dim)Specify velocity boundary condition 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 timetat 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_nameand receives the reference position of a point asSVector{3}and the current timetat every time step. This makes it possible to specify conditions that also depend on the position of a point.
body::AbstractBody:Bodythe 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:
:xor1 - y-direction:
:yor2 - z-direction:
:zor3
- x-direction:
Throws
- Error if the body does not contain a set with
set_name. - Error if the direction is not correctly specified.
- Error 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=2Peridynamics.velocity_ic! — Functionvelocity_ic!(body, set_name, dim, value)
velocity_ic!(fun, body, set_name, dim)Specify velocity initial condition 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:Bodythe 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:
:xor1 - y-direction:
:yor2 - z-direction:
:zor3
- 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 positionpof a point asSVector{3}.
Throws
- Error if the body does not contain a set with
set_name. - Error if the direction is not correctly specified.
- Error 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=1Peridynamics.forcedensity_bc! — Functionforcedensity_bc!(fun, body, set, dim)Specify force density boundary condition 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 timetat 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_nameand receives the reference position of a point asSVector{3}and the current timetat every time step. This makes it possible to specify conditions that also depend on the position of a point.
body::AbstractBody:Bodythe 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:
:xor1 - y-direction:
:yor2 - z-direction:
:zor3
- x-direction:
Throws
- Error if the body does not contain a set with
set_name. - Error if the direction is not correctly specified.
- Error 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=2Peridynamics.precrack! — Functionprecrack!(body, set_a, set_b; update_dmg=true)Create 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
- Error if the body does not contain sets with name
set_aandset_b. - Error 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...)Define 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_aand a point in bodyname_body_bis 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
- Error if
multibody_setupdoes not contain bodies with namename_body_aandname_body_b. - Error if the keyword
radiusis not specified orradius ≤ 0. - Error 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...)Create 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_pointsmatrix 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
8Peridynamics.uniform_sphere — Functionuniform_sphere(diameter, ΔX0; kwargs...)Create 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_pointsmatrix 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
8Peridynamics.uniform_cylinder — Functionuniform_cylinder(diameter::Real, height::Real, ΔX0::Real; kwargs...)Create 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_pointsmatrix 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
8Peridynamics.round_sphere — Functionround_sphere(diameter, ΔX0; kwargs...)Create 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_pointsmatrix 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.311091676163475Peridynamics.round_cylinder — Functionround_cylinder(diameter::Real, height::Real, ΔX0::Real; kwargs...)Create 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_pointsmatrix 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.089969389957473Peridynamics.n_points — Functionn_points(body)Return 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)
1000n_points(multibody_setup)Return 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)
2000Preprocessing & 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!.
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 expressionTime the expression if the mpi rank is zero. Lowers to:
if mpi_isroot()
@time expression
else
expression
endSee also: mpi_isroot.
Peridynamics.@mpiroot — Macro@mpiroot [option] expressionRun the code if the mpi rank is zero. Lowers to something similar as:
if mpi_isroot()
expression
endOptions
: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 keywordstepsis no longer allowed. (optional)steps::Int: Number of calculated time steps. If this keyword is specified, the keywordtimeis 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
- Error if both
timeandstepsare specified as keywords. - Error if neither
timenorstepsare specified as keywords. - Error if
safety_factor < 0orsafety_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.7Peridynamics.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 keywordtimeis 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
- Error if
steps < 0. - Error if
stepsize < 0. - Error 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
Λ 2Peridynamics.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: ABodyorMultibodySetup.time_solver:VelocityVerletorDynamicRelaxation.
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_setupis aBody, thefieldskeyword 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_setupis aMultibodySetup, thefieldskeyword can also be specified for every body separately:fields::Dict{Symbol,T}: A Dictionary containing the fields separately for every body.Tis here every possible type of thefieldskeyword that can be used for a single body.
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=10Peridynamics.Study — TypeStudy(jobcreator::Function, setups::Vector{<:NamedTuple}; root::String)A structure for managing parameter studies with multiple peridynamic simulations. The Study type coordinates the execution of multiple simulation jobs with different parameter configurations, tracks their status, and logs all results to a central logfile.
If a logfile already exists at the study root (from a previous interrupted run), the constructor automatically reads it and initializes the sim_success field based on the recorded completion status. This enables seamless resumption of interrupted studies.
Arguments
jobcreator::Function: A function with signaturejobcreator(setup::NamedTuple, root::String)that creates and returns aJobobject. The function receives:setup: ANamedTuplecontaining the parameters for one simulationroot: The root directory path for the study (to construct individual job paths)
setups::Vector{<:NamedTuple}: A vector of parameter configurations. Each element must be aNamedTuplewith the same field names.
Keywords
root::String: Root directory path where all simulation results and the study logfile will be stored. This directory and all job subdirectories will be created duringsubmit!.
Fields
jobcreator::Function: The job creation functionsetups: Vector of parameter configurationsjobs: Vector of createdJobobjectsjobpaths::Vector{String}: Paths to individual job directories (must be unique)root::String: Root directory for the studylogfile::String: Path to the central study logfilesim_success::Vector{Bool}: Status flags indicating successful completion of each job
Throws
ArgumentError: Ifsetupsis emptyArgumentError: If setups don't have consistent field namesArgumentError: If job paths are not uniqueArgumentError: Ifjobcreatorfails to create a valid job
Example
# some function that creates a Job from a parameter setup
function create_job(setup::NamedTuple, root::String)
body = Body(BBMaterial(), uniform_box(1.0, 1.0, 1.0, 0.1))
material!(body, horizon=0.3, E=setup.E, rho=1000, Gc=100)
velocity_ic!(body, :all_points, :x, setup.velocity)
solver = VelocityVerlet(steps=1000)
# create a unique path for this job based on parameters
path = joinpath(root, "sim_E$(setup.E)_v$(setup.velocity)")
return Job(body, solver; path=path, freq=10)
end
setups = [
(; E=1e9, velocity=1.0),
(; E=2e9, velocity=1.0),
(; E=1e9, velocity=2.0),
]
study = Study(create_job, setups; root="my_parameter_study")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)
Peridynamics.submit! — Functionsubmit!(study::Study; kwargs...)Submit and execute all simulation jobs in a parameter study. Jobs are run sequentially, and each job can utilize MPI or multithreading as configured. If a job fails, the error is logged and execution continues with the remaining jobs.
This function creates the study root directory and logfile, then submits each job using submit. All simulation results and individual job logs are stored in their respective job directories, while a central study logfile tracks the overall progress and status of all simulations.
Resuming Interrupted Runs
If the study logfile already exists (from a previous run), submit! will:
- Read the logfile and refresh
study.sim_successto detect previously completed jobs - Skip jobs that already completed successfully (marked with "skipped" in the logfile)
- Only execute jobs that failed or were not yet started
- Append a "RESUMED" marker with timestamp to the logfile
This allows you to safely resume a study after an interruption (crash, timeout, manual stop) without re-running successful simulations. Simply call submit!(study) again with the same study configuration.
Arguments
study::Study: The study containing all jobs to execute
Keywords
quiet::Bool: Iftrue, suppress terminal output for individual jobs (default:false)- Additional keywords are passed to
submitfor each job
Behavior
- Jobs execute sequentially (one after another)
- Each job inherits MPI/threading configuration from the Julia session
- Failed jobs don't stop execution of remaining jobs
study.sim_successis updated with the status of each job- A study logfile is created at
study.logfilewith:- Study header with timestamp
- Parameters and status for each simulation
- Execution time for successful jobs
- Individual job logs are written to their respective directories
Example
study = Study(create_job, setups; root="my_study")
submit!(study; quiet=true)
# If the process was interrupted, simply call submit! again to resume:
# It will skip completed jobs and only run remaining ones
study = Study(create_job, setups; root="my_study") # Reads existing logfile
submit!(study; quiet=true) # Resumes from where it left off
# Check which simulations succeeded
successful_indices = findall(study.sim_success)
println("Successful simulations: ", successful_indices)
# Read the study logfile
logcontent = read(study.logfile, String)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.
Peridynamics.process_each_job — Functionprocess_each_job(f::Function, study::Study, default_result::NamedTuple)Apply a processing function to each successfully completed job in a parameter study. This function iterates through all jobs in the study and applies the user-defined processing function f to jobs that completed successfully. Failed jobs or jobs where the processing function errors will use the default_result instead.
Before processing, this function automatically refreshes study.sim_success from the logfile if it exists. This ensures that jobs completed in a previous interrupted run are properly detected and processed, even if the current Julia session doesn't reflect their completion status yet.
Arguments
f::Function: A processing function with signaturef(job::Job, setup::NamedTuple)that returns aNamedTuplecontaining the processed results. The function receives:job: TheJobobject for the simulationsetup: The parameter configurationNamedTuplefor this job
study::Study: The study containing the jobs to processdefault_result::NamedTuple: The default result to use when a job failed or when the processing function throws an error
Returns
Vector{<:NamedTuple}: A vector of results with the same length as the number of jobs in the study. Each element is either the result from applyingfor thedefault_resultfor failed/errored cases.
Behavior
- Refreshes job completion status from logfile (if exists) before processing
- Only processes jobs where
study.sim_success[i] == true - Failed jobs automatically receive
default_result(warning logged) - If processing function
fthrows an error, that job receivesdefault_result(error logged) - Processing is sequential (one job at a time)
- All jobs in the study will have a corresponding entry in the results vector
Example
# After running a parameter study (possibly interrupted and resumed)
study = Study(create_job, setups; root="my_study")
# Status is automatically refreshed from logfile in constructor
# Define a function to extract maximum displacement from results
function extract_max_displacement(job::Job, setup::NamedTuple)
# Calculate maximum displacement
max_u = ...
return (; E=setup.E, velocity=setup.velocity, max_displacement=max_u)
end
default = (; E=0.0, velocity=0.0, max_displacement=NaN)
# This will process all successfully completed jobs, even if from a previous run
results = process_each_job(extract_max_displacement, study, default)
# Filter successful results
successful_results = [r for r in results if !isnan(r.max_displacement)]