Function reference

Core functions

CALFEM Core module

Contains all the functions implementing CALFEM standard functionality

calfem.core.assem(edof, K, Ke, f=None, fe=None)[source]

Assemble element matrices Ke ( and fe ) into the global stiffness matrix K ( and the global force vector f ) according to the topology matrix edof.

Parameters:

edof dof topology array K the global stiffness matrix Ke element stiffness matrix f the global force vector fe element force vector

Output parameters:

K the new global stiffness matrix f the new global force vector fe element force vector

calfem.core.bar1e(ex, ep, eq=None)[source]

Ke = bar1e (ex, ep) Ke, fe = bar1e(ex, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a onedimensional bar element.

INPUT: ex = [x1 x2] element node coordinates

ep = [E A] element properties;

E: Young’s modulus A: cross section area ka: axial spring stiffness

eq = [qX] distributed load

OUTPUT: Kebar stiffness matrix [2 x 2]

fe : element load vector [2 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-10-22

O Dahlblom 2022-11-14 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar1s(ex, ep, ed, eq=None, nep=None)[source]

es = bar1s(ex, ep, ed) es = bar1s(ex, ep, ed, eq) es, edi, eci = bar1s(ex, ep, ed, eq, nep) ————————————————————- PURPOSE Compute section forces in one dimensional bar element

INPUT: ex = [x1 x2] element node coordinates

ep = [E A] element properties,

E: Young’s modulus A: cross section area

ed = [u1 u2] element displacement vector

eq = [qX] distributed load

nep : number of evaluation points ( default=2 )

OUTPUT: es = [N1 ; section forces, local directions, in

N2 ; nep points along the beam, dim(es)= nep x 1 …]

edi = [u1 ; element displacements, local directions,

u2 ; in n points along the bar, dim(edi)= nep x 1 …]

eci = [x1; evaluation points on the local x-axis,

x2; (x1=0 and xn=L) …]


LAST MODIFIED: O Dahlblom 2021-02-25

O Dahlblom 2022-11-14 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar1we(ex, ep, eq=None)[source]

Ke = bar1we (ex, ep) Ke, fe = bar1we(ex, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a onedimensional bar element with axial springs.

INPUT: ex = [x1 x2] element node coordinates

ep = [E A kX] element properties;

E: Young’s modulus A: cross section area kX: axial spring stiffness

eq = [qX] distributed load

OUTPUT: Kebar stiffness matrix [2 x 2]

fe : element load vector [2 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-12-17

O Dahlblom 2022-10-19 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar1ws(ex, ep, ed, eq=None, nep=None)[source]

es = bar1ws(ex, ep, ed) es = bar1ws(ex, ep, ed, eq) es, edi, eci = bar1ws(ex, ep, ed, eq, nep) ————————————————————- PURPOSE Compute section forces in one dimensional bar element

INPUT: ex = [x1 x2] element node coordinates

ep = [E A kX] element properties,

E: Young’s modulus A: cross section area kX: axial spring stiffness

ed = [u1 u2] element displacement vector

eq = [qX] distributed load

nep : number of evaluation points ( default=2 )

OUTPUT: es = [N1 ; section forces, local directions, in

N2 ; nep points along the beam, dim(es)= nep x 1 …]

edi = [u1 ; element displacements, local directions,

u2 ; in n points along the bar, dim(edi)= nep x 1 …]

eci = [x1; evaluation points on the local x-axis,

x2; (x1=0 and xn=L) …]


LAST MODIFIED: O Dahlblom 2021-02-25

O Dahlblom 2022-11-14 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar2e(ex, ey, ep, eq=None)[source]

Ke = bar2e(ex, ey, ep) Ke, fe = bar2e(ex, ey, ep, eq) ———————————————————————- PURPOSE Compute the element stiffness matrix for two dimensional bar element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] element node coordinates

ep = [E A] element properties;

E: Young’s modulus A: cross section area

eq = [qX] distributed load

OUTPUT: Kebar stiffness matrix [4 x 4]

fe : element load vector [4 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-10-20

O Dahlblom 2022-11-16 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar2ge(ex, ey, ep, QX)[source]

Ke = bar2ge(ex, ey, ep, QX)

PURPOSE Compute element stiffness matrix for two dimensional geometric nonlinear bar element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] element node coordinates

ep = [E A] element properties;

E: Young’s modulus A: cross section area

QX: axial force in the bar

OUTPUT: Ke : bar stiffness matrix [4 x 4]

LAST MODIFIED: O Dahlblom 2015-12-17

O Dahlblom 2022-11-16 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar2gs(ex, ey, ep, ed, nep=None)[source]

es, QX, edi, eci = bar2s(ex, ey, ep, ed) es, QX, edi, eci = bar2s(ex, ey, ep, ed, nep)

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] element node coordinates

ep = [E A] element properties,

E: Young’s modulus A: cross section area

ed = [u1 … u4] element displacement vector

nep : number of evaluation points ( default=2 )

OUTPUT: es = [N1 ; section forces, local directions, in

N2 ; nep points along the beam, dim(es)= nep x 1 …]

QX: axial force

edi = [u1 ; element displacements, local directions,

u2 ; in n points along the bar, dim(edi)= nep x 1 …]

eci = [x1; evaluation points on the local x-axis,

x2; (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2015-10-20

O Dahlblom 2022-11-16 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.bar2s(ex, ey, ep, ed, eq=None, nep=None)[source]

es = bar2s(ex, ey, ep, ed) es = bar2s(ex, ey, ep, ed, eq) es, edi, eci = bar2s(ex, ey, ep, ed, eq, nep) ————————————————————- PURPOSE Compute normal force in two dimensional bar element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] element node coordinates

ep = [E A] element properties,

E: Young’s modulus A: cross section area

ed = [u1 … u4] element displacement vector

eq = [qX] distributed load

nep : number of evaluation points ( default=2 )

OUTPUT: es = [N1 ; section forces, local directions, in

N2 ; nep points along the beam, dim(es)= nep x 1 …]

edi = [u1 ; element displacements, local directions,

u2 ; in n points along the bar, dim(edi)= nep x 1 …]

eci = [x1; evaluation points on the local x-axis,

x2; (x1=0 and xn=L) …]


LAST MODIFIED: O Dahlblom 2015-12-04

O Dahlblom 2022-11-16 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar3e(ex, ey, ez, ep, eq=None)[source]

Ke = bar2e(ex, ey, ez, ep) Ke, fe = bar2e(ex, ey, ez, ep, eq) ———————————————————————- PURPOSE Compute the element stiffness matrix for three dimensional bar element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] ez = [z1 z2]

ep = [E A] element properties;

E: Young’s modulus A: cross section area

eq = [qX] distributed load

OUTPUT: Kebar stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-10-19

O Dahlblom 2022-11-18 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.bar3s(ex, ey, ez, ep, ed, eq=None, nep=None)[source]

es = bar3s(ex, ey, ez, ep, ed) es = bar3s(ex, ey, ez, ep, ed, eq) es, edi, eci = bar3s(ex, ey, ez, ep, ed, eq, nep) ————————————————————- PURPOSE Compute normal force in three dimensional bar element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2] ez = [z1 z2]

ep = [E A] element properties,

E: Young’s modulus A: cross section area

ed = [u1 … u4] element displacement vector

eq = [qX] distributed load

nep : number of evaluation points ( default=2 )

OUTPUT: es = [N1 ; section forces, local directions, in

N2 ; nep points along the beam, dim(es)= nep x 1 …]

edi = [u1 ; element displacements, local directions,

u2 ; in n points along the bar, dim(edi)= nep x 1 …]

eci = [x1; evaluation points on the local x-axis,

x2; (x1=0 and xn=L) …]


LAST MODIFIED: O Dahlblom 2021-09-01

O Dahlblom 2022-11-18 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam1e(ex, ep, eq=None)[source]

Ke = beam1e(ex, ep) Ke, fe = beam1e(ex, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a one dimensional beam element.

INPUT: ex = [x1 x2] element node coordinates

ep = [E I] element properties;

E: Young’s modulus I: moment of inertia

eq = [qY] distributed load

OUTPUT: Kebeam stiffness matrix [4 x 4]

fe : element load vector [4 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2019-01-09

O Dahlblom 2022-10-25 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam1s(ex, ep, ed, eq=None, nep=None)[source]

es = beam1s(ex, ep, ed) es = beam1s(ex, ep, ed, eq) es, ed, ec = beam1s(ex, ep, ed, eq, nep) ————————————————————- PURPOSE Compute section forces in one dimensional beam element (beam1e).

INPUT ex = [x1 x2] element node coordinates

ep = [E I] element properties,

E: Young’s modulus I: moment of inertia

ed = [u1 … u4] element displacements

eq = [qy] distributed loads, local directions

nep : number of evaluation points ( default=2 )

OUTPUT: es = [V1 M1 ; section forces, local directions, in

V2 M2 ; nep points along the beam, dim(es)= nep x 2 ……]

edi = [v1 ; element displacements, local directions,

v2 ; in nep points along the beam, dim(edi)= nep x 1

….]

eci = [x1; evaluation points on the local x-axis

x2; ..]


LAST MODIFIED: O Dahlblom 2021-09-01

O Dahlblom 2022-10-25 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam1we(ex, ep, eq=None)[source]

Ke = beam1we(ex, ep) Ke, fe = beam1we(ex, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a one dimensional beam element on elastic foundation.

INPUT: ex = [x1 x2] element node coordinates

ep = [E I kY] element properties;

E: Young’s modulus I: moment of inertia kY: transversal found. stiffness

eq = [qY] distributed load

OUTPUT: Ke: beam stiffness matrix [4 x 4]

fe: element load vector [4 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2016-02-17

O Dahlblom 2022-10-18 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam1ws(ex, ep, ed, eq=None, nep=None)[source]

es = beam1ws(ex, ep, ed) es = beam1ws(ex, ep, ed, eq) es, ed, ec = beam1ws(ex, ep, ed, eq, nep) ————————————————————- PURPOSE Compute section forces in one dimensional beam element on elastic foundation (beam1we).

INPUT: ex = [x1 x2] element node coordinates

ep = [E I kY] element properties,

E: Young’s modulus I: moment of inertia kY: transversal foundation stiffness

ed = [u1 … u4] element displacements

eq = [qy] distributed loads, local directions

nep number of evaluation points ( default=2 )

OUTPUT: es = [V1 M1 ; section forces, local directions, in

V2 M2 ; nep points along the beam, dim(es)= n x 2 ……]

edi = [v1 ; element displacements, local directions,

v2 ; in nep points along the beam, dim(edi)= n x 1 …]

eci = [x1 ; evaluation points on the local x-axis

x2 ; …]


LAST MODIFIED: O Dahlblom 2021-09-01

O Dahlblom 2022-10-18 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2crd(ex=None, ey=None, ed=None, mag=None)[source]

number of identical 2D Bernoulli beam elements.

INPUT: ex,ey,

ed, mag

OUTPUT: excd,eycd


LAST MODIFIED: P-E AUSTRELL 1993-10-15 Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2crd_old(ex, ey, ed, mag)[source]

number of identical 2D Bernoulli beam elements.

INPUT: ex,ey,

ed, mag

OUTPUT: excd,eycd


LAST MODIFIED: P-E AUSTRELL 1993-10-15

J Lindemann 2021-12-30 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2de(ex, ey, ep)[source]

Ke, Me = beam2de(ex, ey, ep) Ke, Me, Ce = beam2de(ex, ey, ep) ————————————————————- PURPOSE Calculate the stiffness matrix Ke, the mass matrix Me and the damping matrix Ce for a 2D elastic Bernoulli beam element.

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E,A,I,m,(a,b)] element properties;

E: Young’s modulus A: cross section area I: moment of inertia m: mass per unit length a,b: damping coefficients,

Ce=aMe+bKe

OUTPUT: Ke element stiffness matrix (6 x 6)

Me element mass martix Ce element damping matrix, optional


LAST MODIFIED: K Persson 1995-08-23

O Dahlblom 2022-12-08 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2ds(ex, ey, ep, ed, ev, ea)[source]

es = beam2ds(ex, ey, ep, ed, ev, ea)

PURPOSE Calculate the element forces for a number of identical (nie) 2D Bernoulli beam elements in dynamic analysis.

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E,A,I,m,(a,b)] element properties;

E: Young’s modulus A: cross section area I: moment of inertia m: mass per unit length a,b: damping coefficients,

Ce=aMe+bKe

ed : element displacement matrix

ev : element velocity matrix

ea : element acceleration matrix

OUTPUT: eselement forces in local directions,
= [-N1 -V1 -M1 N2 V2 M2;

……. ……] ; dim(es)= nie x 6


LAST MODIFIED: K Persson 1995-08-23

O Dahlblom 2022-12-08 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2e(ex, ey, ep, eq=None)[source]

Ke = beam2e(ex, ey, ep) Ke, fe = beam2e(ex, ey, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a two dimensional beam element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2]

ep = [E A I] element properties;

E: Young’s modulus A: Cross section area I: moment of inertia

eq = [qX qY] distributed loads, local directions

OUTPUT: Keelement stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-08-17

O Dahlblom 2022-11-21 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2ge(ex, ey, ep, QX, eq=None)[source]

Ke = beam2ge(ex, ey, ep, QX) Ke, fe = beam2ge(ex, ey, ep, QX, eq) ————————————————————- PURPOSE

Compute the element stiffness matrix for a two dimensional beam element with respect to geometric nonlinearity.

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E, A, I] element properties;

E: Young’s modulus A: cross section area I: moment of inertia

QX axial force in the beam

eq = [qY] distributed transverse load

OUTPUT: Keelement stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-12-17

O Dahlblom 2022-12-08 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2gs(ex, ey, ep, ed, QX, eq=None, nep=None)[source]

es, QX = beam2gs(ex, ey, ep, ed, QX) es, QX = beam2gs(ex, ey, ep, ed, QX, eq) es, QX, edi, eci = beam2gs(ex, ey, ep, ed, QX, eq, nep)

beam element (beam2ge).

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E, A, I] element properties;

E: Young’s modulus A: cross section area I: moment of inertia

ed = [u1, … ,u6] element displacement vector

QX axial force

eq = [qy] distributed transverse load

nep number of evaluation points ( default=2 )

OUTPUT: es = [ N1 V1 M1 section forces, local directions, in

N2 V2 M2 n points along the beam, dim(es)= n x 3 ……..]

QX axial force

edi = [ u1 v1 element displacements, local directions,

u2 v2 in n points along the beam, dim(es)= n x 2 …..]

eci = [ x1 local x-coordinates of the evaluation

x2 points, (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2021-09-01

O Dahlblom 2022-12-06 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.beam2gxe(ex, ey, ep, QX, eq=None)[source]

Ke = beam2gxe(ex, ey, ep, QX) Ke, fe = beam2gxe(ex, ey, ep, QX, eq) ————————————————————- PURPOSE

Compute the element stiffness matrix for a two dimensional beam element with respect to geometric nonlinearity with exact solution.

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E, A, I] element properties;

E: Young’s modulus A: cross section area I: moment of inertia

QX axial force in the beam

eq distributed transverse load

OUTPUT: Keelement stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2021-06-21

O Dahlblom 2022-12-06 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2gxs(ex, ey, ep, ed, QX, eq=None, nep=None)[source]

es, QX = beam2gxs(ex, ey, ep, ed, QX) es, QX = beam2gxs(ex, ey, ep, ed, QX, eq) es, QX, edi, eci = beam2gxs(ex, ey, ep, ed, QX, eq, nep)

beam element (beam2gxe).

INPUT: ex = [x1, x2]

ey = [y1, y2] element node coordinates

ep = [E, A, I] element properties;

E: Young’s modulus A: cross section area I: moment of inertia

ed = [u1, … ,u6] element displacement vector

QX axial force

eq = [qy] distributed transverse load

nep number of evaluation points ( default=2 )

OUTPUT: es = [ N1 V1 M1 section forces, local directions, in

N2 V2 M2 n points along the beam, dim(es)= n x 3 ……..]

QX axial force

edi = [ u1 v1 element displacements, local directions,

u2 v2 in n points along the beam, dim(es)= n x 2 …..]

eci = [ x1 local x-coordinates of the evaluation

x2 points, (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2021-09-17

O Dahlblom 2022-12-06 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.beam2s(ex, ey, ep, ed, eq=None, nep=None)[source]

es = beam2s(ex, ey, ep, ed) es = beam2s(ex, ey, ep, ed, eq) es, edi, eci = beam2s(ex, ey, ep, ed, eq, nep)

INPUT: ex = [x1 x2]

ey = [y1 y2] element node coordinates

ep = [E A I] element properties,

E: Young’s modulus A: cross section area I: moment of inertia

ed = [u1 … u6] element displacements

eq = [qx qy] distributed loads, local directions

nep number of evaluation points ( default=2 )

OUTPUT: es = [ N1 V1 M1 section forces, local directions, in

N2 V2 M2 n points along the beam, dim(es)= n x 3 ……..]

edi = [ u1 v1 element displacements, local directions,

u2 v2 in n points along the beam, dim(es)= n x 2 …..]

eci = [ x1 local x-coordinates of the evaluation

x2 points, (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2021-09-08

O Dahlblom 2022-11-21 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.beam2te(ex, ey, ep, eq=None)[source]

Ke = beam2te(ex, ey, ep) Ke, fe = beam2te(ex, ey, ep, eq) ————————————————————- PURPOSE Compute the stiffness matrix for a two dimensional Timoshenko beam element.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2]

ep = [E Gm A I ks] element properties;

E: Young’s modulus G: shear modulus A: Cross section area I: moment of inertia ks: shear correction factor

eq = [qX qY] distributed loads, local directions

OUTPUT: Keelement stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2021-11-05

O Dahlblom 2022-12-08 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2ts(ex, ey, ep, ed, eq=None, nep=None)[source]

es = beam2ts(ex, ey, ep, ed) es = beam2ts(ex, ey, ep, ed, eq) es, edi, eci = beam2s(ex, ey, ep, ed, eq, nep)

element (beam2te).

INPUT: ex = [x1 x2]

ey = [y1 y2] element node coordinates

ep = [E G A I ks] element properties,

E: Young’s modulus G: shear modulus A: cross section area I: moment of inertia ks: shear correction factor

ed = [u1 … u6] element displacements

eq = [qx qy] distributed loads, local directions

nep number of evaluation points ( default=2 )

OUTPUT: es = [ N1 V1 M1 section forces, local directions, in

N2 V2 M2 n points along the beam, dim(es)= n x 3 ……..]

edi = [ u1 v1 teta1 element displacements, local directions,

u2 v2 teta2 in n points along the beam, dim(es)= n x 2 ………..] (Note! For Timoshenko beam element the rotation of the cross

section is not equal to dv/dx)

eci = [ x1 local x-coordinates of the evaluation

x2 points, (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2021-11-05

O Dahlblom 2022-12-08 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.beam2we(ex, ey, ep, eq=None)[source]

Ke = beam2we(ex, ey, ep) Ke, fe = beam2we(ex, ey, ep, eq) ————————————————————- PURPOSE

Compute the stiffness matrix for a two dimensional beam element on elastic foundation.

INPUT: ex = [x1 x2] element node coordinates

ey = [y1 y2]

ep = [E,A,I,kX,kY] element properties;

E: Young’s modulus A: Cross section area I: moment of inertia kX: axial foundation stiffness kY: transversal foundation stiffness

eq = [qX qY] distributed loads, local directions

OUTPUT: Keelement stiffness matrix [6 x 6]

fe : element load vector [6 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-08-07

O Dahlblom 2022-11-21 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam2ws(ex, ey, ep, ed, eq=None, nep=None)[source]

es = beam2ws(ex, ey, ep, ed) es = beam2ws(ex, ey, ep, ed, eq) es, edi, eci = beam2ws(ex, ey, ep, ed, eq, nep)

on elastic foundation.

INPUT: ex = [x1 x2]

ey = [y1 y2] element node coordinates

ep = [E,A,I,kX,kY] element properties,

E: Young’s modulus A: cross section area I: moment of inertia kX: axial foundation stiffness kY: transversal foundation stiffness

ed = [u1 … u6] element displacements

eq = [qx qy] distributed loads, local directions

nep number of evaluation points ( default=2 )

OUTPUT: es = [ N1 V1 M1 section forces, local directions, in

N2 V2 M2 n points along the beam, dim(es)= n x 3 ……..]

edi = [ u1 v1 element displacements, local directions,

u2 v2 in n points along the beam, dim(es)= n x 2 …..]

eci = [ x1 local x-coordinates of the evaluation

x2 points, (x1=0 and xn=L) …]

LAST MODIFIED: O Dahlblom 2022-09-30

O Dahlblom 2022-11-21 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University

calfem.core.beam3e(ex, ey, ez, eo, ep, eq=None)[source]

Ke = beam3e(ex, ey, ez, eo, ep) Ke, fe = beam3e(ex, ey, ez, eo, ep, eq) ————————————————————- PURPOSE Calculate the stiffness matrix for a 3D elastic Bernoulli beam element.

INPUT: ex = [x1 x2]

ey = [y1 y2] ez = [z1 z2] element node coordinates

eo = [xz yz zz] orientation of local z-axis

ep = [E G A Iy Iz Kv] element properties

E: Young’s modulus G: Shear modulus A: Cross section area Iy: Moment of inertia, local y-axis Iz: Moment of inertia, local z-axis Kv: Saint-Venant’s torsion constant

eq = [qX qY qZ qW] distributed loads, local directions

OUTPUT: Ke : element stiffness matrix [12 x 12]

fe : element load vector [12 x 1] (if eq!=None)


LAST MODIFIED: O Dahlblom 2015-10-19

O Dahlblom 2022-11-21 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.beam3s(ex, ey, ez, eo, ep, ed, eq=None, nep=None)[source]

es = beam3s(ex, ey, ez, eo, ep, ed) es = beam3s(ex, ey, ez, eo, ep, ed, eq) es, edi, eci = beam3s(ex, ey, ez, eo, ep, ed, eq, nep) ————————————————————- PURPOSE Calculate the variation of the section forces and displacements along a three-dimensional beam element.

INPUT: ex = [x1 x2]

ey = [y1 y2] ez = [z1 z2] element node coordinates

eo = [xz yz zz] orientation of local z-axis

ep = [E G A Iy Iz Kv] element properties

E: Young’s modulus G: Shear modulus A: Cross section area Iy: Moment of inertia, local y-axis Iz: Moment of inertia, local z-axis Kv: Saint-Venant’s torsion constant

ed = [u1 … u12] element displacements

eq = [qX qY qZ qW] distributed loads, local directions

nep number of evaluation points ( default=2 )

OUTPUT: es = [[N1,Vy1,Vz1,T1,My1,Mz1], section forces in n points along

[N2,Vy2,Vz2,T2,My2,Mz2], the local x-axis [..,…,…,..,…,…], [Nn,Vyn,Vzn,Tn,Myn,Mzn]]

edi = [[u1,v1,w1,fi1], displacements in n points along

[u2,v2,w2,fi2], the local x-axis [..,..,..,…], [un,vn,wn,fin]]

eci = [[x1], local x-coordinates of the evaluation

[x2], points [..], [xn]]

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.coord_extract(edof, coords, dofs, nen=-1)

Create element coordinate matrices ex, ey, ez from edof coord and dofs matrices.

Parameters:

edof [nel x (nen * nnd)], nnd = number of node dofs coords [ncoords x ndims], ndims = node dimensions dofs [ncoords x nnd]

Returns:

ex if ndims = 1 ex, ey if ndims = 2 ex, ey, ez if ndims = 3

calfem.core.coordxtr(edof, coords, dofs, nen=-1)[source]

Create element coordinate matrices ex, ey, ez from edof coord and dofs matrices.

Parameters:

edof [nel x (nen * nnd)], nnd = number of node dofs coords [ncoords x ndims], ndims = node dimensions dofs [ncoords x nnd]

Returns:

ex if ndims = 1 ex, ey if ndims = 2 ex, ey, ez if ndims = 3

calfem.core.create_dofs(nCoords, nDof)[source]

Create dof array [nCoords x nDof]

calfem.core.createdofs(nCoords, nDof)

Create dof array [nCoords x nDof]

calfem.core.effmises(es, ptype)[source]

Calculate effective von mises stresses.

Parameters:

es

ptype= 1: plane stress

2: plane strain 3: axisymmetry 4: three dimensional

es = [[sigx,sigy,[sigz],tauxy] element stress matrix

[ …… ]] one row for each element

Returns:

eseff = [eseff_0 .. eseff_nel-1]

calfem.core.eigen(K, M, b=None)[source]

Solve the generalized eigenvalue problem |K-LM|X = 0, considering boundary conditions

Parameters:

K global stiffness matrix, dim(K) = ndof x ndof M global mass matrix, dim(M) = ndof x ndof b boundary condition vector, dim(b) = nbc x 1

Returns:

L eigenvalue vector, dim(L) = (ndof-nbc) x 1 X eigenvectors, dim(X) = ndof x (ndof-nbc)

calfem.core.error(msg)[source]

Write msg to error log.

calfem.core.exception_logging(exctype, value, tb)[source]

Log exception by using the root logger.

Parameters

exctype : type value : NameError tb : traceback

calfem.core.extractEldisp(edof, a)

Extract element displacements from the global displacement vector according to the topology matrix edof.

Parameters:

a the global displacement vector edof dof topology array

Returns:

ed: element displacement array

calfem.core.extract_ed(edof, a)

Extract element displacements from the global displacement vector according to the topology matrix edof.

Parameters:

a the global displacement vector edof dof topology array

Returns:

ed: element displacement array

calfem.core.extract_eldisp(edof, a)[source]

Extract element displacements from the global displacement vector according to the topology matrix edof.

Parameters:

a the global displacement vector edof dof topology array

Returns:

ed: element displacement array

calfem.core.flw2i4e(ex, ey, ep, D, eq=None)[source]

Compute element stiffness (conductivity) matrix for 4 node isoparametric field element

Parameters:

ex = [x1 x2 x3 x4] element coordinates ey = [y1 y2 y3 y4]

ep = [t ir] thickness and integration rule

D = [[kxx kxy],

[kyx kyy]] constitutive matrix

eq heat supply per unit volume

Returns:

Ke element ‘stiffness’ matrix (4 x 4) fe element load vector (4 x 1)

calfem.core.flw2i4s(ex, ey, ep, D, ed)[source]

Compute flows or corresponding quantities in the 4 node isoparametric element.

Parameters:

ex = [x1 x2 x3 x4] element coordinates ey = [y1 y2 y3 y4]

ep = [t ir] thickness and integration rule

D = [[kxx kxy],

[kyx kyy]] constitutive matrix

ed = [u1, u2, u3, u4] u1,u2,u3,u4: nodal values

Returns:
es = [[qx, qy],

[.., ..]] element flows

et = [[qx, qy],

[… ..]] element gradients

eci=[[ix1, iy1], Gauss point location vector

[… …], nint: number of integration points [ix(nint), iy(nint)]

calfem.core.flw2i8e(ex, ey, ep, D, eq=None)[source]

Compute element stiffness (conductivity) matrix for 8 node isoparametric field element.

Parameters:

ex = [x1, …, x8] element coordinates ey = [y1, …, y8]

ep = [t, ir] thickness and integration rule

D = [[kxx, kxy],

[kyx, kyy]] constitutive matrix

eq heat supply per unit volume

Returns:

Ke element ‘stiffness’ matrix (8 x 8) fe element load vector (8 x 1)

calfem.core.flw2i8s(ex, ey, ep, D, ed)[source]

Compute flows or corresponding quantities in the 8 node isoparametric element.

Parameters:

ex = [x1,x2,x3….,x8] element coordinates ey = [y1,y2,y3….,y8]

ep = [t,ir] thickness and integration rule

D = [[kxx,kxy],

[kyx,kyy]] constitutive matrix

ed = [u1,….,u8] u1,….,u8: nodal values

Returns:
es = [[qx,qy],

[..,..]] element flows

et = [[qx,qy],

[..,..]] element gradients

eci=[[ix1,iy1], Gauss point location vector

[…,…], nint: number of integration points [ix(nint),iy(nint)]]

calfem.core.flw2qe(ex, ey, ep, D, eq=None)[source]

Compute element stiffness (conductivity) matrix for a triangular field element.

Parameters:

ex = [x1, x2, x3, x4] ey = [y1, y2, y3, y4] element coordinates

ep = [t] element thickness

D = [[kxx, kxy],

[kyx, kyy]] constitutive matrix

eq heat supply per unit volume

Returns:

Ke element ‘stiffness’ matrix (4 x 4)

fe element load vector (4 x 1)

calfem.core.flw2qs(ex, ey, ep, D, ed, eq=None)[source]

Compute flows or corresponding quantities in the quadrilateral field element.

Parameters:

ex = [x1, x2, x3, x4] ey = [y1, y2, y3, y4] element coordinates

ep = [t] element thickness

D = [[kxx, kxy],

[kyx, kyy]] constitutive matrix

ed = [[u1, u2, u3, u4],

[.., .., .., ..]] u1,u2,u3,u4: nodal values

eq heat supply per unit volume

Returns:

es = [[qx, qy],

[.., ..]] element flows

et = [[gx, gy],

[.., ..]] element gradients

calfem.core.flw2te(ex, ey, ep, D, eq=None)[source]

Compute element stiffness (conductivity) matrix for a triangular field element.

Parameters:

ex = [x1 x2 x3] ey = [y1 y2 y3] element coordinates

ep = [t] element thickness

D = [kxx kxy;

kyx kyy] constitutive matrix

eq heat supply per unit volume

Returns:

Ke element ‘stiffness’ matrix (3 x 3)

fe element load vector (3 x 1)

calfem.core.flw2ts(ex, ey, D, ed)[source]

Compute flows or corresponding quantities in the triangular field element.

Parameters:

ex = [x1 x2 x3] ey = [y1 y2 y3] element coordinates

D = [kxx kxy

kyx kyy] constitutive matrix

ed =[u1 u2 u3] u1,u2,u3: nodal values

Returns:

es=[ qx qy ]

… ..] element flows

et=[ gx gy ]

… ..] element gradients

calfem.core.flw3i8e(ex, ey, ez, ep, D, eq=None)[source]

Compute element stiffness (conductivity) matrix for 8 node isoparametric field element.

Parameters:

ex = [x1,x2,x3,…,x8] ey = [y1,y2,y3,…,y8] element coordinates ez = [z1,z2,z3,…,z8]

ep = [ir] Ir: Integration rule

D = [[kxx,kxy,kxz],

[kyx,kyy,kyz], [kzx,kzy,kzz]] constitutive matrix

eq heat supply per unit volume

Output:

Ke element ‘stiffness’ matrix (8 x 8) fe element load vector (8 x 1)

calfem.core.flw3i8s(ex, ey, ez, ep, D, ed)[source]

Compute flows or corresponding quantities in the 8 node (3-dim) isoparametric field element.

Parameters:

ex = [x1,x2,x3,…,x8] ey = [y1,y2,y3,…,y8] element coordinates ez = [z1,z2,z3,…,z8]

ep = [ir] Ir: Integration rule

D = [[kxx,kxy,kxz],

[kyx,kyy,kyz], [kzx,kzy,kzz]] constitutive matrix

ed = [[u1,….,u8], element nodal values

[..,….,..]]

Output:

es = [[qx,qy,qz],

[..,..,..]] element flows(s)

et = [[qx,qy,qz], element gradients(s)

[..,..,..]]

eci = [[ix1,ix1,iz1], location vector

[…,…,…], nint: number of integration points [ix(nint),iy(nint),iz(nint)]]

calfem.core.gfunc(G, dt)[source]

Form vector with function values at equally spaced points by linear interpolation

Parameters:

G = [t_i, g_i] t_i: time i, g_i: g(t_i)

dim(G) = np x 2, np = number of points

dt time step

Returns:

t 1-D vector with equally spaced time points g 1-D vector with corresponding function values

calfem.core.hooke(ptype, E, v)[source]

Calculate the material matrix for a linear elastic and isotropic material.

Parameters:

ptype= 1: plane stress

2: plane strain 3: axisymmetry 4: three dimensional

E Young’s modulus v Poissons const.

Returns:

D material matrix

calfem.core.info(msg)[source]

Write msg to info log.

calfem.core.plani4e(ex, ey, ep, D, eq=None)[source]

Calculate the stiffness matrix for a 4 node isoparametric element in plane strain or plane stress.

Parameters:

ex = [x1 … x4] element coordinates. Row array ey = [y1 … y4]

ep =[ptype, t, ir] ptype: analysis type

t : thickness ir: integration rule

D constitutive matrix

eq = [bx; by] bx: body force in x direction
by: body force in y direction

Any array with 2 elements acceptable

Returns:

Ke : element stiffness matrix (8 x 8) fe : equivalent nodal forces (8 x 1)

calfem.core.planqe(ex, ey, ep, D, eq=None)[source]

Calculate the stiffness matrix for a quadrilateral plane stress or plane strain element.

Parameters:

ex=[x1 x2 x3 x4] element coordinates ey=[y1 y2 y3 y4]

ep = [ptype, t] ptype: analysis type

t: element thickness

D constitutive matrix

eq = [bx; bx: body force in x direction

by] by: body force in y direction

OUTPUT: Keelement stiffness matrix (8 x 8)

fe : equivalent nodal forces (row array)

calfem.core.planqs(ex, ey, ep, D, ed, eq=None)[source]

Calculate element normal and shear stress for a quadrilateral plane stress or plane strain element.

Parameters:

ex = [x1 x2 x3 x4] element coordinates ey = [y1 y2 y3 y4]

ep = [ptype, t] ptype: analysis type

t: thickness

D constitutive matrix

ed = [u1 u2 ..u8] element displacement vector

eq = [[bx] bx: body force in x direction

[by]] by: body force in y direction

OUTPUT: es = [ sigx sigy (sigz) tauxy] element stress array

et = [ epsx epsy (epsz) gamxy] element strain array

calfem.core.plante(ex, ey, ep, D, eq=None)[source]

Calculate the stiffness matrix for a triangular plane stress or plane strain element.

Parameters:

ex = [x1,x2,x3] element coordinates ey = [y1,y2,y3]

ep = [ptype,t] ptype: analysis type

t: thickness

D constitutive matrix

eq = [[bx], bx: body force x-dir

[by]] by: body force y-dir

Returns:

Ke element stiffness matrix (6 x 6) fe equivalent nodal forces (6 x 1) (if eq is given)

calfem.core.plantf(ex, ey, ep, es)[source]

Compute internal element force vector in a triangular element in plane stress or plane strain.

Parameters:

ex = [x1,x2,x3] node coordinates ey = [y1,y2,y3]

ep = [ptype,t] ptype: analysis type

t: thickness

es = [[sigx,sigy,[sigz],tauxy] element stress matrix

[ …… ]] one row for each element

OUTPUT:

fe = [[f1],[f2],…,[f8]] internal force vector

calfem.core.plants(ex, ey, ep, D, ed)[source]

Calculate element normal and shear stress for a triangular plane stress or plane strain element.

INPUT: ex = [x1 x2 x3] element coordinates

ey = [y1 y2 y3]

ep = [ptype t ] ptype: analysis type

t: thickness

D constitutive matrix

ed =[u1 u2 …u6 element displacement vector

…… ] one row for each element

OUTPUT: es = [ sigx sigy [sigz] tauxy element stress matrix

…… ] one row for each element

et = [ epsx epsy [epsz] gamxy element strain matrix

…… ] one row for each element

calfem.core.platre(ex, ey, ep, D, eq=None)[source]

Calculate the stiffness matrix for a rectangular plate element. NOTE! Element sides must be parallel to the coordinate axis.

Parameters:

ex = [x1,x2,x3,x4] element coordinates ey = [y1,y2,y3,y4]

ep = [t] thicknes

D constitutive matrix for

plane stress

eq = [qz] load/unit area

Returns:

Ke element stiffness matrix (12 x 12) fe equivalent nodal forces (12 x 1)

calfem.core.soli8e(ex, ey, ez, ep, D, eqp=None)[source]

Ke=soli8e(ex,ey,ez,ep,D) [Ke,fe]=soli8e(ex,ey,ez,ep,D,eq) ————————————————————- PURPOSE Calculate the stiffness matrix for a 8 node (brick) isoparametric element.

INPUT: ex = [x1 x2 x3 … x8]

ey = [y1 y2 y3 … y8] element coordinates ez = [z1 z2 z3 … z8]

ep = [ir] ir integration rule

D constitutive matrix

eq = [bx; by; bz] bx: body force in x direction

by: body force in y direction bz: body force in z direction

OUTPUT: Keelement stiffness matrix

fe : equivalent nodal forces


LAST MODIFIED: M Ristinmaa 1995-10-25

J Lindemann 2022-01-24 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.soli8s(ex, ey, ez, ep, D, ed)[source]

[es,et]=soli8s(ex,ey,ez,ep,D,ed)

8 node (brick) isoparametric element.

INPUT: ex = [x1 x2 x3 … x8]

ey = [y1 y2 y3 … y8] element coordinates ez = [z1 z2 z3 … z8]

ep = [Ir] Ir: integration rule

D constitutive matrix

ed = [u1 u2 ..u24] element displacement vector

OUTPUT: es = [ sigx sigy sigz sigxy sigyz sigxz ;

…… … ]

element stress matrix, one row for each integration point

es = [ eps epsy epsz epsxy epsyz epsxz ;

…… … ]

element strain matrix, one row for each integration point


LAST MODIFIED: M Ristinmaa 1995-10-25

J Lindemann 2022-02-23 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.solveq(K, f, bcPrescr=None, bcVal=None)[source]

Solve static FE-equations considering boundary conditions.

Parameters:

K global stiffness matrix, dim(K)= nd x nd f global load vector, dim(f)= nd x 1

bcPrescr 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values.

If not given all prescribed dofs are assumed 0.

Returns:

a solution including boundary values Q reaction force vector

dim(a)=dim(Q)= nd x 1, nd : number of dof’s

calfem.core.spring1e(ep)[source]

Ke = spring1e(ep)

PURPOSE Compute element stiffness matrix for spring element.

INPUT: ep = [k] spring stiffness or analog quantity

OUTPUT: Ke : spring stiffness matrix, [2 x 2]

LAST MODIFIED: P-E Austrell 1994-11-02

O Dahlblom 2022-11-15 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.spring1s(ep, ed)[source]

es = spring1s(ep, ed)

PURPOSE Compute element force in spring element (spring1e).

INPUT: ep = [k] spring stiffness or analog quantity

ed = [u1 u2] element displacement vector

OUTPUT: es = [N] element force

LAST MODIFIED: P-E AUSTRELL 1994-11-02

O Dahlblom 2022-11-14 (Python version)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.core.spsolveq(K, f, bcPrescr, bcVal=None)[source]

Solve static FE-equations considering boundary conditions.

Parameters:

K global stiffness matrix, dim(K)= nd x nd f global load vector, dim(f)= nd x 1

bcPrescr 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values.

If not given all prescribed dofs are assumed 0.

Returns:

a solution including boundary values Q reaction force vector

dim(a)=dim(Q)= nd x 1, nd : number of dof’s

calfem.core.statcon(K, f, cd)[source]

Condensation of static FE-equations according to the vector cd.

Parameters:

K global stiffness matrix, dim(K) = nd x nd f global load vector, dim(f)= nd x 1

cd vector containing dof’s to be eliminated

dim(cd)= nc x 1, nc: number of condensed dof’s

Returns:

K1 condensed stiffness matrix,

dim(K1)= (nd-nc) x (nd-nc)

f1 condensed load vector, dim(f1)= (nd-nc) x 1

calfem.core.step1(K, C, f, a0, bc, ip, times, dofs)[source]

Algorithm for dynamic solution of first-order FE equations considering boundary conditions.

Parameters:

K conductivity matrix, dim(K) = ndof x ndof C capacity matrix, dim(C) = ndof x ndof f load vector, dim(f) = ndof x (nstep + 1),

If dim(f) = ndof x 1, the values are kept constant during time integration

a0 initial vector a(0), dim(a0) = ndof x 1 bc boundary condition matrix, dim(bc) = nbc x (nstep + 2)

where nbc = number of prescribed degrees of freedom (either constant or time-dependent) The first column contains the numbers of the prescribed degrees of freedom and the subsequent columns contain the time history. If dim(bc) = nbc x 2, the values from the second column are kept constant during time integration

ip array [dt, tottime, alpha], where

dt is the size of the time increment, tottime is the total time, alpha is time integration constant. Frequently used values of alpha are: alpha=0: forward difference; forward Euler, alpha=1/2: trapezoidal rule; Crank-Nicholson alpha=1: backward difference; backward Euler

times array [t(i) …] of times at which output should be written to a and da dofs array [dof(i) …] of degree of freedom numbers for which history output

should be written to ahist and dahist

Returns:

modelhist dictionary containing solution history for the whole model at following keys:
modelhist[‘a’] constains values of a at all timesteps,

alternatively at times specified in ‘times’ dim(modelhist[‘a’]) = ndof x (nstep + 1) or ndof x ntimes

modelhist[‘da’] constains values of da at all timesteps,

alternatively at times specified in ‘times’ dim(modelhist[‘da’]) = ndof x (nstep + 1) or ndof x ntimes

dofhist dictionary containing solution history for the degrees of freedom selected in ‘dofs’:
dofhist[‘a’] constains time history of a at the dofs specified in ‘dofs’

dim(dofhist[‘ahist’]) = ndof x (nstep + 1)

dofhist[‘da’] constains time history of daat the dofs specified in ‘dofs’

dim(dofhist[‘dahist’]) = ndof x (nstep + 1)

calfem.core.step2(K, C, M, f, a0, da0, bc, ip, times, dofs)[source]

Algorithm for dynamic solution of second-order FE equations considering boundary conditions.

Parameters:

K global stiffness matrix, dim(K) = ndof x ndof C global damping matrix, dim(C) = ndof x ndof

If there is no damping in the system, simply set C=[]

M global mass matrix, dim(M) = ndof x ndof f global load vector, dim(f) = ndof x (nstep + 1),

If dim(f) = ndof x 1, the values are kept constant during time integration

a0 initial displacement vector a(0), dim(a0) = ndof x 1 da0 initial velocity vector v(0), dim(da0) = ndof x 1 bc boundary condition matrix, dim(bc) = nbc x (nstep + 2)

where nbc = number of prescribed degrees of freedom (either constant or time-dependent) The first column contains the numbers of the prescribed degrees of freedom and the subsequent columns contain the time history. If dim(bc) = nbc x 2, the values from the second column are kept constant during time integration

ip array [dt, tottime, alpha, delta], where

dt is the size of the time increment, tottime is the total time, alpha and delta are time integration constants for the Newmark family of methods. Frequently used values of alpha and delta are: alpha=1/4, delta=1/2: average acceleration (trapezoidal) rule, alpha=1/6, delta=1/2: linear acceleration alpha=0, delta=1/2: central difference

times array [t(i) …] of times at which output should be written to a, da and d2a dofs array [dof(i) …] of degree of freedom numbers for which history output

should be written to ahist, dahist and d2ahist

Returns:

modelhist dictionary containing solution history for the whole model at following keys:
modelhist[‘a’] constains displacement values at all timesteps,

alternatively at times specified in ‘times’ dim(modelhist[‘a’]) = ndof x (nstep + 1) or ndof x ntimes

modelhist[‘da’] constains velocity values at all timesteps,

alternatively at times specified in ‘times’ dim(modelhist[‘da’]) = ndof x (nstep + 1) or ndof x ntimes

modelhist[‘d2a’] constains acceleration values at all timesteps,

alternatively at times specified in ‘times’ dim(modelhist[‘d2a’]) = ndof x (nstep + 1) or ndof x ntimes

dofhist dictionary containing solution history for the degrees of freedom selected in ‘dofs’:
dofhist[‘a’] constains displacement time history at the dofs specified in ‘dofs’

dim(dofhist[‘ahist’]) = ndof x (nstep + 1)

dofhist[‘da’] constains velocity time history at the dofs specified in ‘dofs’

dim(dofhist[‘dahist’]) = ndof x (nstep + 1)

dofhist[‘d2a’] constains acceleration time history at the dofs specified in ‘dofs’

dim(dofhist[‘d2ahist’]) = ndof x (nstep + 1)

calfem.core.stress2nodal(eseff, edof)[source]

Convert element effective stresses to nodal effective stresses.

Parameters:

eseff = [eseff_0 .. eseff_nel-1] edof = [dof topology array]

Returns:

ev: element value array [[ev_0_0 ev_0_1 ev_0_nen-1 ]

ev_nel-1_0 ev_nel-1_1 ev_nel-1_nen-1]

Geometry functions

CALFEM Geometry module

Contains functions and classes for describing geometry.

class calfem.geometry.Geometry[source]

Instances of GeoData can hold geometric data and be passed to GmshMesher in pycalfem_Mesh to mesh the geometry.

addBSpline(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)[source]

Adds a B-Spline curve

points - List of indices of control points that make a B-spline

[p1, p2, … , pn]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

addCircle(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)[source]

Adds a Circle arc curve.

points - list of 3 indices of point that make a circle arc smaller

than Pi. [startpoint, centerpoint, endpoint]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Marker applied to this curve. Default 0.

elOnCurv - Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

addEllipse(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)[source]

Adds a Ellipse arc curve.

points - List of 4 indices of point that make a ellipse arc smaller

than Pi. [startpoint, centerpoint, mAxisPoint, endpoint] Startpoint is the starting point of the arc. Centerpoint is the point at the center of the ellipse. MAxisPoint is any point on the major axis of the ellipse. Endpoint is the end point of the arc.

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

addPoint(coord, ID=None, marker=0, el_size=1)[source]

Adds a point.

Parameters: coord - [x, y] or [x, y, z].

List, not array.

ID - Positive integer ID of this point. If left unspecified the

point will be assigned the smallest unused point-ID. It is recommended to specify all point-IDs or none.

marker - Marker applied to this point. Default 0.

It is not a good idea to apply non-zero markers to points that are control points on B-splines or center points on circles/ellipses, since this can lead to “loose” nodes that are not part of any elements.

elSize - The size of elements at this point. Default 1. Use to make

a mesh denser or sparser here. Only affects unstructured meshes

addPoints(points, markers=None, ids=None, elSizes=None)[source]

Add points from a numpy-array

addRuledSurface(outer_loop, ID=None, marker=0)[source]

Adds a Ruled Surface (bent surface). Parameters: outer_loop - List of 3 or 4 curve IDs that make up the boundary of

the surface.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

addSpline(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)[source]

Adds a Spline curve

points - List of indices of control points that make a Spline

[p1, p2, … , pn]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

addSplines(points)[source]

Add splines from numpy array

addStructuredSurface(outer_loop, ID=None, marker=0)[source]

Adds a Structured Surface. Parameters: outer_loop - List of 4 curve IDs that make up the boundary of

the surface. The curves must be structured, i.e. their parameter ‘elOnCurv’ must be defined.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

addStructuredVolume(outer_surfaces, ID=None, marker=0)[source]

Adds a Structured Volume Parameters: outer_surfaces - List of surface IDs that make up the outer boundary of

the volume. The surfaces must be Structured Surfaces.

ID - Positive integer ID of this volume. If left unspecified

the volume will be assigned the smallest unused volume-ID. It is recommended to specify all volume-IDs or none.

marker - Integer. Marker applied to this volume. Default 0.

addSurface(outer_loop, holes=[], ID=None, marker=0)[source]

Adds a plane surface (flat). Parameters: outer_loop - List of curve IDs that make up the outer boundary of

the surface. The curves must lie in the same plane.

holes - List of lists of curve IDs that make up the inner

boundaries of the surface. The curves must lie in the same plane.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

addVolume(outer_surfaces, holes=[], ID=None, marker=0)[source]

Adds a Volume Parameters: outer_surfaces - List of surface IDs that make up the outer boundary of

the volume.

holes - List of lists of surface IDs that make up the inner

boundaries of the volume.

ID - Positive integer ID of this volume. If left unspecified

the volume will be assigned the smallest unused volume-ID. It is recommended to specify all volume-IDs or none.

marker - Integer. Marker applied to this volume. Default 0.

bounding_box_2d()[source]

Calculate bounding box geometry

bspline(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)

Adds a B-Spline curve

points - List of indices of control points that make a B-spline

[p1, p2, … , pn]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

circle(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)

Adds a Circle arc curve.

points - list of 3 indices of point that make a circle arc smaller

than Pi. [startpoint, centerpoint, endpoint]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Marker applied to this curve. Default 0.

elOnCurv - Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

ellipse(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)

Adds a Ellipse arc curve.

points - List of 4 indices of point that make a ellipse arc smaller

than Pi. [startpoint, centerpoint, mAxisPoint, endpoint] Startpoint is the starting point of the arc. Centerpoint is the point at the center of the ellipse. MAxisPoint is any point on the major axis of the ellipse. Endpoint is the end point of the arc.

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

getPointCoords(IDs=None)[source]

Returns an N-by-3 list of point coordinates if the parameter is a list of IDs. If the parameter is just a single integer then a single coordinate (simple 3-element list) is returned. If the parameter is undefined (or None) all point coords will be returned

get_point_coords(IDs=None)

Returns an N-by-3 list of point coordinates if the parameter is a list of IDs. If the parameter is just a single integer then a single coordinate (simple 3-element list) is returned. If the parameter is undefined (or None) all point coords will be returned

line(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)

Adds a Spline curve

points - List of indices of control points that make a Spline

[p1, p2, … , pn]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

point(coord, ID=None, marker=0, el_size=1)

Adds a point.

Parameters: coord - [x, y] or [x, y, z].

List, not array.

ID - Positive integer ID of this point. If left unspecified the

point will be assigned the smallest unused point-ID. It is recommended to specify all point-IDs or none.

marker - Marker applied to this point. Default 0.

It is not a good idea to apply non-zero markers to points that are control points on B-splines or center points on circles/ellipses, since this can lead to “loose” nodes that are not part of any elements.

elSize - The size of elements at this point. Default 1. Use to make

a mesh denser or sparser here. Only affects unstructured meshes

pointsOnCurves(IDs)[source]

Returns a list of all geometric points (not nodes) on the curves specified in IDs. IDs may be an integer or a list of integers.

removeCurve(ID)[source]

Removes the curve with this ID

removePoint(ID)[source]

Removes the point with this ID

removeSurface(ID)[source]

Removes the surface with this ID

removeVolume(ID)[source]

Removes the volume with this ID

ruledSurface(outer_loop, ID=None, marker=0)

Adds a Ruled Surface (bent surface). Parameters: outer_loop - List of 3 or 4 curve IDs that make up the boundary of

the surface.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

ruled_surf(outer_loop, ID=None, marker=0)

Adds a Ruled Surface (bent surface). Parameters: outer_loop - List of 3 or 4 curve IDs that make up the boundary of

the surface.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

ruled_surface(outer_loop, ID=None, marker=0)

Adds a Ruled Surface (bent surface). Parameters: outer_loop - List of 3 or 4 curve IDs that make up the boundary of

the surface.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

setCurveMarker(ID, marker)[source]

Sets the marker of the curve with the ID

setPointMarker(ID, marker)[source]

Sets the marker of the point with the ID

setSurfaceMarker(ID, marker)[source]

Sets the marker of the surface with the ID

setVolumeMarker(ID, marker)[source]

Sets the marker of the volume with the ID

spline(points, ID=None, marker=0, el_on_curve=None, el_distrib_type=None, el_distrib_val=None)

Adds a Spline curve

points - List of indices of control points that make a Spline

[p1, p2, … , pn]

ID - Positive integer ID of this curve. If left unspecified the

curve will be assigned the smallest unused curve-ID. It is recommended to specify all curve-IDs or none.

marker - Integer. Marker applied to this curve. Default 0.

elOnCurv - Positive integer. Elements on curve.

The number of element edges that will be distributed along this curve. Only works for structured meshes.

el_distrib_type -

String. Either “bump” or “progression”. Determines how the density of elements vary along the curve for structured meshes. Only works for structured meshes. elOnCurv and el_distrib_val must be be defined if this param is used.

el_distrib_val -

Float. Determines how severe the element distribution is. Only works for structured meshes. elOnCurv and el_distrib_type must be be defined if this param is used.

bump:

Smaller value means elements are bunched up at the edges of the curve, larger means bunched in the middle.

progression:

The edge of each element along this curve (from starting point to end) will be larger than the preceding one by this factor. el_distrib_val = 2 meaning for example that each line element in the series will be twice as long as the preceding one. el_distrib_val < 1 makes each element smaller than the preceeding one.

struct_surf(outer_loop, ID=None, marker=0)

Adds a Structured Surface. Parameters: outer_loop - List of 4 curve IDs that make up the boundary of

the surface. The curves must be structured, i.e. their parameter ‘elOnCurv’ must be defined.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

structuredSurface(outer_loop, ID=None, marker=0)

Adds a Structured Surface. Parameters: outer_loop - List of 4 curve IDs that make up the boundary of

the surface. The curves must be structured, i.e. their parameter ‘elOnCurv’ must be defined.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

structuredVolume(outer_surfaces, ID=None, marker=0)

Adds a Structured Volume Parameters: outer_surfaces - List of surface IDs that make up the outer boundary of

the volume. The surfaces must be Structured Surfaces.

ID - Positive integer ID of this volume. If left unspecified

the volume will be assigned the smallest unused volume-ID. It is recommended to specify all volume-IDs or none.

marker - Integer. Marker applied to this volume. Default 0.

structured_surface(outer_loop, ID=None, marker=0)

Adds a Structured Surface. Parameters: outer_loop - List of 4 curve IDs that make up the boundary of

the surface. The curves must be structured, i.e. their parameter ‘elOnCurv’ must be defined.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

structured_volume(outer_surfaces, ID=None, marker=0)

Adds a Structured Volume Parameters: outer_surfaces - List of surface IDs that make up the outer boundary of

the volume. The surfaces must be Structured Surfaces.

ID - Positive integer ID of this volume. If left unspecified

the volume will be assigned the smallest unused volume-ID. It is recommended to specify all volume-IDs or none.

marker - Integer. Marker applied to this volume. Default 0.

stuffOnSurfaces(IDs)[source]

Returns lists of all geometric points and curves on the surfaces specified in IDs. IDs may be an integer or a list of integers

stuffOnVolumes(IDs)[source]

Returns lists of all geometric points, curves, and surfaces on the volumes specified in IDs. IDs may be an integer or a list of integers

surf(outer_loop, holes=[], ID=None, marker=0)

Adds a plane surface (flat). Parameters: outer_loop - List of curve IDs that make up the outer boundary of

the surface. The curves must lie in the same plane.

holes - List of lists of curve IDs that make up the inner

boundaries of the surface. The curves must lie in the same plane.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

surface(outer_loop, holes=[], ID=None, marker=0)

Adds a plane surface (flat). Parameters: outer_loop - List of curve IDs that make up the outer boundary of

the surface. The curves must lie in the same plane.

holes - List of lists of curve IDs that make up the inner

boundaries of the surface. The curves must lie in the same plane.

ID - Positive integer ID of this surface. If left unspecified

the surface will be assigned the smallest unused surface-ID. It is recommended to specify all surface-IDs or none.

marker - Integer. Marker applied to this surface. Default 0.

volume(outer_surfaces, holes=[], ID=None, marker=0)

Adds a Volume Parameters: outer_surfaces - List of surface IDs that make up the outer boundary of

the volume.

holes - List of lists of surface IDs that make up the inner

boundaries of the volume.

ID - Positive integer ID of this volume. If left unspecified

the volume will be assigned the smallest unused volume-ID. It is recommended to specify all volume-IDs or none.

marker - Integer. Marker applied to this volume. Default 0.

Mesh functions

User interface functions

CALFEM UI module

Routines for interfacing wirt user interface toolkits.

calfem.ui.appInstance(useVisVis=True)[source]

Create a suitable application instance

calfem.ui.app_instance(useVisVis=True)

Create a suitable application instance

calfem.ui.loadUiWidget(uifilename, parent=None)[source]

Load user interface file and return object model

calfem.ui.load_ui_widget(uifilename, parent=None)

Load user interface file and return object model

Utility functions

calfem.utils.applyTractionLinearElement(boundaryElements, coords, dofs, F, marker, q)

Apply traction on part of boundarty with marker. q is added to all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node. The function works with linear line elements. (elm-type 1 in GMSH).

Parameters:

boundaryElements Dictionary with boundary elements, the key is a marker and the values are lists of elements. coords Coordinates matrix dofs Dofs matrix F force matrix. marker Boundary marker to assign boundary condition. q Value to assign boundary condition.

shape = [qx qy] in global coordinates

calfem.utils.apply_bc(boundaryDofs, bcPrescr, bcVal, marker, value=0.0, dimension=0)[source]

Apply boundary condition to bcPresc and bcVal matrices. For 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. bcPresc 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply bc. 0 - all, 1 - x, 2 - y

Returns:

bcPresc Updated 1-dim integer array containing prescribed dofs. bcVal Updated 1-dim float array containing prescribed values.

calfem.utils.apply_bc_3d(boundaryDofs, bcPrescr, bcVal, marker, value=0.0, dimension=0)[source]

Apply boundary condition to bcPresc and bcVal matrices. For 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. bcPresc 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply bc. 0 - all, 1 - x, 2 - y,

3 - z

Returns:

bcPresc Updated 1-dim integer array containing prescribed dofs. bcVal Updated 1-dim float array containing prescribed values.

calfem.utils.apply_force(boundaryDofs, f, marker, value=0.0, dimension=0)[source]

Apply boundary force to f matrix. The value is added to all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y

calfem.utils.apply_force_3d(boundaryDofs, f, marker, value=0.0, dimension=0)[source]

Apply boundary force to f matrix. The value is added to all boundaryDofs defined by marker. Applicable to 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y,

3 - z

calfem.utils.apply_force_total(boundaryDofs, f, marker, value=0.0, dimension=0)[source]

Apply boundary force to f matrix. Total force, value, is distributed over all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Total force value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y

calfem.utils.apply_force_total_3d(boundaryDofs, f, marker, value=0.0, dimension=0)[source]

Apply boundary force to f matrix. Total force, value, is distributed over all boundaryDofs defined by marker. Applicable to 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Total force value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y,

3 - z

calfem.utils.apply_traction_linear_element(boundaryElements, coords, dofs, F, marker, q)[source]

Apply traction on part of boundarty with marker. q is added to all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node. The function works with linear line elements. (elm-type 1 in GMSH).

Parameters:

boundaryElements Dictionary with boundary elements, the key is a marker and the values are lists of elements. coords Coordinates matrix dofs Dofs matrix F force matrix. marker Boundary marker to assign boundary condition. q Value to assign boundary condition.

shape = [qx qy] in global coordinates

calfem.utils.applybc(boundaryDofs, bcPrescr, bcVal, marker, value=0.0, dimension=0)

Apply boundary condition to bcPresc and bcVal matrices. For 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. bcPresc 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply bc. 0 - all, 1 - x, 2 - y

Returns:

bcPresc Updated 1-dim integer array containing prescribed dofs. bcVal Updated 1-dim float array containing prescribed values.

calfem.utils.applybc3D(boundaryDofs, bcPrescr, bcVal, marker, value=0.0, dimension=0)

Apply boundary condition to bcPresc and bcVal matrices. For 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. bcPresc 1-dim integer array containing prescribed dofs. bcVal 1-dim float array containing prescribed values. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply bc. 0 - all, 1 - x, 2 - y,

3 - z

Returns:

bcPresc Updated 1-dim integer array containing prescribed dofs. bcVal Updated 1-dim float array containing prescribed values.

calfem.utils.applyforce(boundaryDofs, f, marker, value=0.0, dimension=0)

Apply boundary force to f matrix. The value is added to all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y

calfem.utils.applyforce3D(boundaryDofs, f, marker, value=0.0, dimension=0)

Apply boundary force to f matrix. The value is added to all boundaryDofs defined by marker. Applicable to 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y,

3 - z

calfem.utils.applyforcetotal(boundaryDofs, f, marker, value=0.0, dimension=0)

Apply boundary force to f matrix. Total force, value, is distributed over all boundaryDofs defined by marker. Applicable to 2D problems with 2 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Total force value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y

calfem.utils.applyforcetotal3D(boundaryDofs, f, marker, value=0.0, dimension=0)

Apply boundary force to f matrix. Total force, value, is distributed over all boundaryDofs defined by marker. Applicable to 3D problems with 3 dofs per node.

Parameters:

boundaryDofs Dictionary with boundary dofs. f force matrix. marker Boundary marker to assign boundary condition. value Total force value to assign boundary condition.

If not given 0.0 is assigned.

dimension dimension to apply force. 0 - all, 1 - x, 2 - y,

3 - z

calfem.utils.calc_bar_displ_limits(a, coords, edof, dofs)[source]

Calc max and min global displacements for bars.

:param array a global displacement array with 3 dofs / node. :param array node coordinates :param array edof beam topology :param array dofs node dofs

calfem.utils.calc_beam_displ_limits(a, coords, edof, dofs)[source]

Calculate max and min displacements for beams.

:param array a global displacement array with 6 dofs / node. :param array node coordinates :param array edof beam topology :param array dofs node dofs

calfem.utils.calc_limits(coords)[source]

Calculate max an min limits of 3d coordinates

calfem.utils.calc_size(coords)[source]

Calculate max and min sizes of 3d coordinates.

calfem.utils.convert_to_node_topo(edof, ex, ey, ez, n_dofs_per_node=3, ignore_first=True)[source]

Routine to convert dof based topology and element coordinates to node based topology required for visualisation with VTK and other visualisation frameworks

Parameters:
  • edof (array) – element topology [nel x (n_dofs_per_node)|(n_dofs_per_node+1)*n_nodes ]

  • ex (array) – element x coordinates [nel x n_nodes]

  • ey (array) – element y coordinates [nel x n_nodes]

  • ez (array) – element z coordinates [nel x n_nodes]

  • n_dofs_per_node (array) – number of dofs per node. (default = 3)

  • ignore_first (boolean) – ignore first column of edof. (default = True)

Return array coords:

Array of node coordinates. [n_nodes x 3]

Return array topo:

Node topology. [nel x n_nodes]

Return array node_dofs:

Dofs for each node. [n_nodes x n_dofs_per_node]

calfem.utils.disp_array(a, headers=[], fmt='.4e', tablefmt='psql', showindex=False)[source]

Print a numpy array in a nice way.

calfem.utils.export_vtk_stress(filename, coords, topo, a=None, el_scalar=None, el_vec1=None, el_vec2=None)[source]

Export mesh and results for a 2D stress problem.

Parameters:

filename Filename of vtk-file coords Element coordinates (np.array) topo Element topology (not dof topology). mesh.topo. (np.array) a Element displacements 2-dof (np.array) el_scalar Scalar values for each element (list) el_vec1 Vector value for each element (list) el_vec2 Vector value for each element (list)

calfem.utils.load_arrays(name)[source]

Load arrays from file.

calfem.utils.load_geometry(name)[source]

Loads a geometry from a file.

calfem.utils.load_mesh(name)[source]

Load a mesh from file.

calfem.utils.readFloat(f)

Read a row from file, f, and return a list of floats.

calfem.utils.readInt(f)

Read a row from file, f, and return a list of integers.

calfem.utils.readSingleFloat(f)

Read a single float from a row in file f. All other values on row are discarded.

calfem.utils.readSingleInt(f)

Read a single integer from a row in file f. All other values on row are discarded.

calfem.utils.read_float(f)[source]

Read a row from file, f, and return a list of floats.

calfem.utils.read_int(f)[source]

Read a row from file, f, and return a list of integers.

calfem.utils.read_single_float(f)[source]

Read a single float from a row in file f. All other values on row are discarded.

calfem.utils.read_single_int(f)[source]

Read a single integer from a row in file f. All other values on row are discarded.

calfem.utils.save_arrays(coords, edof, dofs, bdofs, elementmarkers, boundaryElements, markerDict, name='unnamed_arrays')[source]

Save arrays to file.

calfem.utils.save_geometry(g, name='unnamed_geometry')[source]

Save a geometry to file.

calfem.utils.save_matlab_arrays(coords, edof, dofs, bdofs, elementmarkers, boundaryElements, markerDict, name='Untitled')[source]

Save arrays as MATLAB .mat files.

calfem.utils.save_mesh(mesh, name='Untitled')[source]

Save a mesh to file.

calfem.utils.scalfact2(ex, ey, ed, rat=0.2)[source]

Determine scale factor for drawing computational results, such as displacements, section forces or flux.

Parameters:

ex, ey element node coordinates

ed element displacement matrix or section force matrix

rat relation between illustrated quantity and element size.

If not specified, 0.2 is used.

calfem.utils.which(filename)[source]

Return complete path to executable given by filename.

Visualisation functions (Matplotlib)

CALFEM Visualisation module (matplotlib)

Contains all the functions implementing visualisation routines.

calfem.vis_mpl.axis(*args, **kwargs)[source]

Define axis of figure (Matplotlib passthrough)

calfem.vis_mpl.camera3d()[source]

Get visvis 3D camera.

calfem.vis_mpl.ce2vf(coords, edof, dofs_per_node, el_type)[source]

Duplicate code. Extracts verts, faces and verticesPerFace from input.

calfem.vis_mpl.clf()[source]

Clear visvis figure

calfem.vis_mpl.close(fig=None)[source]

Close visvis figure

calfem.vis_mpl.closeAll()

Close all visvis windows.

calfem.vis_mpl.close_all()[source]

Close all visvis windows.

calfem.vis_mpl.colorbar(**kwargs)[source]

Add a colorbar to current figure

calfem.vis_mpl.create_ordered_polys(geom, N=10)[source]

Creates ordered polygons from the geometry definition

calfem.vis_mpl.dispbeam2(ex, ey, edi, plotpar=[2, 1, 1], sfac=None)[source]

dispbeam2(ex,ey,edi,plotpar,sfac) [sfac]=dispbeam2(ex,ey,edi) [sfac]=dispbeam2(ex,ey,edi,plotpar)

INPUT: ex = [ x1 x2 ]

ey = [ y1 y2 ] element node coordinates.

edi = [ u1 v1;
u2 v2;
…..] matrix containing the displacements

in Nbr evaluation points along the beam.

plotpar=[linetype, linecolour, nodemark]

linetype=1 -> solid linecolour=1 -> black

2 -> dashed 2 -> blue 3 -> dotted 3 -> magenta

4 -> red

nodemark=0 -> no mark

1 -> circle 2 -> star 3 -> point

sfac = [scalar] scale factor for displacements.

Rem. Default if sfac and plotpar is left out is auto magnification

and dashed black lines with circles at nodes -> plotpar=[1 1 1]


LAST MODIFIED: O Dahlblom 2015-11-18

O Dahlblom 2023-01-31 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.drawMesh(coords, edof, dofs_per_node, el_type, title=None, color=(0, 0, 0), face_color=(0.8, 0.8, 0.8), node_color=(0, 0, 0), filled=False, show_nodes=False)

Draws wire mesh of model in 2D or 3D. Returns the Mesh object that represents the mesh. Args:

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_nodes:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

axes:

Matplotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

axes_adjust:

Boolean. True if the view should be changed to show the whole model. Default True.

title:

String. Changes title of the figure. Default “Mesh”.

color:

3-tuple or char. Color of the wire. Defaults to black (0,0,0). Can also be given as a character in ‘rgbycmkw’.

face_color:

3-tuple or char. Color of the faces. Defaults to white (1,1,1). Parameter filled must be True or faces will not be drawn at all.

filled:

Boolean. Faces will be drawn if True. Otherwise only the wire is drawn. Default False.

calfem.vis_mpl.draw_displacements(a, coords, edof, dofs_per_node, el_type, draw_undisplaced_mesh=False, magnfac=-1.0, magscale=0.25, title=None, color=(0, 0, 0), node_color=(0, 0, 0))[source]

Draws scalar element values in 2D or 3D. Returns the world object elementsWobject that represents the mesh.

Args:
ev:

An N-by-1 array or a list of scalars. The Scalar values of the elements. ev[i] should be the value of element i.

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_node:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

displacements:

An N-by-2 or N-by-3 array. Row i contains the x,y,z displacements of node i.

axes:

Matlotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

draw_undisplaced_mesh:

Boolean. True if the wire of the undisplaced mesh should be drawn on top of the displaced mesh. Default False. Use only if displacements != None.

magnfac:

Float. Magnification factor. Displacements are multiplied by this value. Use this to make small displacements more visible.

title:

String. Changes title of the figure. Default “Element Values”.

calfem.vis_mpl.draw_element_values(values, coords, edof, dofs_per_node, el_type, displacements=None, draw_elements=True, draw_undisplaced_mesh=False, magnfac=1.0, title=None, color=(0, 0, 0), node_color=(0, 0, 0))[source]

Draws scalar element values in 2D or 3D.

Args:
values:

An N-by-1 array or a list of scalars. The Scalar values of the elements. ev[i] should be the value of element i.

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_node:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

displacements:

An N-by-2 or N-by-3 array. Row i contains the x,y,z displacements of node i.

draw_elements:

Boolean. True if mesh wire should be drawn. Default True.

draw_undisplaced_mesh:

Boolean. True if the wire of the undisplaced mesh should be drawn on top of the displaced mesh. Default False. Use only if displacements != None.

magnfac:

Float. Magnification factor. Displacements are multiplied by this value. Use this to make small displacements more visible.

title:

String. Changes title of the figure. Default “Element Values”.

calfem.vis_mpl.draw_elements(ex, ey, title='', color=(0, 0, 0), face_color=(0.8, 0.8, 0.8), node_color=(0, 0, 0), line_style='solid', filled=False, closed=True, show_nodes=False)[source]

Draws wire mesh of model in 2D or 3D. Returns the Mesh object that represents the mesh. Args:

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_nodes:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

axes:

Matplotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

axes_adjust:

Boolean. True if the view should be changed to show the whole model. Default True.

title:

String. Changes title of the figure. Default “Mesh”.

color:

3-tuple or char. Color of the wire. Defaults to black (0,0,0). Can also be given as a character in ‘rgbycmkw’.

face_color:

3-tuple or char. Color of the faces. Defaults to white (1,1,1). Parameter filled must be True or faces will not be drawn at all.

filled:

Boolean. Faces will be drawn if True. Otherwise only the wire is drawn. Default False.

calfem.vis_mpl.draw_geometry(geometry, draw_points=True, label_points=True, label_curves=True, title=None, font_size=11, N=20, rel_margin=0.05, draw_axis=False, axes=None)[source]

Draws the geometry (points and curves) in geoData Args:

geoData:

GeoData object. Geodata contains geometric information of the model.

axes:

Matplotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

axes_adjust:

Boolean. If True the view will be changed to show the whole model. Default True.

draw_points:

Boolean. If True points will be drawn.

label_points:

Boolean. If True Points will be labeled. The format is: ID[marker]. If a point has marker==0 only the ID is written.

label_curves:

Boolean. If True Curves will be labeled. The format is: ID(elementsOnCurve)[marker].

font_size:

Integer. Size of the text in the text labels. Default 11.

N:

Integer. The number of discrete points per curve segment. Default 20. Increase for smoother curves. Decrease for better performance.

rel_margin:

Extra spacing between geometry and axis

calfem.vis_mpl.draw_mesh(coords, edof, dofs_per_node, el_type, title=None, color=(0, 0, 0), face_color=(0.8, 0.8, 0.8), node_color=(0, 0, 0), filled=False, show_nodes=False)[source]

Draws wire mesh of model in 2D or 3D. Returns the Mesh object that represents the mesh. Args:

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_nodes:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

axes:

Matplotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

axes_adjust:

Boolean. True if the view should be changed to show the whole model. Default True.

title:

String. Changes title of the figure. Default “Mesh”.

color:

3-tuple or char. Color of the wire. Defaults to black (0,0,0). Can also be given as a character in ‘rgbycmkw’.

face_color:

3-tuple or char. Color of the faces. Defaults to white (1,1,1). Parameter filled must be True or faces will not be drawn at all.

filled:

Boolean. Faces will be drawn if True. Otherwise only the wire is drawn. Default False.

calfem.vis_mpl.draw_nodal_values(values, coords, edof, levels=12, title=None, dofs_per_node=None, el_type=None, draw_elements=False)

Draws element nodal values as filled contours. Element topologies supported are triangles, 4-node quads and 8-node quads.

calfem.vis_mpl.draw_nodal_values_contour(values, coords, edof, levels=12, title=None, dofs_per_node=None, el_type=None, draw_elements=False)[source]

Draws element nodal values as filled contours. Element topologies supported are triangles, 4-node quads and 8-node quads.

calfem.vis_mpl.draw_nodal_values_contourf(values, coords, edof, levels=12, title=None, dofs_per_node=None, el_type=None, draw_elements=False)[source]

Draws element nodal values as filled contours. Element topologies supported are triangles, 4-node quads and 8-node quads.

calfem.vis_mpl.draw_nodal_values_shaded(values, coords, edof, title=None, dofs_per_node=None, el_type=None, draw_elements=False)[source]

Draws element nodal values as shaded triangles. Element topologies supported are triangles, 4-node quads and 8-node quads.

calfem.vis_mpl.draw_node_circles(ex, ey, title='', color=(0, 0, 0), face_color=(0.8, 0.8, 0.8), filled=False, marker_type='o')[source]

Draws wire mesh of model in 2D or 3D. Returns the Mesh object that represents the mesh. Args:

coords:

An N-by-2 or N-by-3 array. Row i contains the x,y,z coordinates of node i.

edof:

An E-by-L array. Element topology. (E is the number of elements and L is the number of dofs per element)

dofs_per_nodes:

Integer. Dofs per node.

el_type:

Integer. Element Type. See Gmsh manual for details. Usually 2 for triangles or 3 for quadrangles.

axes:

Matplotlib Axes. The Axes where the model will be drawn. If unspecified the current Axes will be used, or a new Axes will be created if none exist.

axes_adjust:

Boolean. True if the view should be changed to show the whole model. Default True.

title:

String. Changes title of the figure. Default “Mesh”.

color:

3-tuple or char. Color of the wire. Defaults to black (0,0,0). Can also be given as a character in ‘rgbycmkw’.

face_color:

3-tuple or char. Color of the faces. Defaults to white (1,1,1). Parameter filled must be True or faces will not be drawn at all.

filled:

Boolean. Faces will be drawn if True. Otherwise only the wire is drawn. Default False.

calfem.vis_mpl.eldisp2(ex, ey, ed, plotpar, sfac)[source]

[sfac]=eldisp2(ex,ey,ed,plotpar) [sfac]=eldisp2(ex,ey,ed) ————————————————————-

PURPOSE

Draw the deformed 2D mesh for a number of elements of the same type. Supported elements are:

1) -> bar element 2) -> beam el. 3) -> triangular 3 node el. 4) -> quadrilateral 4 node el. 5) -> 8-node isopar. element

INPUT
ex,ey:………. nen: number of element nodes

nel: number of elements

ed: element displacement matrix

plotpar=[ linetype, linecolor, nodemark]

linetype=1 -> solid linecolor=1 -> black

2 -> dashed 2 -> blue 3 -> dotted 3 -> magenta

4 -> red

nodemark=1 -> circle

2 -> star 0 -> no mark

sfac: scale factor for displacements

Rem. Default if sfac and plotpar is left out is auto magnification

and dashed black lines with circles at nodes -> plotpar=[2 1 1]


LAST MODIFIED: O Dahlblom 2004-10-01

J Lindemann 2021-12-30 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.eldraw2(ex, ey, plotpar, elnum)[source]
calfem.vis_mpl.eldraw2(ex, ey, plotpar) None
calfem.vis_mpl.eldraw2(ex, ey) None
PURPOSE

Draw the undeformed 2D mesh for a number of elements of the same type. Supported elements are:

1) -> bar element 2) -> beam el. 3) -> triangular 3 node el. 4) -> quadrilateral 4 node el. 5) -> 8-node isopar. elemen

