Reference Manual

Simple FE model (volume)

FinEtools.FEMMBaseModule.inspectintegpointsMethod
inspectintegpoints(
    self::FEMM,
    geom::NodalField{GFT},
    u::NodalField{UFT},
    dT::NodalField{TFT},
    felist::AbstractVector{IT},
    inspector::F,
    idat,
    quantity = :Cauchy;
    context...,
) where {
    FEMM<:AbstractFEMMDeforLinear,
    GFT<:Number,
    UFT<:Number,
    TFT<:Number,
    IT<:Integer,
    F<:Function,
}

Inspect integration point quantities.

Arguments

  • geom - reference geometry field
  • u - displacement field
  • dT - temperature difference field
  • felist - indexes of the finite elements that are to be inspected: The fes to be included are: fes[felist].
  • context - structure: see the update!() method of the material.
  • inspector - functionwith the signature idat = inspector(idat, j, conn, x, out, loc); where idat - a structure or an array that the inspector may use to maintain some state, for instance minimum or maximum of stress, j is the element number, conn is the element connectivity, out is the output of the update!() method, loc is the location of the integration point in the reference configuration.

Return

The updated inspector data is returned.

FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.infsup_ghMethod
infsup_gh(
    self::AbstractFEMMDeforLinear,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {A<:AbstractSysmatAssembler,GFT,UFT}

Compute the matrix to produce the norm of the divergence of the displacement.

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.infsup_shMethod
infsup_sh(
    self::AbstractFEMMDeforLinear,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the matrix to produce the seminorm of the displacement (square root of the sum of the squares of the derivatives of the components of displacement).

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.massMethod
mass(
    self::AbstractFEMMDeforLinear,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the consistent mass matrix

This is a general routine for the abstract linear-deformation FEMM.

FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.stiffnessMethod
stiffness(
    self::FEMM,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {FEMM<:AbstractFEMMDeforLinear,A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute and assemble stiffness matrix.

Only for homogeneous materials

The material stiffness matrix is assumed to be the same at all the points of the domain (homogeneous material).

FinEtoolsDeforLinear.FEMMDeforLinearModule.FEMMDeforLinearMethod
FEMMDeforLinear(
    mr::Type{MR},
    integdomain::IntegDomain{S,F},
    material::M,
) where {
    MR<:AbstractDeforModelRed,
    S<:AbstractFESet,
    F<:Function,
    M<:AbstractMatDeforLinearElastic,
}

Constructor of linear deformation finite element modeling machine.

Simple FE models (surface)

FinEtoolsDeforLinear.FEMMDeforWinklerModule.surfacenormalspringstiffnessMethod
surfacenormalspringstiffness(
    self::FEMMDeforWinkler,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
    springconstant::UFT,
    surfacenormal::SurfaceNormal,
) where {A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the stiffness matrix of surface normal spring.

Rationale: consider continuously distributed springs between the surface of the solid body and the 'ground', in the direction normal to the surface. If the spring coefficient becomes large, we have an approximate method of enforcing the normal displacement to the surface.

FinEtoolsDeforLinear.FEMMDeforSurfaceDampingModule.dampingABCMethod
dampingABC(
    self::FEMMDeforSurfaceDamping,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
    impedance::FT,
    surfacenormal::SurfaceNormal,
) where {A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number,FT<:Number}

Compute the damping matrix associated with absorbing boundary conditions (ABC).

Compute the damping matrix associated with absorbing boundary conditions (ABC) representation of the effect of infinite extent of inviscid fluid next to the surface.

Advanced: Mean-strain FEM

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(
    self::FEMMDeforLinearMSH8,
    geom::NodalField{GFT},
) where {GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(
    self::FEMMDeforLinearMST10,
    geom::NodalField{GFT},
) where {GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtools.FEMMBaseModule.inspectintegpointsMethod
inspectintegpoints(
    self::FEMM,
    geom::NodalField{GFT},
    u::NodalField{UFT},
    dT::NodalField{TFT},
    felist::Vector{IT},
    inspector::F,
    idat,
    quantity = :Cauchy;
    context...,
) where {FEMM<:AbstractFEMMDeforLinearMS,GFT<:Number,UFT<:Number,TFT<:Number,IT,F<:Function}

Inspect integration point quantities.

Arguments

  • geom - reference geometry field
  • u - displacement field
  • dT - temperature difference field
  • felist - indexes of the finite elements that are to be inspected: The fes to be included are: fes[felist].
  • context - structure: see the update!() method of the material.
  • inspector - functionwith the signature idat = inspector(idat, j, conn, x, out, loc); where idat - a structure or an array that the inspector may use to maintain some state, for instance minimum or maximum of stress, j is the element number, conn is the element connectivity, out is the output of the update!() method, loc is the location of the integration point in the reference configuration.

Return

The updated inspector data is returned.

FinEtoolsDeforLinear.FEMMDeforLinearMSModule.infsup_ghMethod
infsup_gh(
    self::FEMM,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {FEMM<:AbstractFEMMDeforLinearMS,A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the matrix to produce the norm of the divergence of the displacement.

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearMSModule.infsup_shMethod
infsup_sh(
    self::FEMM,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {FEMM<:AbstractFEMMDeforLinearMS,A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the matrix to produce the seminorm of the displacement (square root of the sum of the squares of the derivatives of the components of displacement).

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMSH8Type
mutable struct FEMMDeforLinearMSH8{
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S,F} where {S<:FESetH8,F<:Function},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
    MS<:MatDeforElastIso,
} <: AbstractFEMMDeforLinearMS

Type for mean-strain linear deformation FEMM based on eight-node hexahedral elements.

FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMST10Type
mutable struct FEMMDeforLinearMST10{
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S,F} where {S<:FESetT10,F<:Function},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
    MS<:MatDeforElastIso,
} <: AbstractFEMMDeforLinearMS

Type for mean-strain linear deformation FEMM based on 10-node tetrahedral elements.

Advanced: Nodal integration

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(self::F,
    geom::NodalField{GFT}) where {F <: AbstractFEMMDeforLinearNICE, GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(
    self::FEMM,
    geom::NodalField{GFT},
) where {FEMM<:FEMMDeforLinearESNICEH8,GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(
    self::FEMM,
    geom::NodalField{GFT};
    stabilization_parameters = _T4_stabilization_parameters,
) where {FEMM<:FEMMDeforLinearESNICET4,GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtools.FEMMBaseModule.inspectintegpointsMethod
inspectintegpoints(
    self::FEMM,
    geom::NodalField{GFT},
    u::NodalField{UFT},
    dT::NodalField{TFT},
    felist::Vector{IT},
    inspector::F,
    idat,
    quantity = :Cauchy;
    context...,
) where {FEMM<:AbstractFEMMDeforLinearESNICE,GFT<:Number,UFT<:Number,TFT<:Number,IT<:Integer,F<:Function}

Inspect integration point quantities.

Arguments

  • geom - reference geometry field
  • u - displacement field
  • dT - temperature difference field
  • felist - indexes of the finite elements that are to be inspected: The fes to be included are: fes[felist].
  • context - structure: see the update!() method of the material.
  • inspector - functionwith the signature idat = inspector(idat, j, conn, x, out, loc); where idat - a structure or an array that the inspector may use to maintain some state, for instance minimum or maximum of stress, j is the element number, conn is the element connectivity, out is the output of the update!() method, loc is the location of the integration point in the reference configuration.

Return

The updated inspector data is returned.

FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.infsup_ghMethod
infsup_gh(
    self::FEMM,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {FEMM<:AbstractFEMMDeforLinearESNICE,A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the matrix to produce the norm of the divergence of the displacement.

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.infsup_shMethod
infsup_sh(
    self::AbstractFEMMDeforLinearESNICE,
    assembler::A,
    geom::NodalField{GFT},
    u::NodalField{UFT},
) where {A<:AbstractSysmatAssembler,GFT<:Number,UFT<:Number}

Compute the matrix to produce the seminorm of the displacement (square root of the sum of the squares of the derivatives of the components of displacement).

This matrix is used in the numerical infsup test (Klaus-Jurgen Bathe, The inf-sup condition and its evaluation for mixed finite element methods, Computers and Structures 79 (2001) 243-252.)

Note

This computation has not been optimized in any way. It can be expected to be inefficient.

FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICEH8Type
mutable struct FEMMDeforLinearESNICEH8{
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S} where {S<:FESetH8},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
    MS<:MatDeforElastIso,
} <: AbstractFEMMDeforLinearESNICE

FEMM type for Energy-sampling Stabilized Nodally Integrated Continuum Elements (ESNICE) based on the 8-node hexahedron.

FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICET4Type
mutable struct FEMMDeforLinearESNICET4{
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S} where {S<:FESetT4},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
    MS<:MatDeforElastIso,
} <: AbstractFEMMDeforLinearESNICE

FEMM type for Energy-sampling Stabilized Nodally Integrated Continuum Elements (ESNICE) based on the 4-node tetrahedron.

Advanced: Incompatible modes

FinEtools.FEMMBaseModule.associategeometry!Method
associategeometry!(self::F,  geom::NodalField{GFT}) where {F<:FEMMDeforLinearIMH8, GFT}

Associate geometry field with the FEMM.

Compute the correction factors to account for the shape of the elements.

FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8Type
mutable struct FEMMDeforLinearIMH8{
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S,F} where {S<:FESetH8,F<:Function},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
} <: AbstractFEMMDeforLinear

Type for mean-strain linear deformation FEMM based on eight-node hexahedral elements with incompatible modes.

Default number of incompatible modes is 12 (Simo formulation). Alternative is 9 incompatible modes (Wilson formulation).

FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8Method
FEMMDeforLinearIMH8(
    mr::Type{MR},
    integdomain::ID,
    mcsys::CS,
    material::M,
) where {
    MR<:AbstractDeforModelRed,
    ID<:IntegDomain{S,F} where {S<:FESetH8,F<:Function},
    CS<:CSys,
    M<:AbstractMatDeforLinearElastic,
}

Constructor. Default number of incompatible modes.

FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8Method
FEMMDeforLinearIMH8(
    mr::Type{MR},
    integdomain::ID,
    material::M,
    nmodes::Int64,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S,F} where {S<:FESetH8,F<:Function}, M<:AbstractMatDeforLinearElastic}

Constructor, with optional configuration of the number of incompatible modes.

Algorithms

Linear deformation

FinEtoolsDeforLinear.AlgoDeforLinearModule.exportdeformationMethod
AlgoDeforLinearModule.exportdeformation(modeldata::FDataDict)

Algorithm for exporting of the deformation for visualization in Paraview.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set

  • "regions" = array of region dictionaries

  • "geom" = geometry field

  • "u" = displacement field, or

  • "us" = array of tuples (name, displacement field)

  • "postprocessing" = dictionary with values for keys

    • "boundary_only" = should only the boundary of the regions be rendered? Default is render the interior.
    • "file" = name of the postprocessing file

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element mmodel machine (mandatory);

Output

modeldata updated with

  • modeldata["postprocessing"]["exported"] = array of data dictionaries, one for each exported file. The data is stored with the keys:

    • "file" - names of exported file
    • "field" - nodal or elemental field
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportmodeMethod
AlgoDeforLinearModule.exportmode(modeldata::FDataDict)

Algorithm for exporting of the mmode shape for visualization in Paraview.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set

  • "regions" = array of region dictionaries

  • "geom" = geometry field

  • "u" = displacement field

  • "W" = Computed free-vibration eigenvectors, neigvs columns

  • "omega" = Computed free-vibration angular frequencies, array of length neigvs

  • "postprocessing" = dictionary with values for keys

    • "boundary_only" = should only the boundary of the regions be rendered? Default is render the interior.
    • "file" = name of the postprocessing file
    • "mode" = which mode should be visualized?
    • "component" = which component of the quantity?
    • "outputcsys" = output coordinate system

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element mmodel machine (mandatory);

Output

modeldata updated with

  • modeldata["postprocessing"]["exported"] = see exportdeformation()
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportstressMethod
AlgoDeforLinearModule.exportstress(modeldata::FDataDict)

Algorithm for exporting of the stress for visualization in Paraview.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set

  • "regions" = array of region dictionaries

  • "geom" = geometry field

  • "u" = displacement field

  • "postprocessing" = dictionary with values for keys

    • "boundary_only" = should only the boundary of the regions be rendered? Default is render the interior.
    • "file" = name of the postprocessing file
    • "quantity" = quantity to be exported (default :Cauchy)
    • "component" = which component of the quantity?
    • "outputcsys" = output coordinate system
    • "inspectormeth" = inspector method to pass to inspectintegpoints()
    • "extrap" = method for extrapolating from the quadrature points to the nodes within one element

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element mmodel machine (mandatory);

Output

modeldata updated with

  • modeldata["postprocessing"]["exported"] = array of data dictionaries, one for each exported file. The data is stored with the keys:

    • "file" - name of exported file
    • "field" - nodal field
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportstresselementwiseMethod
AlgoDeforLinearModule.exportstresselementwise(modeldata::FDataDict)

Algorithm for exporting of the elementwise stress for visualization in Paraview.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set

  • "regions" = array of region dictionaries

  • "geom" = geometry field

  • "u" = displacement field

  • "postprocessing" = dictionary with values for keys

    • "boundary_only" = should only the boundary of the regions be rendered? Default is render the interior.
    • "file" = name of the postprocessing file
    • "quantity" = quantity to be exported (default :Cauchy)
    • "component" = which component of the quantity?
    • "outputcsys" = output coordinate system

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element mmodel machine (mandatory);

Output

modeldata updated with

  • modeldata["postprocessing"]["exported"] = array of data dictionaries, one for each exported file. The data is stored with the keys:

    • "file" - name of exported file
    • "field" - elemental field
FinEtoolsDeforLinear.AlgoDeforLinearModule.linearstaticsMethod
AlgoDeforLinearModule.linearstatics(modeldata::FDataDict)

Algorithm for static linear deformation (stress) analysis.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set
  • "regions" = array of region dictionaries
  • "essential_bcs" = array of essential boundary condition dictionaries
  • "traction_bcs" = array of traction boundary condition dictionaries
  • "temperature_change" = dictionary of data for temperature change

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element model machine (mandatory);

For essential boundary conditions (optional) each dictionary would hold

  • "displacement" = fixed (prescribed) displacement (scalar), or a function with signature function w = f(x). If this value is not given, zero displacement is assumed.
  • "component" = which component is prescribed (1, 2, 3)?
  • "node_list" = list of nodes on the boundary to which the condition applies (mandatory)

For traction boundary conditions (optional) each dictionary would hold key-value pairs

  • "femm" = finite element model machine (mandatory);
  • "traction_vector" = traction vector, either a constant numerical vector, or a function to be used to construct a ForceIntensity object, or it could be the ForceIntensity object itself.

Output

modeldata = the dictionary on input is augmented with the keys

  • "geom" = the nodal field that is the geometry
  • "u" = the nodal field that is the computed displacement
  • "temp" = the nodal field that is the temperature change
  • "work" = work of the applied loads
  • "timing" = dictionary with timing results
FinEtoolsDeforLinear.AlgoDeforLinearModule.modalMethod
AlgoDeforLinearModule.modal(modeldata::FDataDict)

Modal (free-vibration) analysis solver.

Argument

modeldata = dictionary with values for keys

  • "fens" = finite element node set
  • "regions" = array of region dictionaries
  • "essential_bcs" = array of essential boundary condition dictionaries

For each region (connected piece of the domain made of a particular material), mandatory, the region dictionary contains values for keys:

  • "femm" = finite element mmodel machine (mandatory);

For essential boundary conditions (optional) each dictionary would hold

  • "displacement" = fixed (prescribed) displacement (scalar): only zero displacement is allowed for modal analysis.
  • "component" = which component is prescribed (1, 2, 3)?
  • "node_list" = list of nodes on the boundary to which the condition applies (mandatory)

Control parameters:

  • "neigvs" = number of eigenvalues/eigenvectors to compute
  • "omega_shift"= angular frequency shift for mass shifting
  • "use_lumped_mass" = true or false? (Default is false: consistent mass)

Output

modeldata= the dictionary on input is augmented with

  • "geom" = the nodal field that is the geometry
  • "u" = the nodal field that is the computed displacement
  • "neigvs" = Number of computed eigenvectors
  • "W" = Computed eigenvectors, neigvs columns
  • "omega" = Computed angular frequencies, array of length neigvs # For multi point constraints (MPC) (optional):
  • "raw_eigenvalues" = Raw computed eigenvalues # model_data.mpc= cell array of structs, each for one MPC.

Material models

Material for deformation, base functionality

FinEtoolsDeforLinear.MatDeforModule.rotstressvec!Method
rotstressvec!(::Type{DeforModelRed1D},  outstress::Vector{T},  instress::Vector{T},  Rm::_RotationMatrix) where {T}

Rotate the stress vector by the supplied rotation matrix.

Calculate the rotation of the stress vector to the 'bar' coordinate system given by the columns of the rotation matrix Rm.

  • outstress = output stress vector, overwritten inside
  • instress = input stress vector
  • Rm = columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!Method
rotstressvec!(::Type{DeforModelRed2DAxisymm},  outstress::Vector{T},  instress::Vector{T},  Rm::_RotationMatrix) where {T}

Rotate the stress vector by the supplied rotation matrix.

Calculate the rotation of the stress vector to the 'bar' coordinate system given by the columns of the rotation matrix Rm.

  • outstress = output stress vector, overwritten inside
  • instress = input stress vector
  • Rm = columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!Method
rotstressvec!(::Type{DeforModelRed2DStrain},  outstress::Vector{T},  instress::Vector{T},  Rm::_RotationMatrix) where {T}

Rotate the stress vector by the supplied rotation matrix.

Calculate the rotation of the stress vector to the 'bar' coordinate system given by the columns of the rotation matrix Rm.

  • outstress = output stress vector, overwritten inside
  • instress = input stress vector
  • Rm = columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!Method
rotstressvec!(::Type{DeforModelRed2DStress},  outstress::Vector{T},  instress::Vector{T},  Rm::_RotationMatrix) where {T}

Rotate the stress vector by the supplied rotation matrix.

Calculate the rotation of the stress vector to the 'bar' coordinate system given by the columns of the rotation matrix Rm.

  • outstress = output stress vector, overwritten inside
  • instress = input stress vector
  • Rm = columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!Method
rotstressvec!(::Type{DeforModelRed3D},  outstress::Vector{T},  instress::Vector{T}, Rm::_RotationMatrix) where {T}

Rotate the stress vector by the supplied rotation matrix.

Calculate the rotation of the stress vector to the 'bar' coordinate system given by the columns of the rotation matrix Rm.

  • outstress = output stress vector, overwritten inside
  • instress = input stress vector
  • Rm = columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.strainvdetMethod
strainvdet(::Type{DeforModelRed2DStrain},  Cv::Vector{T}) where {T}

Compute the determinant of a symmetric strain-like square matrix represented as a vector. Remember that the shear strain components are twice the entries of the matrix representation.

FinEtoolsDeforLinear.MatDeforModule.strainvdetMethod
strainvdet(::Type{DeforModelRed3D},  Cv::Vector{T}) where {T}

Compute the determinant of a symmetric strain-like square matrix represented as a vector. Remember that the shear strain components are twice the entries of the matrix representation.

FinEtoolsDeforLinear.MatDeforModule.stressvtot!Method
stressvtot!(::Type{DeforModelRed2DAxisymm}, t::Matrix{T}, v::Vector{T}) where {T}

Convert a 4-vector to a matrix of 3x3 stress components (tensor).

Convert a 4-vector to a symmetric matrix of 3x3 stress components (tensor).

The stress vector components need to be ordered as: sigmax, sigmay, sigmaz, tauxy.

FinEtoolsDeforLinear.MatDeforModule.stressvtot!Method
stressvtot!(::Type{DeforModelRed2DStrain}, t::Matrix{T}, v::Vector{T}) where {T}

Convert a vector to a matrix of 2x2 stress components (symmetric tensor).

If v has 4 entries, also the t[3,3] matrix entry is set.

The stress vector components need to be ordered as: sigmax, sigmay, tauxy, sigmaz, which is the ordering used for the plane-strain model reduction.

FinEtoolsDeforLinear.MatDeforModule.tens4deviator!Method
tens4deviator!(t::Array{T, 4}) where {T}

Compute 4th-order deviator tensor.

Double contraction of a second order tensor with this fourth-order tensor produces the deviator part of the second order tensor.

Example

The product of the deviator tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4deviator!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show tr((S - tr(S)/3*I) ), tr(tS)
FinEtoolsDeforLinear.MatDeforModule.tens4dot2!Method
tens4dot2!(R::Array{T, 2}, F::Array{T, 4}, S::Array{T, 2}) where {T}

Compute the double contraction of a 4th-order and a 2nd-order tensors.

Note

The double contraction of two second-order sensors is defined as A:B = tr(A'*B) = A_ij B_ij

The resulting second-order tensor is first zeroed out, and then the result is accumulated.

FinEtoolsDeforLinear.MatDeforModule.tens4identity!Method
tens4identity!(t::Array{T, 4}) where {T}

Compute 4th-order identity tensor.

Example

The product of the identity tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4identity!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show S - tS
FinEtoolsDeforLinear.MatDeforModule.tens4ijkl!Method
tens4ijkl!(t::Array{T, 4}, A::FA, B::FB, op = :+) where {T, FA, FB}

Fill a 4th-order tensor as a dyadic product of two 2nd-order tensors.

The i,j,k,l component is given as t[i,j,k,l]=A(i,j)*B(k,l).

Note

The tensor is accumulated to. It needs to be initialized to zero, if that is desired as the initial state.

Example

t = fill(0.0, 3, 3, 3, 3)
delta = (I, J) -> I == J ? 1.0 : 0.0
tens4ijkl!(t, delta, delta)
S = rand(3, 3)
@show tr(S) * I
tS = fill(0.0, 3, 3)
@show tens4dot2!(tS, t, S)
FinEtoolsDeforLinear.MatDeforModule.tens4ikjl!Method
tens4ikjl!(t::Array{T, 4}, A::FA, B::FB) where {T, FA, FB}

Fill a 4th-order tensor as a dyadic product of two 2nd-order tensors.

The i,j,k,l component is given as t[i,j,k,l]=A(i,k)*B(j,l).

Note

The tensor is accumulated to. It needs to be initialized to zero, if that is desired as the initial state.

Example

t = fill(0.0, 3, 3, 3, 3)
delta = (I, J) -> I == J ? 1.0 : 0.0
tens4ikjl!(t, delta, delta)
S = rand(3, 3)
@show transpose(S) 
tS = fill(0.0, 3, 3)
@show transpose(S) - tens4dot2!(tS, t, S)
FinEtoolsDeforLinear.MatDeforModule.tens4iljk!Method
tens4iljk!(t::Array{T, 4}, A::FA, B::FB) where {T, FA, FB}

Fill a 4th-order tensor as a dyadic product of two 2nd-order tensors.

The i,j,k,l component is given as t[i,j,k,l]=A(i,l)*B(j,k).

Note

The tensor is accumulated to. It needs to be initialized to zero, if that is desired as the initial state.

Example

t = fill(0.0, 3, 3, 3, 3)
delta = (I, J) -> I == J ? 1.0 : 0.0
tens4iljk!(t, delta, delta)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
@show S - tens4dot2!(tS, t, S)
FinEtoolsDeforLinear.MatDeforModule.tens4skewor!Method
tens4skewor!(t::Array{T, 4}) where {T}

Compute 4th-order skewor tensor.

Double contraction of a second order tensor with this fourth-order tensor produces the skew part of the second order tensor.

Example

The product of the skewor tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4skewor!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show (S - S')/2 * I - tS
FinEtoolsDeforLinear.MatDeforModule.tens4symmetrizor!Method
tens4symmetrizor!(t::Array{T, 4}) where {T}

Compute 4th-order symmetrizor tensor.

Double contraction of a second order tensor with this fourth-order tensor produces the symmetric part of the second order tensor.

Example

The product of the symmetrizor tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4symmetrizor!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show (S + S')/2 * I - tS
FinEtoolsDeforLinear.MatDeforModule.tens4symmt6x6tot!Method
tens4symmt6x6tot!(ST::Array{T, 4}, M::Matrix{T}) where {T}

Convert a symmetric 6 x 6 matrix to a symmetric 4th-order tensor.

!!! Note The order corresponds to the arrangement of the components of stress (or strain) tensor, symmetric, three-dimensional, into a 6-component vector.

FinEtoolsDeforLinear.MatDeforModule.tens4symmtto6x6t!Method
tens4symmtto6x6t!(M::Matrix{T}, ST::Array{T, 4}) where {T}

Convert a symmetric 4th-order tensor to a 6 x 6 matrix.

!!! Note The order corresponds to the arrangement of the components of stress (or strain) tensor, symmetric, three-dimensional, into a 6-component vector.

Example

J=tens4_ijkl(eye(3),eye(3))
produces the tracor:
T=rand(3); 
sum(diag(T))*eye(3)
t= tens4_dot_2(J,T)
M= tens4_symm_to_6(ST)
FinEtoolsDeforLinear.MatDeforModule.tens4tracor!Method
tens4tracor!(t::Array{T, 4}) where {T}

Compute 4th-order tracor tensor.

Double contraction of a second order tensor with this fourth-order tensor produces the spherical part of the second order tensor.

Example

The product of the tracor tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4tracor!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show tr(S) * I - tS
FinEtoolsDeforLinear.MatDeforModule.tens4transposor!Method
tens4transposor!(t::Array{T, 4}) where {T}

Compute 4th-order transposor tensor.

Example

The product of the transposor tensor with the second-order tensor S is

t = fill(0.0, 3, 3, 3, 3)
tens4transposor!(t)
S = rand(3, 3)
tS = fill(0.0, 3, 3)
tens4dot2!(tS, t, S)
@show S' - tS

Elasticity

FinEtoolsDeforLinear.MatDeforLinearElasticModule.tangentmoduli!Method
tangentmoduli!(
    self::AbstractMatDeforLinearElastic,
    D::Matrix{FT},
    t::FT,
    dt::FT,
    loc::Matrix{FT},
    label::Int,
) where {FT}

Calculate the material stiffness matrix.

  • D = matrix of tangent moduli, supplied as a buffer and overwritten. Returned

as output.

FinEtoolsDeforLinear.MatDeforLinearElasticModule.update!Method
update!(
    self::AbstractMatDeforLinearElastic,
    stress::Vector{FT},
    output::Vector{FT},
    strain::Vector{FT},
    thstrain::Vector{FT} = zeros(6),
    t::FT = 0.0,
    dt::FT = 0.0,
    loc::Matrix{FT} = zeros(3, 1),
    label::Int = 0,
    quantity = :nothing,
) where {FT}

Update material state.

  • strain = strain vector,
  • thstrain = thermal strain vector,
  • t = current time,
  • dt = current time step,
  • loc = location of the quadrature point in global Cartesian coordinates,
  • label = label of the finite element in which the quadrature point is found.

Output

  • stress = stress vector, allocated by the caller with a size of the number of stress and

strain components, nstressstrain. The components of the stress vector are calculated and stored in the stress vector.

  • output = array which is (if necessary) allocated in an appropriate size, filled with the output quantity, and returned.

Isotropic elasticity

FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIsoMethod
MatDeforElastIso(
    mr::Type{MR},
    mass_density::N,
    E::N,
    nu::N,
    CTE::N,
) where {MR<:AbstractDeforModelRed,N<:Number}

Create an isotropic elastic material providing all material parameters.

Arguments

  • mr::Type{MR}: The type of the deformation model.
  • mass_density::N: The mass density of the material.
  • E::N: The Young's modulus of the material.
  • nu::N: The Poisson's ratio of the material.
  • CTE::N: The coefficient of thermal expansion of the material.

Orthotropic elasticity

FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrthoMethod
MatDeforElastOrtho(
    mr::Type{MR},
    E::N,
    nu::N,
    CTE::N,
) where {MR<:AbstractDeforModelRed,N<:Number}

Create elastic orthotropic material which is really isotropic.

Convenience version with only the specification of the elastic and thermal expansion properties.

FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrthoMethod
MatDeforElastOrtho(
    mr::Type{MR},
    mass_density::N,
    E::N,
    nu::N,
    CTE::N,
) where {MR<:AbstractDeforModelRed,N<:Number}

Create elastic orthotropic material which is really isotropic.

Convenience version with only the specification of the mass density, and the elastic and thermal expansion properties.

FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrthoMethod
MatDeforElastOrtho(
    mr::Type{MR},
    E1::N,
    E2::N,
    E3::N,
    nu12::N,
    nu13::N,
    nu23::N,
    G12::N,
    G13::N,
    G23::N,
) where {MR<:AbstractDeforModelRed,N<:Number}

Create elastic orthotropic material.

Convenience version with only the specification of the elastic properties.

Modules

FinEtoolsDeforLinear.FEMMDeforLinearESNICEModuleModule

Formulation for the small displacement, small strain deformation model for Nodally-Integrated Continuum Elements (NICE).

The approximation is originally from Dohrmann et al IJNME 47 (2000). The formulation was subsequently developed in Krysl, P. and Zhu, B. Locking-free continuum displacement finite elements with nodal integration, International Journal for Numerical Methods in Engineering, 76,7,1020-1043,2008.

The stabilization scheme comes from papers on energy-sampling stabilization for mean-strain elements (Krysl and coauthors).

FinEtoolsDeforLinear.FinEtoolsDeforLinearModule

FinEtoolsDeforLinear (C) 2017-2024, Petr Krysl

Finite Element tools. Julia implementation of the finite element method for continuum mechanics. Package for linear static and dynamic stress analysis problems.

FinEtoolsDeforLinear.FEMMDeforLinearNICEModuleModule

Formulation for the small displacement, small strain deformation model for Nodally-Integrated Continuum Elements (NICE).

The approximation is originally from Dohrmann et al IJNME 47 (2000). The formulation was subsequently developed in Krysl, P. and Zhu, B. Locking-free continuum displacement finite elements with nodal integration, International Journal for Numerical Methods in Engineering, 76,7,1020-1043,2008.

This formulation is at this point obsolete (replaced with ESNICE).

FinEtoolsDeforLinear.FEMMDeforWinklerModuleModule

Module for operations on boundaries of domains to construct system matrices and system vectors for linear deformation models with distributed-spring supports (Winkler foundation model).