DDD.AbstractCrystalStructType
abstract type AbstractCrystalStruct end
struct BCC <: AbstractCrystalStruct end
struct FCC <: AbstractCrystalStruct end
struct HCP <: AbstractCrystalStruct end

Crystal structure types.

source
DDD.AbstractDistributionType
abstract type AbstractDistribution end
struct Zeros <: AbstractDistribution end
struct Rand <: AbstractDistribution end
struct Randn <: AbstractDistribution end
struct Regular <: AbstractDistribution end

Spatial distributions for dislocation sources. These are used to automatically generate networks with a given distribution.

  • Zeros makes the network generation functions place the center of the generated dislocation loops at the origin. This can be used to generate a network and loops can be manually or pseudo-manually distributed in the domain.
  • Rand makes the network generation functions uniformly distribute the dislocations according to the range and buffer values in the dislocation loop structure.
  • Rand makes the network generation functions normally distribute the dislocations according to the range and buffer values in the dislocation loop structure.
  • Rand TBA, will regularly distribute dislocations according to the range, buffer and other args given to the dislocation network generator.
source
DDD.AbstractDlnSegType
abstract type AbstractDlnSeg end
struct segNone <: AbstractDlnSeg end    # Undefined segment
struct segEdge <: AbstractDlnSeg end    # Edge segment
struct segEdgeN <: AbstractDlnSeg end   # Edge segment
struct segScrew <: AbstractDlnSeg end   # Screw segment
struct segMixed <: AbstractDlnSeg end   # Mixed segment

These types are used to automatically generate segments out of Burgers vectors b, slip planes n, and/or line direction l.

  • segEdge: b ⟂ t,
  • segEdgeN: b ⟂ t and b ∥ n ,
  • segScrew: b ∥ t ,
  • segMixed: none of the above.
source
DDD.AbstractDlnStrType
abstract type AbstractDlnStr end
struct loopPrism <: AbstractDlnStr end
struct loopShear <: AbstractDlnStr end
const loopPure = Union{loopPrism,loopShear}
struct loopMixed <: AbstractDlnStr end
struct loopJog <: AbstractDlnStr end
struct loopKink <: AbstractDlnStr end
const loopImpure = Union{loopMixed,loopJog,loopKink}
const loopDefined = Union{loopPure,loopImpure}
struct loopDln <: AbstractDlnStr end

These types are used to automatically generate dislocation loops for simulation initialisation.

  • loopPrism: prismatic loops, their Burgers vectors are perpendicular to the their line direction. They are idealised loops that can be automatically generated as n-gons
  • loopShear: shear loops, their line direction goes through edge, screw and line segments as the loop goes round. They are idealised loops that can be automatically generated as n-gons
  • loopPure: idealised loops
  • loopMixed: loops with prismatic and shear character. They have to be hand-made or require a heuristic to automatically generate
  • loopDln: generic loop used for adding methods to Base functions
  • loopKink and loopJog are structures formed by colliding dislocations. They are not currently used
  • loopImpure: non-idealised loops
  • loopDefined: defined loop types
source
DDD.AbstractElementOrderType
abstract type AbstractElementOrder end
struct LinearElement <: AbstractElementOrder end

Finite element orders for dispatch.

source
DDD.AbstractIntegratorType
abstract type AbstractIntegrator end
struct AdaptiveEulerTrapezoid <: AbstractIntegrator end

Integrator types for dispatch.

source
DDD.AbstractMeshType
abstract type AbstractMesh end
abstract type AbstractRegularCuboidMesh <: AbstractMesh end
struct DispatchRegularCuboidMesh <: AbstractRegularCuboidMesh end

FE mesh types for dispatch.

source
DDD.AbstractMobilityType
abstract type AbstractMobility end
struct mobBCC <: AbstractMobility end
struct mobFCC <: AbstractMobility end
struct mobHCP <: AbstractMobility end

Types to dispatch different mobility functions.

  • mobBCC: used to dispatch the default BCC mobility function.
  • mobFCC: used to dispatch the default FCC mobility function.
  • mobHCP: used to dispatch the default HCP mobility function.
source
DDD.AbstractModelType
abstract type AbstractModel end
abstract type AbstractCantilever <: AbstractModel end
struct CantileverLoad <: AbstractCantilever end

Abstract types for dispatching different models.

source
DDD.AbstractShapeFunctionType
abstract type AbstractShapeFunction end
abstract type AbstractShapeFunction3D <: AbstractShapeFunction end
abstract type AbstractShapeFunction2D <: AbstractShapeFunction end
struct LinearQuadrangle3D <:AbstractShapeFunction3D end
struct LinearQuadrangle2D <:AbstractShapeFunction2D end

Shape function types for dispatch.

