Reference Manual
Simple FE model (volume)
FinEtools.FEMMBaseModule.inspectintegpoints
— Methodinspectintegpoints(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.
geom
- reference geometry fieldu
- displacement fielddT
- temperature difference fieldfelist
- 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); whereidat
- 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_gh
— Methodinfsup_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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.infsup_sh
— Methodinfsup_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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.mass
— Methodmass(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.stiffness
— Methodstiffness(self::FEMM, assembler::A, geom::NodalField{GFT}, u::NodalField{UFT}) where {FEMM<:AbstractFEMMDeforLinear, A<:AbstractSysmatAssembler, GFT<:Number, UFT<:Number}
Compute and assemble stiffness matrix.
The material stiffness matrix is assumed to be the same at all the points of the domain (homogeneous material).
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.thermalstrainloads
— Methodthermalstrainloads(self::AbstractFEMMDeforLinear, assembler::A, geom::NodalField{GFT}, u::NodalField{UFT}, dT::NodalField{TFT}) where {A<:AbstractSysvecAssembler, GFT<:Number, UFT<:Number, TFT<:Number}
Compute the thermal-strain load vector.
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.AbstractFEMMDeforLinear
— TypeAbstractFEMMDeforLinear <: AbstractFEMMBase
Abstract type of FEMM for linear deformation.
FinEtoolsDeforLinear.FEMMDeforLinearModule.FEMMDeforLinear
— TypeFEMMDeforLinear{MR<:AbstractDeforModelRed, S<:AbstractFESet, F<:Function, M<:AbstractMatDeforLinearElastic} <: AbstractFEMMDeforLinear
Class for linear deformation finite element modeling machine.
FinEtoolsDeforLinear.FEMMDeforLinearModule.FEMMDeforLinear
— MethodFEMMDeforLinear(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.surfacenormalspringstiffness
— Methodsurfacenormalspringstiffness(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.FEMMDeforWinklerModule.FEMMDeforWinkler
— TypeFEMMDeforWinkler{S<:AbstractFESet, F<:Function} <: AbstractFEMM
Type for normal spring support (Winkler).
FinEtoolsDeforLinear.FEMMDeforSurfaceDampingModule.dampingABC
— MethoddampingABC(self::FEMMDeforSurfaceDamping, assembler::A, geom::NodalField{GFT}, u::NodalField{UFT}, impedance::UFT, surfacenormal::SurfaceNormal) where {A<:AbstractSysmatAssembler, GFT<:Number, UFT<: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.
FinEtoolsDeforLinear.FEMMDeforSurfaceDampingModule.FEMMDeforSurfaceDamping
— TypeFEMMDeforSurfaceDamping{S<:AbstractFESet, F<:Function} <: AbstractFEMM
Type for surface damping model.
Advanced: Mean-strain FEM
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(self::F, geom::NodalField{GFT}) where {F<:FEMMDeforLinearMSH8, GFT}
Associate geometry field with the FEMM.
Compute the correction factors to account for the shape of the elements.
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(self::F, geom::NodalField{GFT}) where {F<:FEMMDeforLinearMST10, GFT}
Associate geometry field with the FEMM.
Compute the correction factors to account for the shape of the elements.
FinEtools.FEMMBaseModule.inspectintegpoints
— Methodinspectintegpoints(
self::AbstractFEMMDeforLinearMS,
geom::NodalField{GFT},
u::NodalField{UFT},
dT::NodalField{TFT},
felist::Vector{IT},
inspector::F,
idat,
quantity = :Cauchy;
context...,
) where {GFT<:Number,UFT<:Number,TFT<:Number,IT,F<:Function}
Inspect integration point quantities.
Arguments
geom
- reference geometry fieldu
- displacement fielddT
- temperature difference fieldfelist
- 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); whereidat
- 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.stiffness
— Methodstiffness(self::AbstractFEMMDeforLinearMS, assembler::A,
geom::NodalField{GFT},
u::NodalField{UFT}) where {A<:AbstractSysmatAssembler, GFT<:Number, UFT<:Number}
Compute and assemble stiffness matrix.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.infsup_gh
— Methodinfsup_gh(self::AbstractFEMMDeforLinearMS, assembler::A, geom::NodalField{GFT}, u::NodalField{UFT}) where {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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.infsup_sh
— Methodinfsup_sh(self::AbstractFEMMDeforLinearMS, 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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.AbstractFEMMDeforLinearMS
— TypeAbstractFEMMDeforLinearMS <: AbstractFEMMDeforLinear
Abstract type for mean-strain linear deformation FEMM.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMSH8
— Typemutable 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.FEMMDeforLinearMSH8
— MethodFEMMDeforLinearMSH8(
mr::Type{MR},
integdomain::ID,
mcsys::CS,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetH8}, CS<:CSys, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMSH8
— MethodFEMMDeforLinearMSH8(
mr::Type{MR},
integdomain::ID,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetH8}, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMST10
— Typemutable 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.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMST10
— MethodFEMMDeforLinearMST10(
mr::Type{MR},
integdomain::ID,
mcsys::CS,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetT10}, CS<:CSys, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule.FEMMDeforLinearMST10
— MethodFEMMDeforLinearMST10(
mr::Type{MR},
integdomain::ID,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetT10}, M<:AbstractMatDeforLinearElastic}
Constructor.
Advanced: Nodal integration
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(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.
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.stiffness
— Methodstiffness(self::AbstractFEMMDeforLinearNICE,
assembler::A,
geom::NodalField{GFT},
u::NodalField{UFT}) where {A <: AbstractSysmatAssembler, GFT <: Number, UFT <: Number}
Compute and assemble stiffness matrix.
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule.stiffness
— Methodstiffness(self::AbstractFEMMDeforLinearNICE,
geom::NodalField{GFT},
u::NodalField{UFT}) where {GFT <: Number, UFT <: Number}
Compute and assemble stiffness matrix.
FinEtoolsDeforLinear.FEMMDeforLinearNICEModule.AbstractFEMMDeforLinearNICE
— TypeAbstractFEMMDeforLinearNICE <: AbstractFEMMDeforLinear
Abstract FEMM type for Nodally Integrated Continuum Elements (NICE).
FinEtoolsDeforLinear.FEMMDeforLinearNICEModule.FEMMDeforLinearNICEH8
— Typemutable struct FEMMDeforLinearNICEH8{
MR <: AbstractDeforModelRed,
S <: FESetH8,
F <: Function,
M <: AbstractMatDeforLinearElastic,
} <: AbstractFEMMDeforLinearNICE
FEMM type for Nodally Integrated Continuum Elements (NICE) based on the eight-node hexahedron.
FinEtoolsDeforLinear.FEMMDeforLinearNICEModule.FEMMDeforLinearNICET4
— Typemutable struct FEMMDeforLinearNICET4{
MR <: AbstractDeforModelRed,
S <: FESetT4,
F <: Function,
M <: AbstractMatDeforLinearElastic,
} <: AbstractFEMMDeforLinearNICE
FEMM type for Nodally Integrated Continuum Elements (NICE) based on the 4-node tetrahedron.
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(self::F, geom::NodalField{GFT}) where {F<:FEMMDeforLinearESNICEH8, GFT}
Associate geometry field with the FEMM.
Compute the correction factors to account for the shape of the elements.
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(self::F, geom::NodalField{GFT}; stabilization_parameters = _T4_stabilization_parameters) where {F<:FEMMDeforLinearESNICET4, GFT}
Associate geometry field with the FEMM.
Compute the correction factors to account for the shape of the elements.
FinEtools.FEMMBaseModule.inspectintegpoints
— Methodinspectintegpoints(self::AbstractFEMMDeforLinearESNICE, geom::NodalField{GFT}, u::NodalField{UFT}, dT::NodalField{TFT}, felist::Vector{IT}, inspector::F, idat, quantity=:Cauchy; context...) where {GFT<:Number, UFT<:Number, TFT<:Number, IT<:Integer, F<:Function}
Inspect integration point quantities.
Arguments
geom
- reference geometry fieldu
- displacement fielddT
- temperature difference fieldfelist
- 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); whereidat
- 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.stiffness
— Methodstiffness(self::AbstractFEMMDeforLinearESNICE, assembler::A, geom::NodalField{GFT}, u::NodalField{T}) where {A<:AbstractSysmatAssembler, GFT<:Number, T<:Number}
Compute and assemble stiffness matrix.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.infsup_gh
— Methodinfsup_gh(self::AbstractFEMMDeforLinearESNICE, assembler::A, geom::NodalField{GFT}, u::NodalField{UFT}) where {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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.infsup_sh
— Methodinfsup_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.)
This computation has not been optimized in any way. It can be expected to be inefficient.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.AbstractFEMMDeforLinearESNICE
— TypeAbstractFEMMDeforLinearESNICE <: AbstractFEMMDeforLinear
Abstract FEMM type for Nodally Integrated Continuum Elements (NICE) with energy-sampling stabilization.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICEH8
— Typemutable 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 (NICE) based on the p-node hexahedron.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICEH8
— MethodFEMMDeforLinearESNICEH8(
mr::Type{MR},
integdomain::ID,
mcsys::CS,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetH8}, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICEH8
— MethodFEMMDeforLinearESNICEH8(
mr::Type{MR},
integdomain::ID,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetH8}, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICET4
— Typemutable 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 (NICE) based on the 4-node tetrahedron.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICET4
— MethodFEMMDeforLinearESNICET4(
mr::Type{MR},
integdomain::ID,
mcsys::CS,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetT4}, M<:AbstractMatDeforLinearElastic}
Constructor.
FinEtoolsDeforLinear.FEMMDeforLinearESNICEModule.FEMMDeforLinearESNICET4
— MethodFEMMDeforLinearESNICET4(
mr::Type{MR},
integdomain::ID,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S} where {S<:FESetT4}, M<:AbstractMatDeforLinearElastic}
Constructor.
Advanced: Incompatible modes
FinEtools.FEMMBaseModule.associategeometry!
— Methodassociategeometry!(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.FEMMDeforLinearBaseModule.stiffness
— Methodstiffness(self::AbstractFEMMDeforLinearIM, assembler::A, geom::NodalField{FFlt}, u::NodalField{T}) where {A<:AbstractSysmatAssembler, T<:Number}
Compute and assemble stiffness matrix.
FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8
— TypeFEMMDeforLinearIMH8{MR<:AbstractDeforModelRed, S<:FESetH8, F<:Function, M<:AbstractMatDeforLinearElastic}
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.FEMMDeforLinearIMH8
— MethodFEMMDeforLinearIMH8(
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.
FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8
— MethodFEMMDeforLinearIMH8(
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.
FinEtoolsDeforLinear.FEMMDeforLinearIMModule.FEMMDeforLinearIMH8
— MethodFEMMDeforLinearIMH8(
mr::Type{MR},
integdomain::ID,
material::M,
) where {MR<:AbstractDeforModelRed, ID<:IntegDomain{S,F} where {S<:FESetH8,F<:Function}, M<:AbstractMatDeforLinearElastic}
Constructor.
Algorithms
Linear deformation
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportdeformation
— MethodAlgoDeforLinearModule.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 # Defaults"field"
- nodal or elemental field
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportmode
— MethodAlgoDeforLinearModule.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 lengthneigvs
"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"]
= seeexportdeformation()
FinEtoolsDeforLinear.AlgoDeforLinearModule.exportstress
— MethodAlgoDeforLinearModule.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 toinspectintegpoints()
"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.exportstresselementwise
— MethodAlgoDeforLinearModule.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.linearstatics
— MethodAlgoDeforLinearModule.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 functionw = 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 aForceIntensity
object, or it could be theForceIntensity
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 # For body loads (optional):"work"
= work of the applied loads # modeldata.bodyload = cell array of struct,"timing"
= dictionary with timing results # each piece of the domain can have each its own body load
FinEtoolsDeforLinear.AlgoDeforLinearModule.modal
— MethodAlgoDeforLinearModule.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.dett
— Methoddett(::Type{DeforModelRed2DStrain}, C::Matrix{T}) where {T}
Compute the determinant of a general square matrix.
FinEtoolsDeforLinear.MatDeforModule.dett
— Methoddett(::Type{DeforModelRed3D}, C::Matrix{T}) where {T}
Compute the determinant of a general square matrix.
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!
— Methodrotstressvec!(::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 insideinstress
= input stress vectorRm
= columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!
— Methodrotstressvec!(::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 insideinstress
= input stress vectorRm
= columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!
— Methodrotstressvec!(::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 insideinstress
= input stress vectorRm
= columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!
— Methodrotstressvec!(::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 insideinstress
= input stress vectorRm
= columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.rotstressvec!
— Methodrotstressvec!(::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 insideinstress
= input stress vectorRm
= columns are components of 'bar' basis vectors on the 'plain' basis vectors
FinEtoolsDeforLinear.MatDeforModule.strainttov!
— Methodstrainttov!(::Type{DeforModelRed2DStrain}, v::Vector{T}, t::Matrix{T}) where {T}
Convert a symmetric matrix of 2x2 strain components into a 3-component vector.
FinEtoolsDeforLinear.MatDeforModule.strainttov!
— Methodstrainttov!(::Type{DeforModelRed3D}, v::Vector{T}, t::Matrix{T}) where {T}
Convert a symmetric matrix of 3x3 strain components into a 6-component vector.
FinEtoolsDeforLinear.MatDeforModule.strainvdet
— Methodstrainvdet(::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.strainvdet
— Methodstrainvdet(::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.strainvtot!
— Methodstrainvtot!(::Type{DeforModelRed2DStrain}, t::Matrix{T}, v::Vector{T}) where {T}
Convert a strain 3-vector to a matrix of 2x2 strain components (symmetric tensor).
FinEtoolsDeforLinear.MatDeforModule.strainvtot!
— Methodstrainvtot!(::Type{DeforModelRed3D}, t::Matrix{T}, v::Vector{T}) where {T}
Convert a strain 3-vector to a matrix of 2x2 strain components (symmetric tensor).
FinEtoolsDeforLinear.MatDeforModule.strainvtr
— Methodstrainvtr(::Type{DeforModelRed2DStrain}, Cv::Vector{T}) where {T}
Compute the trace of a symmetric strain-like square matrix represented as a vector.
FinEtoolsDeforLinear.MatDeforModule.strainvtr
— Methodstrainvtr(::Type{DeforModelRed3D}, Cv::Vector{T}) where {T}
Compute the trace of a symmetric strain-like square matrix represented as a vector.
FinEtoolsDeforLinear.MatDeforModule.stressttov!
— Methodstressttov!(::Type{DeforModelRed2DStrain}, v::Vector{T}, t::Matrix{T}) where {T}
Convert a symmetric matrix of 2x2 stress components to a 3-component vector.
FinEtoolsDeforLinear.MatDeforModule.stressttov!
— Methodstressttov!(::Type{DeforModelRed2DStress}, v::Vector{T}, t::Matrix{T}) where {T}
Convert a symmetric matrix of 2x2 stress components to a 3-component vector.
FinEtoolsDeforLinear.MatDeforModule.stressttov!
— Methodstressttov!(::Type{DeforModelRed3D}, v::Vector{T}, t::Matrix{T}) where {T}
Convert a symmetric matrix of 3x3 stress components to a 6-component vector.
FinEtoolsDeforLinear.MatDeforModule.stressvtot!
— Methodstressvtot!(::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!
— Methodstressvtot!(::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.stressvtot!
— Methodstressvtot!(::Type{DeforModelRed2DStress}, t::Matrix{T}, v::Vector{T}) where {T}
Convert a 3-vector to a matrix of 2x2 stress components (symmetric tensor).
FinEtoolsDeforLinear.MatDeforModule.stressvtot!
— Methodstressvtot!(::Type{DeforModelRed3D}, t::Matrix{T}, v::Vector{T}) where {T}
Convert a 6-vector to a matrix of 3x3 stress components (symmetric tensor).
FinEtoolsDeforLinear.MatDeforModule.tens4checksymmetry
— Methodtens4checksymmetry(C4th)
If the fourth-order tensor of material elasticity has the full set of symmetries, return true; otherwise false.
FinEtoolsDeforLinear.MatDeforModule.tens4deviator!
— Methodtens4deviator!(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!
— Methodtens4dot2!(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.
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!
— Methodtens4identity!(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!
— Methodtens4ijkl!(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)
.
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!
— Methodtens4ikjl!(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)
.
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!
— Methodtens4iljk!(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)
.
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!
— Methodtens4skewor!(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!
— Methodtens4symmetrizor!(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!
— Methodtens4symmt6x6tot!(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!
— Methodtens4symmtto6x6t!(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!
— Methodtens4tracor!(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!
— Methodtens4transposor!(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
FinEtoolsDeforLinear.MatDeforModule.AbstractMatDefor
— TypeAbstractMatDefor
Abstract type that represents deformable materials.
Elasticity
FinEtoolsDeforLinear.MatDeforLinearElasticModule.tangentmoduli!
— Methodtangentmoduli!(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.thermalstrain!
— Methodthermalstrain!(self::AbstractMatDeforLinearElastic, thstrain::Vector{FT}, dT= 0.0) where {FT}
Compute thermal strain from the supplied temperature increment.
thstrain
= thermal strain vector, supplied as buffer, returned as output.
FinEtoolsDeforLinear.MatDeforLinearElasticModule.update!
— Methodupdate!(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.
FinEtoolsDeforLinear.MatDeforLinearElasticModule.AbstractMatDeforLinearElastic
— TypeAbstractMatDeforLinearElastic <: AbstractMatDefor
Abstract Linear Elasticity material.
Isotropic elasticity
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— TypeMatDeforElastIso{MR<:AbstractDeforModelRed, MTAN<:Function, MUPD<:Function, MTHS<:Function} <: AbstractMatDeforLinearElastic
Linear isotropic elasticity material.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed1DStrain}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 1D models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed1DStress}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 1D models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed2DAxisymm}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 2D axially symmetric models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed2DStrain}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 2D plane strain models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed2DStress}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 2D plane stress models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{DeforModelRed3D}, args::NTuple{4, FT}) where FT
Create elastic isotropic material for 3D stress models.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{MR}, E, nu) where {MR}
Create isotropic elastic material with default mass density and thermal expansion.
FinEtoolsDeforLinear.MatDeforElastIsoModule.MatDeforElastIso
— MethodMatDeforElastIso(mr::Type{MR}, mass_density, E, nu, CTE) where {MR<:AbstractDeforModelRed}
Create an isotropic elastic material providing all material parameters.
Orthotropic elasticity
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— TypeMatDeforElastOrtho{MR<:AbstractDeforModelRed, MTAN<:Function, MUPD<:Function, MTHS<:Function} <: AbstractMatDeforLinearElastic
Linear orthotropic elasticity material.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(mr::Type{DeforModelRed2DAxisymm}, args::NTuple{13, FT}) where FT
Create elastic orthotropic material for 2D axially symmetric models.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(mr::Type{DeforModelRed2DStrain}, args::NTuple{13, FT}) where FT
Create elastic orthotropic material for 2D plane strain models.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(mr::Type{DeforModelRed2DStress}, args::NTuple{13, FT}) where FT
Create elastic orthotropic material for 2D plane stress models.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(mr::Type{DeforModelRed3D}, args::NTuple{13, FT}) where FT
Create elastic orthotropic material for 3D models.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(
mr::Type{MR},
E::N,
nu::N,
) where {MR<:AbstractDeforModelRed,N<:Number}
Create elastic orthotropic material which is really isotropic.
Convenience version with only the specification of the elastic properties.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(
mr::Type{MR},
mass_density::N,
E1::N,
E2::N,
E3::N,
nu12::N,
nu13::N,
nu23::N,
G12::N,
G13::N,
G23::N,
CTE1::N,
CTE2::N,
CTE3::N,
) where {MR<:AbstractDeforModelRed,N<:Number}
Create elastic orthotropic material.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(
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 elastic and thermal expansion properties.
FinEtoolsDeforLinear.MatDeforElastOrthoModule.MatDeforElastOrtho
— MethodMatDeforElastOrtho(
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.FEMMDeforLinearESNICEModule
— ModuleFEMMDeforLinearESNICEModule
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).
FinEtoolsDeforLinear.FEMMDeforLinearBaseModule
— ModuleAbstractFEMMDeforLinearBaseModule
Base module for operations on interiors of domains to construct system matrices and system vectors for linear deformation models.
FinEtoolsDeforLinear.FinEtoolsDeforLinear
— ModuleFinEtoolsDeforLinear (C) 2017-2024, Petr Krysl
Finite Element tools. Julia implementation of the finite element method for continuum mechanics. Package for linear stress analysis problems.
FinEtoolsDeforLinear.FEMMDeforLinearMSModule
— ModuleFEMMDeforLinearMSModule
Module for operations on interiors of domains to construct system matrices and system vectors for linear deformation models: mean-strain formulation.
FinEtoolsDeforLinear.MatDeforModule
— ModuleMatDeforModule
Module to support general operations for deformation material models.
FinEtoolsDeforLinear.FEMMDeforLinearModule
— ModuleFEMMDeforLinearModule
Module for operations on interiors of domains to construct system matrices and system vectors for linear deformation models.
FinEtoolsDeforLinear.AlgoDeforLinearModule
— ModuleAlgoDeforLinearModule
Module for algorithms used in linear deformation models.
FinEtoolsDeforLinear.FEMMDeforLinearNICEModule
— ModuleFEMMDeforLinearNICEModule
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.
FinEtoolsDeforLinear.FEMMDeforSurfaceDampingModule
— ModuleFEMMDeforSurfaceDampingModule
Module for operations on the damping associated with absorbing boundary conditions (ABC) representation of the effect of infinite extent of inviscid fluid next to the surface.
Missing docstring for FinEtoolsDeforLinear.MatDeforElastIsoModule
. Check Documenter's build log for details.
Missing docstring for FinEtoolsDeforLinear.MatDeforLinearElasticModule
. Check Documenter's build log for details.
FinEtoolsDeforLinear.FEMMDeforWinklerModule
— ModuleFEMMDeforWinklerModule
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).
FinEtoolsDeforLinear.MatDeforElastOrthoModule
— ModuleMatDeforElastOrthoModule
Module for orthotropic elastic material.
FinEtoolsDeforLinear.FEMMDeforLinearIMModule
— ModuleFEMMDeforLinearIMModule
Module for operations on interiors of domains to construct system matrices and system vectors for linear deformation models: incompatible-mode formulation.