TopOptProblems

This sub-module of TopOpt defines a number of standard topology optimization problems for the convenient testing of algorithms.

Problem types

Abstract type

StiffnessTopOptProblem is an abstract type that a number of linear elasticity, quasi-static, topology optimization problems subtype.

TopOpt.TopOptProblems.StiffnessTopOptProblemType
abstract type StiffnessTopOptProblem{dim, T} <: AbstractTopOptProblem end

An abstract stiffness topology optimization problem. All subtypes must have the following fields:

  • ch: a Ferrite.ConstraintHandler struct
  • metadata: Metadata having various cell-node-dof relationships
  • black: a BitVector of length equal to the number of elements where black[e] is 1 iff the e^th element must be part of the final design
  • white: a BitVector of length equal to the number of elements where white[e] is 1 iff the e^th element must not be part of the final design
  • varind: an AbstractVector{Int} of length equal to the number of elements where varind[e] gives the index of the decision variable corresponding to element e. Because some elements can be fixed to be black or white, not every element has a decision variable associated.
source

Test problems

The following types are all concrete subtypes of StiffnessTopOptProblem. PointLoadCantilever is a cantilever beam problem with a point load as shown below. HalfMBB is the half Messerschmitt-Bölkow-Blohm (MBB) beam problem commonly used in topology optimization literature. LBeam and TieBeam are the common L-beam and tie-beam test problem used in topology optimization literature. The PointLoadCantilever and HalfMBB problems can be either 2D or 3D depending on the type of the inputs to the constructor. If the number of elements and sizes of elements are 2-tuples, the problem constructed will be 2D. And if they are 3-tuples, the problem constructed will be 3D. For the 3D versions, the point loads are applied at approximately the mid-depth point. The TieBeam and LBeam problems are always 2D.

Missing docstring.

Missing docstring for PointLoadCantilever. Check Documenter's build log for details.

TopOpt.TopOptProblems.PointLoadCantileverMethod
PointLoadCantilever(::Type{Val{CellType}}, nels::NTuple{dim,Int}, sizes::NTuple{dim}, E, ν, force) where {dim, CellType}
  • dim: dimension of the problem
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the center right of the cantilever beam (positive is downward)
  • nels: number of elements in each direction, a 2-tuple for 2D problems and a 3-tuple for 3D problems
  • sizes: the size of each element in each direction, a 2-tuple for 2D problems and a 3-tuple for 3D problems
  • CellType: can be either :Linear or :Quadratic to determine the order of the geometric and field basis functions and element type. Only isoparametric elements are supported for now.

Example:

nels = (60,20);
sizes = (1.0,1.0);
E = 1.0;
ν = 0.3;
force = 1.0;

# Linear elements and linear basis functions
celltype = :Linear

# Quadratic elements and quadratic basis functions
#celltype = :Quadratic

problem = PointLoadCantilever(Val{celltype}, nels, sizes, E, ν, force)
source
TopOpt.TopOptProblems.HalfMBBType
 |
 |
 v
O*********************************
O*                               *
O*                               *
O*                               *
O*********************************
                                 O

struct HalfMBB{dim, T, N, M} <: StiffnessTopOptProblem{dim, T}
    rect_grid::RectilinearGrid{dim, T, N, M}
    E::T
    ν::T
    ch::ConstraintHandler{<:DofHandler{dim, <:Cell{dim,N,M}, T}, T}
    force::T
    force_dof::Integer
    black::AbstractVector
    white::AbstractVector
    varind::AbstractVector{Int}
    metadata::Metadata
end
  • dim: dimension of the problem
  • T: number type for computations and coordinates
  • N: number of nodes in a cell of the grid
  • M: number of faces in a cell of the grid
  • rect_grid: a RectilinearGrid struct
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the top left of half the MBB (positive is downward)
  • force_dof: dof number at which the force is applied
  • ch: a Ferrite.ConstraintHandler struct
  • metadata: Metadata having various cell-node-dof relationships
  • black: a BitVector of length equal to the number of elements where black[e] is 1 iff the e^th element must be part of the final design
  • white: a BitVector of length equal to the number of elements where white[e] is 1 iff the e^th element must not be part of the final design
  • varind: an AbstractVector{Int} of length equal to the number of elements where varind[e] gives the index of the decision variable corresponding to element e. Because some elements can be fixed to be black or white, not every element has a decision variable associated.