INPUT
ex,ey:………. nen: number of element nodes

nel: number of elements

plotpar=[ linetype, linecolor, nodemark]

linetype=1 -> solid linecolor=1 -> black

2 -> dashed 2 -> blue 3 -> dotted 3 -> magenta

4 -> red

nodemark=1 -> circle

2 -> star 0 -> no mark

elnum=edof(:,1) ; i.e. the first column in the topology matrix

Rem. Default is solid white lines with circles at nodes.

calfem.vis_mpl.error(msg)[source]

Log error message

calfem.vis_mpl.figure(figure=None, show=True, fig_size=(6, 5.33))[source]

Create a visvis figure with extras.

calfem.vis_mpl.figureClass()

Return visvis Figure class.

calfem.vis_mpl.figure_class()[source]

Return visvis Figure class.

calfem.vis_mpl.gca()[source]

Get current axis of the current visvis figure.

calfem.vis_mpl.info(msg)[source]

Log information message

calfem.vis_mpl.pltstyle(plotpar)[source]
INPUT

plotpar=[ linetype, linecolor, nodemark ]

linetype=1 -> solid linecolor=1 -> black

2 -> dashed 2 -> blue 3 -> dotted 3 -> magenta

4 -> red

nodemark=1 -> circle

2 -> star 0 -> no mark