source
DDD.BoundariesType
struct Boundaries{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12}
    model::T1
    noExit::T2
    uGammaDln::T3
    tGammaDln::T4
    uDofsDln::T5
    tDofsDln::T6
    uGamma::T7
    tGamma::T8
    mGamma::T9
    uDofs::T10
    tDofs::T11
    mDofs::T12
    tK::T13
end

Stores the nodes and degrees of freedom upon which the different boundary conditions are applied.

Fields

  • noExit: Faces that are impenetrable to dislocations
  • uGammaDln: Nodes on which dislocation displacements are calculated
  • tGammaDln: Nodes on which dislocation tractions are calculated
  • uDofsDln: Degrees of freedom on which dislocation displacements are calculated
  • tDofsDln: Degrees of freedom on which dislocation tractions are calculated
  • uGamma: Nodes with displacement boundaries
  • tGamma: Nodes with traction boundaries
  • mGamma: Nodes with displacement and traction boundaries
  • uDofs: Degrees of freedom with specified displacements
  • tDofs: Degrees of feedom with specified tractions
  • mDofs0: Degrees of feedom with specified displacements and tractions
  • tK1: Stiffness matrix of traction degrees of freedom
source
DDD.BoundariesMethod
Boundaries(
    ::FEMParameters{T1,T2,T3,T4,T5} where {T1,T2,T3<:CantileverLoad,T4,T5},
    femMesh::RegularCuboidMesh; 
    kw...
)

Creates Boundaries for loading a hexahedral cantilever.

source
DDD.BoundaryNodeType
struct BoundaryNode{T1,T2}
    index::T1
    node::T2
end

Stores corresponding type, indices and node number of boundary nodes.

source
DDD.DislocationLoopType
struct DislocationLoop{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14}
    loopType::T1        # Loop type.
    numSides::T2        # Number of sides in the loop.
    nodeSide::T3        # Nodes per side of the loop.
    numLoops::T4        # Number of loops to generate when making the network.
    segLen::T5          # Segment lengths.
    slipSystemIdx::T6   # Slip system.
    label::T7           # Node labels.
    links::T8           # Links.
    slipPlane::T9       # Slip planes.
    bVec::T10           # Burgers vectors.
    coord::T11          # Coordinates.
    buffer::T12         # Buffer for distributions.
    range::T13          # Range for distributions.
    dist::T14           # Distribution.
end

Stores a dislocation loop.

source
DDD.DislocationLoopMethod
DislocationLoop(
    loopType::AbstractDlnStr,
    numSides,
    nodeSide,
    numLoops,
    segLen,
    slipSystem,
    _slipPlane,
    _bVec,
    label,
    buffer,
    range,
    dist,
)

Fallback for creating a generic DislocationLoop.

source
DDD.DislocationLoopMethod
DislocationLoop(
    loopType::loopImpure,
    numSides,
    nodeSide,
    numLoops,
    segLen,
    slipSystem,
    _slipPlane::AbstractArray,
    _bVec::AbstractArray,
    label::AbstractVector{nodeTypeDln},
    buffer,
    range,
    dist::AbstractDistribution,
)

Fallback DislocationLoop constructor for other as of yet unimplemented loopImpure.

source
DDD.DislocationLoopMethod
DislocationLoop(
    loopType::loopPure,
    numSides,
    nodeSide,
    numLoops,
    segLen,
    slipSystem,
    _slipPlane::AbstractArray,
    _bVec::AbstractArray,
    label::AbstractVector{nodeTypeDln},
    buffer,
    range,
    dist::AbstractDistribution,
)

Constructor for loopPure DislocationLoops.

Inputs

  • loopType::loopPure: can be either loopPrism() or loopShear() to make prismatic or shear loops.
  • numSides: number of sides in the loop
  • nodeSide: nodes per side of the loop
  • numLoops: number of loops to generate when making the dislocation network
  • segLen: length of each initial dislocation segment
  • slipSystem: slip system from SlipSystem the loop belongs to
  • _slipPlane: slip plane vector
  • _bVec: Burgers vector
  • label: node labels of type nodeTypeDln
  • buffer: buffer for increasing the spread of the generated dislocation loops in the network
  • range: range on which the loops will be distributed in the network
  • dist: distribution used to generate the network
source
DDD.DislocationLoopMethod
DislocationLoop(;
    loopType::AbstractDlnStr,
    numSides,
    nodeSide,
    numLoops,
    segLen,
    slipSystem,
    _slipPlane,
    _bVec,
    label,
    buffer,
    range,
    dist,
)

Create a DislocationLoop.

source
DDD.DislocationNetworkType
struct DislocationNetwork{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14}
    numNode::T1
    numSeg::T2
    maxConnect::T3
    label::T4
    links::T5
    connectivity::T6
    linksConnect::T7
    slipPlane::T8
    segIdx::T9
    bVec::T10
    coord::T11
    nodeVel::T12
    nodeForce::T13
    segForce::T14