source
TopOpt.TopOptProblems.HalfMBBMethod
HalfMBB(::Type{Val{CellType}}, nels::NTuple{dim,Int}, sizes::NTuple{dim}, E, ν, force) where {dim, CellType}
  • dim: dimension of the problem
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the top left of half the MBB (positive is downward)
  • nels: number of elements in each direction, a 2-tuple for 2D problems and a 3-tuple for 3D problems
  • sizes: the size of each element in each direction, a 2-tuple for 2D problems and a 3-tuple for 3D problems
  • CellType: can be either :Linear or :Quadratic to determine the order of the geometric and field basis functions and element type. Only isoparametric elements are supported for now.

Example:

nels = (60,20);
sizes = (1.0,1.0);
E = 1.0;
ν = 0.3;
force = -1.0;

# Linear elements and linear basis functions
celltype = :Linear

# Quadratic elements and quadratic basis functions
#celltype = :Quadratic

problem = HalfMBB(Val{celltype}, nels, sizes, E, ν, force)
source
TopOpt.TopOptProblems.LBeamType
////////////
............
.          .
.          .
.          . 
.          .                    
.          ......................
.                               .
.                               . 
.                               . |
................................. v
                                force

struct LBeam{T, N, M} <: StiffnessTopOptProblem{2, T}
    E::T
    ν::T
    ch::ConstraintHandler{<:DofHandler{2, <:Cell{2,N,M}, T}, T}
    force::T
    force_dof::Integer
    black::AbstractVector
    white::AbstractVector
    varind::AbstractVector{Int}
    metadata::Metadata
end
  • T: number type for computations and coordinates
  • N: number of nodes in a cell of the grid
  • M: number of faces in a cell of the grid
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the center right of the cantilever beam (positive is downward)
  • force_dof: dof number at which the force is applied
  • ch: a Ferrite.ConstraintHandler struct
  • metadata: Metadata having various cell-node-dof relationships
  • black: a BitVector of length equal to the number of elements where black[e] is 1 iff the e^th element must be part of the final design
  • white: a BitVector of length equal to the number of elements where white[e] is 1 iff the e^th element must not be part of the final design
  • varind: an AbstractVector{Int} of length equal to the number of elements where varind[e] gives the index of the decision variable corresponding to element e. Because some elements can be fixed to be black or white, not every element has a decision variable associated.
source
TopOpt.TopOptProblems.LBeamMethod
LBeam(::Type{Val{CellType}}, ::Type{T}=Float64; length = 100, height = 100, upperslab = 50, lowerslab = 50, E = 1.0, ν = 0.3, force = 1.0) where {T, CellType}
  • T: number type for computations and coordinates
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the center right of the cantilever beam (positive is downward)
  • length, height, upperslab and lowerslab are explained in LGrid.
  • CellType: can be either :Linear or :Quadratic to determine the order of the geometric and field basis functions and element type. Only isoparametric elements are supported for now.

Example:

E = 1.0;
ν = 0.3;
force = 1.0;

# Linear elements and linear basis functions
celltype = :Linear

# Quadratic elements and quadratic basis functions
#celltype = :Quadratic

problem = LBeam(Val{celltype}, E = E, ν = ν, force = force)
source
TopOpt.TopOptProblems.TieBeamType
                                                               1
                                                               
                                                              OOO
                                                              ...
                                                              . .
                                                           4  . . 
                                30                            . .   
/ . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . <-
/ .                                                                 . <- 2 f 
/ .    3                                                            . <- 
/ . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . <-
                                                              ^^^
                                                              |||
                                                              1 f

struct TieBeam{T, N, M} <: StiffnessTopOptProblem{2, T}
    E::T
    ν::T
    force::T
    ch::ConstraintHandler{<:DofHandler{2, N, T, M}, T}
    black::AbstractVector
    white::AbstractVector
    varind::AbstractVector{Int}
    metadata::Metadata
end
  • T: number type for computations and coordinates
  • N: number of nodes in a cell of the grid
  • M: number of faces in a cell of the grid
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the center right of the cantilever beam (positive is downward)
  • ch: a Ferrite.ConstraintHandler struct
  • metadata: Metadata having various cell-node-dof relationships
  • black: a BitVector of length equal to the number of elements where black[e] is 1 iff the e^th element must be part of the final design
  • white: a BitVector of length equal to the number of elements where white[e] is 1 iff the e^th element must not be part of the final design
  • varind: an AbstractVector{Int} of length equal to the number of elements where varind[e] gives the index of the decision variable corresponding to element e. Because some elements can be fixed to be black or white, not every element has a decision variable associated.
