7. Interface and Gap Elements - Inter2D, Inter3D, Gap

Interface and gap elements are designed to model surface interactions between two surfaces termed the master and slave surfaces. In structural analysis these interactions include various types of surface contact in which the relative displacement of the interacting faces remains small and the nodes on opposite faces are closely positioned. In thermal analysis these interactions include thermal contact resistance, surface convection and radiation effects. In thermal analysis the nodes on opposite faces are not restrained to be closely positioned.

When used to model contacting surfaces, the surface interaction can consist of a normal pressure and two tangential tractions which would typically arise from contact and sliding frictional effects. The contact pressures are functions of the relative displacement of the master and slave surfaces. The exact nature of the pressure, relative displacement relationship is determined by the underlying interface material model (usually an InterProp object). The relative displacements, dx, dy, dz along the element local axes are computed given the total displacement of the master surface Dxm,Dym,Dzm and the slave surface Dxs,Dys,Dzs in the element local system.

dx = Dxm - Dxs - ci
dy = Dym - Dys
dz = Dzm - Dzs

where ci is the initial clearance. The initial clearance is the sum of the specified element node clearance and the actual geometric distance between the master and slave surfaces along the master surface normal. The inclusion of the actual geometric distance in the initial clearance calculation is optional.

When used to model thermal interactions, in all cases the geometry of the surface is determined solely from the coordinates and other element node properties associated with the nodes on the master surface. The type of thermal interaction selected determines the expression used to compute the generalized temperature difference, dt, given the temperature Tm of the master surface and the temperature Ts of the slave surface. The flux q across the interface is given by:

q = -h*dt

where h is the film coefficient supplied from the underlying interface material (usually an InterProp object) and as such may be temperature dependent. The type of interface is specified by selecting an exchange law and optional exponent, a. Currently supported thermal interactions are:

  • VFE_HEATEXCH_NONE, dt = 0.

  • VFE_HEATEXCH_LINEAR, linear law, dt = Ts - Tm

  • VFE_HEATEXCH_EXP, exponential law, dt = Ts**a - Tm**a

  • VFE_HEATEXCH_POW, power law, dt = (Ts - Tm)**a

  • VFE_HEATEXCH_RAD, black body radiation law, dt = (Ts**4 - Tm**4)

7.1. 2D Interface Elements - Inter2D

The Inter2D module is used to simulate surface interactions in 2D space. The interface element may be used in structural and thermal analysis.

The methods associated with a Inter2D object are the following.

Instance a Inter2D object using vfe_Inter2DBegin(). Once a Inter2D is instanced, set the material function attribute object MatlFun using vfe_Inter2DSetObject(). The current topology of the element is specified using vfe_Inter2DSetTopology(). Certain geometric details of the interface such as initial gap separation are specified using vfe_Inter2DSetPropPtr(). Query the element degree of freedom map using vfe_Inter2DDofMap() and vfe_Inter2DNumDof().

If the element is to support a material model, such as a surface friction material model, which requires a material history then the user must manage the material history information using vfe_Inter2DSetHistPtr() and vfe_Inter2DInitHist().

7.2. Element Geometry

The interface element consists of two matched surfaces termed the master surface and the slave surface. The basic interface element geometry is defined by the node coordinates of the the master surface and possibly the coordinates of the slave surface. The separation of the two surfaces is defined by normals and clearances specified at the nodes of the master surface. These properties are set using vfe_Inter2DSetPropPtr(). Other element node properties which may be set are temperature, depth and area factor.

The master and slave surfaces may be either explicitly defined as surfaces in terms of the nodal geometry as shown in Figure 7-1a or implicitly as lines in terms of the nodal geometry and associated nodal area factors as shown in Figure 7-1b. For 2D planar analysis the surface area is multiplied by the element node depth. If the interface element is being used to model surface contact in structural analysis, the separation of the master and slave surfaces is always measured along the normal to the master surface. The master surface normal may be determined from the nodal geometry or specified at the master surface nodes using vfe_Inter2DSetPropPtr(). The element local coordinate system is constructed so that the x’ direction is along the specified normal to the master surface. The y’ direction is perpendicular to x’.

In thermal analysis the geometry of the slave surface is ignored and only the temperatures of the slave surface nodes are used. For thermal analysis the element coordinate system is unused and is undefined.