OUTPUT

s1: linetype and color for mesh lines s2: type and color for node markers


LAST MODIFIED: Ola Dahlblom 2004-09-15 Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.pltstyle2(plotpar)[source]
INPUT

plotpar=[ linetype, linecolor, nodemark ]

linetype=1 -> solid linecolor=1 -> black

2 -> dashed 2 -> blue 3 -> dotted 3 -> magenta

4 -> red

nodemark=1 -> circle

2 -> star 0 -> no mark

OUTPUT

s1: linetype and color for mesh lines s2: type and color for node markers


LAST MODIFIED: Ola Dahlblom 2004-09-15 Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.scalfact2(ex, ey, ed, rat=0.2)[source]

[sfac]=scalfact2(ex,ey,ed,rat) [sfac]=scalfact2(ex,ey,ed) ————————————————————- PURPOSE Determine scale factor for drawing computational results, such as displacements, section forces or flux.

INPUT

ex,ey: element node coordinates

ed: element displacement matrix or section force matrix

rat: relation between illustrated quantity and element size. If not specified, 0.2 is used.


LAST MODIFIED: O Dahlblom 2004-09-15

J Lindemann 2021-12-29 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.scalgraph2(sfac, magnitude, plotpar)[source]
calfem.vis_mpl.scalgraph2(sfac, magnitude) None