end

Stores a dislocation network.

source
DDD.DislocationNetworkType
DislocationNetwork(
    sources::DislocationLoopCollection,
    maxConnect = 4,
    args...;
    memBuffer = nothing,
    checkConsistency = true,
    kw...,
)

Creates a DislocationNetwork out of a DislocationLoopCollection.

Inputs

  • args... are optional arguments that will be passed on to the loopDistribution function which distributes the loops in sources according to the type of their dist variable.
  • kw... are optional keyword arguments that will also be passed to loopDistribution.
  • memBuffer is the numerical value for allocating memory in advance. The quantity, memBuffer × N, where N is the total number of nodes in sources, will be the initial number of entries allocated in the matrices that keep the network's data. If no memBuffer is provided, the number of entries allocated will be `round(N*log2(N)).
source
DDD.DislocationNetworkMethod
DislocationNetwork(;
    links::AbstractArray,
    slipPlane::AbstractArray,
    bVec::AbstractArray,
    coord::AbstractArray,
    label::AbstractVector{nodeTypeDln},
    nodeVel::AbstractArray,
    nodeForce::AbstractArray,
    numNode = length(label),
    numSeg = size(links, 2),
    maxConnect = 4,
    connectivity::AbstractArray = zeros(Int, 1 + 2 * maxConnect, length(label)),
    linksConnect::AbstractArray = zeros(Int, 2, size(links, 2)),
    segIdx::AbstractArray = zeros(Int, size(links, 2), 3),
    segForce::AbstractArray = zeros(3, 2, size(links, 2)),
)

Create a DislocationNetwork. We recommend generating networks from DislocationLoop unless you want a special case.

source
DDD.DislocationParametersType
struct DislocationParameters{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14,T15,T16,T17,T18,T19,T20,T21,T22,T23,T24}
    mobility::T1
    dragCoeffs::T2
    coreRad::T3
    coreRadSq::T4
    coreRadMag::T5
    coreEnergy::T6
    minSegLen::T7
    minSegLenSq::T8
    maxSegLen::T9
    twoMinSegLen::T10
    minArea::T11
    maxArea::T12
    minAreaSq::T13
    maxAreaSq::T14
    collisionDist::T15
    collisionDistSq::T16
    slipStepCritLen::T17
    slipStepCritArea::T18
    remesh::T19
    collision::T20
    separation::T21
    virtualRemesh::T22
    parCPU::T23
    parGPU::T24
end

Stores the dislocation parameters.

Fieldnames

  • mobility: mobility law
  • dragCoeffs: drag coefficients, use of a named tuple is recommended
  • coreRad: dislocation core radius
  • coreRadSq: square of the dislocation core radius
  • coreRadMag: magnitude of the dislocation core radius
  • coreEnergy: core energy
  • minSegLen: minimum segment length
  • maxSegLen: maximum segment length
  • twoMinSegLen: two times minimum segment length
  • minSegLenSq: square of the minimum segment length
  • minArea: minimum area
  • maxArea: maximum area
  • minAreaSq: square of the minimum area
  • maxAreaSq: sqare of the maximum area
  • collisionDist: collision distance
  • slipStepCritLen: critical length for slip step tracking
  • slipStepCritArea: critical area for slip step tracking
  • remesh: remeshing flag
  • collision: collision flag
  • separation: separation flag
  • virtualRemesh: virtual remeshing flag
  • parCPU: parallelise over CPU flag
  • parGPU: parallelise over GPU flag
source
DDD.DislocationParametersMethod
DislocationParameters(;
    mobility::AbstractMobility,
    dragCoeffs = (edge = 1.0, screw = 2.0, climb = 1e9),
    coreRad = 1.0,
    coreRadMag = 1.0,
    coreEnergy = 1 / (4 * π) * log(coreRad / 0.1),
    minSegLen = 2 * coreRad,
    maxSegLen = 20 * coreRad,
    minArea = minSegLen^2 / sqrt(2),
    maxArea = 100 * minArea,
    collisionDist = minSegLen / 2,
    slipStepCritLen = maxSegLen / 2,
    slipStepCritArea = 0.5 * (slipStepCritLen^2) * sind(1),
    remesh = true,
    collision = true,
    separation = true,
    virtualRemesh = true,
    parCPU = false,
    parGPU = false,
)

Creates DislocationParameters. Automatically calculates derived values, asserts values are reasonable and provides sensible default values.

source
DDD.FEMParametersType
struct FEMParameters{T1,T2,T3,T4,T5,T6,T7,T8,T9}
    type::T1
    order::T2
    model::T3
    dx::T4
    dy::T5
    dz::T6
    mx::T7
    my::T8
    mz::T9
end

Stores the finite element parameters.

