6. Discrete Elements - ConElas, ConMass, Bulk, Spring
Discrete elements are a general class of point or line elements which have applicability in both 2D and 3D space and are used in practice to idealize certain mass, stiffness or heat capacitance modelling situations. The ConElas and ConMass modules are useful for representing various types of concentrated or degree of freedom springs and masses. The Bulk module is used to represent a volumetric mass or heat capacitance. The Spring module is used to model various type of specialized elastic connections such as spot welds and flexible links.
6.1. Concentrated and Degree of Freedom Springs - ConElas
The ConElas module is used to simulate abstract ground and line elastic springs and damper. Spring and damper stiffness coefficients may be entered explicitly at specified degrees of freedom at each of two points. This feature is sometimes referred to as a degree of freedom or “scalar” spring/damper. Intrinsic extensional and/or torsional spring/damper stiffnesses may be specified and the ConElas module will generate either a rank-one or a rank-two matrix with an extensional and/or torsional component introduced in the direction defined by the end point coordinates.
The methods associated with a ConElas object are the following.
Begin and end an instance of an object, generic object functions
vfe_ConElasBegin()
- create an instance of a ConElas objectvfe_ConElasEnd()
- destroy an instance of a ConElas objectvfe_ConElasError()
- return ConElas object error flag
Attributes and Parameters
vfe_ConElasDef()
- define spring type and dimensionvfe_ConElasInq()
- inquire spring type and dimensionvfe_ConElasSetDofScalar()
- set dof spring scalar propertiesvfe_ConElasSetDofVector()
- set dof spring vector propertiesvfe_ConElasSetIntrinsic()
- set intrinsic spring stiffnessesvfe_ConElasSetLocalSystem()
- set element local system directionvfe_ConElasSetPropPtr()
- set property pointer
Degree of Freedom Information
vfe_ConElasDofMap()
- query element degree of freedom mapvfe_ConElasNumDof()
- query number of element degrees of freedom
Structural analysis computations
vfe_ConElasReact()
- reaction vectorvfe_ConElasReactStiff()
- reaction vector, stiffness matrixvfe_ConElasReactStiffDamp()
- reaction vector, stiffness and damping matricesvfe_ConElasStiff()
- linear stiffness matrixvfe_ConElasStiffDamp()
- linear stiffness and damping matrices
Instance a ConElas object using vfe_ConElasBegin()
.
Once a ConElas is instanced,
define the spring/damper type and dimension using vfe_ConElasDef()
.
If the spring/damper is defined to exist in 3D space, then all 3 translations
and 3 rotations may exist in the element degree of freedom map.
If the spring/damper is defined to exist in 2D space, then
only the x and y translations and z rotation may appear in the
element degree of freedom map.
There are two basic types of spring/dampers available in the ConElas module.
The first type are termed degree-of-freedom or “dof” spring/dampers. These
may contain either one or two nodes in which each node is connected
to a single degree of freedom. A dof spring/damper may be either a scalar
or a vector spring/damper. The scalar spring is
similar to NASTRAN-like CELASx elements, while the scalar damper is similar to
NASTRAN-like CDAMPx elements. The directions of the
scalar spring degree of freedoms are assumed to be predetermined.
Use vfe_ConElasSetDofScalar()
to define scalar springs
A vector spring is similar to ABAQUS like SPRING1 and SPRING2 elements.
The directions of the spring degrees of freedom are determined by the
spring element local system which is evaluated at the spring nodes using
the prescription specified by vfe_ConElasSetLocalSystem()
.
Use vfe_ConElasSetDofVector()
to define vector spring/dampers.
ConElas supports damping properties. For reaction and damping matrix
calculations the nodal velocities may be defined using vfe_ConElasSetPropPtr()
.
If left undefined they are assumed zero.
The second basic spring/damper type is the intrinsic stiffness spring/damper
with two nodes.
An extension and/or torsional component may be specified for the intrinsic
spring. An extensional component activates all translational freedoms
at the spring nodes and a torsional component activates all rotational
freedoms. The extensional and torsional components act along
the line between the element endpoints. Therefore the element endpoints
must not be coincident. If spring stiffness must be modelled between
two coincident nodes a combination of vector dof springs may be used.
Use vfe_ConElasSetIntrinsic()
to define
the intrinsic extensional and/or torsional stiffness for an intrinsic spring.
A purely extensional intrinsic spring is similar to ABAQUS SPRINGA elements.
The degree of freedom configuration for the element may be returned
using vfe_ConElasNumDof()
and vfe_ConElasDofMap()
.
Unlike other, more general, finite elements, these elements do not require
a material function interface using the MatlFun object.
Element reactions, stiffness, and damping matrices are computed using
vfe_ConElasStiff()
, vfe_ConElasReactStiff()
, vfe_ConElasReact()
,
vfe_ConElasStiffDamp()
, and vfe_ConElasReactStiffDamp()
.
6.2. Function Descriptions
The currently available ConElas functions are described in detail in this section.
-
vfe_ConElas *vfe_ConElasBegin(void)
create an instance of a ConElas object
Create an instance of a ConElas object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a ConElas object using
void vfe_ConElasEnd (vfe_ConElas *conelas)
Return the current value of a ConElas object error flag using
Vint vfe_ConElasError (vfe_ConElas *conelas)
- Returns:
The function returns a pointer to the newly created ConElas object. If the object creation fails, NULL is returned.
-
void vfe_ConElasEnd(vfe_ConElas *p)
destroy an instance of a ConElas object
-
Vint vfe_ConElasError(vfe_ConElas *p)
return the current value of a ConElas object error flag
-
void vfe_ConElasDef(vfe_ConElas *p, Vint type, Vint dimension)
define spring type and dimension
Define the type and dimension of the spring to be specified. The type
CONELAS_DOFSCALAR
requires a call tovfe_ConElasSetDofScalar()
,CONELAS_DOFVECTOR
requires a call tovfe_ConElasSetDofVector()
andCONELAS_INTRINSIC
requires a call tovfe_ConElasSetIntrinsic()
.Inquire of a defined type and dimension as output arguments using
void vfe_ConElasInq (vfe_ConElas *conelas, Vint *type, Vint *dime)
- Errors
SYS_ERROR_ENUM
is generated if an improper type or dimension is specified.
- Parameters:
p – Pointer to ConElas object.
type – Spring type
x=CONELAS_DOFSCALAR Degree of freedom scalar spring =CONELAS_DOFVECTOR Degree of freedom vector spring =CONELAS_INTRINSIC Intrinsic spring stiffness
dimension – Spring dimension
=VFE_2D 2D analysis =VFE_3D 3D analysis
-
void vfe_ConElasInq(vfe_ConElas *p, Vint *type, Vint *dim)
inquire of a defined type and dime as output arguments
See
vfe_ConElasDef()
-
void vfe_ConElasSetDofScalar(vfe_ConElas *p, Vint tag1, Vint tag2, Vdouble kl[3], Vdouble dl[3])
set dof spring scalar properties
Define a scalar spring/damper symmetric matrix connecting a degree of freedom at each of two points. If the dimension of the element is 2D, then tag1 and tag2 must be one of
SYS_DOF_TX
,SYS_DOF_TY
, orSYS_DOF_RZ
, otherwise any translational or rotational degree of freedom type is allowed. The element stiffness and reaction calculation are independent of node coordinates.- Errors
SYS_ERROR_ENUM
is generated if an improper tag1 or tag2 is specified.
- Parameters:
p – Pointer to ConElas object.
tag1 – First degree of freedom type
tag2 – Second degree of freedom type
kl – Lower triangle of symmetric stiffness matrix
dl – Lower triangle of damping matrix
-
void vfe_ConElasSetDofVector(vfe_ConElas *p, Vint tag1, Vint tag2, Vdouble kl[3], Vdouble dl[3])
set dof spring vector properties
Define a vector spring/damper symmetric matrix connecting a degree of freedom at each of two points in the element local system at each point. If the dimension of the element is 2D, then tag1 and tag2 must be one of
SYS_DOF_TX
,SYS_DOF_TY
, orSYS_DOF_RZ
, otherwise any translational or rotational degree of freedom type is allowed. The element local system at each point is specified byvfe_ConElasSetLocalSystem()
.- Errors
SYS_ERROR_ENUM
is generated if an improper tag1 or tag2 is specified.
- Parameters:
p – Pointer to ConElas object.
tag1 – First degree of freedom type
tag2 – Second degree of freedom type
kl – Lower triangle of element local symmetric stiffness matrix
dl – Lower triangle of element local damping matrix
-
void vfe_ConElasSetIntrinsic(vfe_ConElas *p, Vint type, Vdouble kext, Vdouble ktor, Vdouble dext, Vdouble dtor)
set intrinsic spring stiffnesses
Define the type and magnitude of intrinsic stiffness. The type
CONELAS_EXT
ignores ktor and dtor. The typeCONELAS_TOR
ignores kext and dext. If the dimension of the element is 2D, then only a type ofCONELAS_EXT
is allowed.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_OPERATION
is generated if an unallowed type is specified for 2D analysis.
- Parameters:
p – Pointer to ConElas object.
type – Spring stiffness type
x=CONELAS_EXT Extensional only =CONELAS_TOR Torsional only =CONELAS_EXTTOR Extensional and torsional
kext – Intrinsic extensional stiffness
ktor – Intrinsic torsional stiffness
dext – Intrinsic extensional damping coefficient
dtor – Intrinsic torsional damping coefficient
-
void vfe_ConElasNumDof(vfe_ConElas *p, Vint analysistype, Vint *nedofs)
query number of element degrees of freedom
Query for number of element degree of freedom nedofs. The number of degrees of freedom will generally be equal to the number of nodal degrees of freedom per node times the number of nodes. Use
vfe_ConElasDofMap()
to return the location and type of each degree of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to ConElas object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
nedofs – [out] Number of element degrees of freedom
-
void vfe_ConElasDofMap(vfe_ConElas *p, Vint analysistype, Vint loc[], Vint tag[])
query element degree of freedom map
Query for element degree of freedom map. The degree of freedom map consists of a location index, loc and type, tag for each degree of freedom used by the element.
The location index is a positive node index into the element connectivity indicating a nodal freedom. The tag indicates the type of the degree of freedom. Tag values are one of a set of enumerated types which indicate whether the degree of freedom is a translation or rotation. The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_ConElasNumDof()
to return the number of element degrees of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to ConElas object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
loc – [out] Vector of degree of freedom locations
tag – [out] Vector of degree of freedom types
-
void vfe_ConElasReact(vfe_ConElas *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
reaction vector
Compute the reaction vector, r, given the node coordinates, x, and the degree of freedom displacement vector, u. If damping is defined the nodal velocities are set with
vfe_ConElasSetPropPtr()
.- Errors
SYS_ERROR_COMPUTE
is generated if element geometry can not be defined.
- Parameters:
p – Pointer to ConElas object.
x – Array of node locations.
u – Degree of freedom vector of displacements
r – [out] Degree of freedom reaction vector
-
void vfe_ConElasStiff(vfe_ConElas *p, Vdouble x[][3], Vdouble kl[])
linear stiffness matrix
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_COMPUTE
is generated if element geometry can not be defined.
- Parameters:
p – Pointer to ConElas object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
-
void vfe_ConElasReactStiff(vfe_ConElas *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
reaction vector, stiffness matrix
Compute the reaction vector, r, and optionally the stiffness matrix, k, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the stiffness matrix is returned. If damping is defined the nodal velocities are set with
vfe_ConElasSetPropPtr()
.- Errors
SYS_ERROR_COMPUTE
is generated if element geometry can not be defined.
- Parameters:
p – Pointer to ConElas object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kflag – Flag to compute stiffness matrix, k
=SYS_OFF Do not compute stiffness matrix =SYS_ON Compute and return stiffness matrix
r – [out] Degree of freedom reaction vector
k – [out] Degree of freedom stiffness matrix
-
void vfe_ConElasStiffDamp(vfe_ConElas *p, Vdouble x[][3], Vdouble kl[], Vdouble dl[])
linear stiffness and damping matrices
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_COMPUTE
is generated if element geometry can not be defined.
- Parameters:
p – Pointer to ConElas object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
dl – [out] Degree of freedom damping matrix
-
void vfe_ConElasReactStiffDamp(vfe_ConElas *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vint dflag, Vdouble r[], Vdouble k[], Vdouble d[])
reaction vector, stiffness and damping matrices
Compute the reaction vector, r, and optionally the stiffness matrix, k, and the damping matrix, d, given the node coordinates, x, and the degree of freedom displacement vector, u. If damping is defined the nodal velocities are set with
vfe_ConElasSetPropPtr()
. The lower triangle of the stiffness matrix is returned.- Errors
SYS_ERROR_COMPUTE
is generated if element geometry can not be defined.
- Parameters:
p – Pointer to ConElas object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kflag – Flag to compute stiffness matrix, k
=SYS_OFF Do not compute stiffness matrix =SYS_ON Compute and return stiffness matrix
dflag – Flag to compute damping matrix, d
=SYS_OFF Do not compute damping matrix =SYS_ON Compute and return damping matrix
r – [out] Degree of freedom reaction vector
k – [out] Degree of freedom stiffness matrix
d – [out] Degree of freedom damping matrix
-
void vfe_ConElasSetLocalSystem(vfe_ConElas *p, Vint type, Vdouble vec[], Vdouble angle)
set element local system direction
Specify the convention to be used to construct the orientation of the element local x’,y’,z’ coordinate system. This local system is computed for vector spring elements only at each node location to determine the local directions of the node degree of freedom.
For a description of element coordinate systems, type, and associated orientation vector data, please see Element Coordinate Systems
- Errors
SYS_ERROR_ENUM
is generated if an improper type is input.
- Parameters:
p – Pointer to ConElas object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate y’,z’ axes about the element x’ axis in degrees.
-
void vfe_ConElasSetPropPtr(vfe_ConElas *p, Vint type, Vdouble *propptr)
set pointer to element nodal properties
Set a pointer to the start of a specified type of element properties. Note that the properties are not copied by this function, only the pointer itself is copied. If a property pointer is not set the element assumes a default value. By default the velocity used in damping functions is 0.0.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to ConElas object.
type – Type of element property
x=VFE_PROP_VELOCITY Velocity
propptr – Pointer to start of element nodal properties
6.3. Concentrated and Degree of Freedom Masses - ConMass
The ConMass module is used to simulate concentrated mass and degree of freedom masses. Mass matrix coefficients may be entered explicitly at specified degrees of freedom at each of two points. This feature is sometimes referred to as a degree of freedom or “scalar” mass. A variety of other concentrated mass types are available.
The methods associated with a ConMass object are the following.
Begin and end an instance of an object, generic object functions
vfe_ConMassBegin()
- create an instance of a ConMass objectvfe_ConMassEnd()
- destroy an instance of a ConMass objectvfe_ConMassError()
- return ConMass object error flag
Attributes and Parameters
vfe_ConMassDef()
- define mass type and dimensionvfe_ConMassInq()
- inquire mass type and dimensionvfe_ConMassSetParami()
- set integer parametersvfe_ConMassSetDofLink()
- set degree of freedom link propertiesvfe_ConMassSetLocalSystem()
- set local inertia tensor systemvfe_ConMassSetMatrix()
- set symmetric mass matrix at a pointvfe_ConMassSetOffset()
- set offset of center of massvfe_ConMassSetPrinc()
- set translational mass and rotary inertiavfe_ConMassSetTrans()
- set single translational mass
Degree of Freedom Information
vfe_ConMassDofMap()
- query element degree of freedom mapvfe_ConMassNumDof()
- query number of element degrees of freedom
Structural analysis computations
vfe_ConMassElemLoad()
- body force vectorvfe_ConMassMass()
- consistent mass matrixvfe_ConMassMassDiag()
- diagonal mass matrix
Instance a ConMass object using vfe_ConMassBegin()
.
Once a ConMass is instanced,
define the mass type and dimension using vfe_ConMassDef()
.
If the mass is defined to exist in 3D space, then all 3 translations
and 3 rotations may exist in the element degree of freedom map.
If the mass is defined to exist in 2D space, then
only the x and y translations and z rotation may appear in the
element degree of freedom map.
Use vfe_ConMassSetDofLink()
to define scalar masses similar to
NASTRAN like CMASSx elements.
Use vfe_ConMassSetMatrix()
to define
the components of a complete symmetric mass matrix at a point similar
to NASTRAN like CONM1 elements.
Use vfe_ConMassSetPrinc()
to define
a concentrated mass as a translational mass and a rotary inertia tensor similar
to NASTRAN like CONM2 elements.
Use vfe_ConMassSetTrans()
to define
a single concentrated translational mass at a point.
The degree of freedom configuration for the element may be returned
using vfe_ConMassNumDof()
and vfe_ConMassDofMap()
.
Unlike other, more general, finite elements, these elements do not require
a material function interface using the MatlFun object.
Element consistent and diagonal mass are computed using
vfe_ConMassMass()
and vfe_ConMassMassDiag()
.
6.4. Function Descriptions
The currently available ConMass functions are described in detail in this section.
-
vfe_ConMass *vfe_ConMassBegin(void)
create an instance of a ConMass object
Create an instance of a ConMass object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a ConMass object using
void vfe_ConMassEnd (vfe_ConMass *conmass)
Return the current value of a ConMass object error flag using
Vint vfe_ConMassError (vfe_ConMass *conmass)
- Returns:
The function returns a pointer to the newly created ConMass object. If the object creation fails, NULL is returned.
-
void vfe_ConMassEnd(vfe_ConMass *p)
destroy an instance of a ConMass object
-
Vint vfe_ConMassError(vfe_ConMass *p)
return the current value of a ConMass object error flag
-
void vfe_ConMassDef(vfe_ConMass *p, Vint type, Vint dimension)
define mass type and dimension
Define the type and dimension of the mass to be specified. The type
CONMASS_DOFLINK
requires a call tovfe_ConMassSetDofLink()
,CONMASS_TRANS
requires a call tovfe_ConMassSetTrans()
,CONMASS_PRINC
requires a call tovfe_ConMassSetPrinc()
andCONMASS_MATRIX
requires a call tovfe_ConMassSetMatrix()
, By default type isCONMASS_DOFLINK
and dime isVFE_3D
.Inquire of a defined type and dime as output arguments using
void vfe_ConMassInq (vfe_ConMass *conmass, Vint *type, Vint *dime)
- Errors
SYS_ERROR_ENUM
is generated if an improper type or dime is specified.
- Parameters:
p – Pointer to ConMass object.
type – Mass type
x=CONMASS_DOFLINK Degree of freedom link =CONMASS_TRANS Translational mass =CONMASS_PRINC Translational mass, rotary inertia =CONMASS_MATRIX Symmetric mass matrix
dimension – Mass dimension
=VFE_2D 2D analysis =VFE_3D 3D analysis
-
void vfe_ConMassInq(vfe_ConMass *p, Vint *type, Vint *dim)
inquire of a defined type and dime as output arguments
See
vfe_ConMassDef()
-
void vfe_ConMassSetDofLink(vfe_ConMass *p, Vint tag1, Vint tag2, Vdouble kl[3])
set degree of freedom link properties
Define a general mass symmetric matrix connecting a degree of freedom at each of two points. If the dimension of the element is 2D, then tag1 and tag2 must be one of
SYS_DOF_TX
,SYS_DOF_TY
, orSYS_DOF_RZ
, otherwise any translational or rotational degree of freedom type is allowed.- Errors
SYS_ERROR_ENUM
is generated if an improper tag1 or tag2 is specified.
- Parameters:
p – Pointer to ConMass object.
tag1 – First degree of freedom type
tag2 – Second degree of freedom type
kl – Lower triangle of symmetric mass matrix
-
void vfe_ConMassSetTrans(vfe_ConMass *p, Vdouble mass)
set single translational mass
Define a single translational mass. This mass value is inserted as the diagonal translational mass for all appropriate translational degrees of freedom.
- Parameters:
p – Pointer to ConMass object.
mass – Single translational mass
-
void vfe_ConMassSetPrinc(vfe_ConMass *p, Vdouble mass, Vdouble rotary[3], Vdouble tm[3][3], Vdouble offset[3])
set translational mass and rotary inertia
Define a concentrated mass matrix at a point as a single translational mass and the principal moments of rotary inertia, rotary. The orientation of the principal axes is specified with the direction cosine matrix, tm. The center of gravity of the mass may be offset from its attachment point.
- Parameters:
p – Pointer to ConMass object.
mass – Single translational mass
rotary – Principal values of rotary inertia.
tm – Orientation of principal axes of rotary inertia tensor
offset – Offset of the center of gravity of the mass.
-
void vfe_ConMassSetMatrix(vfe_ConMass *p, Vdouble a[])
set symmetric mass matrix at a point
Define a concentrated mass matrix at a point as a complete symmetric mass matrix. The mass matrix entries are dependent upon the dimension of the element. In 3D the entries of a are as follows in which the order of the degrees of freedom at the point is
SYS_DOF_TX
,SYS_DOF_TY
,SYS_DOF_TZ
,SYS_DOF_RX
,SYS_DOF_RY
,SYS_DOF_RZ
.In 2D the entries of a are as follows in which the order of the degrees of freedom at the point is- - | a[ 0] a[ 1] a[ 3] a[ 6] a[10] a[15] | | a[ 1] a[ 2] a[ 4] a[ 7] a[11] a[16] | M = | a[ 3] a[ 4] a[ 5] a[ 8] a[12] a[17] | | a[ 6] a[ 7] a[ 8] a[ 9] a[13] a[18] | | a[10] a[11] a[12] a[13] a[14] a[19] | | a[15] a[16] a[17] a[18] a[19] a[20] | - -
SYS_DOF_TX
,SYS_DOF_TY
,SYS_DOF_RZ
.- - | a[ 0] a[ 1] a[ 3] | M = | a[ 1] a[ 2] a[ 4] | | a[ 3] a[ 4] a[ 5] | - -
- Parameters:
p – Pointer to ConMass object.
a – Complete lower triangle of symmetric mass matrix at a node
-
void vfe_ConMassNumDof(vfe_ConMass *p, Vint analysistype, Vint *nedofs)
query number of element degrees of freedom
Query for number of element degree of freedom nedofs. The number of degrees of freedom will generally be equal to the number of nodal degrees of freedom per node times the number of nodes. Use
vfe_ConMassDofMap()
to return the location and type of each degree of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to ConMass object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
nedofs – [out] Number of element degrees of freedom
-
void vfe_ConMassDofMap(vfe_ConMass *p, Vint analysistype, Vint loc[], Vint tag[])
query element degree of freedom map
Query for element degree of freedom map. The degree of freedom map consists of a location index, loc and type, tag for each degree of freedom used by the element.
The location index is a positive node index into the element connectivity indicating a nodal freedom. The tag indicates the type of the degree of freedom. Tag values are one of a set of enumerated types which indicate whether the degree of freedom is a translation or rotation. The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_ConMassNumDof()
to return the number of element degrees of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to ConMass object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
loc – [out] Vector of degree of freedom locations
tag – [out] Vector of degree of freedom types
-
void vfe_ConMassSetParami(vfe_ConMass *p, Vint type, Vint iparam)
set integer parameters
Set element formulation parameters. The inclusion of rotary inertias to the element mass matrix is set with the
VFE_ROTINERTIA
parameter. By defaultVFE_ROTINERTIA
is set toSYS_OFF
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to ConMass object.
type – Type of formulation parameter to set
x=VFE_ROTINERTIA Enable rotary inertias
iparam – Integer parameter value.
x=SYS_OFF Disable =SYS_ON Enable
-
void vfe_ConMassMass(vfe_ConMass *p, Vdouble x[][3], Vdouble m[])
consistent mass matrix
Compute the consistent mass matrix, m. The lower triangle of the consistent mass is returned. Use
vfe_ConMassMassDiag()
to compute diagonal mass.- Parameters:
p – Pointer to ConMass object.
x – Array of node locations.
m – [out] Consistent mass matrix
-
void vfe_ConMassMassDiag(vfe_ConMass *p, Vdouble x[][3], Vdouble md[])
diagonal mass matrix
Compute the diagonal mass matrix, md. Use
vfe_ConMassMass()
to compute consistent mass.- Parameters:
p – Pointer to ConMass object.
x – Array of node locations.
md – [out] Diagonal mass matrix
-
void vfe_ConMassElemLoad(vfe_ConMass *p, Vdouble x[][3], Vdouble q[][3], Vdouble f[])
body force vector
Compute the consistent degree of freedom body loads given acceleration load vector, q on an element. The vector q contains an acceleration vector for for each node in the element. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The input element loads are in the units of force per unit mass.
- Parameters:
p – Pointer to ConMass object.
x – Array of node locations.
q – Array of node accelerations
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_ConMassSetLocalSystem(vfe_ConMass *p, Vint type, Vdouble vec[], Vdouble angle)
set local inertia tensor system
Specify the local coordinate system of the inertia tensor. The default is
SYS_ELEMSYS_GLOBAL
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is input.
- Parameters:
p – Pointer to ConMass object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate y’,z’ axes about the element x’ axis in degrees.
6.5. Bulk Elements - Bulk
The Bulk module is used to represent a volumetric mass or heat capacitance at a single point.
The methods associated with a Bulk object are the following.
Begin and end an instance of an object, generic object functions
vfe_BulkBegin()
- create an instance of a Bulk objectvfe_BulkEnd()
- destroy an instance of a Bulk objectvfe_BulkError()
- return Bulk object error flag
Attributes and Parameters
vfe_BulkDef()
- define bulk dimensionvfe_BulkInq()
- inquire bulk dimensionvfe_BulkSetObject()
- set attribute objectvfe_BulkSetPropPtr()
- set pointer to element nodal properties
Degree of Freedom and Integration Information
vfe_BulkDofMap()
- query element degree of freedom mapvfe_BulkNumDof()
- query number of element degrees of freedom
Structural analysis computations
vfe_BulkElemLoad()
- body force vectorvfe_BulkMass()
- consistent mass matrixvfe_BulkMassDiag()
- diagonal mass matrix
Thermal analysis computations
vfe_BulkElemHeat()
- body heat generationvfe_BulkCap()
- consistent capacitance matrixvfe_BulkCapDiag()
- diagonal capacitance matrix
Instance a Bulk object using vfe_BulkBegin()
.
Once a Bulk is instanced,
define the bulk dimension using vfe_BulkDef()
.
In structural analysis, if the bulk mass is defined in 2D space, then
only x and y translations appear in the element degree of freedom map.
If the bulk mass is defined in 3D space, then
x, y and z translations appear in the element degree of freedom map.
The bulk mass does not contain any rotary inertia.
In thermal analysis, the bulk capacitance will always have the
temperature appear in the element degree of freedom map.
The Bulk module accesses material density and specific heat through
a MatlFun interface to a primitive material. Use vfe_BulkSetObject()
to set the MatlFun object.
The bulk element temperature and volume are set using vfe_BulkSetPropPtr()
.
Query the element degree of freedom map using vfe_BulkDofMap()
and
vfe_BulkNumDof()
.
For structural analysis the bulk element supports the computation of body
forces and consistent and diagonal mass using vfe_BulkElemLoad()
,
vfe_BulkMass()
and vfe_BulkMassDiag()
.
For thermal analysis the bulk element supports the computation of body
heat generation and consistent and diagonal capacitance
using vfe_BulkElemHeat()
,
vfe_BulkCap()
and vfe_BulkCapDiag()
. Note that the coordinate location of the
bulk element is immaterial. The functions for the computation of body loads,
mass and capacitance have the coordinates in their input argument lists
only for consistency with the argument lists of the analogous functions
for all other element types such as Solid3D, etc.
6.6. Function Descriptions
The currently available Bulk functions are described in detail in this section.
-
vfe_Bulk *vfe_BulkBegin(void)
create an instance of a Bulk object
Create an instance of a MatlFun object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a MatlFun object using
void vfe_BulkEnd (vfe_Bulk *bulk)
Return the current value of a MatlFun object error flag using
Vint vfe_BulkError (vfe_Bulk *bulk)
- Returns:
The function returns a pointer to the newly created MatlFun object. If the object creation fails, NULL is returned.
-
void vfe_BulkEnd(vfe_Bulk *p)
destroy an instance of a Bulk object
See
vfe_BulkBegin()
-
Vint vfe_BulkError(vfe_Bulk *p)
return the current value of a Bulk object error flag
See
vfe_BulkBegin()
-
void vfe_BulkSetObject(vfe_Bulk *p, Vint objecttype, Vobject *object)
set attribute object
Set a pointer to an attribute object. Users must supply a MatlFun object prior to computing any quantity that requires a material model definition.
- Errors
SYS_ERROR_OBJECTTYPE
is generated if an improper objecttype is specified.
- Parameters:
p – Pointer to Bulk object.
objecttype – The object type identifier
x=VFE_MATLFUN MatlFun object
object – Pointer to the object to be set.
-
void vfe_BulkDef(vfe_Bulk *p, Vint dimension)
define bulk dimension
Define the dimension of the bulk element.
Inquire of a defined dimension as an output argument using
void vfe_BulkInq (vfe_Bulk *bulk, Vint *dime)
- Errors
SYS_ERROR_ENUM
is generated if an improper dimension is specified.
- Parameters:
p – Pointer to Bulk object.
dimension – Bulk dimension
=VFE_2D 2D analysis =VFE_3D 3D analysis
-
void vfe_BulkSetPropPtr(vfe_Bulk *p, Vint type, Vdouble *propptr)
set pointer to element nodal properties
Set a pointer to the start of a specified type of element properties. Note that the properties are not copied by this function, only the pointer itself is copied. If a property pointer is not set the element assumes a default value. By default the temperature is 0. and the volume is 1.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to Bulk object.
type – Type of element property
x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_VOLUME Volume
propptr – Pointer to start of element nodal properties
-
void vfe_BulkNumDof(vfe_Bulk *p, Vint analysistype, Vint *nedofs)
query number of element degrees of freedom
Query for number of element degrees of freedom nedofs. For structural analysis, nedofs is 2 for 2D space and 3 for 3D space. For thermal analysis, nedofs is 1. Use
vfe_BulkDofMap()
to return the location and type of each degree of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to Bulk object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
nedofs – [out] Number of element degrees of freedom
-
void vfe_BulkDofMap(vfe_Bulk *p, Vint analysistype, Vint loc[], Vint tag[])
query element degree of freedom map
Query for element degree of freedom map. The degree of freedom map consists of a location index, loc and type, tag for each degree of freedom used by the element.
The location index is a positive node index into the element connectivity indicating a nodal freedom. The tag indicates the type of the degree of freedom. Tag values are one of a set of enumerated types which indicate whether the degree of freedom is a translation, temperature, etc.
The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_BulkNumDof()
to return the number of element degrees of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to Bulk object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
loc – [out] Vector of degree of freedom locations
tag – [out] Vector of degree of freedom types
-
void vfe_BulkMass(vfe_Bulk *p, Vdouble x[][3], Vdouble m[])
consistent mass matrix
Compute the consistent mass matrix, m, given the node coordinates, x. The lower triangle of the consistent mass is returned. Use
vfe_BulkMassDiag()
to compute a diagonal mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
m – [out] Degree of freedom consistent mass matrix
-
void vfe_BulkMassDiag(vfe_Bulk *p, Vdouble x[][3], Vdouble md[])
diagonal mass matrix
Compute the diagonal mass matrix, md, given the node coordinates, x. The diagonal mass is returned as a degree of freedom length vector. Use
vfe_BulkMass()
to compute a consistent mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
md – [out] Degree of freedom diagonal mass vector
-
void vfe_BulkElemLoad(vfe_Bulk *p, Vdouble x[][3], Vdouble q[][3], Vdouble f[])
body force vector
Compute the consistent degree of freedom body loads given acceleration load vector, q on an element. The vector q contains an acceleration vector for for each node in the element. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The input element loads are in the units of force per unit mass. Note that the computation of consistent loads uses the material density.
- Errors
SYS_ERROR_COMPUTE
is generated if a non-positive Jacobian is computed.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
q – Array of node accelerations
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_BulkElemHeat(vfe_Bulk *p, Vdouble x[][3], Vdouble q[], Vdouble f[])
body heat generation
Compute the consistent degree of freedom body heat generation given nodal heat generation vector, q on an element. The vector q contains heat generation per volume for each node in the element. The output array of consistent degree of freedom loads, f, contains the heat generation in the element. The input element loads are in the units of power per unit volume.
- Errors
SYS_ERROR_COMPUTE
is generated if a non-positive Jacobian is computed.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
q – Array of node heat fluxes
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_BulkCap(vfe_Bulk *p, Vdouble x[][3], Vdouble u[], Vdouble c[])
consistent capacitance matrix
Compute the consistent capacitance matrix, c, given the node coordinates, x, and temperatures, u. The lower triangle of the consistent capacitance is returned. Use
vfe_BulkCapDiag()
to compute a diagonal capacitance matrix. This calculation requires material density and specific heat.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
u – Degree of freedom vector of temperatures
c – [out] Degree of freedom consistent capacitance matrix
-
void vfe_BulkCapDiag(vfe_Bulk *p, Vdouble x[][3], Vdouble u[], Vdouble cd[])
diagonal capacitance matrix
Compute the diagonal capacitance matrix, cd, given the node coordinates, x, and temperatures, u. The diagonal capacitance is returned as a degree of freedom length vector. Use
vfe_BulkCap()
to compute a consistent capacitance matrix. This calculation requires material density and specific heat.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.
- Parameters:
p – Pointer to Bulk object.
x – Array of node locations.
u – Degree of freedom vector of temperatures
cd – [out] Degree of freedom diagonal capacitance matrix
6.7. Springs - Spring
The Spring module is used to model specialized elastic connections such as spot welds and flexible joints. In their simplest form, springs act like connectors between two end points. The Spring module generalizes this concept to allow for each end point connection to be distributed over an element surface.
The methods associated with a Spring object are the following.
Begin and end an instance of an object, generic object functions
vfe_SpringBegin()
- create an instance of a Spring objectvfe_SpringEnd()
- destroy an instance of a Spring objectvfe_SpringError()
- return Spring object error flag
Attributes and Parameters
vfe_SpringDef()
- specify spring typevfe_SpringInq()
- inquire spring typevfe_SpringSetHistPtr()
- set pointers to material historyvfe_SpringSetLocalSystem()
- set stress axes directionvfe_SpringSetObject()
- set attribute objectvfe_SpringSetParamd()
- set element formulation parametersvfe_SpringSetParami()
- set element formulation parametersvfe_SpringSetBush()
- set generalized spring propertiesvfe_SpringSetLink()
- set link propertiesvfe_SpringSetWeld()
- set spot weld propertiesvfe_SpringSetPntPnt()
- set point to point connectionvfe_SpringSetPntSurf()
- set point to surface connectionvfe_SpringSetSurfSurf()
- set surface to surface connectionvfe_SpringDirCos()
- compute spring local direction cosines
Degree of Freedom and Integration Information
vfe_SpringDofMap()
- query element degree of freedom mapvfe_SpringNumDof()
- query number of element degrees of freedomvfe_SpringNumIntPnt()
- query number of element integration points
Structural analysis computations
vfe_SpringGeomStiff()
- geometric stiffness matrixvfe_SpringInitHist()
- initialize material historyvfe_SpringMass()
- consistent mass matrixvfe_SpringMassDiag()
- diagonal mass matrixvfe_SpringReact()
- reaction vectorvfe_SpringReactStiff()
- reaction vector, stiffness matrixvfe_SpringStiff()
- linear stiffness matrixvfe_SpringStrsStrn()
- stress and strain
Instance a Spring object using vfe_SpringBegin()
.
Once a Spring object is instanced,
define the spring type using vfe_SpringDef()
.
Currently options exist for spot welds, and generalized springs.
The spot weld is similar to the NASTRAN CWELD element and is specified
using vfe_SpringSetWeld()
. The generalized spring is similar to the
NASTRAN CBUSH element and is specified using vfe_SpringSetBush()
.
Set the material function attribute object MatlFun using
vfe_SpringSetObject()
.
Query the element degree of freedom map using vfe_SpringDofMap()
and
vfe_SpringNumDof()
.
If the element is to support a material model,
such as a plastic material model,
which requires a material history then the user must manage the material
history information using vfe_SpringSetHistPtr()
and vfe_SpringInitHist()
.
Since all element geometry is available directly to the element
the material function should point to a primitive material.
6.8. Element Geometry
The basic spring element connection is between two end points, A and B.
The Spring module allows for the motion of an end point to be represented
by the motion of a single node or the motion of an element face at a
point. The three possibilities become node to node, node to element face
and element face to element face.
Therefore the nodal connections of a spring may range from two nodes
to the set of nodes on two element faces.
If the connection is node to node the end point coordinates A and B are
entered using vfe_SpringSetPntPnt()
.
If the connection is node to element face the end point coordinates A and B
and the topology of the element face to be connected to point B are
entered using vfe_SpringSetPntSurf()
.
If the connection is element face
to element face the end point coordinates A and B
and the topologies of the associated element faces are
entered using vfe_SpringSetSurfSurf()
. The topologies of the two element
faces need not be the same.
Options exist for the Spring module to alter the end point coordinates under
certain conditions. One option specifies that the end point coordinates will
be projected to the associated node or element face.
This option is
set using vfe_SpringSetParami()
with type VIS_ENDPROJECT.
An additional end point alteration may occur in the case of the weld type
spring if the length/diameter ratio exceeds certain bounds.
These limits are specified using vfe_SpringSetParamd()
with the types
VFE_MINRATIO and VFE_MAXRATIO.
The element coordinate system x’,y’,z’ is constructed so that the x’
coordinate system is directed from end point A to B
in their final adjusted locations. The orientation of the y’ and z’ axes
is determined using vfe_SpringSetLocalSystem()
.
The link type spring is the simplest: it connects two nodes with a simple spring constant. The force along the direction defined by the two nodes depends on the relative change in distance between the nodes.
The weld type spring is characterized by a specialized beam connection
characterized by a radius and a length from end point A to B.
The computed stress and strain are the same as an equivalent
beam element connected to points A and B.
A MatlFun of a primitive material must be set using vfe_SpringSetObject()
.
The weld type spring geometry is
illustrated below in Figure 6-1.
The generalized bush type spring is characterized by a spring location S and a set of 6 stiffness coefficients and 6 stress and 6 strain recovery coefficients. Each of the stiffness coefficients relates the relative displacement of the 3 translations and 3 rotations of the end points of the spring in the spring element local system to the element forces. The stress is reported as the element force times the associated stress recovery coefficient. The strain is reported as the element relative displacement times the associated strain recovery coefficient. A MatlFun is not required. The generalized (bush) type spring geometry is illustrated below in Figure 6-2.
The nodal freedoms used by spring elements depend upon the spring element type and the end point connection. If the motion of an end point is connected to an element face, only the translational freedoms of the element face nodes are used to compute the motion of the connected end point.
6.9. Function Descriptions
The currently available Spring functions are described in detail in this section.
-
vfe_Spring *vfe_SpringBegin(void)
create an instance of a Spring object
Create an instance of a Spring object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a Spring object using
void vfe_SpringEnd (vfe_Spring *spring)
Return the current value of a Spring object error flag using
Vint vfe_SpringError (vfe_Spring *spring)
- Returns:
The function returns a pointer to the newly created Spring object. If the object creation fails, NULL is returned.
-
void vfe_SpringEnd(vfe_Spring *p)
destroy an instance of a Spring object
-
Vint vfe_SpringError(vfe_Spring *p)
return the current value of a Spring object error flag
-
void vfe_SpringSetObject(vfe_Spring *p, Vint objecttype, Vobject *object)
set attribute object
Set a pointer to an attribute object. Users must supply a MatlFun object prior to computing any quantity that requires a material model definition. The bush type spring does not require a MatlFun object.
- Errors
SYS_ERROR_OBJECTTYPE
is generated if an improper objecttype is specified.
- Parameters:
p – Pointer to Spring object.
objecttype – The object type identifier
x=VFE_MATLFUN MatlFun object
object – Pointer to the object to be set.
-
void vfe_SpringSetParami(vfe_Spring *p, Vint type, Vint iparam)
set element formulation parameters
Set element formulation parameters.
Use
VFE_ENDPROJECT
to specify that an end point coordinate is to be projected to its connected node location or to the surface of its connected element face. By defaultVFE_ENDPROJECT
is set toSYS_ON
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper iparam is specified.
- Parameters:
p – Pointer to Spring object.
type – Type of formulation parameter to set
x=VFE_ENDPROJECT Project end point
iparam – Integer parameter value.
x=SYS_OFF Disable =SYS_ON Enable
-
void vfe_SpringSetParamd(vfe_Spring *p, Vint type, Vdouble dparam)
set element formulation parameters
Set element formulation parameters.
Use
VFE_MINRATIO
andVFE_MAXRATIO
to set the minimum and maximum length/diameter ratios for the weld type spring. By defaultVFE_MINRATIO
is .2 andVFE_MAXRATIO
is 5.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper dparam is specified.
- Parameters:
p – Pointer to Spring object.
type – Type of formulation parameter to set
x=VFE_MINRATIO Minimum length/diameter ratio =VFE_MAXRATIO Maximum length/diameter ratio
dparam – Double precision parameter value.
-
void vfe_SpringDef(vfe_Spring *p, Vint type, Vint contype)
specify spring type
Define the type of spring to be specified.
SPRING_WELD
requires a call tovfe_SpringSetWeld()
,SPRING_BUSH
requires a call tovfe_SpringSetBush()
,SPRING_LINK
requires a call tovfe_SpringSetLink()
, Specify the connection type to be made at end point A and B.SPRING_PNTPNT
requires a call tovfe_SpringSetPntPnt()
,SPRING_PNTSURF
requires a call tovfe_SpringSetPntSurf()
,SPRING_SURFSURF
requires a call tovfe_SpringSetSurfSurf()
.SPRING_NONE
does not required any additional call, as the direction is specified by the coordinates passed to the element routines. This is only valid for type=SPRING_LINK
.Inquire of the defined type and contype as output arguments using
void vfe_SpringInq (vfe_Spring *spring, Vint *type, Vint *contype)
- Errors
SYS_ERROR_ENUM
is generated if an improper type or contype are specified.
- Parameters:
p – Pointer to Spring object.
type – Spring type
x=SPRING_WELD Weld type spring =SPRING_BUSH Bush type spring =SPRING_LINK Link type spring
contype – End point A and B connection type
x=SPRING_PNTPNT Point to point connection =SPRING_PNTSURF Point to surface connection =SPRING_SURFSURF Surface to surface connection =SPRING_NONE Coordinate connection
-
void vfe_SpringInq(vfe_Spring *p, Vint *type, Vint *contype)
inquire of the defined type and contype as output arguments
See
vfe_SpringDef()
-
void vfe_SpringSetWeld(vfe_Spring *p, Vdouble r)
set spot weld properties
Set the radius, r, of the spot weld.
- Errors
SYS_ERROR_OPERATION
is generated if the spring type is notSPRING_WELD
.
- Parameters:
p – Pointer to Spring object.
r – Radius of spot weld
-
void vfe_SpringSetBush(vfe_Spring *p, Vdouble xs[3], Vdouble ks[6], Vdouble sc[6], Vdouble ec[6])
set generalized spring properties
Set the location, xs, of the discrete spring, the spring stiffnesses, ks, and stress and strain recovery coefficients sc and ec.
- Errors
SYS_ERROR_OPERATION
is generated if the spring type is notSPRING_BUSH
.
- Parameters:
p – Pointer to Spring object.
xs – Location of spring
ks – Stiffness coefficients
sc – Stress recovery coefficients
ec – Strain recovery coefficients
-
void vfe_SpringSetLink(vfe_Spring *p, Vdouble stiff)
set link properties
Set the stiffness, stiff, of the spring.
- Errors
SYS_ERROR_OPERATION
is generated if the spring type is notSPRING_LINK
.
- Parameters:
p – Pointer to Spring object.
stiff – Spring stiffness
-
void vfe_SpringSetPntPnt(vfe_Spring *p, Vdouble xa[3], Vdouble xb[3])
set point to point connection
Set the locations of the end points, xa and xb.
- Errors
SYS_ERROR_OPERATION
is generated if the connection type is notSPRING_PNTPNT
.
- Parameters:
p – Pointer to Spring object.
xa – Location of end point A
xb – Location of end point B
-
void vfe_SpringSetPntSurf(vfe_Spring *p, Vdouble xa[3], Vdouble xb[3], Vint shapeb, Vint maxib, Vint maxjb)
set point to surface connection
Set the locations of the end points, xa and xb and the topology of the element face to be connected to end point B.
- Errors
SYS_ERROR_OPERATION
is generated if the connection type is notSPRING_PNTPNT
.
- Parameters:
p – Pointer to Spring object.
xa – Location of end point A
xb – Location of end point B
shapeb – The shape for the element face connected to B
x=SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral
maxib – The number of points along the i direction.
maxjb – The number of points along the j direction.
-
void vfe_SpringSetSurfSurf(vfe_Spring *p, Vdouble xa[3], Vint shapea, Vint maxia, Vint maxja, Vdouble xb[3], Vint shapeb, Vint maxib, Vint maxjb)
set surface to surface connection
Set the locations of the end points, xa and xb and the topology of the element faces to be connected to end points A and B.
- Errors
SYS_ERROR_OPERATION
is generated if the connection type is notSPRING_PNTPNT
.
- Parameters:
p – Pointer to Spring object.
xa – Location of end point A
shapea – The shape for the element face connected to A
x=SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral
maxia – The number of points along the i direction.
maxja – The number of points along the j direction.
xb – Location of end point B
shapeb – The shape for the element face connected to B
x=SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral
maxib – The number of points along the i direction.
maxjb – The number of points along the j direction.
-
void vfe_SpringReactStiff(vfe_Spring *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
reaction vector, stiffness matrix
Compute the reaction vector, r, and optionally the stiffness matrix, k, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kflag – Flag to compute stiffness matrix, k
=SYS_OFF Do not compute stiffness matrix =SYS_ON Compute and return stiffness matrix
r – [out] Degree of freedom reaction vector
k – [out] Degree of freedom stiffness matrix
-
void vfe_SpringStiff(vfe_Spring *p, Vdouble x[][3], Vdouble kl[])
linear stiffness matrix
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
-
void vfe_SpringReact(vfe_Spring *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
reaction vector
Compute the reaction vector, r, given the node coordinates, x, and the degree of freedom displacement vector, u.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
u – Degree of freedom vector of displacements
r – [out] Degree of freedom reaction vector
-
void vfe_SpringNumDof(vfe_Spring *p, Vint analysistype, Vint *nedofs)
query number of element degrees of freedom
Query for number of element degree of freedom nedofs. The number of degrees of freedom will generally be equal to the number of nodal degrees of freedom per node times the number of nodes plus the number of elemental degrees of freedom. Use
vfe_SpringDofMap()
to return the location and type of each degree of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to Spring object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
nedofs – [out] Number of element degrees of freedom
-
void vfe_SpringDofMap(vfe_Spring *p, Vint analysistype, Vint loc[], Vint tag[])
query element degree of freedom map
Query for element degree of freedom map. The degree of freedom map consists of a location index, loc and type, tag for each degree of freedom used by the element.
The location index is either a positive node index into the element connectivity indicating a nodal freedom or a zero value indicating an elemental degree of freedom. The tag indicates the type of the degree of freedom. Tag values are one of a set of enumerated types which indicate whether the degree of freedom is a translation or rotation, etc.
The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_SpringNumDof()
to return the number of element degrees of freedom.- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to Spring object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
loc – [out] Vector of degree of freedom locations
tag – [out] Vector of degree of freedom types
-
void vfe_SpringNumIntPnt(vfe_Spring *p, Vint analysistype, Vint *nepnts)
query number of element integration points
Query for number of element integration points nepnts.
- Errors
SYS_ERROR_ENUM
is generated if an improper analysistype is specified.
- Parameters:
p – Pointer to Spring object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
nepnts – [out] Number of element integration points
-
void vfe_SpringSetHistPtr(vfe_Spring *p, Vdouble *oldhist, Vdouble *newhist)
set pointers to material history
Set pointers to the start of the material history data at the previous step, oldhist and the current step newhist. This function is required when an underlying material model such as plasticity is used. Note that the material history data is not copied by this function, only the pointers themselves are copied. The number of double precision values required for the material history at a step is the number of history variables at an integration point times the number of element integration points. The number of history variables depends on the underlying material model and may be queried using
vfe_MatlFunNumHist()
. The number of element integration points is returned usingvfe_SpringNumIntPnt()
. By default the pointers to the material history are NULL. If the number of history variables is zero, this function need not be called.- Parameters:
p – Pointer to Spring object.
oldhist – Pointer to start of material history at previous step
newhist – Pointer to start of material history at current step
-
void vfe_SpringInitHist(vfe_Spring *p)
initialize material history
Initialize the values of the history variables used in the underlying element or primitive material model for the element. This operation should be performed once for each element (at the first load or time step for example) to initialize the old history variables to reflect the initial configuration condition. If the number of history variables is zero, this function need not be called.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.
- Parameters:
p – Pointer to Spring object.
-
void vfe_SpringMass(vfe_Spring *p, Vdouble x[][3], Vdouble m[])
consistent mass matrix
Compute the consistent mass matrix, m, given the node coordinates, x. The lower triangle of the consistent mass is returned. Use
vfe_SpringMassDiag()
to compute a diagonal mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
m – [out] Degree of freedom consistent mass matrix
-
void vfe_SpringMassDiag(vfe_Spring *p, Vdouble x[][3], Vdouble md[])
diagonal mass matrix
Compute the diagonal mass matrix, md, given the node coordinates, x. The diagonal mass is returned as a degree of freedom length vector. Use
vfe_SpringMass()
to compute a consistent mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
md – [out] Degree of freedom diagonal mass vector
-
void vfe_SpringGeomStiff(vfe_Spring *p, Vdouble x[][3], Vdouble u[], Vdouble kg[])
geometric stiffness matrix
Compute the geometric stiffness matrix, kg, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the geometric stiffness is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kg – [out] Degree of freedom geometric stiffness matrix
-
void vfe_SpringSetLocalSystem(vfe_Spring *p, Vint type, Vdouble vec[], Vdouble angle)
set local coordinate system convention
Specify the convention to be used to construct the orientation of the element local x’,y’,z’ coordinate system with respect to the spring axis. This local system is computed at each integration point location on the element axis and is assumed to be the coordinate system in which the material properties of the element at the integration point are expressed. For stress and strain computation for output using
vfe_SpringStrsStrn()
, the local coordinate system is evaluated at each output location and is the coordinate system in which the output stresses and strains at the output location are expressed.The x’ axis is always constructed to be directed from the first to the second spring end point. The orientation of the y’ and z’ axes perpendicular to the element x’ axis is determined by type. The vec array is only used if the specified type requires position or direction vectors. An additional rotation of the y’,z’ axes about the x’ axis can be specified with angle. By default the local system convention is
SYS_ELEMSYS_STANDARD
with angle set to 0.For a description of element coordinate systems, type, and associated orientation vector data, please see Element Coordinate Systems
- Errors
SYS_ERROR_ENUM
is generated if an improper type is input.
- Parameters:
p – Pointer to Spring object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate y’,z’ axes about the element x’ axis in degrees.
-
void vfe_SpringDirCos(vfe_Spring *p, Vdouble x[][3], Vdouble tm[][3][3])
compute spring local direction cosines
Compute the direction cosine matrices of the element local coordinate system. For stress and strain computation the local coordinate system at each stress output location is the coordinate system in which the output stresses and strains at the location using
vfe_SpringStrsStrn()
are expressed. Given that X’,Y’ and Z’ are three orthonormal vectors indicating the direction of the local coordinate axes in the global coordinate system (x,y,z), then the direction cosine matrix, tm for this local coordinate system is defined as:The local coordinate system is determined by the local system convention set usingtm[0][0] = X'x tm[0][1] = X'y tm[0][2] = X'z tm[1][0] = Y'x tm[1][1] = Y'y tm[1][2] = Y'z tm[2][0] = Z'x tm[2][1] = Z'y tm[2][2] = Z'z
vfe_SpringSetLocalSystem()
.- Parameters:
p – Pointer to Spring object.
x – Array of point locations defining spring axis
tm – [out] Array of direction cosine matrices at the element nodes.
-
void vfe_SpringStrsStrn(vfe_Spring *p, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[])
stress and strain
Compute nodal stresses and strains, strs and strn, given the node coordinates, x, and the degree of freedom displacement vector, u. The stresses and strains are computed in the element local system at end points A and B. The values at end A are then output at all nodes connected to end point A. Similarly for end point B, the values at end B are then output at all nodes connected to end point B. The stress and strain values are ordered first by the 6 components followed by the the number of element nodes.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Spring object.
x – Array of node locations.
u – Degree of freedom vector of displacements
strs – [out] Array of nodal stresses
strn – [out] Array of nodal strains