source
TopOpt.TopOptProblems.TieBeamMethod
TieBeam(::Type{Val{CellType}}, ::Type{T} = Float64, refine = 1, force = T(1); E = T(1), ν = T(0.3)) where {T, CellType}
  • T: number type for computations and coordinates
  • E: Young's modulus
  • ν: Poisson's ration
  • force: force at the center right of the cantilever beam (positive is downward)
  • refine: an integer value of 1 or greater that specifies the mesh refinement extent. A value of 1 gives the standard tie-beam problem in literature.
  • CellType: can be either :Linear or :Quadratic to determine the order of the geometric and field basis functions and element type. Only isoparametric elements are supported for now.
source

Reading INP Files

In TopOpt.jl, you can import a .inp file to an instance of the problem struct InpStiffness. This can be used to construct problems with arbitrary unstructured ground meshes, complex boundary condition domains and load specifications. The .inp file can be exported from a number of common finite element software such as: FreeCAD or ABAQUS.

TopOpt.TopOptProblems.InputOutput.INP.InpStiffnessType
struct InpStiffness{dim, N, TF, TI, TBool, Tch <: ConstraintHandler, GO, TInds <: AbstractVector{TI}, TMeta<:Metadata} <: StiffnessTopOptProblem{dim, TF}
    inp_content::InpContent{dim, TF, N, TI}
    geom_order::Type{Val{GO}}
    ch::Tch
    black::TBool
    white::TBool
    varind::TInds
    metadata::TMeta