Fields

  • type: Mesh type
  • order: Element order
  • model: Experimental model
  • dx: Dimension in x
  • dy: Dimension in y
  • dz: Dimension in z
  • mx: Elements in x
  • my: Elements in y
  • mz: Elements in z
source
DDD.ForceDisplacementType
struct ForceDisplacement{T1,T2,T3,T4,T5,T6}
    uTilde::T1
    uHat::T2
    u::T3
    fTilde::T4
    fHat::T5
    f::T6
end

Stores displacements and forces applied on the FE nodes.

Fields

  • uTilde: Dislocation displacements
  • uHat: Corrective displacements
  • u: Displacements
  • fTilde Dislocation forces
  • fHat: Corrective forces
  • f: Forces
source
DDD.ForceDisplacementDotType
struct ForceDisplacementDot{T1,T2,T3,T4}
    uDotDofs::T1
    uDot::T2
    fDotDofs::T3
    fDot::T4
end

Loading and displacement rate. Stores the degrees of freedom on which the loading is applied as well as the loading values.

Fields

  • uDotDofs: degrees of freedom on which a displacement is applied
  • uDot: displacement rate
  • fDotDofs: degrees of freedom on which a load is applied
  • fDot: loading rate
source
DDD.IntegrationParametersType
struct IntegrationParameters{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10}
    method::T1
    tmin::T2
    tmax::T3
    dtmin::T4
    dtmax::T5
    abstol::T6
    reltol::T7
    maxchange::T8
    exponent::T9
    maxiter::T10
end

Stores integration parameters.

source
DDD.IntegrationTimeType
struct IntegrationTime{T1,T2,T3}
    dt::T1      # Current time step.
    time::T2    # Current simulation time.
    step::T3    # Current simulation step.
end

Store integration time and steps.

source
DDD.MaterialParametersType
struct MaterialParameters{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14,T15}
    crystalStruct::T1   # Crystal structure.
    μ::T2               # Shear modulus.
    μMag::3             # Magnitude of shear modulus.
    ν::T4               # Poisson ratio.
    E::T5               # Young's modulus.
    omνInv::T6          # 1 / (1 - ν)
    opνInv::T7          # 1 / (1 + ν)
    νomνInv::T8         # ν / (1 - ν)
    νopνInv::T9         # v / (1 + ν)
    μ4π::T10            # μ / (4π)
    μ8π::T11            # μ / (8π)
    μ4πν::T12           # μ / (4π (1 - ν))
    omνInv8π::T13       # 1 / (8π (1 - ν))
    om2νomνInv8π::T14   # (1 - 2 * ν) / (8π (1 - ν))
    σPN::T15            # Peierls-Nabarro stress.
end

Store material parameters.

source
DDD.RegularCuboidMeshType
struct RegularCuboidMesh{
    T1,
    T2,
    T3,
    T4,
    T5,
    T6,
    T7,
    T8,
    T9,
    T10,
    T11,
    T12,
    T13,
    T14,
    T15,
    T16,
    T17,
    T18,
    T19,
    T20,
    T21,
    T22,
    T23,
    T24,
    T25,
    T26,
    T27,
    T28,
} <: AbstractRegularCuboidMesh
    order::T1
    dx::T2
    dy::T3
    dz::T4
    mx::T5
    my::T6
    mz::T7
    w::T8
    h::T9
    d::T10
    scale::T11
    numElem::T12
    numNode::T13
    C::T14
    vertices::T15
    faces::T16
    faceNorm::T17
    faceMidPt::T18
    cornerNode::T19
    edgeNode::T20
    faceNode::T21
    surfNode::T22
    surfNodeArea::T23
    surfNodeNorm::T24
    surfElemNode::T25
    coord::T26
    connectivity::T27
    K::T28
end

Stores regular a cuboid mesh.

Fields

  • order: Element order
  • dx: Size in x
  • dy: Size in y
  • dz: Size in z
  • mx: Elements in x
  • my: Elements in y
  • mz: Elements in z
  • w: Width
  • h: Height
  • d: Depth
  • scale: Mesh scale
  • numElem: Number of elements
  • numNode: Number of nodes
  • C: Stiffness tensor
  • vertices: Vertices
  • faces: Faces
  • faceNorm: Face normals
  • faceMidPt: Face mid-points
  • cornerNode: Corner nodes set
  • edgeNode: Edge nodes set
  • faceNode: Face node set
  • coord: Node coordinates
  • connectivity: Node connectivity
  • K: Stiffness matrix
source
DDD.RegularCuboidMeshMethod
RegularCuboidMesh(
    matParams::MaterialParameters,
    femParams::FEMParameters{F1,F2,F3,F4,F5} where {F1<:DispatchRegularCuboidMesh,F2<:LinearElement,F3,F4,F5}
)

