Functions

Helpers

FinEtools.AssemblyModule.assemble!Method
assemble!(self::SysvecAssemblerOpt{T}, vec::MV,
  dofnums::D) where {T<:Number, MV<:AbstractArray{T}, D<:AbstractArray{FInt}}

Assemble an elementwise vector.

The method assembles a column element vector using the vector of degree of freedom numbers for the rows.

source
FinEtools.AssemblyModule.startassembly!Method
startassembly!(self::SysvecAssemblerOpt{T},
  ndofs_row::FInt) where {T<:Number}

Start assembly.

The method makes the buffer for the vector assembly. It must be called before the first call to the method assemble.

ndofs_row= Total number of degrees of freedom.

source

FEM machines

Nonlinear deformation

Simple FE models

FinEtoolsDeforNonlinear.FEMMDeforNonlinearBaseModule.geostiffnessMethod
geostiffness(self::AbstractFEMMDeforNonlinear, assembler::A, geom::NodalField{FFlt}, un1::NodalField{T}, un::NodalField{T}, tn::FFlt, dtn::FFlt) where {A<:AbstractSysmatAssembler, T<:Number}

Compute and assemble geometric stiffness matrix.

Arguments

  • assembler = matrix assembler,
  • geom = geometry field: reference coordinates of the nodes,
  • un1 = displacement field at time tn1 = tn + dtn,
  • un = displacement field at time tn,
  • tn = time in step n
  • dtn = increment of time
source
FinEtoolsDeforNonlinear.FEMMDeforNonlinearBaseModule.nzebcloadsMethod
nzebcloads(self::AbstractFEMMDeforNonlinear, assembler::A, geom::NodalField{FFlt}, un1::NodalField{T}, un::NodalField{T}, du::NodalField{T}, tn::FFlt, dtn::FFlt) where {A<:SysvecAssembler, T<:Number}

Compute and assemble load vector due to prescribed increment of displacements.

Arguments

  • assembler = matrix assembler,
  • geom = geometry field: reference coordinates of the nodes,
  • un1 = displacement field at time tn1 = tn + dtn,
  • un = displacement field at time tn,
  • du = field of displacement increment prescribed at time tn1 = tn + dtn, The increment is stored in the fixed values of this field.
  • tn = time in step n
  • dtn = increment of time
source
FinEtoolsDeforNonlinear.FEMMDeforNonlinearBaseModule.restoringforceMethod
restoringforce(self::AbstractFEMMDeforNonlinear, assembler::A, geom::NodalField{FFlt}, un1::NodalField{T}, un::NodalField{T}, tn::FFlt, dtn::FFlt, savestate = false) where {A<:AbstractSysvecAssembler, T<:Number}

Compute the restoring force vector.

Note: This method UPDATES the state of the FEMM object. In particular, the material state gets updated.

Arguments

  • assembler = vector assembler,
  • geom = geometry field: reference coordinates of the nodes,
  • un1 = displacement field at time tn1 = tn + dtn,
  • un = displacement field at time tn,
  • tn = time in step n
  • dtn = increment of time
  • savestate = bool flag: should we modify the material states (savestate = true)? Otherwise work with a copy of the material state.
source
FinEtoolsDeforNonlinear.FEMMDeforNonlinearBaseModule.stiffnessMethod
stiffness(self::AbstractFEMMDeforNonlinear, assembler::A, geom::NodalField{FFlt}, un1::NodalField{T}, un::NodalField{T}, tn::FFlt, dtn::FFlt) where {A<:AbstractSysmatAssembler, T<:Number}

Compute and assemble stiffness matrix.

Arguments

  • assembler = matrix assembler,
  • geom = geometry field: reference coordinates of the nodes,
  • un1 = displacement field at time tn1 = tn + dtn,
  • un = displacement field at time tn,
  • tn = time in step n
  • dtn = increment of time
source

Algorithms

Nonlinear deformation

FinEtoolsDeforNonlinear.AlgoDeforNonlinearModule.nonlinearstaticsMethod
AlgoDeforNonlinearModule.nonlinearstatics(modeldata::FDataDict)

Algorithm for static nonlinear deformation (stress) analysis.

The algorithm chooses steps from the array of load multipliers: the step takes it precisely from the preceding step to the next step in one go.

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" = when this key is not present, the assumption is that the displacement is fixed at zero (0). Otherwise, this needs to be set to a function with signature f(x, lambda), where x is the location of the node in the reference configuration, and lambda is the load factor. In other words, whenever this quantity is supplied, it is implied that the displacement depends on the load factor.
  • "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

  • "femm" = finite element model machine (mandatory);
  • "traction_vector" = traction vector, a force-intensity (ForceIntensity) object.

