5. Truss Elements - Truss3D

Truss elements are a specialization of a beam element in which only axial forces are supported. Truss elements are one-dimensional structural members in 3D space, that are usually long and slender, which only transmit axial force. The Truss3D module is used to simulate 3D truss elements. The 3 node version of the truss is useful in modelling curved cable structures.

The methods associated with a Truss3D object are the following.

Instance a Truss3D object using vfe_Truss3DBegin(). Once a Truss3D is instanced, set the material function attribute object MatlFun using vfe_Truss3DSetObject(). The current topology of the element is specified using vfe_Truss3DSetTopology(). The area of the truss element is specified using vfe_Truss3DSetPropPtr(). Query the element degree of freedom map using vfe_Truss3DDofMap() and vfe_Truss3DNumDof().

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_Truss3DSetHistPtr() and vfe_Truss3DInitHist(). Since all element geometry is available directly to the element in terms of node coordinates and areas, the material function should point to a primitive material.

5.1. Element Geometry

The basic truss element geometry is defined by the node coordinates of the truss axis and the cross section properties of the truss. The truss section is defined by a cross sectional area. The cross section is always assumed to be perpendicular to the truss axis. The truss section properties which are set using vfe_Truss3DSetPropPtr() are illustrated below in Figure 5-2. Other element quantities which may be set using this function are node temperature and reference temperature. The x’ axis of the truss is always tangent to the axis of the truss. The components of stress y’ and z’ in the plane perpendicular to x’ are always zero as well as all shear stresses. The element natural coordinate, r, is tangent to the the truss centroidal axis.

../../_images/vfetools-truss2.gif

Figure 5-2, Geometry for 3 Node 3D Truss Element

Use the function vfe_Truss3DStrsAdapt() to aid in computing element strain energy, strain energy error and other useful quantities to aid in solution error estimation and mesh adaptation. The function vfe_Truss3DHFlxAdapt() performs a similar computation for heat transfer analysis.

5.2. Function Descriptions

The currently available Truss3D functions are described in detail in this section.

vfe_Truss3D *vfe_Truss3DBegin(void)

create an instance of a Truss3D object

Create an instance of a Truss3D object. Memory is allocated for the object private data and the pointer to the object is returned. Default topology is the 2-noded truss.

Destroy an instance of a Truss3D object using

void vfe_Truss3DEnd (vfe_Truss3D *truss3d)

Return the current value of a Truss3D object error flag using

Vint vfe_Truss3DError (vfe_Truss3D *truss3d)

Returns:

The function returns a pointer to the newly created Truss3D object. If the object creation fails, NULL is returned.

void vfe_Truss3DEnd(vfe_Truss3D *p)

destroy an instance of a Truss3D object

See vfe_Truss3DBegin()

Vint vfe_Truss3DError(vfe_Truss3D *p)

return the current value of a Truss3D object error

See vfe_Truss3DBegin()