FE domain with linear hexahedral elements with 2 × 2 × 2 Gauss nodes. Uses the canonical node, face and edge ordering of https://www.mcs.anl.gov/uploads/cels/papers/P1573A.pdf.

Corners:

      7----------8
     /|         /|
    / |        / |
   /  |       /  |
  /   4------/---3
 /   /      /   /
5----------6   /
|  /       |  /
| /        | /
|/         |/
1----------2

Edges:

      .----11-----.
     /|          /|
    / 8         / 7
  12  |       10  |
  /   .-----3-/---.
 /   /       /   /
.-----9-----.   /
|  4        |  2
5 /         6 /
|/          |/
.-----1-----.

Faces:

      .----------.
     /|         /|
    / |    3   / |
   /  | 6     /  |
  /   .------/---.
 / 4 /      / 2 /
.----------.   /
|  /     5 |  /
| /  1     | /
|/         |/
.----------.


Z    
^   y
|  ^  
| / 
.-----> x
source
DDD.SlipSystemType
struct SlipSystem{T1,T2,T3}
    crystalStruct::T1
    slipPlane::T2
    bVec::T3
end

Stores slip systems.

source
DDD.SlipSystemMethod
SlipSystem(;
    crystalStruct::AbstractCrystalStruct,
    slipPlane::AbstractArray,
    bVec::AbstractArray
)

Creates a SlipSystem.

Ensures slipPlane ⟂ bVec.

source
DDD.nodeTypeDlnType
@enum nodeTypeDln begin
    noneDln = 0    # Undefined node, value at initialisation
    intMobDln = 1  # Internal mobile node
    intFixDln = 2  # Internal fixed node
    srfMobDln = 3  # Mobile surface node
    srfFixDln = 4  # Fixed surface node
    extDln = 5     # External node
    tmpDln = 6     # Temporary flag, used during topological operations
end

Different types of nodes behave differently. There are only a finite number of them so an enumerated type provides safety and efficiency. Each value represents a different type of node and therefore its behaviour.

Instances

  • noneDln: uninitialised nodes.
  • intMobDln: mobile nodes internal to the convex hull of the domain. They take part in tractions, displacements and dislocation interactions.
  • intFixDln: fixed nodes internal to the convex hull of the domain. They participate in the same way as intMobDln nodes except for the fact that their velocities is fixed are zero.
  • srfMobDln: mobile nodes that live on the surface of the convex hull of the domain, they are used to track slip steps and therefore participate in the same things as internal nodes but their velocities are restricted to the convex hull surface.
  • srfFixDln: fixed surface nodes and have the same characteristics as mobile surface nodes except for having zero velocity.
  • extDln: external nodes that do not participate in dislocation interactions or forces but are used to calculate displacements and track slip steps.
  • tmpDln: nodes that are temporarily flagged before they are assigned another type.
source
DDD.nodeTypeFEType
@enum nodeTypeFE begin
    noneFE = 0
    corner = 1  # Corner node
    edge = 2    # Edge node
    face = 3    # Face node
    intFE = 4   # Internal node
end

Finite element node type.

Instances

  • noneFE = 0: Uninitialised
  • corner = 1: Corner
  • edge = 2: Edge
  • face = 3: Face
  • intFE = 4: Internal
source
DDD.:⊗Method
⊗(x::AbstractVector, y::AbstractVector)

Tensor product.

source
DDD.buildMeshMethod
buildMesh(
    matParams::MaterialParameters, 
    femParams::FEMParameters{F1,F2,F3,F4,F5} where {F1<:DispatchRegularCuboidMesh,F2,F3,F4,F5}
)

Creates a RegularCuboidMesh.

source
DDD.calcDisplacementDislocationTriangle!Method
calcDisplacementDislocationTriangle!(uTilde, uDofs, matParams, A, B, C, b, P)

Written by F.Hofmann 9/7/09 Routine to compute the displacement field for a triangular dislocation loop ABC at point P.

Modified by F.Hofmann 5/11/18

Inputs

  • A, B, C: three column vectors defining the nodes of the dislocation loop.
  • P: column vector with coordinates of the point at which the displacement is evaluated in dimension 1. Then in dimension 2 this is a list of points at which to evaluate the field.
  • b: column vector with 3 burgers vector components.

Translated by Daniel Celis Garza.

source
DDD.calcPKForceFunction
calcPKForce(
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
    idx = nothing,
)

calcPKForce!(
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
    idx = nothing,
)

Compute the Peach-Koehler force on segments by using calc_σHat.

$f = (\hat{\mathbb{\sigma}} \cdot \overrightarrow{b}) \times \overrightarrow{t}$

source
DDD.calcSegForceFunction
calcSegForce(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
    idx = nothing,
)

calcSegForce!(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
    idx = nothing,
)

Compute total force on dislocation segments.

source
DDD.calcSelfForceFunction
calcSelfForce(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    network::DislocationNetwork,
    idx = nothing,
)