Control parameters The following attributes may be supplied:

  • "load_multipliers" = For what load multipliers should the solution be calculated? Array of monotonically increasing numbers.
  • "line_search" = Should we use line search? Boolean. Default = true.
  • "maxdu_tol" = Tolerance on the magnitude of the largest incremental displacement component.
  • "maxbal_tol" = Tolerance on the magnitude of the largest out-of-balance force component.
  • "iteration_observer" = observer function to be called after each iteration. Default is to do nothing. The observer function has a signature iteration_observer(lambda,iter,du,modeldata) where lambda is the current load factor, iter is the iteration number, du is the nodal field of current displacement increments.
  • "increment_observer" = observer function to be called after convergence is reached in each step (optional) The observer function has a signature output(lambda, modeldata) where lambda is the current load factor. Default is to do nothing. The increment observer can refer to the following key-value pairs in modeldata: - "un1": converged displacement in current step - "un": converged displacement in the last step - "t": current value of load factor - "dt": increment of the load factor

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
  • "reactions" = computed reaction nodal field
  • "timing" = dictionary with timing results
source

Material models

Material models for nonlinear deformation

FinEtoolsDeforNonlinear.MatDeforNonlinearModule.tangentmoduli!Method
tangentmoduli!(self::M, D::FFltMat, statev::FFltVec, Fn1::FFltMat,
    Fn::FFltMat, tn::FFlt, dtn::FFlt, loc::FFltMat, label::FInt)
    where {M<:AbstractMatDeforNonlinear}

Calculate the material stiffness matrix.

Arguments

  • self = material
  • D = matrix of tangent moduli, supplied as a buffer and overwritten. Returned as output.
  • statev = material state vector, the content of this vector must not change inside this function.
  • Fn1 = deformation gradient at time tn1 = tn + dtn,
  • Fn = deformation gradient at time tn,
  • tn = time in step n
  • dtn = increment of time
  • loc = location of the integration point in the reference coordinates (time t0),
  • label = label of the element containing the integration point
Note

The deformation gradients and the matrix of the tangent moduli are expressed with respect to the local material coordinate system.

source
FinEtoolsDeforNonlinear.MatDeforNonlinearModule.totlag2curr!Method
totlag2curr!(c, C, F)

Convert a total Lagrangean constitutive matrix to a current Lagrangean one (sometimes known as "Eulerian").

  • C = Lagrangean constitutive matrix, 6x6, symmetric
  • F = current deformation gradient, FiJ = partial xi / partial X_J

The transformation is cijkl = 1/J CIJKL FiI FjJ FkK FlL. In the present case the fourth-order tensor is represented with a 6 x 6 matrix.

source
FinEtoolsDeforNonlinear.MatDeforNonlinearModule.totlag2curr4th!Method
totlag2curr4th!(c, C, F)

Convert a total Lagrangean constitutive matrix to a current Lagrangean one (sometimes known as "Eulerian").

  • C = Lagrangean constitutive matrix, fourth-order tensor
  • F = current deformation gradient, FiJ = partial xi / partial X_J

The transformation is cijkl = 1/J CIJKL FiI FjJ FkK FlL. Both the input and the output are fourth-order tensors.

source
FinEtoolsDeforNonlinear.MatDeforNonlinearModule.totlag2currsymm!Method
totlag2currsymm!(c, C, F)

Convert a total Lagrangean constitutive matrix to a current Lagrangean one (sometimes known as "Eulerian").

  • C = Lagrangean constitutive matrix, 6x6, symmetric
  • F = current deformation gradient, FiJ = partial xi / partial X_J

The transformation is cijkl = 1/J CIJKL FiI FjJ FkK FlL. In the present case the fourth-order tensor is represented with a 6 x 6 matrix.

Note

The Lagrangean material stiffness matrices, both input and output, are presumed symmetric.

source
FinEtoolsDeforNonlinear.MatDeforNonlinearModule.update!Method
update!(self::M, statev::FFltVec, stress::FFltVec, output::FFltVec,
    Fn1::FFltMat, Fn::FFltMat, tn::FFlt, dtn::FFlt,
    loc::FFltMat=zeros(3,1), label::FInt=0, quantity=:nothing)
    where {M<:AbstractMatDeforNonlinear}

Update material state.

Arguments

  • self = material
  • statev = state variables: array which is (if necessary) allocated in an appropriate size, filled with the state variables, and returned. The contents of this vector may change as the state of the material is updated by the logic inside this function. If this change is to be saved, it must happen outside of this function.
  • cauchy = Cauchy 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.
  • Fn1 = deformation gradient at time tn1 = tn + dtn,
  • Fn = deformation gradient at time tn,
  • tn = time in step n
  • dtn = increment of time
  • loc = location of the integration point in the reference coordinates (time t0),
  • label = label of the element containing the integration point

Output

  • cauchy = Cauchy stress vector
  • output = output array
Note

The deformation gradients and the stress vector are expressed with respect to the local material coordinate system.

source

Material models for neohookean hyperelasticity

Material models for St Venant-Kirchhoff hyperelasticity

Material models for Mooney-Rivlin hyperelasticity