void vfe_Truss3DSetObject(vfe_Truss3D *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 Truss3D object.

  • objecttype – The object type identifier

    x=VFE_MATLFUN            MatlFun object
    

  • object – Pointer to the object to be set.

void vfe_Truss3DSetParami(vfe_Truss3D *p, Vint type, Vint iparam)

set element formulation parameters

Set element formulation parameters. Set element strain type using VFE_STRAINTYPE with a value of either VFE_LARGESTRAIN to enable large strain or VFE_SMALLSTRAIN to enable small strains. By default VFE_STRAINTYPE is set to VFE_SMALLSTRAIN. If a small strain material is used under VFE_LARGESTRAIN then it is assumed that the material law relates the Green-Lagrange strain tensor with the second second Piola-Kirchhoff stress tensor. Otherwise, the material law relates the deformation gradient with the Cauchy stress.

The parameter VFE_TEMPMATLAVG toggles the method for computing the temperature used for evaluating temperature dependent material properties. If enabled, the temperature used for temperature dependent material properties is the average of the element node point temperatures. If disabled, the temperature is isoparametrically interpolated from the node point temperatures at each element integration point. By default VFE_TEMPMATLAVG is set to SYS_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 Truss3D object.

  • type – Type of formulation parameter to set

    x=VFE_STRAINTYPE         Element strain type
     =VFE_TEMPMATLAVE        Average material temperature flag
    

  • iparam – Integer parameter value.

    x=VFE_SMALLSTRAIN        Small strain
     =VFE_LARGESTRAIN        Large strain
    

void vfe_Truss3DSetParamd(vfe_Truss3D *p, Vint type, Vdouble dparam)

set element formulation parameters

Set element formulation parameters. Use VFE_MAXPROJANG to set the maximum angle in degrees between the truss tangent at an element integration point and the truss tangent at a node. By default VFE_MAXPROJANG is set to 90. degrees.

Use VFE_ADDEDMASS to add non-structural mass/length to the truss element. By default VFE_ADDEDMASS is set to 0..

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 Truss3D object.

  • type – Type of formulation parameter to set

    x=VFE_MAXPROJANG         Maximum nodal projection angle
     =VFE_ADDEDMASS          Additional mass/length
    

  • dparam – Double precision parameter value.

void vfe_Truss3DSetTopology(vfe_Truss3D *p, Vint maxi)

set element topology

Specify the topology of a 3D truss element. If maxi is set to 3 then a quadratic element form is assumed. The default topology is maxi = 0.

Errors

SYS_ERROR_VALUE is generated if an improper maxi is specified.

Parameters:
  • p – Pointer to Truss3D object.

  • maxi – The number of points along the i direction. If maxi = 0 then the linear Serendipity element form is assumed.

void vfe_Truss3DSetPropPtr(vfe_Truss3D *p, Vint type, Vdouble *propptr)

set pointer to element nodal properties

Set a pointer to the start of a specified type of element nodal properties. Note that the properties are not copied by this function, only the pointer itself is copied. Element nodal properties must be defined as an array of Vdouble, as they should be defined for all nodes. If a property pointer is not set the element assumes a default value. for the associated property. By default the area is 1., the temperature is 0. and the reference temperature is 0. If a pre-stress is set it is added to the computed stress resulting from the material model.

Errors

SYS_ERROR_ENUM is generated if an improper type is specified.

Parameters:
  • p – Pointer to Truss3D object.

  • type – Type of element property

    x=VFE_PROP_AREA             Areas
     =VFE_PROP_TEMPERATURE      Temperatures
     =VFE_PROP_TEMPREF          Reference temperatures
     =VFE_PROP_PRESTRESS        Prestresses
    

  • propptr – Pointer to start of element nodal properties

void vfe_Truss3DStiff(vfe_Truss3D *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 Truss3D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom stiffness matrix

void vfe_Truss3DStrsStrn(vfe_Truss3D *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 truss local system. The x’ axis of the system is along the axis of the element. The convention used to generate local coordinate systems is specified using vfe_Truss3DSetLocalSystem(). The orientation of the y’ and z’ axes is in some sense arbitrary since there is only non-zero stress along the x’ axis. The stress and strain values are ordered first by the 6 tensor 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 Truss3D 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

void vfe_Truss3DReact(vfe_Truss3D *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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • r[out] Degree of freedom reaction vector

void vfe_Truss3DReactStiff(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])

reaction vector, stiffness 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 Truss3D 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_Truss3DInitHist(vfe_Truss3D *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 Truss3D object.

void vfe_Truss3DNumDof(vfe_Truss3D *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_Truss3DDofMap() 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 Truss3D 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_Truss3DNumIntPnt(vfe_Truss3D *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 Truss3D object.

  • analysistype – The type of analysis

    x=VFE_ANALYSIS_STRUCTURAL Structural analysis
     =VFE_ANALYSIS_THERMAL    Thermal analysis
    

  • nepnts[out] Number of element integration points

void vfe_Truss3DDofMap(vfe_Truss3D *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, temperature, etc.

The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use vfe_Truss3DNumDof() 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 Truss3D 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_Truss3DSetHistPtr(vfe_Truss3D *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 using vfe_Truss3DNumIntPnt(). 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 Truss3D object.

  • oldhist – Pointer to start of material history at previous step

  • newhist – Pointer to start of material history at current step

void vfe_Truss3DMass(vfe_Truss3D *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_Truss3DMassDiag() 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 Truss3D object.

  • x – Array of node locations.

  • m[out] Degree of freedom consistent mass matrix

void vfe_Truss3DMassDiag(vfe_Truss3D *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_Truss3DMass() 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 Truss3D object.

  • x – Array of node locations.

  • md[out] Degree of freedom diagonal mass vector

void vfe_Truss3DGeomStiff(vfe_Truss3D *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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • kg[out] Degree of freedom geometric stiffness matrix

void vfe_Truss3DDistLoad(vfe_Truss3D *p, Vdouble x[][3], Vint enttype, Vint no, Vint loadtype, Vdouble q[], Vdouble f[])

distributed load vector

Compute the consistent degree of freedom loads given a distributed load, q on a 3D truss element edge. The vector q contains values for the load type for each node in the element. If the loadtype is VFE_DISTLOAD_TRAC then q contains a vector traction at each element node. If the traction is applied to an edge the units are force/length. If the loadtype is VFE_DISTLOAD_TANGFORCE then q contains a scalar force/length at each element node.

Errors

  • SYS_ERROR_ENUM is generated if an improper enttype or loadtype is specified.

  • SYS_ERROR_OPERATION is generated if an invalid combination of enttype and loadtype is specified.

  • SYS_ERROR_VALUE is generated if an improper no is specified.

  • SYS_ERROR_COMPUTE is generated if a zero edge Jacobian is computed.

Parameters:
  • p – Pointer to Truss3D object.

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_EDGE     Element edge
    

  • no – Element edge number

  • loadtype – Distributed load type

    x=VFE_DISTLOAD_TRAC      Load directed along vector
     =VFE_DISTLOAD_TANGFORCE Load directed along element edge
    

  • q – Vector of distributed load values

  • f[out] Degree of freedom vector of consistent loads.

void vfe_Truss3DElemLoad(vfe_Truss3D *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 Truss3D object.

  • x – Array of node locations.

  • q – Array of node accelerations

  • f[out] Degree of freedom vector of consistent loads.

void vfe_Truss3DDistHeat(vfe_Truss3D *p, Vdouble x[][3], Vint enttype, Vint no, Vdouble q[], Vdouble f[])

distributed heat loads

Compute the consistent degree of freedom loads given a distributed heat load, q along a 3D truss element edge. The vector q contains values for the heat flux for each node in the element. The distributed loads are in units of heat flux per unit length.

Errors

  • SYS_ERROR_ENUM is generated if an improper enttype is specified.

  • SYS_ERROR_VALUE is generated if an improper no is specified.

  • SYS_ERROR_COMPUTE is generated if a zero edge Jacobian is computed.

Parameters:
  • p – Pointer to Truss3D object.

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_EDGE     Element edge
    

  • no – Element edge number

  • q – Vector of distributed load values

  • f[out] Degree of freedom vector of consistent loads.

void vfe_Truss3DElemHeat(vfe_Truss3D *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 Truss3D object.

  • x – Array of node locations.

  • q – Array of node heat fluxes

  • f[out] Degree of freedom vector of consistent loads.

void vfe_Truss3DPower(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vdouble r[])

thermal power vector

Compute the power vector, r, given the node coordinates, x, and the degree of freedom temperature 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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • r[out] Degree of freedom power vector

void vfe_Truss3DPowerCond(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])

thermal power, conductance matrix

Compute the power vector, r, and optionally the conductance matrix, k, given the node coordinates, x, and the degree of freedom temperature vector, u. The lower triangle of the conductance 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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • kflag – Flag to compute conductance matrix, k

    =SYS_OFF      Do not compute conductance matrix
    =SYS_ON       Compute and return conductance matrix
    

  • r[out] Degree of freedom power vector

  • k[out] Degree of freedom conductance matrix

void vfe_Truss3DCond(vfe_Truss3D *p, Vdouble x[][3], Vdouble kl[])

thermal conductance matrix

Compute the linear conductance matrix, kl, given the node coordinates, x. The lower triangle of the conductance 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 Truss3D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom conductance matrix

void vfe_Truss3DHFlxTGrd(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vdouble hflx[], Vdouble tgrd[])

heat flux and thermal gradient

Compute nodal heat fluxes and temperature gradients, hflx and tgrd, given the node coordinates, x, and the degree of freedom temperature vector, u. The flux and gradients are computed in the truss local system. The x’ axis of the system is along the axis of the element. The convention used to generate local coordinate systems is specified using vfe_Truss3DSetLocalSystem(). The orientation of the y’ and z’ axes is in some sense arbitrary since there is only non-zero flux along the x’ axis.

The flux and gradient values are ordered first by the 3 vectoral components followed by 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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • hflx[out] Array of nodal heat fluxes

  • tgrd[out] Array of nodal temperature gradients

void vfe_Truss3DCap(vfe_Truss3D *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_Truss3DCapDiag() 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 or the maximum nodal projection angle is exceeded.

Parameters:
  • p – Pointer to Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • c[out] Degree of freedom consistent capacitance matrix

void vfe_Truss3DCapDiag(vfe_Truss3D *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_Truss3DCap() 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.

  • SYS_ERROR_COMPUTE is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.

Parameters:
  • p – Pointer to Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • cd[out] Degree of freedom diagonal capacitance matrix

void vfe_Truss3DDirCos(vfe_Truss3D *p, Vdouble x[][3], Vdouble tm[][3][3])

compute truss 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_Truss3DStrsStrn() 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:

 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_Truss3DSetLocalSystem().

Parameters:
  • p – Pointer to Truss3D object.

  • x – Array of point locations defining truss axis

  • tm[out] Array of direction cosine matrices at the element nodes.

void vfe_Truss3DSetLocalSystem(vfe_Truss3D *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 truss axis. This local system is computed at each integration point location on the truss 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_Truss3DStrsStrn(), 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 tangent to the truss axis. The orientation of the y’ and z’ axes perpendicular to the truss 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 Truss3D object.

  • type – Local system convention

  • vec – Orientation vector data

  • angle – Angle to rotate truss y’,z’ axes about the truss x’ axis in degrees.

void vfe_Truss3DStrsAdapt(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vdouble strss[], Vdouble *setot, Vdouble *seerr, Vdouble *hsize, Vdouble *order, Vdouble *d)

stress based error analysis

Compute the element total strain energy, setot, and strain energy error, seerr, given the element displacements, u, and an estimate of the exact nodal stresses, strss. In addition useful quantities such as the characteristic length, effective polynomial order and dimension of the element are returned. The element dimension, d, is 1. These quantities are useful for computing new characteristic element length for mesh adaptation.

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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • strss – Array of recovered nodal stresses

  • setot[out] Total strain energy

  • seerr[out] Strain energy error

  • hsize[out] Characteristic length

  • order[out] Effective polynomial order

  • d[out] Dimension

void vfe_Truss3DHFlxAdapt(vfe_Truss3D *p, Vdouble x[][3], Vdouble u[], Vdouble hflxs[], Vdouble *hetot, Vdouble *heerr, Vdouble *hsize, Vdouble *order, Vdouble *d)

heat flux based error analysis

Compute the element total heat energy, hetot, and heat energy error, heerr, given the element temperatures, u, and an estimate of the exact nodal heat flux, hflxs. In addition useful quantities such as the characteristic length, effective polynomial order and dimension of the element are returned. The element dimension, d, is 1. These quantities are useful for computing new characteristic element length for mesh adaptation.

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 Truss3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • hflxs – Array of recovered nodal heat flux

  • hetot[out] Total heat energy

  • heerr[out] Heat energy error

  • hsize[out] Characteristic length

  • order[out] Effective polynomial order

  • d[out] Dimension