calcSelfForce!(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    network::DislocationNetwork,
    idx = nothing,
)

Compute the self-interaction force on dislocation segments. calcSelfForce! for its mutating form.

source
DDD.calc_uTilde!Method

Calculates displacements from dislocations. Bruce Bromage, bruce.bromage@materials.ox.ac.uk Github: @brucebromage

Calculating dislocation displacements on the surface of a volume B Bromage and E Tarleton Published 29 October 2018 • © 2018 IOP Publishing Ltd Modelling and Simulation in Materials Science and Engineering, Volume 26, Number 8

@article{bromage2018calculating,
  title={Calculating dislocation displacements on the surface of a volume},
  author={Bromage, B and Tarleton, E},
  journal={Modelling and Simulation in Materials Science and Engineering},
  volume={26},
  number={8},
pages={085007},
  year={2018},
  publisher={IOP Publishing}
}
source
DDD.calc_σHatMethod
calc_σHat(
    mesh::RegularCuboidMesh{T1,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14} where {T1 <: LinearElement,T2,T3,T4,T5,T6,T7,T8,T9,T10,T11,T12,T13,T14},
    forceDisplacement::ForceDisplacement,
    x0,
)

Compute the stress, ̂σ, on a dislocation segment x0 as a result of body forces on a RegularCuboidMesh composed of LinearElement()(@ref). Used by calcPKForce.

Returns

σ = [
        σxx σxy σxz
        σxy σyy σyz
        σxz σyz σzz
    ]
source
DDD.calc_σTildeFunction
calc_σTilde(
    x0,
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    network::DislocationNetwork,
    idx = nothing,
)

Computes the stress tensor, ̃σ, induced by dislocations on points x0.

Returns

σxx = σ[1, :]
σyy = σ[2, :]
σzz = σ[3, :]
σxy = σ[4, :]
σxz = σ[5, :]
σyz = σ[6, :]
source
DDD.calc_σTilde!Function
calc_σTilde!(
    σ,
    x0,
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    network::DislocationNetwork,
idx = nothing,
)

In-place computation of the stress tensor,̃σ, induced by dislocations on points x0.

Returns

σxx = σ[1, :]
σyy = σ[2, :]
σzz = σ[3, :]
σxy = σ[4, :]
σxz = σ[5, :]
σyz = σ[6, :]
source
DDD.checkNetworkMethod
checkNetwork(network::DislocationNetwork)

Checks the validity of the dislocation network. It ensures the following conditions are met by the member variables of network:

  1. connectivity and links have the same number of non-zero entries;
  2. all entries in bVec are non-zero;
  3. only the trailing columns of connectivity are zeros;
  4. consistency between connectivity and links;
  5. bVec is conserved among connected nodes;
  6. entries in links are unique;
  7. consistency betwen connectivity and linksConnect
source
DDD.coarsenNetwork!Method
coarsenNetwork!(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
)

Coarsens network such that no links are smaller than the minimum allowable length and so that no two links form triangles with area under the minimum allowed.

source
DDD.coarsenVirtualNetwork!Method
coarsenVirtualNetwork!(dlnParams::DislocationParameters, network::DislocationNetwork)

Check whether virtual nodes can be eliminated based on:

  1. If they are not connected to any surface nodes
  2. If they are not due to an angle change in the simulated volume surface

Bruce Bromage, Github @brucebromage Michromechanical Testing Group Department of Materials, University of Oxford bruce.bromage@materials.ox.ac.uk May 2017

Adapted Jan 2021 Daniel Celis Garza, Github @dcelisgarza

source
DDD.compStructMethod
compStruct(arg1, arg2; verbose::Bool = false)

Compares values of the fields of two variables arg1 and arg2 with the same structure. If verbose = true, it will print which fields are different from each other.

Examples

julia> struct MyStruct1; x; end
julia> test1 = MyStruct1(1)
MyStruct1(1)
julia> test2 = MyStruct1(5)
MyStruct1(5)
julia> compStruct(test1, test2; verbose = true)
Structures differ in field: x.
false
julia> compStruct(1, 1; verbose = true)
true
julia> compStruct(1, [1]; verbose = true)
false
source
DDD.deriv!Method
deriv!(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
)

Computes the nodal velocities of a network.

source
DDD.dlnMobilityFunction
dlnMobility(
    dlnParams::DislocationParameters{T1,T2,T3,T4} where {T1 <: mobBCC,T2,T3,T4},
    matParams::MaterialParameters,
    network::DislocationNetwork,
    idx = nothing,
)

Compute the nodal force and velocities for a BCC material.

Original by Bruce Bromage at the Department of Materials of the University of Oxford, @brucebromage on github.

This is outdated, new capabilities include rotating the frame of reference and better handling of cross-slip.

