Programs
4 Static Equilibrium Programs
PtFEM.p41 — Function.Method p41
One dimensional analysis of an axially loaded elastic Rod using 2-node Line elements.
Constructors
p41(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Type of structural fin_el
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Float64} : x-coordinate vectorOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1
* eq_nodal_forces_and_moments : Contribution of distributed loads to loaded_nodesReturn values
* (jfem, dis_dt, fm_dt) : Tuple of jFem, dis_dt and fm_dt
where:
jfem::jFem : Computational result type
dis_dt : Displacement data table
fm_dt : Forces and moments data tableRelated help
?StructuralElement : List of available structural element types
?Rod : Help on a Rod structural element
?FiniteElement : List finite element types
?Line : Help on Line finite elementMethod p41
One dimensional analysis of an axially loaded elastic Rod using 2-node Line elements.
Constructors
p41(m, data) # Re-use factored global stiffness matrixArguments
* `m::jFEM` : Previously created jFEM model
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Type of structural fin_el
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Float64} : x-coordinate vectorOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1
* eq_nodal_forces_and_moments : Contribution of distributed loads to loaded_nodesReturn values
* (jfem, dis_dt, fm_dt) : Tuple of jFem, dis_dt and fm_dt
where:
jfem::jFem : Computational result type
dis_dt : Displacement data table
fm_dt : Forces and moments data tableRelated help
?StructuralElement : List of available structural element types
?Rod : Help on a Rod structural element
?FiniteElement : List finite element types
?Line : Help on Line finite elementPtFEM.p42 — Method.Method p42
Analysis of elastic pin-jointed frames using 2-node rod elements in 2- or 3-dimensions.
Constructors
p42(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataDictionary keys
* struc_el::StructuralElement : Type of structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::Vector{Float64} : x coordinate vector
* y_coords::Vector{Float64} : y coordinate vector
* g_num::Array{Int,2} : Element node connectionsOptional additional dictionary keys
* penalty::Float64 : Penalty for fixed freedoms
* etype::Vector{Int} : Element material vector
* z_coords::Vector{Float64} : z coordinate vector (3D)
* eq_nodal_forces_and_moments : Contribution of distributed loads to loaded_nodesReturn values
* (jfem, dis_dt, fm_dt) : Tuple of jFem, dis_dt and fm_dt
where:
jfem::jFem : Computational result type
dis_dt : Displacement data table
fm_dt : Forces and moments data tableRelated help
?StructuralElement : List structural element types
?Frame : Help on a Rod structural fin_el
?FiniteElement : List finite element types
?Line : Help on Line finite elementPtFEM.p43 — Method.p43
Analysis of elastic beams using 2-node Beam structural elements and Line finite elements. Elastic foundation is optional.
Constructors
p43(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataDictionary keys
* struc_el::StructuralElement : Type of structural fin_el
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::LinSpace{Float64} : x coordinate vector
* g_num::Array{Int,2} : Element node connections
* fixed_freedoms::Array{Tuple{Vector{Int}} : Fixed freedomsOptional additional dictionary keys
* etype::Vector{Int} : Element material vector
* penalty::Float64 : Penalty for fixed freedoms
* eq_nodal_forces_and_moments : Equivalent nodal loadsReturn values
* (jfem, dis_dt, fm_dt) : Tuple of jFem, dis_dt and fm_dt
where:
jfem::jFem : Computational result type
dis_dt : Displacement data table
fm_dt : Forces and moments data tableRelated help
?StructuralElement : Help on structural elements
?Rod : Help on a Rod structural fin_el
?FiniteElement : Help on finite element typesPtFEM.p44 — Method.p44
Analysis of elastic rigid-jointed frames using a 2-node Frame structural element and Line finite elements in 2 or 3 dimensions.
Constructors
p44(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataDictionary keys
* struc_el::StructuralElement : Type of structural fin_el
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Float64} : x coordinate vector
* y_coords::FloatRange{Float64} : y coordinate vector
* g_num::Array{Int,2} : Element node connections
* fixed_freedoms::Array{Tuple{Vector{Int}} : Fixed freedomsOptional additional dictionary keys
* etype::Vector{Int} : Element material vector
* penalty::Float64 : Penalty for fixed freedoms
* z_coords::FloatRange{Float64} : z coordinate vector
* eq_nodal_forces_and_moments : Equivalent nodal loadsReturn values
* (jfem, dis_dt, fm_dt) : Tuple of jFem, dis_dt and fm_dt
where:
jfem::jFem : Computational result type
dis_dt : Displacement data table
fm_dt : Forces and moments data tableRelated help
?StructuralElement : Help on structural elements
?Beam : Help on a Beam structural fin_el
?FiniteElement : Help on finite element typesPtFEM.p45 — Method.Method p45
Analysis of elasto-plastic beams or rigid-jointed frames using a 2-node Frame structural element in 1, 2 or 3 dimensions.
Constructors
p45(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Type of structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Float64} : x-coordinate vector
* dload::FloatRange{Float64} : load stepsOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1
* y_coords::FloatRange{Float64} : y-coordinate vector (2D)
* z_coords::FloatRange{Float64} : x-coordinate vector (3D)
* limit = 250 : Iteration limit
* tol = 0.0001 : Tolerance for iteration convergenceRelated help
?StructuralElement : List of available structural element types
?Frame : Help on a Frame structural element
?FiniteElement : List finite element types
?Line : Help on Line finite elementPtFEM.p46 — Method.Method p46
Stability (buckling) analysis of elastic beams using a 2-node Beam structural element and Line finite elements. Elastic foundation is optional.
Constructors
p46(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Type of structural fin_el
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Float64} : x-coordinate vectorOptional additional data dictionary keys
* etype::Vector{Int} : Element material vector if np_types > 1
* limit = 250 : Iteration limit
* tol = 0.0001 : Tolerance for iteration convergenceRelated help
?StructuralElement : List of available structural element types
?Beam : Help on a Beam structural element
?FiniteElement : List finite element types
?Line : Help on Line finite elementPtFEM.p47 — Method.Method p47
Analysis of plates (Plane structural element) using 4-node Quadrilateral finite elements. Homogeneous material with identical elements. Mesh numbered in x or y direction.
Constructors
p47(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Floalt64} : x-coordinate vector
* y_coords::FloatRange{Floalt64} : y-coordinate vector
* thickness:: Float64 : Thickness of plateOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1Return values
* (fm_dt, sigma_dt) : Tuple of jFem, dis_dt and fm_dt
where:
fm_dt : Forces and moments data table
sigma_dt : Stresses data tableRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite element5 Elastic Solids Programs
PtFEM.p51 — Method.Method p51
Plane or axisymmetric strain analysis of an elastic solid (Plane structural element) using 3-, 6-, 10- or 15-node right-angled triangles (Triangle finite elements) or 4-, 8- or 9-node rectangular quadrilaterals (Quadrilateral finite elements). Mesh numbered in x(r)- or y(z)- direction.
Constructors
p51(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Floalt64} : x-coordinate vector
* y_coords::FloatRange{Floalt64} : y-coordinate vector
* thickness:: Float64 : Thickness of plateOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1Return values
* (fem, fm_dt, sigma_dt) : Tuple of jFem, dis_dt and fm_dt
where:
fm_dt : Forces and moments data table
sigma_dt : Stresses data tableRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite elementPtFEM.p52 — Method.Method p52
Non-axisymmetric analysis of an axisymmetric elastic solid using 8-node rectangular quadrilaterals. Mesh numbered in r- or z- direction.
Constructors
p52(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material propertiesOptional additional data dictionary keys
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1
* lth::Int :
* iflag::Int :
* chi::Float64 :Return values
* (fem, fm_dt, sigma_dt) : Tuple of jFem, dis_dt and fm_dt
where:
fm_dt : Forces and moments data table
sigma_dt : Stresses data tableRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite element6 Material Nonlinearity Programs
PtFEM.p61 — Method.Method p61
Plane strain bearing capacity analysis of an elastic-plastic (von Mises) material using 8-node rectangular quadrilaterals. Viscoplastic strain method.
Constructors
p61(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Floalt64} : x-coordinate vector
* y_coords::FloatRange{Floalt64} : y-coordinate vector
* thickness:: Float64 : Thickness of plate
* tol::Float64 : Convergence tolerance
* qincs::Vector{Float64} : Incremental load stepsOptional additional data dictionary keys
* limit = 250 : Iteration limit
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1Return values
* (g_coord, g_num, disp) : where:
g_coord : Coordinates
g_num : Node numbering
disp : Matrix of displacementsRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite elementPtFEM.p62 — Method.Method p62
Plane strain bearing capacity analysis of an elastic-plastic (von Mises) material using 8-node rectangular quadrilaterals.
Viscoplastic strain method.
No global stiffness matrix assembly.
Diagonally preconditioned conjugate gradient solver.
Constructors
p62(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* support::Array{Tuple{Int,Array{Int,1}},1} : Fixed-displacements vector
* loaded_nodes::Array{Tuple{Int,Array{Float64,1}},1} : Node load vector
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Floalt64} : x-coordinate vector
* y_coords::FloatRange{Floalt64} : y-coordinate vector
* thickness:: Float64 : Thickness of plate
* tol::Float64 : Convergence tolerance
* qincs::Vector{Float64} : Incremental load stepsOptional additional data dictionary keys
* limit = 250 : Iteration limit
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1Return values
* (g_coord, g_num, disp) : where:
g_coord : Coordinates
g_num : Node numbering
disp : Matrix of displacementsRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite elementPtFEM.p63 — Method.Method p63
Plane strain bearing capacity analysis of an elastic-plastic (Mohr-Coulomb) material using 8-node rectangular quadrilaterals. Rigid smooth footing. Displacement control. Viscoplastic strain method.
Constructors
p63(data)Arguments
* `data::Dict{Symbol, Any}` : Dictionary containing all input dataRequired data dictionary keys
* struc_el::StructuralElement : Structural element
* properties::Vector{Float64} : Material properties
* x_coords::FloatRange{Floalt64} : x-coordinate vector
* y_coords::FloatRange{Floalt64} : y-coordinate vectorOptional additional data dictionary keys
* tol::Float64 : Convergence tolerance
* limit = 250 : Iteration limit
* incs::Int : Incremental load steps
* presc::Float64 : Wall displacement increment
* penalty = 1e20 : Penalty used for fixed degrees of freedoms
* etype::Vector{Int} : Element material vector if np_types > 1Return values
* (g_coord, g_num, disp) : where:
g_coord : Coordinates
g_num : Node numbering
disp : Matrix of displacementsRelated help
?StructuralElement : List of available structural element types
?Plane : Help on a Plane structural element
?FiniteElement : List finite element types
?Quadrilateral : Help on Quadrilateral finite elementStructural Element Types
PtFEM.StructuralElement — Type.StructuralElement
Abstract structural element type.
Type
abstract StructuralElementSubtypes
* Rod::StructuralElement : Rod(nxe, np_types, nip, fin_el)
* Beam::StructuralElement : Beam(nod, nodof)
* Frame::StructuralElement : Frame(nod, nodof)
* Plane::StructuralElement : Plane(nod, nodof)
* Solid::StructuralElement : Solid(nod, nodof)
* GenericSolid::StructuralElement : GenericSolid(nod, nodof)Related help
?FiniteElement : Show all finite elements
?Rod : Help on Rod structural element
?Beam : Help on Beam structural element
?Frame : Help on Frame structural element
?Plane : Help on Plane structural element
?Solid : Help on Solid structural element
?GenericSolid : Help on GenericSolid structural elementPtFEM.Rod — Type.Rod
Concrete 1D structural element with only axial stresses.
Constructor
Rod(nels, np_types, nip, fin_el)Arguments
* nels::Int : Number of fin_els (stored in field nxe)
* np_types::Int : Number of different property types
* nip::Int : Number of integration points
* fin_el::FiniteElement : Line(nod, nodof)Related help
?StructuralElement : Help on structural elements
?FiniteElement : Help on finite element types
?Line : Help on a Line finite elementPtFEM.Beam — Type.Beam
Concrete structural element with transverse and moment loading.
Constructor
Beam(ndim, nip, fin_el)Arguments
* ndim::Int : Number of dimensions
* nst::Int : Number of stress terms
* nxe::Int : Number of different property types
* nip::Int : Number of integration points
* direction::Symbol : Number of integration points
* fin_el::FiniteElement : Line(nod, nodof)
* axisymmetric::Bool : Axisymmetric if trueRelated help
?StructuralElement : Help on structural elements
?FiniteElement : Help on finite element types
?Line : Help on a Line finite elementPtFEM.Frame — Type.Frame
Pin- or rigid-jointed structural element.
Constructor
Frame(nels, nn, ndim, finite_element(nod, nodof))Arguments
* nels::Int : Number of elements
* nn:Int : Number of nodes
* ndim::Int : Number of dimensions
* nst::Int : Number of stress terms
* nip::Int : Number of integration points
* fin_el::FiniteElement : Line(nod, nodof)Related help
?StructuralElement : List structural elements
?FiniteElement : List finite element types
?Line : Help on a Line finite elementPtFEM.Plane — Type.Plane
Plate structural element.
Constructor
Plane(ndim, nst, nxe, nye, nip, dir, finite_element(nod, nodof), axisymmetric)Arguments
* ndim::Int : Number of dimensions
* nst::Int : Number of stress terms
* nxe::Int : Number of elements in x direction
* nye::Int : Number of elements in y direction
* nip::Int : Number of integration points
* dir::Symbol : Direction of node numbering
* fin_el::FiniteElement : Line(nod, nodof)
* axisymmetric::Bool : AxisymmetricRelated help
?StructuralElement : List structural elements
?FiniteElement : List finite element types
?Line : Help on a Line finite elementPtFEM.Solid — Type.Solid
Solid structural element.
Constructor
Solid(ndim, nst, nxe, nye, nze, nip, finite_element(nod, nodof))Arguments
* ndim::Int : Number of dimensions
* nst::Int : Number of stress terms
* nxe::Int : Number of elements in x direction
* nye::Int : Number of elements in y direction
* nze::Int : Number of elements in z direction
* nip::Int : Number of integration points
* fin_el::FiniteElement : Line(nod, nodof)Related help
?StructuralElement : List structural elements
?FiniteElement : List finite element typesPtFEM.GenericSolid — Type.GenericSolid
Solid structural element.
Constructor
GenericSolid(ndim, nst, nels, nn, nip, finite_element(nod, nodof), axisymmetric)Arguments
* ndim::Int : Number of dimensions
* nst::Int : Number of stress terms
* nels::Int : Number of finite elements
* nn::Int : Number of nodes
* nip::Int : Number of integration points
* fin_el::FiniteElement : Finite element type used
* axisymmetric::Bool : AxisymmetricRelated help
?StructuralElement : List structural elements
?FiniteElement : List finite element typesFinite Element Types
PtFEM.FiniteElement — Type.FiniteElement
Abstract finite element type.
Type
abstract FiniteElementSubtypes
* Line::FiniteElement : 1D Line(nod, nodof)
* Triangle::FiniteElement : 2D Triangle(nod, nodof)
* Quadrilateral::FiniteElement : 2D Quadrilateral(nod, nodof)
* Hexahedron::FiniteElement : 3D Hexahedron(nod, nodof)
* Tetrahedron::FiniteElement : 3D Tetrahedron(nod, nodof)PtFEM.Line — Type.Line (Interval)
1D type finite element
Constructor
Line(nod, nodof)
Line(nodof)Arguments
* nod::Int : Number of nodes for finite element, defaults to 2
* nodof::Int : Number of degrees of freedom per nodeRelated help
?FiniteElement : Help on finite element typesPtFEM.Triangle — Type.Triangle
2D type finite element
Constructor
Triangle(nod, nodof)Arguments
* nod::Int : Number of nodes for finite element (3, 6, 10, 15)
* nodof::Int : Number of degrees of freedom per nodeRelated help
?FiniteElement : Help on finite element typesPtFEM.Quadrilateral — Type.Quadrilateral
2D type finite element
Constructor
Quadrilateral(nod, nodof)Arguments
* nod::Int : Number of nodes for finite element (4, 8, 9)
* nodof::Int : Number of degrees of freedom per nodeRelated help
?FiniteElement : Help on finite element typesPtFEM.Hexahedron — Type.hexahedron
3D type finite element
Constructor
Hexahedron(nod, nodof)Arguments
* nod::Int : Number of nodes for finite element (8, 14, 20)
* nodof::Int : Number of degrees of freedom per nodeRelated help
?FiniteElement : Help on finite element typesPtFEM.Tetrahedron — Type.Tetrahedron
3D type finite element
Constructor
Tetrahedron(nod, nodof)
Tetrahedron(nodof)Arguments
* nod::Int : Number of nodes for finite element (defaults to 4)
* nodof::Int : Number of degrees of freedom per nodeRelated help
?FiniteElement : Help on finite element typesOther Julia Types
PtFEM.FEM — Type.FEM
Computational structure used in chapter 5 (Skyline format used)
PtFEM.jFEM — Type.jFEM
Computational structure used in chapter 4 (Julia Sparse matrices used)
PtFEM - Main
PtFEM.beam_gm — Method.beam_gm
This subroutine forms the beam geometric matrix for stability analysis.
Method
beam_gm(ell::Float64)Arguments
* ell::Float64 : Element lengthReturn value
* gm::::Matrix{Float64}(4,4) : Geometric matrix for beam elementPtFEM.beam_km — Method.beam_km
This subroutine forms the stiffness matrix of a beam element (bending only).
Method
beam_km(ei, ell)Arguments
* ei::Float64 : Element stiffness
* ell::Float64 : Element lengthReturn value
* km::::Matrix{Float64} : Stiiness matrix for beam element (Updated)PtFEM.beam_mm — Method.beam_mm
This subroutine forms the stiffness matrix of a beam element (bending only).
Method
beam_mm(ei, ell)Arguments
* fs::Float64 : Element density
* ell::Float64 : Element lengthReturn value
* mm::::Matrix{Float64} : Mass matrix for beam elembeam_mmated)PtFEM.beemat! — Method.beemat!
This subroutine forms the strain-displacement matrix for axisymmetric solids subjected to non-axisymmetric loading.
Method
beemat!(bee, deriv)Arguments
* bee::Matrix{Float64} : Bee matrix (Updated)
* deriv::Matrix{Float64} : DerivativePtFEM.bmat_nonaxi! — Method.bmat_nonaxi!
This subroutine forms the strain-displacement matrix for axisymmetric solids subjected to non-axisymmetric loading.
Method
bmat_nonaxi!(bee, radius, coord, deriv, fun, iflag, lth)Arguments
* bee::Matrix{Float64} : Bee matrix (Updated)
* radius::Float64 : r coordinate of the Gauss point
* coord::Matrix{Float64} : Nodal coordinate matrix
* deriv::Matrix{Float64} : Derivative
* fun::Vector{Float64} : Shape function
* iflag::Int : 1 = symmetric, -1 = anti-symmetric
* lth::Int : HarmonicPtFEM.checon — Method.checon
This subroutine sets converged to .FALSE. if relative change in loads and oldlds is greater than tol and updates oldlds.
Method
checon(loads, oldlds, tol)Arguments
* loads::Vector{Float64} : Displacements vector/OffsetArray
* oldlds::Vector{Float64} : Previous displacement vector/OffsetArray
* tol::Float64 : Convergence toleranceReturn value
* ::Bool : Convergence achievedPtFEM.deemat! — Method.deemat!
This subroutine returns the elastic dee matrix for ih=3 (plane strain), ih=4 (axisymmetry or plane strain elastoplasticity) or ih=6 (three dimensions).
Method
deemat!(dee, e, v)Arguments
* dee::Matrix{Float64} : Dee matrix (Updated)
* e::Float64 : Young's modulus
* v::Float64 : Poisson's ratioPtFEM.fkdiag! — Method.fkdiag!
This subroutine returns the elastic dee matrix for ih=3 (plane strain), ih=4 (axisymmetry or plane strain elastoplasticity) or ih=6 (three dimensions).
Method
fkdiag!(kdiag, g)Arguments
* kdiag::Vector{Int} : Bandwidth vector (Updated)
* g::Vector{Int} : Element steering vectorPtFEM.fmplat! — Method.fmplat!
This subroutine forms the 2nd derivatives for rectangular plate bending fin_els.
Method
fmplat!(d2x, d2y, d2xy, points, aa, bb, i)Arguments
* d2x::Vector{Float64} : x derivative term (Updated)
* d2y::Vector{Float64} : y derivative term (Updated)
* d2xy::Vector{Float64} : x,y derivative term (Updated)
* points::Matrix{Float64} : Location of Gauss points
* aa::Float64 : Dimension of plate
* bb::Float64 : Dimension of plate
* i::Int : Gauss point to usePtFEM.formm! — Method.formm!
This subroutine forms the derivatives of the invariants with respect to stress in 2- or 3-d. See equation 6.25.
Function
formm!(stress, m1, m2, m3)Arguments
* stress::Vector{Float64} : Stress vector, see eq 6.25
* m1::Matrix{Float64} : m1 matrix
* m2::Matrix{Float64} : m2 matrix
* m3::Matrix{Float64} : m3 matrixReturn values
* m1::Matrix{Float64} : m1 matrix
* m2::Matrix{Float64} : m2 matrix
* m3::Matrix{Float64} : m3 matrixPtFEM.formnf! — Method.formnf!
Returns nodal freedom numbering array nf
Function
formnf!(nodof, nn, nf)Arguments
* nodof::Int : Number of degrees of freedom for each node
* nn::Int : Number of nodes in mesh
* nf::Array{Int,2} : Nodal freedom matrix (updated)PtFEM.fsparm! — Method.fsparm!
Function fsparm assembles fin_el matrices into a Julia sparse global stiffness matrix.
Method
fsparm!(gsm, g, km)Arguments
* gsm::SparseArrays{Float64, Float64} : Sparse stiffnes matrix (Updated)
* g::Vector{Int} : Global coordinate vector.
* km::Matrix{Float64} : Stiffness matrix.PtFEM.glob_to_loc! — Method.glob_to_loc!
This subroutine transforms the global end reactions and moments into the local system (2- or 3-d). Called from hinge!().
Function
glob_to_loc!(loc, glob, gamma, coord)Arguments
* loc::Vector{Float64} : Local force and momemts (Updated)
* glob::Vector{Float64} : Globale forces and moments
* gamma::Float64 : Element orientation angle (3D)
* coord::Matrix{Float64} : Nodal coordinatesPtFEM.glob_to_axial — Method.glob_to_axial
This subroutine transforms the global end reactions into an axial force for rod fin_els (2- or 3-d).
Function
glob_to_axial(glob, coord)Arguments
* glob::Vector{Float64} : Globale forces and moments
* coord::Matrix{Float64} : Nodal coordinatesREturn value
* ::Float64 : Axial forcePtFEM.hinge! — Method.hinge!
This subroutine forms the end forces and moments to be applied to a member if a joint has gone plastic.
Function
hinge!(coord, holdr, action, react, prop, iel, etype, gamma)Arguments
* coord::Matrix{Float64} : Nodal coordinates
* holdr::Matrix{Float64} : Existing reactions
* action::Vector{Float64} : Incremental reactions
* react::Vector{Float64} : Correction to reactions (Updated)
* prop::Matrix{Float64} : Beam properties
* iel::Int : Element number
* etype::Vector{Int} : Element type
* gamma::Vector{Float64} : Element orientation (3D)PtFEM.invar — Method.invar
This subroutine forms the stress invariants in 2- or 3-d. See equations 6.3 and 6.4
Function
invar(stress, sigm, dsbar, theta)Arguments
* stress::Vector{Float64} : Stress vector
* sigm::Float64 : Invariant, eq 6.4 (Updated)
* dsbar::Float64 : Invariant, eq 6.4 (Updated)
* theta::Float64 : Invariant, eq 6.3 (Updated)REturn values
* stress::Vector{Float64} : Stress vector
* sigm::Float64 : Invariant, eq 6.4 (Updated)
* dsbar::Float64 : Invariant, eq 6.4 (Updated)PtFEM.linmul_sky! — Method.linmul_sky!
This subroutine forms the product of symmetric matrix stored as a skyline and a vector.
Method
linmul_sky!(kv, disps, loads, kdiag)Arguments
* kv::Vector{Float64} : Sparse stiffnes matrix (Skyline format)
* disps::Vector{Float64} : Displacements
* loads::Vector{Float64} : Loads vector (Updated)
* kdiag::Vector{Int} : Bandwidth vectorPtFEM.loc_to_glob! — Method.loc_to_glob!
This subroutine transforms the local end reactions and moments into the global system (3-d).
Function
loc_to_glob!(loc, glob, gamma, coord)Arguments
* loc::Vector{Float64} : Local force and momemts (Updated)
* glob::Vector{Float64} : Globale forces and moments
* gamma::Float64 : Element orientation angle (3D)
* coord::Matrix{Float64} : Nodal coordinatesPtFEM.mocouf — Method.mocouf
This subroutine calculates the value of the yield function for a Mohr-Coulomb material (phi in degrees).
Function
mocouf(phi, c, sigm, dsbar, theta)Arguments
* psi::Float64 : Local force and momemts (Updated)
* c::Float64 : Globale forces and moments
* sigm::Float64 : Element orientation angle (3D)
* dsbar::Float64 : Globale forces and moments
* theta::Float64 : Element orientation angle (3D)Return value
* ::Float64 : Value of yield functionPtFEM.mocouq — Method.mocouq
This subroutine forms the derivatives of a Mohr-Coulomb potential function with respect to the three invariants (psi in degrees).
Function
(dq1,dq2,dq3) = mocouq(psi,dsbar,theta)Arguments
* psi::Float64 : Local force and momemts (Updated)
* dsbar::Float64 : Globale forces and moments
* theta::Float64 : Element orientation angle (3D)Return values
* dq1::Float64 : Local force and momemts (Updated)
* dq2::Float64 : Globale forces and moments
* dq3::Float64 : Element orientation angle (3D)PtFEM.num_to_g! — Method.num_to_g!
Returns the element steering vector g from the element node numbering num and the nodal freedom array nf.
Function
num_to_g!(num, nf, g)Arguments
* num::Vector{Int} : Node numbering vector
* nf::Matrix{Int} : Nodal freedom array
* g::Vector{Int} : Element steering vector (Updated)PtFEM.pin_jointed! — Method.pin_jointed!
This subroutine forms the global stiffness matrix of a general pin-joionted structural element (1-, 2- or 3-d).
Function
pin_jointed!(km, ea, coord)Arguments
* km::Matrix{Float64} : Element stiffness matrix (Updated)
* ea::Float64 : Element stiffness
* coord::Matrix{Float64}} : Element nodal coordinatesPtFEM.rigid_jointed! — Method.rigid_jointed!
This subroutine forms the global stiffness matrix of a general pin-joionted structural element (1-, 2- or 3-d).
Function
rigid_jointed!(km, prop, gamma, etype, iel, coord)Arguments
* km::Matrix{Float64} : Element stiffness matrix (Updated)
* prop::Matrix{Float64} : Element properties
* gamma::Vector{Float64} : Element orientations (3D)
* etype::Vector{Int} : Element type vector
* iel::Int : Element number
* coord::Matrix{Float64}} : Element nodal coordinatesPtFEM.rod_km! — Method.rod_km!
This subroutine forms the stiffness matrix of a 1-d "rod" fin_el.
Function
rod_km!(km, ea, length)Arguments
* km::Matrix{Float64} : Element stiffness matrix (Updated)
* ea::Float64 : Element stiffness
* ell::Float64 : Element lengthPtFEM.rod_mm! — Method.rod_mm!
This subroutine forms the consistent mass matrix of a 1-d "rod" fin_el.
Function
rod_mm!(km, ell)Arguments
* mm::Matrix{Float64} : Element mass matrix (Updated)
* ell::Float64 : Element lengthPtFEM.sample! — Function.sample!
This subroutine returns the local coordinates and weighting coefficients of the integrating points.
Function
sample!(fin_el, s, wt)Arguments
* fin_el::FiniteElement : Finite element type
* s::Matrix{Float64} : Local coordinates (Updated)
* wt::Vector{Float64} : Weighting coefficients (Updated)PtFEM.shape_der! — Method.shape_der!
This subroutine produces derivatives of shape functions with respect to local coordinates.
Function
shape_der!(der, point, i)Arguments
* der::Matrix{Float64} : Function derivative (Updated)
* points::Matrix{Float64} : Local coordinates of integration points
* i::Int : Integration pointPtFEM.shape_fun! — Method.shape_fun!
This subroutine produces derivatives of shape functions with respect to local coordinates.
Function
shape_fun!(fun, point, i)Arguments
* fun::Vector{Float64} : Shape function (Updated)
* points::Matrix{Float64} : Local coordinates of integration points
* i::Int : Integration pointPtFEM.stability — Method.stability
Function spabac! performs Cholesky forward and back-substitution on a symmetric skyline global matrix. The loads vector is updated.
###Arguments
stability(gsm, ggm, tol, limit)Arguments
* gsm::SparseMatrixCSC{Float64,Int} : Factored global stiffness matrix
* ggm::SparseMatrixCSC{Float64,Int} : Factored geometric matrix
* tol::Float64 : Convergence tolerance
* limit::Int : Iteration limitPtFEM - Geom
PtFEM.geom_rect! — Function.geom_rect!
This subroutine forms the coordinates and connectivity for a rectangular mesh of right angled triangular elements (3, 6, 10 or 15-node) or quadrilateral elements (4, 8 or 9-node) counting in the x- or y-dir.
Method
geom_rect!(fin_el, iel, x_coords, y_coords, coord, num, dir)Arguments
* fin_el::FiniteElement : Shape of finite element
(Trangle or Quadrilateral)
* iel::Int : Element number
* x_coords::FloatRange{Float64} : x coordinates
* y_coords::FloatRange{Float64} : y coordinates
* coord::Matrix{Float64} : Nodal coordinates (Updated)
* num::Vector{Int} : Node numbers (Updated)
* dir::Symbol : Node numbering directionPtFEM.hexahedron_xz! — Function.hexahedron_xz!
This subroutine generates nodal coordinates and numbering for 8, 14 or 20-node "bricks" counting x-z planes in the y-direction.
Method
hexahedron_xz!(iel, x_coords, y_coords, z_coords, coord, num)Arguments
* iel::Int : Element number
* x_coords::FloatRange{Float64} : x coordinates
* y_coords::FloatRange{Float64} : y coordinates
* z_coords::FloatRange{Float64} : y coordinates
* coord::Matrix{Float64} : Nodal coordinates (Updated)
* num::Vector{Int} : Node numbers (Updated)
* dir::Symbol : Node numbering directionPtFEM.mesh_size — Function.mesh_size
Function mesh_size returns the number of fin_els (nels) and the number of nodes (nn) in a 1, 2 or 3-d geometry-created mesh.
Method
(nels, nn) = mesh_size(fin_el, nxe, [nye[, nze]])Arguments
* fin_el::FiniteElement : Shape of finite element
1D: Line
2D: Trangle or Quadrilateral
3D: Hexahedron
* nxe::Int : Number of fin_els in x direction
* nye::Int : Number of fin_els in y direction (for 2D and 3D)
* nze::Int : Number of fin_els in z direction (3D only)mesh_size
mesh_size: The function mesh_size returns the number of fin_els (nels) and the number of nodes (nn) in a 2-d geometry-created mesh.
Method
(nels, nn) = mesh_size(fin_el, nxe)Arguments
* `fin_el` : Shape of 2D finite element (Triangle)
* `nxe` : Number of fin_els in x direction
* `nxe` : Number of fin_els in y directionmesh_size
mesh_size: The function mesh_size returns the number of fin_els (nels) and the number of nodes (nn) in a 2-d geometry-created mesh.
Method
(nels, nn) = mesh_size(fin_el, nxe, nye)Arguments
* `fin_el` : Shape of 2D finite element (Quadrilateral)
* `nxe` : Number of fin_els in x direction
* `nye` : Number of fin_els in y directionmesh_size
mesh_size: The function mesh_size returns the number of fin_els (nels) and the number of nodes (nn) in a 3-d geometry-created mesh.
Method
(nels, nn) = mesh_size(fin_el, nxe, nye, nze)Arguments
* `fin_el` : Shape of 2D finite element (Hexahedron)
* `nxe` : Number of fin_els in x direction
* `nye` : Number of fin_els in y direction
* `nxe` : Number of fin_els in x directionPtFEM - Plot methods
PtFEM.mesh(data::Dict, g_coord::Array{Float64,2}, g_num::Array{Int, 2}, disp, ampl, pdir)PtFEM - VTK methods
PtFEM.vtk — Method.vtk
Plots displacements and directions
Function
mesh(data, fm_dt, sigma_dt, dir, fname)Arguments
* data::Dict : Input dictionary
* fm_dt::DataTable : Forces and moments DataTable
* sigma_dt::DataTable : Stresses DataTable
* dir::AbstractString : Project directory
* fname::AbstractString : Output VTK file namePtFEM - Julia functions & operators
Base.LinAlg.cholfact — Function.cholfact(A, [uplo::Symbol,] Val{false}) -> CholeskyCompute the Cholesky factorization of a dense symmetric positive definite matrix A and return a Cholesky factorization. The matrix A can either be a Symmetric or Hermitian StridedMatrix or a perfectly symmetric or Hermitian StridedMatrix. In the latter case, the optional argument uplo may be :L for using the lower part or :U for the upper part of A. The default is to use :U. The triangular Cholesky factor can be obtained from the factorization F with: F[:L] and F[:U]. The following functions are available for Cholesky objects: size, \, inv, and det. A PosDefException exception is thrown in case the matrix is not positive definite.
Example
julia> A = [4. 12. -16.; 12. 37. -43.; -16. -43. 98.]
3×3 Array{Float64,2}:
4.0 12.0 -16.0
12.0 37.0 -43.0
-16.0 -43.0 98.0
julia> C = cholfact(A)
Base.LinAlg.Cholesky{Float64,Array{Float64,2}} with factor:
[2.0 6.0 -8.0; 0.0 1.0 5.0; 0.0 0.0 3.0]
julia> C[:U]
3×3 UpperTriangular{Float64,Array{Float64,2}}:
2.0 6.0 -8.0
⋅ 1.0 5.0
⋅ ⋅ 3.0
julia> C[:L]
3×3 LowerTriangular{Float64,Array{Float64,2}}:
2.0 ⋅ ⋅
6.0 1.0 ⋅
-8.0 5.0 3.0
julia> C[:L] * C[:U] == A
truecholfact(A, [uplo::Symbol,] Val{true}; tol = 0.0) -> CholeskyPivotedCompute the pivoted Cholesky factorization of a dense symmetric positive semi-definite matrix A and return a CholeskyPivoted factorization. The matrix A can either be a Symmetric or Hermitian StridedMatrix or a perfectly symmetric or Hermitian StridedMatrix. In the latter case, the optional argument uplo may be :L for using the lower part or :U for the upper part of A. The default is to use :U. The triangular Cholesky factor can be obtained from the factorization F with: F[:L] and F[:U]. The following functions are available for PivotedCholesky objects: size, \, inv, det, and rank. The argument tol determines the tolerance for determining the rank. For negative values, the tolerance is the machine precision.
cholfact(A; shift = 0.0, perm = Int[]) -> CHOLMOD.FactorCompute the Cholesky factorization of a sparse positive definite matrix A. A must be a SparseMatrixCSC, Symmetric{SparseMatrixCSC}, or Hermitian{SparseMatrixCSC}. Note that even if A doesn't have the type tag, it must still be symmetric or Hermitian. A fill-reducing permutation is used. F = cholfact(A) is most frequently used to solve systems of equations with F\b, but also the methods diag, det, and logdet are defined for F. You can also extract individual factors from F, using F[:L]. However, since pivoting is on by default, the factorization is internally represented as A == P'*L*L'*P with a permutation matrix P; using just L without accounting for P will give incorrect answers. To include the effects of permutation, it's typically preferable to extract "combined" factors like PtL = F[:PtL] (the equivalent of P'*L) and LtP = F[:UP] (the equivalent of L'*P).
Setting the optional shift keyword argument computes the factorization of A+shift*I instead of A. If the perm argument is nonempty, it should be a permutation of 1:size(A,1) giving the ordering to use (instead of CHOLMOD's default AMD ordering).
This method uses the CHOLMOD library from SuiteSparse, which only supports doubles or complex doubles. Input matrices not of those element types will be converted to SparseMatrixCSC{Float64} or SparseMatrixCSC{Complex128} as appropriate.
Many other functions from CHOLMOD are wrapped but not exported from the Base.SparseArrays.CHOLMOD module.
Base.:\ — Function.\(x, y)Left division operator: multiplication of y by the inverse of x on the left. Gives floating-point results for integer arguments.
julia> 3 \ 6
2.0
julia> inv(3) * 6
2.0
julia> A = [1 2; 3 4]; x = [5, 6];
julia> A \ x
2-element Array{Float64,1}:
-4.0
4.5
julia> inv(A) * x
2-element Array{Float64,1}:
-4.0
4.5\(A, B)Matrix division using a polyalgorithm. For input matrices A and B, the result X is such that A*X == B when A is square. The solver that is used depends upon the structure of A. If A is upper or lower triangular (or diagonal), no factorization of A is required and the system is solved with either forward or backward substitution. For non-triangular square matrices, an LU factorization is used.
For rectangular A the result is the minimum-norm least squares solution computed by a pivoted QR factorization of A and a rank estimate of A based on the R factor.
When A is sparse, a similar polyalgorithm is used. For indefinite matrices, the LDLt factorization does not use pivoting during the numerical factorization and therefore the procedure can fail even for invertible matrices.
Example
julia> A = [1 0; 1 -2]; B = [32; -4];
julia> X = A \ B
2-element Array{Float64,1}:
32.0
18.0
julia> A * X == B
truePtFEM - Parallel processing
Base.Distributed.pmap — Function.pmap([::AbstractWorkerPool], f, c...; distributed=true, batch_size=1, on_error=nothing, retry_delays=[]), retry_check=nothing) -> collectionTransform collection c by applying f to each element using available workers and tasks.
For multiple collection arguments, apply f elementwise.
Note that f must be made available to all worker processes; see Code Availability and Loading Packages for details.
If a worker pool is not specified, all available workers, i.e., the default worker pool is used.
By default, pmap distributes the computation over all specified workers. To use only the local process and distribute over tasks, specify distributed=false. This is equivalent to using asyncmap. For example, pmap(f, c; distributed=false) is equivalent to asyncmap(f,c; ntasks=()->nworkers())
pmap can also use a mix of processes and tasks via the batch_size argument. For batch sizes greater than 1, the collection is processed in multiple batches, each of length batch_size or less. A batch is sent as a single request to a free worker, where a local asyncmap processes elements from the batch using multiple concurrent tasks.
Any error stops pmap from processing the remainder of the collection. To override this behavior you can specify an error handling function via argument on_error which takes in a single argument, i.e., the exception. The function can stop the processing by rethrowing the error, or, to continue, return any value which is then returned inline with the results to the caller.
Consider the following two examples. The first one returns the exception object inline, the second a 0 in place of any exception:
julia> pmap(x->iseven(x) ? error("foo") : x, 1:4; on_error=identity)
4-element Array{Any,1}:
1
ErrorException("foo")
3
ErrorException("foo")
julia> pmap(x->iseven(x) ? error("foo") : x, 1:4; on_error=ex->0)
4-element Array{Int64,1}:
1
0
3
0Errors can also be handled by retrying failed computations. Keyword arguments retry_delays and retry_check are passed through to retry as keyword arguments delays and check respectively. If batching is specified, and an entire batch fails, all items in the batch are retried.
Note that if both on_error and retry_delays are specified, the on_error hook is called before retrying. If on_error does not throw (or rethrow) an exception, the element will not be retried.
Example: On errors, retry f on an element a maximum of 3 times without any delay between retries.
pmap(f, c; retry_delays = zeros(3))Example: Retry f only if the exception is not of type InexactError, with exponentially increasing delays up to 3 times. Return a NaN in place for all InexactError occurrences.
pmap(f, c; on_error = e->(isa(e, InexactError) ? NaN : rethrow(e)), retry_delays = ExponentialBackOff(n = 3))PtFEM - Deprecated
PtFEM.fkdiag! — Function.fkdiag!
This subroutine returns the elastic dee matrix for ih=3 (plane strain), ih=4 (axisymmetry or plane strain elastoplasticity) or ih=6 (three dimensions).
Method
fkdiag!(kdiag, g)Arguments
* kdiag::Vector{Int} : Bandwidth vector (Updated)
* g::Vector{Int} : Element steering vectorPtFEM.fromSkyline — Method.fromSkyline
Helper function to convert a Skyline vector to a full matrix.
Type
fromSkyline(skyline::Vector{Float64}, kdiag::Vector{Int})Arguments
* skyline::Vector{Float64} : 1D Line(nod, nodof)
* kdiag::Vector{Int} : 2D Triangle(nod, nodof)PtFEM.fsparv! — Method.fsparv!
Function fsparv! assembles fin_el matrices into a symmetric skyline global matrix. The Skyline vector kv is updated.
Method
fsparv!(kv, km, g, km)Arguments
* kv::Vector{Float64} : Sparse stiffnes matrix (Updated)
* km::Matrix{Float64} : Symmetric element stiffnes matrix
* g::Vector{Int} : Global steering vector.
* kdiag::Vector{Int} : Location of diagoinal termsPtFEM.linmul_sky! — Function.linmul_sky!
This subroutine forms the product of symmetric matrix stored as a skyline and a vector.
Method
linmul_sky!(kv, disps, loads, kdiag)Arguments
* kv::Vector{Float64} : Sparse stiffnes matrix (Skyline format)
* disps::Vector{Float64} : Displacements
* loads::Vector{Float64} : Loads vector (Updated)
* kdiag::Vector{Int} : Bandwidth vectorPtFEM.skyline2sparse — Method.skyline2sparse
Converts a Skyline matrix to a Julia Sparse matrix
Function
skyline2sparse(skyline, kdiag)Arguments
* skyline::Vector{Float64} : Skyline matrix
* kdiag::Vector{Int} : Element diagonal index vectorPtFEM.spabac! — Method.spabac!
Function spabac! performs Cholesky forward and back-substitution on a symmetric skyline global matrix. The loads vector is updated.
###Arguments
spabac!(kv, loads, kdiag)Arguments
* kv::Vector{Float64} : Skyline vector of global stiffness matrix
* loads::Vector{Float64} : Load vector (Updated)
* kdiag::Vector{Int} : Diagonal elemnt index vectorPtFEM.sparin! — Method.sparin!
Function sparin! performs Cholesky factorisation on a symmetric skyline global matrix. The vector kv is updated.
###Arguments
sparin!(kv, kdiag)Arguments
* kv::Vector{Float64} : Global stiffness matrix (Updated)
* kdiag::Vector{Int} : Diagonal elemnt index vectorIndex
PtFEM.BeamPtFEM.FEMPtFEM.FiniteElementPtFEM.FramePtFEM.GenericSolidPtFEM.HexahedronPtFEM.LinePtFEM.PlanePtFEM.QuadrilateralPtFEM.RodPtFEM.SolidPtFEM.StructuralElementPtFEM.TetrahedronPtFEM.TrianglePtFEM.jFEMBase.:\Base.Distributed.pmapBase.LinAlg.cholfactPtFEM.beam_gmPtFEM.beam_kmPtFEM.beam_mmPtFEM.beemat!PtFEM.bmat_nonaxi!PtFEM.checonPtFEM.deemat!PtFEM.fkdiag!PtFEM.fkdiag!PtFEM.fmplat!PtFEM.formm!PtFEM.formnf!PtFEM.fromSkylinePtFEM.fsparm!PtFEM.fsparv!PtFEM.geom_rect!PtFEM.glob_to_axialPtFEM.glob_to_loc!PtFEM.hexahedron_xz!PtFEM.hinge!PtFEM.invarPtFEM.linmul_sky!PtFEM.linmul_sky!PtFEM.loc_to_glob!PtFEM.mesh_sizePtFEM.mocoufPtFEM.mocouqPtFEM.num_to_g!PtFEM.p41PtFEM.p42PtFEM.p43PtFEM.p44PtFEM.p45PtFEM.p46PtFEM.p47PtFEM.p51PtFEM.p52PtFEM.p61PtFEM.p62PtFEM.p63PtFEM.pin_jointed!PtFEM.rigid_jointed!PtFEM.rod_km!PtFEM.rod_mm!PtFEM.sample!PtFEM.shape_der!PtFEM.shape_fun!PtFEM.skyline2sparsePtFEM.spabac!PtFEM.sparin!PtFEM.stabilityPtFEM.vtk