The methods associated with a ConElas object are the following.
*vfe_ConElasBegin - create an instance of a ConElas object vfe_ConElasEnd - destroy an instance of a ConElas object vfe_ConElasError - return ConElas object error flag
vfe_ConElasDef - define spring type and dimension Inq - inquire spring type and dimension vfe_ConElasSetDofScalar - set dof spring scalar properties vfe_ConElasSetDofVector - set dof spring vector properties vfe_ConElasSetIntrinsic - set intrinsic spring stiffnesses vfe_ConElasSetLocalSystem - set element local system direction vfe_ConElasSetPropPtr - set property pointer
vfe_ConElasDofMap - query element degree of freedom map vfe_ConElasNumDof - query number of element degrees of freedom
vfe_ConElasReact - reaction vector vfe_ConElasReactStiff - reaction vector, stiffness matrix vfe_ConElasReactStiffDamp - reaction vector, stiffness and damping matrices vfe_ConElasStiff - linear stiffness matrix vfe_ConElasStiffDamp - linear stiffness and damping matrices
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.
*vfe_ConElasBegin - create an instance of a ConElas object
vfe_ConElas *vfe_ConElasBegin ()
None
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)
vfe_ConElasDef - define spring type and dimension
void vfe_ConElasDef (vfe_ConElas *conelas, Vint type, Vint dime)
conelas Pointer to ConElas object. type Spring type =CONELAS_DOFSCALAR Degree of freedom scalar spring =CONELAS_DOFVECTOR Degree of freedom vector spring =CONELAS_INTRINSIC Intrinsic spring stiffness dime Spring dimension =VFE_2D 2D analysis =VFE_3D 3D analysis
None
Inquire of a defined type and dime as output arguments using
void vfe_ConElasInq (vfe_ConElas *conelas, Vint *type, Vint *dime)
vfe_ConElasSetDofScalar - set dof spring scalar properties
void vfe_ConElasSetDofScalar (vfe_ConElas *conelas, Vint tag1, Vint tag2, Vdouble kl[3], Vdouble dl[3])
conelas 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
None
vfe_ConElasSetDofVector - set dof spring vector properties
void vfe_ConElasSetDofVector (vfe_ConElas *conelas, Vint tag1, Vint tag2, Vdouble kl[3], Vdouble dl[3])
conelas 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
None
vfe_ConElasSetIntrinsic - set intrinsic spring stiffnesses
void vfe_ConElasSetIntrinsic (vfe_ConElas *conelas, Vint type, Vdouble kext, Vdouble ktor, Vdouble dext, Vdouble dtor)
conelas Pointer to ConElas object. type Spring stiffness type =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
None
vfe_ConElasSetPropPtr - set pointer to element nodal properties
void vfe_ConElasSetPropPtr (vfe_ConElas *conelas, Vint type, Vdouble *propptr)
conelas Pointer to ConElas object. type Type of element property =VFE_PROP_VELOCITY Velocity propptr Pointer to start of element nodal properties
None
vfe_ConElasSetLocalSystem - set element local system direction
void vfe_ConElasSetLocalSystem (vfe_ConElas *conelas, Vint type, Vdouble vec[], Vdouble angle)
conelas 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.
None
For a description of element coordinate systems, type, and associated
orientation vector data, please see
1.4 Element Coordinate Systems
vfe_ConElasDofMap - query element degree of freedom map
void vfe_ConElasDofMap (vfe_ConElas *conelas, Vint analysistype, Vint loc[], Vint tag[])
conelas Pointer to ConElas object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis
loc Vector of degree of freedom locations tag Vector of degree of freedom types
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.
vfe_ConElasNumDof - query number of element degrees of freedom
void vfe_ConElasNumDof (vfe_ConElas *conelas, Vint analysistype, Vint *nedofs)
conelas Pointer to ConElas object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis
nedofs Number of element degrees of freedom
vfe_ConElasReact - reaction vector
void vfe_ConElasReact (vfe_ConElas *conelas, Vdouble x[][3], Vdouble u[], Vdouble r[])
conelas Pointer to ConElas object. x Array of node locations. u Degree of freedom vector of displacements
r Degree of freedom reaction vector
vfe_ConElasReactStiff - reaction vector, stiffness matrix
void vfe_ConElasReactStiff (vfe_ConElas *conelas, Vdouble x[][3], Vdouble u[], Vdouble kflag, Vdouble r[], Vdouble k[])
conelas 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 Degree of freedom reaction vector k Degree of freedom stiffness matrix
vfe_ConElasReactStiffDamp - reaction vector, stiffness and damping matrices
void vfe_ConElasReactStiffDamp (vfe_ConElas *conelas, Vdouble x[][3], Vdouble u[], Vdouble kflag, Vdouble dflag, Vdouble r[], Vdouble k[], Vdouble d[])
conelas 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 Degree of freedom reaction vector k Degree of freedom stiffness matrix d Degree of freedom damping matrix
vfe_ConElasStiff - linear stiffness matrix
void vfe_ConElasStiff (vfe_ConElas *conelas, Vdouble x[][3], Vdouble kl[])
conelas Pointer to ConElas object. x Array of node locations.
kl Degree of freedom stiffness matrix
vfe_ConElasStiffDamp - linear stiffness and damping matrices
void vfe_ConElasStiffDamp (vfe_ConElas *conelas, Vdouble x[][3], Vdouble kl[], Vdouble dl[])
conelas Pointer to ConElas object. x Array of node locations.
kl Degree of freedom stiffness matrix dl Degree of freedom damping matrix
The methods associated with a ConMass object are the following.
*vfe_ConMassBegin - create an instance of a ConMass object vfe_ConMassEnd - destroy an instance of a ConMass object vfe_ConMassError - return ConMass object error flag
vfe_ConMassDef - define mass type and dimension Inq - inquire mass type and dimension vfe_ConMassSetParami - set integer parameters vfe_ConMassSetDofLink - set degree of freedom link properties vfe_ConMassSetLocalSystem - set local inertia tensor system vfe_ConMassSetMatrix - set symmetric mass matrix at a point vfe_ConMassSetOffset - set offset of center of mass vfe_ConMassSetPrinc - set translational mass and rotary inertia vfe_ConMassSetTrans - set single translational mass
vfe_ConMassDofMap - query element degree of freedom map vfe_ConMassNumDof - query number of element degrees of freedom
vfe_ConMassElemLoad - body force vector vfe_ConMassMass - consistent mass matrix vfe_ConMassMassDiag - diagonal mass matrix
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.
*vfe_ConMassBegin - create an instance of a ConMass object
vfe_ConMass *vfe_ConMassBegin ()
None
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)
vfe_ConMassDef - define mass type and dimension
void vfe_ConMassDef (vfe_ConMass *conmass, Vint type, Vint dime)
conmass Pointer to ConMass object. type Mass type =CONMASS_DOFLINK Degree of freedom link =CONMASS_TRANS Translational mass =CONMASS_PRINC Translational mass, rotary inertia =CONMASS_MATRIX Symmetric mass matrix dime Mass dimension =VFE_2D 2D analysis =VFE_3D 3D analysis
None
Inquire of a defined type and dime as output arguments using
void vfe_ConMassInq (vfe_ConMass *conmass, Vint *type, Vint *dime)
vfe_ConMassSetParami - set integer parameters
void vfe_ConMassSetParami (vfe_ConMass *conmass, Vint type, Vint iparam)
conmass Pointer to ConMass object. type Type of formulation parameter to set =VFE_ROTINERTIA Enable rotary inertias iparam Integer parameter value. =SYS_OFF Disable =SYS_ON Enable
None
vfe_ConMassSetDofLink - set degree of freedom link properties
void vfe_ConMassSetDofLink (vfe_ConMass *conmass, Vint tag1, Vint tag2, Vdouble kl[3])
conmass Pointer to ConMass object. tag1 First degree of freedom type tag2 Second degree of freedom type kl Lower triangle of symmetric mass matrix
None
vfe_ConMassSetLocalSystem - set local inertia tensor system
void vfe_ConMassSetLocalSystem (vfe_ConMass *conmass, Vint type, Vdouble vec[], Vdouble angle)
conmass 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.
None
For a description of element coordinate systems, type, and associated
orientation vector data, please see
1.4 Element Coordinate Systems
vfe_ConMassSetTrans - set single translational mass
void vfe_ConMassSetTrans (vfe_ConMass *conmass, Vdouble mass)
conmass Pointer to ConMass object. mass Single translational mass
None
vfe_ConMassSetPrinc - set translational mass and rotary inertia
void vfe_ConMassSetPrinc (vfe_ConMass *conmass, Vdouble mass, Vdouble rotary[3], Vdouble tm[3][3], Vdouble offset[3]);
conmass 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.
None
vfe_ConMassSetMatrix - set symmetric mass matrix at a point
void vfe_ConMassSetMatrix (vfe_ConMass *conmass, Vdouble a[21])
conmass Pointer to ConMass object. a Complete lower triangle of symmetric mass matrix at a node
None
- - | 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] | - -
In 2D 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_RZ.
- - | a[ 0] a[ 1] a[ 3] | M = | a[ 1] a[ 2] a[ 4] | | a[ 3] a[ 4] a[ 5] | - -
vfe_ConMassSetOffset - set offset of center of mass
void vfe_ConMassSetOffset (vfe_ConMass *conmass, Vdouble offset[3])
conmass Pointer to ConMass object. offset Offset of center of mass in global coordinates
None
vfe_ConMassDofMap - query element degree of freedom map
void vfe_ConMassDofMap (vfe_ConMass *conmass, Vint analysistype, Vint loc[], Vint tag[])
conmass Pointer to ConMass object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
loc Vector of degree of freedom locations tag Vector of degree of freedom types
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.
vfe_ConMassNumDof - query number of element degrees of freedom
void vfe_ConMassNumDof (vfe_ConMass *conmass, Vint analysistype, Vint *nedofs)
conmass Pointer to ConMass object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
nedofs Number of element degrees of freedom
vfe_ConMassElemLoad - body force vector
void vfe_ConMassElemLoad (vfe_ConMass *conmass, Vdouble x[][3], Vdouble q[][3], Vdouble f[],
conmass Pointer to ConMass object. x Array of node locations. q Array of node accelerations
f Degree of freedom vector of consistent loads.
vfe_ConMassMass - consistent mass matrix
void vfe_ConMassMass (vfe_ConMass *conmass, Vdouble x[][3], Vdouble m[])
conmass Pointer to ConMass object. x Array of node locations.
m Consistent mass matrix
vfe_ConMassMassDiag - diagonal mass matrix
void vfe_ConMassMassDiag (vfe_ConMass *conmass, Vdouble x[][3], Vdouble md[])
conmass Pointer to ConMass object. x Array of node locations.
md Diagonal mass matrix
The methods associated with a Bulk object are the following.
*vfe_BulkBegin - create an instance of a Bulk object vfe_BulkEnd - destroy an instance of a Bulk object vfe_BulkError - return Bulk object error flag
vfe_BulkDef - define bulk dimension Inq - inquire bulk dimension vfe_BulkSetObject - set attribute object vfe_BulkSetPropPtr - set pointer to element nodal properties
vfe_BulkDofMap - query element degree of freedom map vfe_BulkNumDof - query number of element degrees of freedom
vfe_BulkElemLoad - body force vector vfe_BulkMass - consistent mass matrix vfe_BulkMassDiag - diagonal mass matrix
vfe_BulkElemHeat - body heat generation vfe_BulkCap - consistent capacitance matrix vfe_BulkCapDiag - diagonal capacitance matrix
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.
*vfe_BulkBegin - create an instance of a Bulk object
vfe_Bulk *vfe_BulkBegin ()
None
Destroy an instance of a Bulk object using
void vfe_BulkEnd (vfe_Bulk *bulk)
Return the current value of a Bulk object error flag using
Vint vfe_BulkError (vfe_Bulk *bulk)
vfe_BulkDef - define bulk dimension
void vfe_BulkDef (vfe_Bulk *bulk, Vint dime)
bulk Pointer to Bulk object. dime Bulk dimension =VFE_2D 2D analysis =VFE_3D 3D analysis
None
Inquire of a defined dime as an output argument using
void vfe_BulkInq (vfe_Bulk *bulk, Vint *dime)
vfe_BulkSetObject - set attribute object
void vfe_BulkSetObject (vfe_Bulk *bulk, Vint objecttype, Vobject *object)
bulk Pointer to Bulk object. objecttype The object type identifier =VFE_MATLFUN MatlFun object object Pointer to the object to be set.
None
vfe_BulkSetPropPtr - set pointer to element nodal properties
void vfe_BulkSetPropPtr (vfe_Bulk *bulk, Vint type, Vdouble *propptr)
bulk Pointer to Bulk object. type Type of element property =VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_VOLUME Volume propptr Pointer to start of element nodal properties
None
vfe_BulkDofMap - query element degree of freedom map
void vfe_BulkDofMap (vfe_Bulk *bulk, Vint analysistype, Vint loc[], Vint tag[])
bulk Pointer to Bulk object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
loc Vector of degree of freedom locations tag Vector of degree of freedom types
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.
vfe_BulkNumDof - query number of element degrees of freedom
void vfe_BulkNumDof (vfe_Bulk *bulk, Vint analysistype, Vint *nedofs)
bulk Pointer to Bulk object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis =VFE_ANALYSIS_THERMAL Thermal analysis
nedofs Number of element degrees of freedom
vfe_BulkElemLoad - body force vector
void vfe_BulkElemLoad (vfe_Bulk *bulk, Vdouble x[][3], Vdouble q[][3], Vdouble f[],
bulk Pointer to Bulk object. x Array of node locations. q Array of node accelerations
f Degree of freedom vector of consistent loads.
vfe_BulkMass - consistent mass matrix
void vfe_BulkMass (vfe_Bulk *bulk, Vdouble x[][3], Vdouble m[])
bulk Pointer to Bulk object. x Array of node locations.
m Degree of freedom consistent mass matrix
vfe_BulkMassDiag - diagonal mass matrix
void vfe_BulkMassDiag (vfe_Bulk *bulk, Vdouble x[][3], Vdouble md[])
bulk Pointer to Bulk object. x Array of node locations.
md Degree of freedom diagonal mass vector
vfe_BulkElemHeat - body heat generation
void vfe_BulkElemHeat (vfe_Bulk *bulk, Vdouble x[][3], Vdouble q[], Vdouble f[],
bulk Pointer to Bulk object. x Array of node locations. q Array of node heat fluxes
f Degree of freedom vector of consistent loads.
vfe_BulkCap - consistent capacitance matrix
void vfe_BulkCap (vfe_Bulk *bulk, Vdouble x[][3], Vdouble u[], Vdouble c[])
bulk Pointer to Bulk object. x Array of node locations. u Degree of freedom vector of temperatures
c Degree of freedom consistent capacitance matrix
vfe_BulkCapDiag - diagonal capacitance matrix
void vfe_BulkCapDiag (vfe_Bulk *bulk, Vdouble x[][3], Vdouble u[], Vdouble cd[])
bulk Pointer to Bulk object. x Array of node locations. u Degree of freedom vector of temperatures
cd Degree of freedom diagonal capacitance matrix
The methods associated with a Spring object are the following.
*vfe_SpringBegin - create an instance of a Spring object vfe_SpringEnd - destroy an instance of a Spring object vfe_SpringError - return Spring object error flag
vfe_SpringDef - specify spring type Inq - inquire spring type vfe_SpringSetHistPtr - set pointers to material history vfe_SpringSetLocalSystem - set stress axes direction vfe_SpringSetObject - set attribute object vfe_SpringSetParamd - set element formulation parameters vfe_SpringSetParami - set element formulation parameters vfe_SpringSetBush - set generalized spring properties vfe_SpringSetLink - set link properties vfe_SpringSetWeld - set spot weld properties vfe_SpringSetPntPnt - set point to point connection vfe_SpringSetPntSurf - set point to surface connection vfe_SpringSetSurfSurf - set surface to surface connection vfe_SpringDirCos - compute spring local direction cosines
vfe_SpringDofMap - query element degree of freedom map vfe_SpringNumDof - query number of element degrees of freedom vfe_SpringNumIntPnt - query number of element integration points
vfe_SpringGeomStiff - geometric stiffness matrix vfe_SpringInitHist - initialize material history vfe_SpringMass - consistent mass matrix vfe_SpringMassDiag - diagonal mass matrix vfe_SpringReact - reaction vector vfe_SpringReactStiff - reaction vector, stiffness matrix vfe_SpringStiff - linear stiffness matrix vfe_SpringStrsStrn - stress and strain
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.
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.
*vfe_SpringBegin - create an instance of a Spring object
vfe_Spring *vfe_SpringBegin ()
None
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)
vfe_SpringDef - specify spring type
void vfe_SpringDef (vfe_Spring *spring, Vint type, Vint contype)
spring Pointer to Spring object. type Spring type =SPRING_WELD Weld type spring =SPRING_BUSH Bush type spring =SPRING_LINK Link type spring contype End point A and B connection type =SPRING_PNTPNT Point to point connection =SPRING_PNTSURF Point to surface connection =SPRING_SURFSURF Surface to surface connection =SPRING_NONE Coordinate connection
None
Inquire of the defined type and contype as output arguments using
void vfe_SpringInq (vfe_Spring *spring, Vint *type, Vint *contype)
vfe_SpringSetHistPtr - set pointers to material history
void vfe_SpringSetHistPtr (vfe_Spring *spring, Vdouble *oldhist, Vdouble *newhist)
spring Pointer to Spring object. oldhist Pointer to start of material history at previous step newhist Pointer to start of material history at current step
None
vfe_SpringSetLocalSystem - set local coordinate system convention
void vfe_SpringSetLocalSystem (vfe_Spring *spring, Vint type, Vdouble vec[], Vdouble angle)
spring 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.
None
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
1.4 Element Coordinate Systems
vfe_SpringSetObject - set attribute object
void vfe_SpringSetObject (vfe_Spring *spring, Vint objecttype, Vobject *object)
spring Pointer to Spring object. objecttype The object type identifier =VFE_MATLFUN MatlFun object object Pointer to the object to be set.
None
vfe_SpringSetParamd - set element formulation parameters
void vfe_SpringSetParamd (vfe_Spring *spring, Vint type, Vdouble dparam)
spring Pointer to Spring object. type Type of formulation parameter to set =VFE_MINRATIO Minimum length/diameter ratio =VFE_MAXRATIO Maximum length/diameter ratio dparam Double precision parameter value.
None
Use VFE_MINRATIO and VFE_MAXRATIO to set the minimum and maximum length/diameter ratios for the weld type spring. By default VFE_MINRATIO is .2 and VFE_MAXRATIO is 5.
vfe_SpringSetParami - set element formulation parameters
void vfe_SpringSetParami (vfe_Spring *spring, Vint type, Vint iparam)
spring Pointer to Spring object. type Type of formulation parameter to set =VFE_ENDPROJECT Project end point iparam Integer parameter value. =SYS_OFF Disable =SYS_ON Enable
None
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 default VFE_ENDPROJECT is set to SYS_ON.
vfe_SpringSetBush - set generalized spring properties
void vfe_SpringSetBush (vfe_Spring *spring, Vdouble xs[3], Vdouble ks[6], Vdouble sc[6], Vdouble ec[6])
spring Pointer to Spring object. xs Location of spring ks Stiffness coefficients sc Stress recovery coefficients ec Strain recovery coefficients
None
vfe_SpringSetLink - set link properties
void vfe_SpringSetLink (vfe_Spring *spring, Vdouble stiff)
spring Pointer to Spring object. stiff Spring stiffness
None
vfe_SpringSetWeld - set spot weld properties
void vfe_SpringSetWeld (vfe_Spring *spring, Vdouble r)
spring Pointer to Spring object. r Radius of spot weld
None
vfe_SpringSetPntPnt - set point to point connection
void vfe_SpringSetPntPnt (vfe_Spring *spring, Vdouble xa[3], Vdouble xb[3])
spring Pointer to Spring object. xa Location of end point A xb Location of end point B
None
vfe_SpringSetPntSurf - set point to surface connection
void vfe_SpringSetPntSurf (vfe_Spring *spring, Vdouble xa[3], Vdouble xb[3], Vint shapeb, Vint maxib, Vint maxjb)
spring 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 =SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral maxib The number of points along the i direction. maxjb The number of points along the j direction.
None
vfe_SpringSetSurfSurf - set surface to surface connection
void vfe_SpringSetSurfSurf (vfe_Spring *spring, Vdouble xa[3], Vint shapea, Vint maxia, Vint maxja, Vdouble xb[3], Vint shapeb, Vint maxib, Vint maxjb)
spring Pointer to Spring object. xa Location of end point A shapea The shape for the element face connected to A =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 =SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral maxib The number of points along the i direction. maxjb The number of points along the j direction.
None
vfe_SpringDofMap - query element degree of freedom map
void vfe_SpringDofMap (vfe_Spring *spring, Vint analysistype, Vint loc[], Vint tag[])
spring Pointer to Spring object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis
loc Vector of degree of freedom locations tag Vector of degree of freedom types
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.
vfe_SpringNumDof - query number of element degrees of freedom
void vfe_SpringNumDof (vfe_Spring *spring, Vint analysistype, Vint *nedofs)
spring Pointer to Spring object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis
nedofs Number of element degrees of freedom
vfe_SpringNumIntPnt - query number of element integration points
void vfe_SpringNumIntPnt (vfe_Spring *spring, Vint analysistype, Vint *nepnts)
spring Pointer to Spring object. analysistype The type of analysis =VFE_ANALYSIS_STRUCTURAL Structural analysis
nepnts Number of element integration points
vfe_SpringDirCos - compute spring local direction cosines
void vfe_SpringDirCos (vfe_Spring *spring, Vdouble x[][3], Vdouble tm[][3][3])
spring Pointer to Spring object. x Array of point locations defining spring axis
tm Array of direction cosine matrices at the element nodes.
tm[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
The local coordinate system is determined by the local system convention set using vfe_SpringSetLocalSystem.
vfe_SpringGeomStiff - geometric stiffness matrix
void vfe_SpringGeomStiff (vfe_Spring *spring, Vdouble x[][3], Vdouble u[], Vdouble kg[])
spring Pointer to Spring object. x Array of node locations. u Degree of freedom vector of displacements
kg Degree of freedom geometric stiffness matrix
vfe_SpringInitHist - initialize material history
void vfe_SpringInitHist (vfe_Spring *spring)
spring Pointer to Spring object.
None
vfe_SpringMass - consistent mass matrix
void vfe_SpringMass (vfe_Spring *spring, Vdouble x[][3], Vdouble m[])
spring Pointer to Spring object. x Array of node locations.
m Degree of freedom consistent mass matrix
vfe_SpringMassDiag - diagonal mass matrix
void vfe_SpringMassDiag (vfe_Spring *spring, Vdouble x[][3], Vdouble md[])
spring Pointer to Spring object. x Array of node locations.
md Degree of freedom diagonal mass vector
vfe_SpringReact - reaction vector
void vfe_SpringReact (vfe_Spring *spring, Vdouble x[][3], Vdouble u[], Vdouble r[])
spring Pointer to Spring object. x Array of node locations. u Degree of freedom vector of displacements
r Degree of freedom reaction vector
vfe_SpringReactStiff - reaction vector, stiffness matrix
void vfe_SpringReactStiff (vfe_Spring *spring, Vdouble x[][3], Vdouble u[], Vdouble kflag, Vdouble r[], Vdouble k[])
spring 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 Degree of freedom reaction vector k Degree of freedom stiffness matrix
vfe_SpringStiff - linear stiffness matrix
void vfe_SpringStiff (vfe_Spring *spring, Vdouble x[][3], Vdouble kl[])
spring Pointer to Spring object. x Array of node locations.
kl Degree of freedom stiffness matrix
vfe_SpringStrsStrn - stress and strain
void vfe_SpringStrsStrn (vfe_Spring *spring, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[])
spring Pointer to Spring object. x Array of node locations. u Degree of freedom vector of displacements
strs Array of nodal stresses strn Array of nodal strains