source
DDD.dlnMobility!Function
dlnMobility!(
    dlnParams::DislocationParameters{T1,T2,T3,T4} where {T1 <: mobBCC,T2,T3,T4},
    matParams::MaterialParameters,
    network::DislocationNetwork,
    idx = nothing,
)

In-place computation of the nodal force and velocities for a BCC material.

Original by Bruce Bromage at the Department of Materials of the University of Oxford, @brucebromage on github.

This is outdated, new capabilities include rotating the frame of reference and better handling of cross-slip.

source
DDD.externalAngleMethod
externalAngle(n::Int)

Compute the exterior angle of a regular polygon with n sides.

source
DDD.findConnectedNodeMethod
findConnectedNode(network::DislocationNetwork, node, connection)

Use a node's connection to find the other node in the link. Returns a tuple with the link corresponding to the connection of node, the row of links the connected node appears, and the connected node.

Examples

We want to find the node connected to node 5 via connection 3.

julia> link, rowColLink, connectedNode = findConnectedNode(network, 5, 10)
julia> link
10
julia> oppRowLink
1
julia> connectedNode
25

Means that the link corresponding to connection 3 of node 5 is network.links[:, 10]; the node connected to node 5 is found on network.links[oppRowLink, 10] and is node 25, so the node we're looking for is network.links[oppRowLink, 10] == 25 and the node whose connection we are looking for is network.links[3-oppRowLink, 10] == 5. In this case oppRowLink == 1, so through its 3rd connection, node 5 is the second node on link 10, where it is connected to node 25.

source
DDD.getSegmentIdxMethod
getSegmentIdx(links, label)

Finds the indices of a link and corresponding nodes.

Returns

n × 3 matrix is of the form [i, node1, node2].

  • i can be used to find the Burgers vector, slip plane and segment forces of segment i, eg bVec[:, i].
  • node1 and node2 can be used to find the coordinate and velocity of the nodes, eg l = coord[:, node2] - coord[:, node1].
source
DDD.integrate!Method
integrate!(
    intParams::IntegrationParameters{T1,T2,T3} where {T1 <: AdaptiveEulerTrapezoid,T2,T3},
    intVars::IntegrationTime,
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
)

Integrates nodal velocities using a time-adaptive Euler-Trapezoid method.

source
DDD.internalAngleMethod
internalAngle(n::Int)

Compute the interior angle of a regular polygon with n sides.

source
DDD.linePlaneIntersectMethod
linePlaneIntersect(n::T, p0::T, l::T, l0::T) where {T <: AbstractVector}

Finds the intersect between a line and a plane. n is the plane normal, p0 is a point on the plane, l is the line vector, l0 is a point on a line.

source
DDD.loadBoundariesMethod
loadBoundaries(dict::Dict{T1,T2}) where {T1,T2}

Constructs Boundaries out of a dictionary.

Note

tK may be null if it was factorised when the variable was saved.

source
DDD.loadJSONMethod
loadJSON(filename::AbstractString)

Wrapper for JSON.parsefile(filename). Loads a JSON file as a dictionary.

source
DDD.makeConnectMethod
makeConnect(links, maxConnect)
makeConnect!(network::DislocationNetwork)

Creates connectivity and linksConnect matrices for DislocationNetwork.

  • connectivity contains the number of other other nodes each node is connected to, up to maxConnect other nodes. Every (2i, j) entry contains a node connected to node j. Every (2i+1, j) coordinate contains whether that node is the first or second node in the link.
  • linksConnect relates connections enabled by a link. Analogous to the connectivity of a link.
source
DDD.makeNetwork!Method
makeNetwork!(
    links,
    slipPlane,
    bVec,
    coord,
    label,
    sources,
    lims,
    initIdx,
    args...;
    kw...,
)

Internal function called by DislocationNetwork to fill the arrays that define the network.

source
DDD.mergeNode!Method
mergeNode!(network::DislocationNetwork, nodeKept, nodeGone)

Merges nodeGone into nodeKept.

source
DDD.minimumDistanceMethod
minimumDistance(x0, x1, y0, y1, vx0, vx1, vy0, vy1)

Calculates the minimum distance between two segments seg1 = x0 → x1, seg2 = y0 → y1.

Returns

  • distSq: minimum distance squared between L1 and L2.
  • dDistSqDt: rate of change with respect to time of the minimum distance squared between L1 and L2.
  • L1: normalised position on seg1 that is closest to seg2.
  • L2: normalised position on seg2 that is closest to seg1.
source
DDD.refineNetwork!Method
refineNetwork!(
    dlnParams::DislocationParameters,
    matParams::MaterialParameters,
    mesh::AbstractMesh,
    forceDisplacement::ForceDisplacement,
    network::DislocationNetwork,
)

Refines a dislocation network to ensure all the segments are shorter than the maximum allowable length and so no two links form a triangle with an area over the maximum allowed.