INPUT: sfac = [scalar] scale factor.

magnitude = [Ref x y] The graphic scale has a length equivalent

to Ref and starts at coordinates (x,y). If no coordinates are given the starting

point will be (0,-0.5).

plotpar=[linecolor]

linecolor=1 -> black 2 -> blue 3 -> magenta 4 -> red


LAST MODIFIED: O Dahlblom 2015-12-02

O Dahlblom 2023-01-23 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.secforce2(ex, ey, es, plotpar, sfac)[source]
calfem.vis_mpl.secforce2(ex, ey, es, plotpar, sfac, eci) None

[sfac]=secforce2(ex,ey,es) [sfac]=secforce2(ex,ey,es,plotpar) ————————————————————————– PURPOSE: Draw section force diagram for a two dimensional bar or beam element.

INPUT: ex = [ x1 x2 ]

ey = [ y1 y2 ] element node coordinates.

es = [ S1;
S2;
… ] vector containing the section force

in Nbr evaluation points along the element.

plotpar=[linecolour, elementcolour]

linecolour=1 -> black elementcolour=1 -> black

2 -> blue 2 -> blue 3 -> magenta 3 -> magenta 4 -> red 4 -> red

sfac = [scalar] scale factor for section force diagrams.

eci = [ x1;

x2;

… ] local x-coordinates of the evaluation points (Nbr).

If not given, the evaluation points are assumed to be uniformly distributed


LAST MODIFIED: O Dahlblom 2019-12-16

O Dahlblom 2023-01-31 (Python)

Copyright (c) Division of Structural Mechanics and

Division of Solid Mechanics. Lund University


calfem.vis_mpl.show()[source]

Use in Qt applications

calfem.vis_mpl.showAndWait()

Wait for plot to show

calfem.vis_mpl.showAndWaitMpl()

Wait for plot to show

calfem.vis_mpl.show_and_wait()[source]

Wait for plot to show

calfem.vis_mpl.show_and_wait_mpl()[source]

Wait for plot to show

calfem.vis_mpl.subplot(*args)[source]

Create a visvis subplot.

calfem.vis_mpl.title(*args, **kwargs)[source]

Define title of figure (Matplotlib passthrough)

calfem.vis_mpl.topo_to_tri(edof)[source]

Converts 2d element topology to triangle topology to be used with the matplotlib functions tricontour and tripcolor.

Visualisation functions (VisVis)

Interactive Geometry Editor