../../_images/vfetools-inter1a.gif

Figure 7-1a, Geometry for 6 Node Surface 2D Interface Element

The implicit surface interface element is only available for thermal analysis. In this case, since the master surface is defined as a single node, the element node area factor is interpreted as an effective thickness at the node for the calculation of surface area.

../../_images/vfetools-inter1b.gif

Figure 7-1b, Geometry for 2 Node Line 2D Interface Element

7.3. Function Descriptions

vfe_Inter2D *vfe_Inter2DBegin(void)

create an instance of an Inter2D object

Create an instance of a Inter2D object. Memory is allocated for the object private data and the pointer to the object is returned. Default topology is the 4-node per surface quadrilateral interface

Destroy an instance of a Inter2D object using

void vfe_Inter2DEnd (vfe_Inter2D *inter2d)

Return the current value of a Inter2D object error flag using

Vint vfe_Inter2DError (vfe_Inter2D *inter2d)

Returns

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

void vfe_Inter2DEnd(vfe_Inter2D *p)

destroy an instance of an Inter2D object

See vfe_Inter2DBegin()

Vint vfe_Inter2DError(vfe_Inter2D *p)

return the current value of an Inter2D object error flag

See vfe_Inter2DBegin()

void vfe_Inter2DSetObject(vfe_Inter2D *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 InterProp module provides a suitable material model.

Errors

SYS_ERROR_OBJECTTYPE is generated if an improper objecttype is specified.

Parameters
  • p – Pointer to Inter2D object.

  • objecttype – The object type identifier

    x=VFE_MATLFUN            MatlFun object
    

  • object – Pointer to the object to be set.

void vfe_Inter2DSetTopology(vfe_Inter2D *p, Vint shape, Vint maxi, Vint maxj)

set element topology

Specify the topology of a 2D interface element. If maxi is set to 3 then a quadratic element form is assumed. The default topology is SYS_SHAPEQUAD with maxi = 0 and maxj = 2.

Errors

  • SYS_ERROR_VALUE is generated if an improper maxi or maxj is specified.

  • SYS_ERROR_ENUM is generated if an improper shape is input.

Parameters
  • p – Pointer to Inter2D object.

  • shape – The topological shape for the element

    x=VIS_SHAPEQUAD          Quadrilateral
     =VIS_SHAPELINE          Line
    

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

  • maxj – The number of points along the j direction. For quadrilateral shapes maxj = 2. For line shapes maxj = 0.

void vfe_Inter2DSetParami(vfe_Inter2D *p, Vint type, Vint iparam)

set element formulation parameters

Set element formulation parameters. Set element heat exchange law using VFE_HEATEXCH. By default VFE_HEATEXCH is set to VFE_HEATEXCH_LINEAR.

The parameter VFE_GEOMCLEARANCE is used to toggle the inclusion of the actual geometric location of the slave nodes in the calculation of the clearance between the master and slave surfaces. If VFE_GEOMCLEARANCE is disabled then the geometry of the slave surface is ignored and the initial clearance between the master and slave surfaces is determined by the element node clearances set using vfe_Inter2DSetPropPtr(). If VFE_GEOMCLEARANCE is enabled then the actual geometric separation of the master and slave surfaces is added to the specified element node clearances. By default VFE_GEOMCLEARANCE is set to SYS_OFF.

The parameter VFE_2D sets the particular 2D approximation for the element. The choices include plane stress, VFE_PLANESTRESS, plane strain, VFE_PLANESTRAIN and axisymmetric, VFE_AXISYMMETRIC. By default VFE_2D is set to VFE_PLANESTRESS.

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

  • type – Type of formulation parameter to set

    x=VFE_HEATEXCH           Element exchange law
     =VFE_GEOMCLEARANCE      Clearance determined by geometry
     =VFE_2D                 2D Approximation
    

  • iparam – Integer parameter value.

    x=VFE_HEATEXCH_NONE      No exchange
     =VFE_HEATEXCH_LINEAR    Linear exchange law
     =VFE_HEATEXCH_EXP       Exponential exchange law
     =VFE_HEATEXCH_POW       Power exchange law
     =VFE_HEATEXCH_RAD       Black body radiation law
     =VFE_AXISYMMETRIC       2D axisymmetric
     =VFE_PLANESTRAIN        2D plane strain
     =VFE_PLANESTRESS        2D plane stress
    

void vfe_Inter2DSetParamd(vfe_Inter2D *p, Vint type, Vdouble dparam)

set double formulation parameters

Set double formulation parameters. Set element heat exchange exponent law using VFE_HEATEXP. This value is used if VFE_HEATEXCH is set to VFE_HEATEXCH_POW. Defaults to 1.0.

Set the absolute zero temperature for radiation exchange. Defaults 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 Inter2D object.

  • type – Type of formulation parameter to set

    x=VFE_HEATEXP            Heat exchange exponent
     =VFE_ABSZEROTEMP        Absolute zero temperature
    

  • dparam – Double parameter value.

void vfe_Inter2DStiff(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom stiffness matrix

void vfe_Inter2DStrsStrn(vfe_Inter2D *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 interface local system.

The strs is composed of 3 contact pressures and the strn is composed of 3 associated relative displacements in the element local system. The z’ values are returned as 0.

The strs and strn values are ordered first by the 3 vectoral 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.

Parameters
  • p – Pointer to Inter2D 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_Inter2DReact(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • r[out] Degree of freedom reaction vector

void vfe_Inter2DReactStiff(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D 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_Inter2DNumDof(vfe_Inter2D *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_Inter2DDofMap() 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 Inter2D 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_Inter2DDofMap(vfe_Inter2D *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_Inter2DNumDof() 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 Inter2D 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_Inter2DInitHist(vfe_Inter2D *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 Inter2D object.

void vfe_Inter2DNumIntPnt(vfe_Inter2D *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 Inter2D 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_Inter2DDirCos(vfe_Inter2D *p, Vdouble x[][3], Vdouble tm[][3][3])

compute interface 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_Inter2DStrsStrn() are expressed. Given that X’ and Y’ are two 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] = 0.
 tm[1][0] = Y'x  tm[1][1] = Y'y  tm[1][2] = 0.
 tm[2][0] = 0.   tm[2][1] = 0.   tm[2][2] = 1.
The local coordinate system is determined by the local system convention and the element node normals which may be optionally set using vfe_Inter2DSetPropPtr().

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

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

void vfe_Inter2DSetPropPtr(vfe_Inter2D *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 for the associated property. By default the temperature is 0., the depth is 1., the area factor is 1., the initial clearance is 0. and the normal is assumed to be computed from the master surface geometry.

Errors

SYS_ERROR_ENUM is generated if an improper type is specified.

Parameters
  • p – Pointer to Inter2D object.

  • type – Type of element property

    x=VFE_PROP_AREAFACT         Area Factor
     =VFE_PROP_DEPTH            Depth
     =VFE_PROP_TEMPERATURE      Temperatures
     =VFE_PROP_NORMAL           Normals
     =VFE_PROP_CLEARANCE        Initial clearance
    

  • propptr – Pointer to start of element nodal properties

void vfe_Inter2DSetHistPtr(vfe_Inter2D *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_Inter2DNumIntPnt(). 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 Inter2D object.

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

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

void vfe_Inter2DPowerCond(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D 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_Inter2DPower(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • r[out] Degree of freedom power vector

void vfe_Inter2DCond(vfe_Inter2D *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.

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom conductance matrix

void vfe_Inter2DHFlxTGrd(vfe_Inter2D *p, Vdouble x[][3], Vdouble u[], Vdouble hflx[], Vdouble tgrd[])

Heat flux and temperature difference.

Compute nodal heat fluxes and temperature differences, hflx and tgrd, given the node coordinates, x, and the degree of freedom temperature vector, u.

The hflx is a scalar flux and the tgrd is the associated temperature difference across the interface at each element node.

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 Inter2D 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 differences

void vfe_Inter2DDistHeat(vfe_Inter2D *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 on a 2D interface element edge or node. The vector q contains values for the heat flux for each node in the element. The distributed flux for both edges and nodes are in units of flux per unit area.

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

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_NODE     Element node
    =SYS_EDGE     Element edge
    

  • no – Element node or edge number

  • q – Vector of distributed load values

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

void vfe_Inter2DDistLoad(vfe_Inter2D *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 2D interface element face. 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 loadtype is VFE_DISTLOAD_PRES then q contains a scalar pressure 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 face Jacobian is computed.

Parameters
  • p – Pointer to Inter2D object.

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_EDGE     Element edge
    

  • no – Element edge number

    =1            Master face
    =3            Slave face
    

  • loadtype – Distributed load type

    x=VFE_DISTLOAD_PRES      Pressure normal to face
     =VFE_DISTLOAD_TRAC      Traction directed along vector
    

  • q – Vector of distributed load values

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

7.4. 3D Interface Elements - Inter3D

The Inter3D module is used to simulate surface interactions in 3D space. The interface element may be used in structural and thermal analysis.

The methods associated with a Inter3D object are the following.

Instance a Inter3D object using vfe_Inter3DBegin(). Once a Inter3D is instanced, set the material function attribute object MatlFun using vfe_Inter3DSetObject(). The current topology of the element is specified using vfe_Inter3DSetTopology(). Certain geometric details of the interface such as initial gap separation are specified using vfe_Inter3DSetPropPtr(). Query the element degree of freedom map using vfe_Inter3DDofMap() and vfe_Inter3DNumDof().

If the element is to support a material model, such as a surface friction material model, which requires a material history then the user must manage the material history information using vfe_Inter3DSetHistPtr() and vfe_Inter3DInitHist().

7.5. Element Geometry

The interface element consists of two matched surfaces termed the master surface and the slave surface. The basic interface element geometry is defined by the node coordinates of the the master surface and possibly the coordinates of the slave surface. The separation of the two surfaces is defined by normals and clearances specified at the nodes of the master surface. These properties are set using vfe_Inter3DSetPropPtr(). Other element node properties which may be set are temperature and area factor.

The master and slave surfaces may be either explicitly defined as surfaces in terms of the nodal geometry as shown in Figure 7-2a or implicitly as lines in terms of the nodal geometry and associated nodal area factors as shown in Figure 7-2b. If the interface element is being used to model surface contact in structural analysis, the separation of the master and slave surfaces is always measured along the normal to the master surface. The master surface normal may be determined from the nodal geometry or specified at the master surface nodes using vfe_Inter3DSetPropPtr(). The element local coordinate system is constructed so that the x’ direction is along the specified normal to the master surface. The y’ and z’ directions are constructed using the prescription specified using vfe_Inter3DSetLocalSystem(). If the normal, x’ direction, is specified using vfe_Inter3DSetPropPtr(), the y’ and z’ tangential directions are adjusted to remain orthogonal to x’.

In thermal analysis the geometry of the slave surface is ignored and only the temperatures of the slave surface nodes are used. For thermal analysis the element coordinate system is unused and is undefined.

../../_images/vfetools-inter2a.gif

Figure 7-2a, Geometry for 16 Node Surface 3D Interface Element

The implicit surface interface element is only available for thermal analysis. In this case, since the master surface is defined as a line of nodes, the element node area factor is interpreted as an effective thickness at each node for the calculation of surface area.

../../_images/vfetools-inter2b.gif

Figure 7-2b, Geometry for 6 Node Line 3D Interface Element

7.6. Function Descriptions

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

vfe_Inter3D *vfe_Inter3DBegin(void)

create an instance of an Inter3D object

Create an instance of a Inter3D object. Memory is allocated for the object private data and the pointer to the object is returned. Default topology is the 4-node per surface quadrilateral interface

Destroy an instance of a Inter3D object using

void vfe_Inter3DEnd (vfe_Inter3D *inter3d)

Return the current value of a Inter3D object error flag using

Vint vfe_Inter3DError (vfe_Inter3D *inter3d)

Returns

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

void vfe_Inter3DEnd(vfe_Inter3D *p)

destroy an instance of an Inter3D object

See vfe_Inter3DBegin()

Vint vfe_Inter3DError(vfe_Inter3D *p)

return the current value of an Inter3D object error flag

See vfe_Inter3DBegin()

void vfe_Inter3DSetObject(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D object.

  • objecttype – The object type identifier

    x=VFE_MATLFUN            MatlFun object
    

  • object – Pointer to the object to be set.

void vfe_Inter3DSetTopology(vfe_Inter3D *p, Vint shape, Vint maxi, Vint maxj, Vint maxk)

set element topology

Specify the topology of a 3D interface element. If maxi is set to 3 then a quadratic element form is assumed. The default topology is SYS_SHAPEHEX with maxi = maxj = 0 and maxk = 2.

Errors

  • SYS_ERROR_VALUE is generated if an improper maxi, maxj or maxk is specified.

  • SYS_ERROR_ENUM is generated if an improper shape is input.

Parameters
  • p – Pointer to Inter3D object.

  • shape – The topological shape for the element

    x=VIS_SHAPEWED           Wedge
     =VIS_SHAPEHEX           Hexahedron
     =VIS_SHAPEQUAD          Quadrilateral
    

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

  • maxj – The number of points along the j direction. If maxj = 0 then a Serendipity finite element is assumed. If 2 <= maxj <= 3 and 2 <= maxi <= 3, then a Lagrange finite element is assumed. If maxi != 0 then maxj = 0 or maxj = maxi. For quadrilateral shapes maxj = 2.

  • maxk – The number of points along the k direction. For wedge and hexahedron shapes maxk = 2. For quadrilateral shapes maxk is ignored.

void vfe_Inter3DSetParami(vfe_Inter3D *p, Vint type, Vint iparam)

set element formulation parameters

Set element formulation parameters. Set element heat exchange law using VFE_HEATEXCH. By default VFE_HEATEXCH is set to VFE_HEATEXCH_LINEAR.

The parameter VFE_GEOMCLEARANCE is used to toggle the inclusion of the actual geometric location of the slave nodes in the calculation of the clearance between the master and slave surfaces. If VFE_GEOMCLEARANCE is disabled then the geometry of the slave surface is ignored and the initial clearance between the master and slave surfaces is determined by the element node clearances set using vfe_Inter3DSetPropPtr(). If VFE_GEOMCLEARANCE is enabled then the actual geometric separation of the master and slave surfaces is added to the specified element node clearances. By default VFE_GEOMCLEARANCE is set to SYS_OFF.

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

  • type – Type of formulation parameter to set

    x=VFE_HEATEXCH           Element exchange law
     =VFE_GEOMCLEARANCE      Clearance determined by geometry
    

  • iparam – Integer parameter value.

    x=VFE_HEATEXCH_NONE      No exchange
     =VFE_HEATEXCH_LINEAR    Linear exchange law
     =VFE_HEATEXCH_EXP       Exponential exchange law
     =VFE_HEATEXCH_POW       Power exchange law
     =VFE_HEATEXCH_RAD       Black body radiation law
    

void vfe_Inter3DSetParamd(vfe_Inter3D *p, Vint type, Vdouble dparam)

set double formulation parameters

Set double formulation parameters. Set element heat exchange exponent law using VFE_HEATEXP. This value is used if VFE_HEATEXCH is set to VFE_HEATEXCH_POW. Defaults to 1.0.

Set the tolerance for clamping the natural coordinates when performing projections. Defaults to 0.1.

Set the absolute zero temperature for radiation exchange. Defaults 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 Inter3D object.

  • type – Type of formulation parameter to set

    x=VFE_HEATEXP       Heat exchange exponent
     =VFE_EXTENSION     Tolerance for natural coordinates clamping
     =VFE_ABSZEROTEMP   Absolute zero temperature
    

  • dparam[in] Double parameter value.

void vfe_Inter3DStiff(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom stiffness matrix

void vfe_Inter3DStrsStrn(vfe_Inter3D *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 interface local system. The convention used to generate local coordinate systems is specified using vfe_Inter3DSetLocalSystem() with an optional surface normal specified using vfe_Inter3DSetPropPtr().

The strs is composed of 3 contact pressures and the strn is composed of 3 associated relative displacements in the element local system.

The strs and strn values are ordered first by the 3 vectoral 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.

Parameters
  • p – Pointer to Inter3D 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_Inter3DReact(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • r[out] Degree of freedom reaction vector

void vfe_Inter3DReactStiff(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D 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_Inter3DNumDof(vfe_Inter3D *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_Inter3DDofMap() 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 Inter3D 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_Inter3DDofMap(vfe_Inter3D *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_Inter3DNumDof() 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 Inter3D 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_Inter3DInitHist(vfe_Inter3D *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 Inter3D object.

void vfe_Inter3DNumIntPnt(vfe_Inter3D *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 Inter3D 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_Inter3DSetLocalSystem(vfe_Inter3D *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 interface surface.

For stress and strain computation for output using vfe_Inter3DStrsStrn(), 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 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 Inter3D object.

  • type – Local system convention

  • vec – Orientation vector data

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

void vfe_Inter3DSetMatlSystem(vfe_Inter3D *p, Vint type, Vdouble vec[], Vdouble angle)

set material coordinate system convention

Specify the convention to be used to construct the orientation of the element material x’,y’,z’ coordinate system with respect to the interface surface. This material system is computed at each integration point location on the interface and is assumed to be the coordinate system in which the material properties of the element at the integration point are expressed.

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

  • type – Local system convention

  • vec – Orientation vector data

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

void vfe_Inter3DDirCos(vfe_Inter3D *p, Vdouble x[][3], Vdouble tm[][3][3])

compute interface 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_Inter3DStrsStrn() 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_Inter3DSetLocalSystem() and the element node normals which may be optionally set using vfe_Inter3DSetPropPtr().

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

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

void vfe_Inter3DSetPropPtr(vfe_Inter3D *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 for the associated property. By default the temperature is 0., the area factor is 1., the initial clearance is 0. and the normal is assumed to be computed from the master surface geometry.

Errors

SYS_ERROR_ENUM is generated if an improper type is specified.

Parameters
  • p – Pointer to Inter3D object.

  • type – Type of element property

    x=VFE_PROP_AREAFACT         Area Factor
     =VFE_PROP_TEMPERATURE      Temperatures
     =VFE_PROP_NORMAL           Normals
     =VFE_PROP_CLEARANCE        Initial clearance
    

  • propptr – Pointer to start of element nodal properties

void vfe_Inter3DSetHistPtr(vfe_Inter3D *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_Inter3DNumIntPnt(). 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 Inter3D object.

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

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

void vfe_Inter3DPowerCond(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D 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_Inter3DPower(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • r[out] Degree of freedom power vector

void vfe_Inter3DCond(vfe_Inter3D *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.

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

  • kl[out] Degree of freedom conductance matrix

void vfe_Inter3DHFlxTGrd(vfe_Inter3D *p, Vdouble x[][3], Vdouble u[], Vdouble hflx[], Vdouble tgrd[])

Heat flux and temperature difference.

Compute nodal heat fluxes and temperature differences, hflx and tgrd, given the node coordinates, x, and the degree of freedom temperature vector, u.

The hflx is a scalar flux and the tgrd is the associated temperature difference across the interface at each element node.

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 Inter3D 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 differences

void vfe_Inter3DDistLoad(vfe_Inter3D *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 interface element face. 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 loadtype is VFE_DISTLOAD_PRES then q contains a scalar pressure 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 face Jacobian is computed.

Parameters
  • p – Pointer to Inter3D object.

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_FACE     Element face
    

  • no – Element face number

    =1            Master face
    

  • loadtype – Distributed load type

    x=VFE_DISTLOAD_PRES      Pressure normal to face
     =VFE_DISTLOAD_TRAC      Traction directed along vector
    

  • q – Vector of distributed load values

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

void vfe_Inter3DDistHeat(vfe_Inter3D *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 on a 3D interface element face or edge. The vector q contains values for the heat flux for each node in the element. The distributed flux for faces and edges are in units of flux per unit area.

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

  • x – Array of node locations.

  • enttype – Entity type on which load is applied

    =SYS_EDGE     Element edge
    =SYS_FACE     Element face
    

  • no – Element edge or face number

    =1            Master edge or face
    

  • q – Vector of distributed load values

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

7.7. Gap Elements - Gap

The Gap module is used to simulate simple surface interactions in 2D or 3D space. The master and slave surfaces are assumed to be infinite, passing through the gap element end points and and perpendicular to the specified normal at the first end point. The gap element may be used in structural and thermal analysis.

The methods associated with a Gap object are the following.

Instance a Gap object using vfe_GapBegin(). Once a Gap is instanced, set the material function attribute object MatlFun using vfe_GapSetObject(). The function vfe_GapDef() is used to defined the gap as either being used in 2D or 3D analysis.

Certain geometric details of the gap such as initial gap separation are specified using vfe_GapSetPropPtr(). Query the element degree of freedom map using vfe_GapDofMap() and vfe_GapNumDof().

If the element is to support a material model, such as a surface friction material model, which requires a material history then the user must manage the material history information using vfe_GapSetHistPtr() and vfe_GapInitHist().

7.8. Element Geometry

The gap element consists of two matched surfaces termed the master surface and the slave surface. The basic gap element geometry is defined by the node coordinates of the the first node and possibly the coordinates of the second node. The separation of the two surfaces is defined by normals and clearances specified at the first of the master surface. These properties are set using vfe_GapSetPropPtr(). Other element node properties which may be set are temperature and area.

The master surface passes through the first node and is perpendicular to the normal at the first node. The normal may be computed to be along the line connecting the first and second node or may be explicitly specified using vfe_GapSetPropPtr(). The slave surface passes through the second node and is parallel to the master surface. Both surfaces are assumed to be of infinite extent for contact purposes. The 3D gap element geometry is shown in Figure 7-3a and the 2D gap geometry is shown in Figure 7-3b.

If the gap element is being used to model surface contact in structural analysis, the separation of the master and slave surfaces is always measured along the normal to the master surface. The element local coordinate system is constructed so that the x’ direction is along the specified normal to the master surface. The y’ and z’ directions are constructed using the prescription specified using vfe_GapSetLocalSystem(). If the normal, x’ direction, is specified using vfe_GapSetPropPtr(), the y’ and z’ tangential directions are adjusted to remain orthogonal to x’.

In thermal analysis the geometry of the slave surface is ignored and only the temperature of the slave surface node is used. The element node area is used at the master node for the calculation of surface area. For thermal analysis the element coordinate system is unused and is undefined.

../../_images/vfetools-gap3a.gif

Figure 7-3a, Geometry for 3D Gap Element

In 2D, the gap element may be used for either 2D axisymmetric or 2D planar analysis.

../../_images/vfetools-gap3b.gif

Figure 7-3b, Geometry for 2D Gap Element

7.9. Function Descriptions

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

vfe_Gap *vfe_GapBegin(void)

create an instance of a Gap object

Create an instance of a Gap object. Memory is allocated for the object private data and the pointer to the object is returned.

Destroy an instance of a Gap object using

void vfe_GapEnd (vfe_Gap *gap)

Return the current value of a Gap object error flag using

Vint vfe_GapError (vfe_Gap *gap)

Returns

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

void vfe_GapEnd(vfe_Gap *p)

destroy an instance of a Gap object

See vfe_GapBegin()

Vint vfe_GapError(vfe_Gap *p)

return the current value of a Gap object error flag

See vfe_GapBegin()

void vfe_GapSetObject(vfe_Gap *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 InterProp module provides a suitable material model.

Errors

SYS_ERROR_OBJECTTYPE is generated if an improper objecttype is specified.

Parameters
  • p – Pointer to Gap object.

  • objecttype – The object type identifier

    x=VFE_MATLFUN            MatlFun object
    

  • object[out] Pointer to the object to be set.

void vfe_GapDef(vfe_Gap *p, Vint dimension)

define gap dimensionality

Define the dimension of the gap element.

Errors

SYS_ERROR_ENUM is generated if an improper dime is specified.

Parameters
  • p – Pointer to Gap object.

  • dimension – Gap dimension

    =VFE_2D       2D analysis
    =VFE_3D       3D analysis
    

void vfe_GapSetParami(vfe_Gap *p, Vint type, Vint iparam)

set element formulation parameters

Set element formulation parameters. Set element heat exchange law using VFE_HEATEXCH. By default VFE_HEATEXCH is set to VFE_HEATEXCH_LINEAR.

The parameter VFE_GEOMCLEARANCE is used to toggle the inclusion of the actual geometric location of the slave nodes in the calculation of the clearance between the master and slave surfaces. If VFE_GEOMCLEARANCE is disabled then the geometry of the slave surface is ignored and the initial clearance between the master and slave surfaces is determined by the element node clearances set using vfe_GapSetPropPtr(). If VFE_GEOMCLEARANCE is enabled then the actual geometric separation of the master and slave surfaces is added to the specified element node clearances. By default VFE_GEOMCLEARANCE is set to SYS_OFF.

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

  • type – Type of formulation parameter to set

    x=VFE_HEATEXCH           Element exchange law
     =VFE_GEOMCLEARANCE      Clearance determined by geometry
    

  • iparam – Integer parameter value.

    x=VFE_HEATEXCH_NONE      No exchange
     =VFE_HEATEXCH_LINEAR    Linear exchange law
     =VFE_HEATEXCH_EXP       Exponential exchange law
     =VFE_HEATEXCH_POW       Power exchange law
    

void vfe_GapSetParamd(vfe_Gap *p, Vint type, Vdouble dparam)

set double formulation parameters

Set element heat exchange exponent law using VFE_HEATEXP. This value is used if VFE_HEATEXCH is set to VFE_HEATEXCH_POW. Defaults to 1.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 Gap object.

  • type – Type of formulation parameter to set

    x=VFE_HEATEXP            Heat exchange exponent
    

  • dparam[out] Double parameter value.

void vfe_GapStiff(vfe_Gap *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.

Parameters
  • p – Pointer to Gap object.

  • x – Array of node locations.

  • kl[out] Degree of freedom stiffness matrix

void vfe_GapStrsStrn(vfe_Gap *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 gap local system. The convention used to generate local coordinate systems is specified using vfe_GapSetLocalSystem() with an optional surface normal specified using vfe_GapSetPropPtr().

The strs is composed of 3 contact pressures and the strn is composed of 3 associated relative displacements in the element local system.

The strs and strn values are ordered first by the 3 vectoral 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.

Parameters
  • p – Pointer to Gap 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_GapReact(vfe_Gap *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

Parameters
  • p – Pointer to Gap object.

  • x – Array of node locations.

  • u – Degree of freedom vector of displacements

  • r[out] Degree of freedom reaction vector

void vfe_GapReactStiff(vfe_Gap *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.

Parameters
  • p – Pointer to Gap 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_GapNumDof(vfe_Gap *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_GapDofMap() 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 Gap 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_GapDofMap(vfe_Gap *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_GapNumDof() 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 Gap 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_GapInitHist(vfe_Gap *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 Gap object.

void vfe_GapNumIntPnt(vfe_Gap *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 Gap 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_GapDirCos(vfe_Gap *p, Vdouble x[][3], Vdouble tm[][3][3])

compute gap 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_GapStrsStrn() 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_GapSetLocalSystem() and the element node normals which may be optionally set using vfe_GapSetPropPtr().

Parameters
  • p – Pointer to Gap object.

  • x – Array of node locations.

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

void vfe_GapSetPropPtr(vfe_Gap *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., the area is 1., the initial clearance is 0. and the normal is assumed to be computed from the gap endpoint locations.

Errors

SYS_ERROR_ENUM is generated if an improper type is specified.

Parameters
  • p – Pointer to Gap object.

  • type – Type of element property

    x=VFE_PROP_AREA             Area
     =VFE_PROP_TEMPERATURE      Temperatures
     =VFE_PROP_NORMAL           Normals
     =VFE_PROP_CLEARANCE        Initial clearance
    

  • propptr – Pointer to start of element nodal properties

void vfe_GapSetHistPtr(vfe_Gap *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_GapNumIntPnt(). 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 Gap object.

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

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

void vfe_GapPowerCond(vfe_Gap *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.

Parameters
  • p – Pointer to Gap 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_GapPower(vfe_Gap *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.

Parameters
  • p – Pointer to Gap object.

  • x – Array of node locations.

  • u – Degree of freedom vector of temperatures

  • r[out] Degree of freedom power vector

void vfe_GapCond(vfe_Gap *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.

Parameters
  • p – Pointer to Gap object.

  • x – Array of node locations.

  • kl[out] Degree of freedom conductance matrix

void vfe_GapHFlxTGrd(vfe_Gap *p, Vdouble x[][3], Vdouble u[], Vdouble hflx[], Vdouble tgrd[])

Heat flux and temperature difference.

Compute nodal heat fluxes and temperature differences, hflx and tgrd, given the node coordinates, x, and the degree of freedom temperature vector, u.

The hflx is a scalar flux and the tgrd is the associated temperature difference across the gap at each element node.

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 Gap 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 difference

void vfe_GapSetLocalSystem(vfe_Gap *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.

For stress and strain computation for output using vfe_GapStrsStrn(), 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 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 Gap object.

  • type – Local system convention

  • vec – Orientation vector data

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