source
DDD.rot3DMethod
rot3D(xyz, uvw, abc, θ)

Rotate point xyz about the vector uvw that crosses point abc by the angle θ. Further details found here.

Examples

julia> rot3D([1;1;1],[1;0;0],[0;0;0],π/2)
3-element Array{Float64,1}:
  1.0
 -0.9999999999999999
  1.0
julia> rot3D([1;1;1],[1;0;0],[0;0;0],-π/2)
3-element Array{Float64,1}:
1.0
1.0
-0.9999999999999999
julia> rot3D([1;1;1],[1;0;0],[0;0;0],π)
3-element Array{Float64,1}:
  1.0
 -1.0000000000000002
 -0.9999999999999999
source
DDD.saveJSONMethod
saveJSON(filename::AbstractString, args...; mode::AbstractString = "w")

Wrapper for JSON.print.

source
DDD.shapeFunctionMethod
shapeFunction(::LinearQuadrangle3D, x::AbstractVector, y::AbstractVector, z::AbstractVector)

Computes the linear shape functions N[1:8][p] for (x, y, z) point p on a 3D linear quadrangle.

source
DDD.shapeFunctionMethod
shapeFunction(::LinearQuadrangle3D, x, y, z)

Computes the linear shape functions N[1:8] for an (x, y, z) point on a 3D linear quadrangle.

source
DDD.shapeFunctionDerivMethod
shapeFunctionDeriv(shape<:AbstractShapeFunction, x::AbstractVector, y::AbstractVector, z::AbstractVector)

Computes the derivatives of the linear shape functions N[1:3, 1:8] for an (x, y, z) point on a 3D linear quadrangle.

Returns

dNdS[x, n, p](x,y,z) := x'th derivative of shape function n for point p.
source
DDD.shapeFunctionDerivMethod
shapeFunctionDeriv(shape<:AbstractShapeFunction, x, y, z)

Computes the derivatives of the linear shape functions N[1:3, 1:8] for an (x, y, z) point on a 3D linear quadrangle.

Returns

dNdS[x, n](x,y,z) := x'th derivative of shape function n.
source
DDD.splitNode!Method
splitNode!(network::DislocationNetwork, splitNode, splitConnect, midCoord, midVel)

Splits node splitNode along connection splitConnect, puts it at coordinate midCoord with velocity midVel. If it is called from within refineNetwork!, midCoord and midVel are the coordinate between splitNode and the node connected to it via splitConnect and gives it the mean velocity of the two nodes.

source
DDD.inclusiveComparisonMethod
inclusiveComparison(data, args...)::Bool

Compare data to a tuple, return true if it is equal to any arg, false if it is not equal to any.

Examples

julia> inclusiveComparison("f", 1,4,5,"f")
true
julia> inclusiveComparison(23.246, 1.5, 4, 5, "f")
false
source
DDD.makeInstanceDictMethod
makeInstanceDict(valType::DataType)

Make a dictionary of enumerated variable instances. Helps in translating JSON files.

source
DDD.makeTypeDictMethod
makeTypeDict(valType::DataType)

Inputs contain strings that correspond to DDD data types. This function atuomatically creates a dictionary for all concrete subtypes of a given valType.

Examples

julia> abstract type MyAbstractType end
julia> struct MyStruct1 <: MyAbstractType end
julia> struct MyStruct2 <: MyAbstractType end
julia> makeTypeDict(MyAbstractType)
Dict{String,Any} with 4 entries:
  "DDD.MyStruct1()" => MyStruct1()
  "DDD.MyStruct2()" => MyStruct2()
  "MyStruct1()"     => MyStruct1()
  "MyStruct2()"     => MyStruct2()
source
DDD.removeConnection!Method
removeConnection!(network::DislocationNetwork, nodeKept, connectGone)

Removes connection connectGone from nodeKept.

source
DDD.removeLink!Function
removeLink!(network::DislocationNetwork, linkGone, lastLink = nothing)

Removes link linkGone and uses lastLink to reoganise the network.

source
DDD.removeNode!Function
removeNode!(network::DislocationNetwork, nodeGone, lastNode = nothing)

Removes node nodeGone from network.

source
DDD.subTypeTreeMethod
subTypeTree(t; dict = Dict(), level = 1, cutoff = 0)

Create subtype dictionary. Adapted from https://github.com/JuliaLang/julia/issues/24741

source
JSON.Writer.lowerMethod
JSON.lower(
    t::T,
) where {
    T <: Union{
        AbstractCrystalStruct,
        AbstractMobility,
        AbstractIntegrator,
        AbstractDlnSeg,
        AbstractDlnStr,
        AbstractDistribution,
    },
}

JSON.lower(t::nodeTypeDln)

Extensions to JSON.lower for custom types. Allows these variables to be serialised properly.

source