end
  • dim: dimension of the problem
  • TF: number type for computations and coordinates
  • N: number of nodes in a cell of the grid
  • inp_content: an instance of InpContent which stores all the information from the `.inp file.
  • geom_order: a field equal to Val{GO} where GO is an integer representing the order of the finite elements. Linear elements have a geom_order of Val{1} and quadratic elements have a geom_order of Val{2}.
  • metadata: Metadata having various cell-node-dof relationships
  • black: a BitVector of length equal to the number of elements where black[e] is 1 iff the e^th element must be part of the final design
  • white: a BitVector of length equal to the number of elements where white[e] is 1 iff the e^th element must not be part of the final design
  • varind: an AbstractVector{Int} of length equal to the number of elements where varind[e] gives the index of the decision variable corresponding to element e. Because some elements can be fixed to be black or white, not every element has a decision variable associated.
source
TopOpt.TopOptProblems.InputOutput.INP.InpStiffnessMethod
InpStiffness(filepath::AbstractString; keep_load_cells = false)

Imports stiffness problem from a .inp file, e.g. InpStiffness("example.inp"). The keep_load_cells keyword argument will enforce that any cell with a load applied on it is to be part of the final optimal design generated by topology optimization algorithms.

source
Missing docstring.

Missing docstring for IO.INP.Parser.InpContent. Check Documenter's build log for details.

Grids

Grid types are defined in TopOptProblems because a number of topology optimization problems share the same underlying grid but apply the loads and boundary conditions at different locations. For example, the PointLoadCantilever and HalfMBB problems use the same rectilinear grid type, RectilinearGrid, under the hood. The LBeam problem uses the LGrid function under the hood to construct an L-shaped Ferrite.Grid. New problem types can be defined using the same grids but different loads or boundary conditions.

TopOpt.TopOptProblems.RectilinearGridType
struct RectilinearGrid{dim, T, N, M, TG<:Ferrite.Grid{dim, <:Ferrite.Cell{dim,N,M}, T}} <: AbstractGrid{dim, T}
    grid::TG
    nels::NTuple{dim, Int}
    sizes::NTuple{dim, T}
    corners::NTuple{2, Vec{dim, T}}
    white_cells::BitVector
    black_cells::BitVector
    constant_cells::BitVector
end

A type that represents a rectilinear grid with corner points corners.

  • dim: dimension of the problem
  • T: number type for computations and coordinates
  • N: number of nodes in a cell of the grid
  • M: number of faces in a cell of the grid
  • grid: a Ferrite.Grid struct
  • nels: number of elements in every dimension
  • sizes: dimensions of each rectilinear cell
  • corners: 2 corner points of the rectilinear grid
  • white_cells: cells fixed to be void during optimization
  • black_cells: cells fixed to have material during optimization
  • constant_cells: cells fixed to be either void or have material during optimization
source
TopOpt.TopOptProblems.RectilinearGridMethod
RectilinearGrid(::Type{Val{CellType}}, nels::NTuple{dim,Int}, sizes::NTuple{dim,T}) where {dim, T, CellType}

Constructs an instance of RectilinearGrid.

  • dim: dimension of the problem
  • T: number type for coordinates
  • nels: number of elements in every dimension
  • sizes: dimensions of each rectilinear cell

Example:

rectgrid = RectilinearGrid((60,20), (1.0,1.0))
source
TopOpt.TopOptProblems.LGridFunction
LGrid(::Type{Val{CellType}}, ::Type{T}; length = 100, height = 100, upperslab = 50, lowerslab = 50) where {T, CellType}
LGrid(::Type{Val{CellType}}, nel1::NTuple{2,Int}, nel2::NTuple{2,Int}, LL::Vec{2,T}, UR::Vec{2,T}, MR::Vec{2,T}) where {CellType, T}

Constructs a Ferrite.Grid that represents the following L-shaped grid.

        upperslab   UR
       ............
       .          .
       .          .
       .          . 
height .          .                     MR
       .          ......................
       .                               .
       .                               . lowerslab
       .                               .
       .................................
     LL             length

Examples:

LGrid(upperslab = 30, lowerslab = 70)
LGrid(Val{:Linear}, (2, 4), (2, 2), Vec{2,Float64}((0.0,0.0)), Vec{2,Float64}((2.0, 4.0)), Vec{2,Float64}((4.0, 2.0)))
source

Finite element backend

Currently, TopOpt uses Ferrite.jl for FEA-related modeling. This means that all the problems above are described in the language and types of Ferrite.

Matrices and vectors

ElementFEAInfo

TopOpt.TopOptProblems.ElementFEAInfoType
struct ElementFEAInfo{dim, T}
    Kes::AbstractVector{<:AbstractMatrix{T}}
    fes::AbstractVector{<:AbstractVector{T}}
    fixedload::AbstractVector{T}
    cellvolumes::AbstractVector{T}
    cellvalues::CellValues{dim, T}
    facevalues::FaceValues{<:Any, T}
    metadata::Metadata
    black::AbstractVector
    white::AbstractVector
    varind::AbstractVector{Int}
    cells
end

An instance of the ElementFEAInfo type stores element information such as:

  • Kes: the element stiffness matrices,
  • fes: the element load vectors,
  • cellvolumes: the element volumes,
  • cellvalues and facevalues: two Ferrite types that facilitate cell and face iteration and queries.
  • metadata: that stores degree of freedom (dof) to node mapping, dof to cell mapping, etc.
  • black: a BitVector such that black[i] is 1 iff element i must be part of any feasible design.
  • white: a BitVector such that white[i] is 1 iff element i must never be part of any feasible design.
  • varind: a vector such that varind[i] gives the decision variable index of element i.
  • cells: the cell connectivities.
source
TopOpt.TopOptProblems.ElementFEAInfoMethod
ElementFEAInfo(sp, quad_order=2, ::Type{Val{mat_type}}=Val{:Static}) where {mat_type}

Constructs an instance of ElementFEAInfo from a stiffness problem sp using a Gaussian quadrature order of quad_order. The element matrix and vector types will be:

  1. SMatrix and SVector if mat_type is :SMatrix or :Static, the default,
  2. MMatrix and MVector if mat_type is :MMatrix, or
  3. Matrix and Vector otherwise.

The static matrices and vectors are more performant and GPU-compatible therefore they are used by default.

source

GlobalFEAInfo

TopOpt.TopOptProblems.GlobalFEAInfoType
struct GlobalFEAInfo{T, TK<:AbstractMatrix{T}, Tf<:AbstractVector{T}, Tchol}
    K::TK
    f::Tf
    cholK::Tchol
end

An instance of GlobalFEAInfo hosts the global stiffness matrix K, the load vector f and the cholesky decomposition of the K, cholK.

source
TopOpt.TopOptProblems.GlobalFEAInfoMethod
GlobalFEAInfo(::Type{T}=Float64) where {T}

Constructs an empty instance of GlobalFEAInfo where the field K is an empty sparse matrix of element type T and the field f is an empty dense vector of element type T.

source
TopOpt.TopOptProblems.GlobalFEAInfoMethod
GlobalFEAInfo(sp::StiffnessTopOptProblem)

Constructs an instance of GlobalFEAInfo where the field K is a sparse matrix with the correct size and sparsity pattern for the problem instance sp. The field f is a dense vector of the appropriate size. The values in K and f are meaningless though and require calling the function assemble! to update.

source