ObjectFFRF
This object is used to represent equations modelled by the FFRF. It contains a RigidBodyNode (always node 0) and a list of other nodes representing the finite element nodes used in the FFRF. Note that temporary matrices and vectors are subject of change in future. Usually you SHOULD NOT USE THIS OBJECT - use the much more efficient ObjectFFRFreducedOrder object with modal reduction instead.
Authors: Gerstmayr Johannes, Zwölfer Andreas
Additional information for ObjectFFRF:
- This
Object
has/provides the following types =Body
,MultiNoded
,SuperElement
- Requested
Node
type: read detailed information of item
The item ObjectFFRF with type = ‘FFRF’ has the following parameters:
- name [type = String, default = ‘’]:objects’s unique name
- nodeNumbers [\(\mathbf{n}\indf = [n_0,\,\ldots,\,n_{n_\mathrm{nf}}]\tp\), type = ArrayNodeIndex, default = []]:node numbers which provide the coordinates for the object (consecutively as provided in this list); the \((n_\mathrm{nf}+1)\) nodes represent the nodes of the FE mesh (except for node 0); the global nodal position needs to be reconstructed from the rigid-body motion of the reference frame
- massMatrixFF [\(\LU{b}{{\mathbf{M}}} \in \Rcal^{n\indf \times n\indf}\), type = PyMatrixContainer, default = PyMatrixContainer[]]:body-fixed and ONLY flexible coordinates part of mass matrix of object given in Python numpy format (sparse (CSR) or dense, converted to sparse matrix); internally data is stored in triplet format
- stiffnessMatrixFF [\(\LU{b}{{\mathbf{K}}} \in \Rcal^{n\indf \times n\indf}\), type = PyMatrixContainer, default = PyMatrixContainer[]]:body-fixed and ONLY flexible coordinates part of stiffness matrix of object in Python numpy format (sparse (CSR) or dense, converted to sparse matrix); internally data is stored in triplet format
- dampingMatrixFF [\(\LU{b}{{\mathbf{D}}} \in \Rcal^{n\indf \times n\indf}\), type = PyMatrixContainer, default = PyMatrixContainer[]]:body-fixed and ONLY flexible coordinates part of damping matrix of object in Python numpy format (sparse (CSR) or dense, converted to sparse matrix); internally data is stored in triplet format
- forceVector [\(\LU{0}{{\mathbf{f}}} = [\LU{0}{{\mathbf{f}}\indr},\; \LU{0}{{\mathbf{f}}\indf}]\tp \in \Rcal^{n_c}\), type = NumpyVector, default = []]:generalized, force vector added to RHS; the rigid body part \({\mathbf{f}}_r\) is directly applied to rigid body coordinates while the flexible part \({\mathbf{f}}\indf\) is transformed from global to local coordinates; note that this force vector only allows to add gravity forces for bodies with COM at the origin of the reference frame
- forceUserFunction [\({\mathbf{f}}_{user} = [\LU{0}{{\mathbf{f}}_{\mathrm{r},user}},\; \LU{b}{{\mathbf{f}}_{\mathrm{f},user}}]\tp \in \Rcal^{n_c}\), type = PyFunctionVectorMbsScalarIndex2Vector, default = 0]:A Python user function which computes the generalized user force vector for the ODE2 equations; note the different coordinate systems for rigid body and flexible part; The function args are mbs, time, objectNumber, coordinates q (without reference values) and coordinate velocities q_t; see description below
- massMatrixUserFunction [\({\mathbf{M}}_{user} \in \Rcal^{n_c\times n_c}\), type = PyFunctionMatrixMbsScalarIndex2Vector, default = 0]:A Python user function which computes the TOTAL mass matrix (including reference node) and adds the local constant mass matrix; note the different coordinate systems as described in the FFRF mass matrix; see description below
- coordinateIndexPerNode [type = ArrayIndex, default = []]:this list contains the local coordinate index for every node, which is needed, e.g., for markers; the list is generated automatically every time parameters have been changed
- physicsMass [\(m\), type = UReal, default = 0.]:total mass [SI:kg] of FFRF object, auto-computed from mass matrix \(\LU{b}{{\mathbf{M}}}\)
- physicsInertia [\(J_r \in \Rcal^{3 \times 3}\), type = Matrix3D, default = [[1,0,0], [0,1,0], [0,0,1]]]:inertia tensor [SI:kgm\(^2\)] of rigid body w.r.t. to the reference point of the body, auto-computed from the mass matrix \(\LU{b}{{\mathbf{M}}}\)
- physicsCenterOfMass [\(\LU{b}{{\mathbf{b}}}_{COM}\), type = Vector3D, size = 3, default = [0.,0.,0.]]:local position of center of mass (COM); auto-computed from mass matrix \(\LU{b}{{\mathbf{M}}}\)
- PHItTM [\(\tPhi\indt\tp \in \Rcal^{n\indf \times 3}\), type = NumpyMatrix, default = Matrix[]]:projector matrix; may be removed in future
- referencePositions [\({\mathbf{x}}\cRef \in \Rcal^{n\indf}\), type = NumpyVector, default = []]:vector containing the reference positions of all flexible nodes
- tempVector [\({\mathbf{v}}_{temp} \in \Rcal^{n\indf}\), type = NumpyVector, default = []]:temporary vector
- tempCoordinates [\({\mathbf{c}}_{temp} \in \Rcal^{n\indf}\), type = NumpyVector, default = []]:temporary vector containing coordinates
- tempCoordinates_t [\(\dot {\mathbf{c}}_{temp} \in \Rcal^{n\indf}\), type = NumpyVector, default = []]:temporary vector containing velocity coordinates
- tempRefPosSkew [\(\tilde{\mathbf{p}}\indf \in \Rcal^{n\indf \times 3}\), type = NumpyMatrix, default = Matrix[]]:temporary matrix with skew symmetric local (deformed) node positions
- tempVelSkew [\(\dot{\tilde{\mathbf{c}}}\indf \in \Rcal^{n\indf \times 3}\), type = NumpyMatrix, default = Matrix[]]:temporary matrix with skew symmetric local node velocities
- visualization [type = VObjectFFRF]:parameters for visualization of item
The item VObjectFFRF has the following parameters:
- show [type = Bool, default = True]:set true, if item is shown in visualization and false if it is not shown; use visualizationSettings.bodies.deformationScaleFactor to draw scaled (local) deformations; the reference frame node is shown with additional letters RF
- color [type = Float4, size = 4, default = [-1.,-1.,-1.,-1.]]:RGBA color for object; 4th value is alpha-transparency; R=-1.f means, that default color is used
- triangleMesh [type = NumpyMatrixI, default = MatrixI[]]:a matrix, containg node number triples in every row, referring to the node numbers of the GenericODE2 object; the mesh uses the nodes to visualize the underlying object; contour plot colors are still computed in the local frame!
- showNodes [type = Bool, default = False]:set true, nodes are drawn uniquely via the mesh, eventually using the floating reference frame, even in the visualization of the node is show=False; node numbers are shown with indicator ‘NF’
DESCRIPTION of ObjectFFRF
The following output variables are available as OutputVariableType in sensors, Get…Output() and other functions:
Additional output variables for superelement node access
Functions like GetObjectOutputSuperElement(...)
, see Section MainSystem: Object,
or SensorSuperElement
, see Section MainSystem: Sensor, directly access special output variables
(OutputVariableType
) of the mesh nodes \(n_i\) of the superelement.
Additionally, the contour drawing of the object can make use the OutputVariableType
of the meshnodes.
Super element output variables
super element output variables
|
symbol
|
description
|
---|---|---|
Position
|
\(\LU{0}{{\mathbf{p}}}\cConfig(n_i) = \LU{0}{\pRef\cConfig} + \LU{0b}{\Rot}\cConfig \LU{b}{{\mathbf{p}}}\cConfig(n_i)\)
|
global position of mesh node \(n_i\) including rigid body motion and flexible deformation
|
Displacement
|
\(\LU{0}{{\mathbf{c}}}\cConfig(n_i) = \LU{0}{{\mathbf{p}}\cConfig(n_i)} - \LU{0}{{\mathbf{p}}\cRef(n_i)}\)
|
global displacement of mesh node \(n_i\) including rigid body motion and flexible deformation
|
Velocity
|
\(\LU{0}{{\mathbf{v}}}\cConfig(n_i) = \LU{0}{\dot \pRef\cConfig} + \LU{0b}{\Rot}\cConfig (\LU{b}{\dot {\mathbf{q}}\indf}\cConfig(n_i) + \LU{b}{\tomega}\cConfig \times \LU{b}{{\mathbf{p}}}\cConfig(n_i))\)
|
global velocity of mesh node \(n_i\) including rigid body motion and flexible deformation
|
Acceleration
|
\(\begin{array}{l} \LU{0}{{\mathbf{a}}}\cConfig(n_i) = \LU{0}{\ddot \pRef\cConfig}\cConfig \\
+ \LU{0b}{\Rot}\cConfig \LU{b}{\ddot {\mathbf{q}}\indf}\cConfig(n_i) \\
+ 2\LU{0}{\tomega}\cConfig \times \LU{0b}{\Rot}\cConfig \LU{b}{\dot {\mathbf{q}}\indf}\cConfig(n_i) \\
+ \LU{0}{\talpha}\cConfig \times \LU{0}{{\mathbf{p}}}\cConfig(n_i) \\
+ \LU{0}{\tomega}\cConfig \times (\LU{0}{\tomega}\cConfig \times \LU{0}{{\mathbf{p}}}\cConfig(n_i)) \end{array}\)
|
global acceleration of mesh
node \(n_i\) including rigid body motion and flexible deformation; note that \(\LU{0}{{\mathbf{p}}}\cConfig(n_i) = \LU{0b}{\Rot} \LU{b}{{\mathbf{p}}}\cConfig(n_i)\)
|
DisplacementLocal
|
\(\LU{b}{{\mathbf{d}}}\cConfig(n_i) = \LU{b}{{\mathbf{p}}}\cConfig(n_i) - \LU{b}{{\mathbf{x}}}\cRef(n_i)\)
|
local displacement of mesh node \(n_i\), representing the flexible deformation within the body frame; note that \(\LU{0}{{\mathbf{u}}}\cConfig \neq \LU{0b}{\Rot}\LU{b}{{\mathbf{d}}}\cConfig\) !
|
VelocityLocal
|
\(\LU{b}{\dot {\mathbf{q}}\indf}\cConfig(n_i)\)
|
local velocity of mesh node \(n_i\), representing the rate of flexible deformation within the body frame
|
Definition of quantities
intermediate variables
|
symbol
|
description
|
---|---|---|
object coordinates
|
\({\mathbf{q}} = [{\mathbf{q}}\indt\tp,\;{\mathbf{q}}\indr\tp,\;{\mathbf{q}}\indf\tp]\tp\)
|
object coordinates
|
rigid body coordinates
|
\({\mathbf{q}}\indrigid = [{\mathbf{q}}\indt\tp,\;{\mathbf{q}}\indr\tp]\tp = [q_0,\,q_1,\,q_2,\,\psi_0,\,\psi_1,\,\psi_2,\,\psi_3]\tp\)
|
rigid body coordinates in case of Euler parameters
|
reference frame (rigid body) position
|
\(\LU{0}{\pRef\cConfig} = \LU{0}{{\mathbf{q}}_\mathrm{t,config}}+\LU{0}{{\mathbf{q}}_\mathrm{t,ref}}\)
|
global position of underlying rigid body node \(n_0\) which defines the reference frame origin
|
reference frame (rigid body) orientation
|
\(\LU{0b}{\Rot(\ttheta)}\cConfig\)
|
transformation matrix for transformation of local (reference frame) to global coordinates, given by underlying rigid body node \(n_0\)
|
local nodal position
|
\(\LU{b}{{\mathbf{p}}^{(i)}} = \LU{b}{{\mathbf{x}}^{(i)}}\cRef + \LU{b}{{\mathbf{q}}\indf^{(i)}}\)
|
vector of body-fixed (local) position of node \((i)\), including flexible part
|
local nodal positions
|
\(\LU{b}{{\mathbf{p}}} = \LU{b}{{\mathbf{x}}}\cRef + \LU{b}{{\mathbf{q}}\indf}\)
|
vector of all body-fixed (local) nodal positions including flexible part
|
rotation coordinates
|
\(\ttheta\cCur = [\psi_0,\,\psi_1,\,\psi_2,\,\psi_3]\tp\cRef + [\psi_0,\,\psi_1,\,\psi_2,\,\psi_3]\cCur\tp\)
|
rigid body coordinates in case of Euler parameters
|
flexible coordinates
|
\(\LU{b}{{\mathbf{q}}\indf}\)
|
flexible, body-fixed coordinates
|
transformation of flexible coordinates
|
\(\LU{0b}{{\mathbf{A}}_{bd}} = \mathrm{diag}([\LU{0b}{{\mathbf{A}}},\;\ldots,\;\LU{0b}{{\mathbf{A}}})\)
|
block diagonal transformation matrix, which transforms all flexible coordinates from local to global coordinates
|
The derivations follow Zwölfer and Gerstmayr with only small modifications in the notation.
Nodal coordinates
Consider an object with \(n = 1 + n_\mathrm{nf}\) nodes, \(n_\mathrm{nf}\) being the number of ‘flexible’ nodes and one additional node is the rigid body node for the reference frame. The list if node numbers is \([n_0,\,\ldots,\,n_{n_\mathrm{nf}}]\) and the according numbers of nodal coordinates are \([n_{c_0},\,\ldots,\,n_{c_n}]\), where \(n_0\) denotes the rigid body node. This gives \(n_c\) total nodal coordinates,
whereof the number of flexible coordinates is
The total number of equations (=coordinates) of the object is \(n_c\). The first node \(n_0\) represents the rigid body motion of the underlying reference frame with \(n_{c\indr} = n_{c_0}\) coordinates (e.g., \(n_{c\indr}=6\) coordinates for Euler angles and \(n_{c\indr}=7\) coordinates in case of Euler parameters; currently only the Euler parameter case is implemented.).
Kinematics
We assume a finite element mesh with The kinematics of the FFRF is based on a splitting of translational (\({\mathbf{c}}_t \in \Rcal^{n\indf}\)), rotational (\({\mathbf{c}}\indr \in \Rcal^{n\indf}\)) and flexible (\({\mathbf{c}}\indf \in \Rcal^{n\indf}\)) nodal displacements,
which are written in global coordinates in Eq. (58) but will be transformed to other coordinates later on.
In the present formulation of ObjectFFRF
, we use the following set of object coordinates (unknowns)
with \(\LU{0}{{\mathbf{q}}}\indt \in \Rcal^{3}\), \(\ttheta \in \Rcal^{4}\) and \(\LU{b}{{\mathbf{q}}}\indf \in \Rcal^{n\indf}\). Note that parts of the coordinates \({\mathbf{q}}\) can be already interpreted in specific coordinate systems, which is therefore added.
With the relations
we obtain the total relation of (global) nodal displacements to the object coordinates
On velocity level, we have
with the matrix \({\mathbf{L}} \in \Rcal^{n\indf \times n_c}\)
with the rotation parameters specific matrix \(\LU{b}{{\mathbf{G}}}\), implicitly defined in the rigid body node by the relation \(\LU{b}{\tomega} = \LU{b}{{\mathbf{G}}} \dot \ttheta\)and the body-fixed nodal position vector (for node \(i\))
and the special tilde matrix for vectors \({\mathbf{p}} \in \Rcal^{3 {n_\mathrm{nf}}}\),
with the tilde operator for a \({\mathbf{p}}^{(i)} \in \Rcal^{3}\) defined in the common notations section.
Equations of motion
We use the Lagrange equations extended for constraint \({\mathbf{g}}\),
with the quantities
Note that \(\LU{b}{{\mathbf{M}}}\) and \(\LU{b}{{\mathbf{K}}}\) are the conventional finite element mass an stiffness matrices defined in the body frame.
Elementary differentiation rules of the Lagrange equations lead to
with \({\mathbf{M}} = \LU{b}{{\mathbf{M}}}\) and \(\hat {\mathbf{K}}\) becoming obvious in Eq. (62). Note that Eq. (61) is given in global coordinates for the translational part, in terms of rotation parameters for the rotation part and in body-fixed coordinates for the flexible part of the equations.
In case that computeFFRFterms = True
, the equations eq-objectffrf-leq can be transformed into the equations of motion,
in which iN
represents the itemNumber (=objectNumber of ObjectFFRF in mbs) in the user function.
The mass terms are given as
In case that computeFFRFterms = False
, the mass terms \({\mathbf{M}}\indtt, {\mathbf{M}}\indtr, {\mathbf{M}}\indtf, {\mathbf{M}}\indrr,
{\mathbf{M}}\indrf, \LU{b}{{\mathbf{M}}}\) in Eq. (62) are set to zero (and not computed) and
the quadratic velocity vector \({\mathbf{f}}_{v} = \Null\).
Note that the user functions \({\mathbf{f}}_{user}(mbs, t, i_N, {\mathbf{q}},\dot {\mathbf{q}})\) and \({\mathbf{M}}_{user}(mbs, t, i_N, {\mathbf{q}},\dot {\mathbf{q}})\) may be empty (=0).
The detailed equations of motion for this element can be found in .
The quadratic velocity vector follows as
with the special matrix
CoordinateLoads are added for each ODE2 coordinate on the RHS of the latter equation.
If the rigid body node is using Euler parameters \(\ttheta = [\theta_0,\,\theta_1,\,\theta_2,\,\theta_3]\tp\), an additional constraint (constraint nr.0) is added automatically for the Euler parameter norm, reading
In order to suppress the rigid body motion of the mesh nodes, you should apply a ObjectConnectorCoordinateVector object with the following constraint equations which impose constraints of a so-called Tisserand frame, giving 3 constraints for the position of the center of mass
and 3 constraints for the rotation,
Userfunction: forceUserFunction(mbs, t, itemNumber, q, q_t)
A user function, which computes a force vector depending on current time and states of object. Can be used to create any kind of mechanical system by using the object states.
arguments / return
|
type or size
|
description
|
---|---|---|
mbs |
MainSystem
|
provides MainSystem mbs to which object belongs
|
t |
Real
|
current time in mbs
|
itemNumber |
Index
|
integer number of the object in mbs, allowing easy access to all object data via mbs.GetObjectParameter(itemNumber, …)
|
q |
Vector \(\in \Rcal^n_c\)
|
object coordinates (nodal displacement coordinates of rigid body and mesh nodes) in current configuration, without reference values
|
q_t |
Vector \(\in \Rcal^n_c\)
|
object velocity coordinates (time derivative of
q ) in current configuration |
returnValue
|
Vector \(\in \Rcal^{n_c}\)
|
returns force vector for object
|
Userfunction: massMatrixUserFunction(mbs, t, itemNumber, q, q_t)
A user function, which computes a mass matrix depending on current time and states of object. Can be used to create any kind of mechanical system by using the object states.
arguments / return
|
type or size
|
description
|
---|---|---|
mbs |
MainSystem
|
provides MainSystem mbs to which object belongs
|
t |
Real
|
current time in mbs
|
itemNumber |
Index
|
integer number of the object in mbs, allowing easy access to all object data via mbs.GetObjectParameter(itemNumber, …)
|
q |
Vector \(\in \Rcal^n_c\)
|
object coordinates (nodal displacement coordinates of rigid body and mesh nodes) in current configuration, without reference values
|
q_t |
Vector \(\in \Rcal^n_c\)
|
object velocity coordinates (time derivative of
q ) in current configuration |
returnValue
|
NumpyMatrix \(\in \Rcal^{n_c \times n_c}\)
|
returns mass matrix for object
|
Relevant Examples and TestModels with weblink:
objectFFRFTest.py (TestModels/), objectFFRFTest2.py (TestModels/)
The web version may not be complete. For details, consider also the Exudyn PDF documentation : theDoc.pdf