2. Solid Elements - Solid2D, Solid3D
Solid elements are the basic volumetric elements in VfeTools. There are two modules for the simulation of volumetric, continuum based finite elements: Solid2D for plane stress, plane strain, and axisymmetric 2D analysis and Solid3D for full 3D analysis.
2.1. 2D Solid Elements - Solid2D
The Solid2D module is used to simulate 2D continuum. It can be used for any planar or axisymmetric simulation. In the latter case, the coordinate system utilized is r-z-theta, where r is the radial direction, z is the axial direction and theta is the circumferential direction. All coordinates are input to Solid2D functions as a three component position vector, in all cases the third component is ignored.
In plane stress and plane strain structural analysis, all symmetric tensors are given in the order xx, yy, zz, xy, yz, zx. Note that in plane strain, the zz, yz, and zx components of strain are identically zero, while the yz and the zx components of stress are zero. In plane stress, yz and zx components of strain are zero, while the zz, yz, and zx components of stress are zero. In the axisymmetric case, the tensor order becomes rr, zz, theta-theta, rz, z-theta, theta-r.
In planar thermal analysis, all vectors are given in the order x, y, z where the z component of the temperature gradient is zero. In the axisymmetric case, the vectoral order becomes r, z, theta where the theta component of the temperature gradient is zero.
The methods associated with a Solid2D object are the following.
- Begin and end an instance of an object, generic object functions - vfe_Solid2DBegin()- create an instance of a Solid2D object
- vfe_Solid2DEnd()- destroy an instance of a Solid2D object
- vfe_Solid2DError()- return Solid2D object error flag
 
- Attributes and Parameters - vfe_Solid2DSetHistPtr()- set pointers to material history
- vfe_Solid2DSetMatlSystem()- set material axes direction
- vfe_Solid2DSetObject()- set attribute object
- vfe_Solid2DSetParami()- set element formulation parameters
- vfe_Solid2DSetPropPtr()- set pointer to element nodal properties
- vfe_Solid2DSetTopology()- set input element topology
 
- Degree of Freedom and Integration Information - vfe_Solid2DDofMap()- query element degree of freedom map
- vfe_Solid2DNumDof()- query number of element degrees of freedom
- vfe_Solid2DNumIntPnt()- query number of element integration points
 
- General - vfe_Solid2DShapeGrad()- shape functions and gradients
 
- Structural analysis computations - vfe_Solid2DConcLoad()- concentrated load vector
- vfe_Solid2DDistLoad()- distributed load vector.
- vfe_Solid2DElemLoad()- body force vector
- vfe_Solid2DGeomStiff()- geometric stiffness matrix
- vfe_Solid2DInitHist()- initialize material history
- vfe_Solid2DMass()- consistent mass matrix
- vfe_Solid2DMassDiag()- diagonal mass matrix
- vfe_Solid2DReact()- reaction vector
- vfe_Solid2DReactStiff()- reaction vector, stiffness matrix
- vfe_Solid2DStiff()- linear stiffness matrix
- vfe_Solid2DStrsAdapt()- stress based error analysis
- vfe_Solid2DStrsStrn()- stress and strain
 
- Thermal analysis computations - vfe_Solid2DDistHeat()- distributed heat loads.
- vfe_Solid2DElemHeat()- body heat generation
- vfe_Solid2DCap()- consistent capacitance matrix
- vfe_Solid2DCapDiag()- diagonal capacitance matrix
- vfe_Solid2DPower()- thermal power
- vfe_Solid2DPowerCond()- thermal power, conductance matrix
- vfe_Solid2DCond()- conductance matrix
- vfe_Solid2DHFlxAdapt()- heat flux based error analysis
- vfe_Solid2DHFlxTGrd()- heat flux and temperature gradient
 
Instance a Solid2D object using vfe_Solid2DBegin().
Once a Solid2D is instanced,
set the material function attribute object MatlFun using
vfe_Solid2DSetObject().
The current topology of the element is specified using
vfe_Solid2DSetTopology().
The user may select an element formulation using vfe_Solid2DSetParami().
If the solid element is used in plane stress or plane strain analysis,
the depth of the element may be specified using vfe_Solid2DSetPropPtr().
Query the element degree of freedom map using vfe_Solid2DDofMap() and
vfe_Solid2DNumDof().
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_Solid2DSetHistPtr() and vfe_Solid2DInitHist().
Use the function
vfe_Solid2DStrsAdapt() 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_Solid2DHFlxAdapt() performs a similar
computation for heat transfer analysis.
It is useful to use reduced integration rules for certain parabolic 2D and 3D solids for performance reasons. In general reduced integration rules can introduce zero energy modes within a single element. These zero energy modes are suppressed in most practical circumstances by attachments to adjacent elements or constraints across the element edges or faces. The number of adjacent elements required to suppress the zero energy modes for each Serendipity element type is listed below.
The element technology to be used is set using the
function vfe_Solid2DSetParami() with type VFE_TECH.
Standard fully integrated, VFE_TECH_ISOP,
and underintegrated, VFE_TECH_URED exist for all elements.
The default technology is VFE_TECH_ISOP which will work in all cases.
However for some element types, alternative technologies should be
considered for improved performance.
The recommended technology varies with element topology, order and
2D approximation and is listed below.
The number of integration points used for fully
integrated technologies appears in paraenthesis.  Underintegrated
technologies use a single integration point for the linear element
forms and for the higher order element forms use the fully integrated
rule for the one-less order.
- Linear triangle axisymmertic VFE_TECH_ISOP (3), planar VFE_TECH_ISOP (1). 
- Parabolic triangle axisymmertic VFE_TECH_ISOP (6), planar VFE_TECH_ISOP (3). 
- Cubic triangle axisymmertic VFE_TECH_ISOP (7), planar VFE_TECH_ISOP (6). 
 
- Linear quadrilateral, VFE_TECH_MIXED or VFE_TECH_ENHANCED. Use VFE_TECH_MIXED unless bending deformations are particularly important with a coarse mesh. The fully integrated technology, VFE_TECH_ISOP (4), is overly stiff. 
- Parabolic quadrilateral axisymetric, VFE_TECH_URED The fully integrated technology, VFE_TECH_ISOP (9), is overly stiff. planar, if Lagrangian VFE_TECH_ISOP, if Serendipity VFE_TECH_URED if adjacent elements exist across 3 or more edges else VFE_TECH_ISOP. 
- Cubic quadrilateral axisymetric, VFE_TECH_ISOP (16) 
2.2. Function Descriptions
The currently available Solid2D functions are described in detail in this section.
- 
vfe_Solid2D *vfe_Solid2DBegin(void)
- create an instance of a Solid2D object - Create an instance of a Solid2D object. Memory is allocated for the object private data and the pointer to the object is returned. By default the element assumes plane stress approximations, with a unit depth. Default topology is the 4-noded Serendipity quadrilateral with isoparametric technology. - Destroy an instance of a Solid2D object using - void vfe_Solid2DEnd (vfe_Solid2D *solid2d) - Return the current value of a Solid2D object error flag using - Vint vfe_Solid2DError (vfe_Solid2D *solid2d) - Returns
- The function returns a pointer to the newly created Solid2D object. If the object creation fails, NULL is returned. 
 
- 
void vfe_Solid2DEnd(vfe_Solid2D *p)
- destroy an instance of a Solid2D object 
- 
Vint vfe_Solid2DError(vfe_Solid2D *p)
- return the current value of a Solid2D object error flag 
- 
void vfe_Solid2DSetObject(vfe_Solid2D *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_OBJECTTYPEis generated if an improper objecttype is specified.
 - Parameters
- p – Pointer to Solid2d object. 
- objecttype – The object type identifier - x=VFE_MATLFUN MatlFun object 
- object – Pointer to the object to be set. 
 
 
- 
void vfe_Solid2DSetParami(vfe_Solid2D *p, Vint type, Vint iparam)
- set element formulation parameters - Set element formulation technology parameters. The basic element technology is set with the - VFE_TECHparameter. By default- VFE_TECHis set to- VFE_TECH_ISOP.- The parameter - VFE_CHECKTOTALVOLtoggles of Jacobian checking performed during element volume integrations. If this parameter is not enabled, then an error is generated if the Jacobian contribution of any integration point is non-positive. If this parameter is enabled then an error is generated only if the total Jacobian contribution over all integration points is non-positive. By default- VFE_CHECKTOTALVOLis set to- SYS_OFF.- Set element strain type using - VFE_STRAINTYPEwith a value of either- VFE_LARGESTRAINto enable large strain or- VFE_SMALLSTRAINto enable small strains. By default- VFE_STRAINTYPEis set to- VFE_SMALLSTRAIN.- The parameter - VFE_2Dsets the particular 2D approximation for the element. The choices include plane stress,- VFE_PLANESTRESS, plane strain,- VFE_PLANESTRAINand axisymmetric,- VFE_AXISYMMETRIC. By default- VFE_2Dis set to- VFE_PLANESTRESS.- The parameter - VFE_TEMPMATLAVGtoggles 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_TEMPMATLAVGis set to- SYS_ON.- Errors
- SYS_ERROR_ENUMis generated if an improper type is specified.
- SYS_ERROR_VALUEis generated if an improper iparam is specified.
 
 - Parameters
- p – Pointer to Solid2d object. 
- type – Type of formulation parameter to set - x=VFE_TECH Element technology =VFE_CHECKTOTALVOL Total Jacobian checking =VFE_STRAINTYPE Element strain type =VFE_TEMPMATLAVE Average material temperature flag =VFE_2D 2D Approximation 
- iparam – Integer parameter value. - x=SYS_ON Enable =SYS_OFF Disable =VFE_TECH_ISOP Standard isoparametric technology =VFE_TECH_URED Uniform reduced technology =VFE_TECH_MIXED Mixed technology =VFE_TECH_ENHANCED Enhanced technology =VFE_AXISYMMETRIC 2D axisymmetric =VFE_PLANESTRAIN 2D plane strain =VFE_PLANESTRESS 2D plane stress =VFE_SMALLSTRAIN Small strain =VFE_LARGESTRAIN Large strain 
 
 
- 
void vfe_Solid2DSetTopology(vfe_Solid2D *p, Vint shape, Vint maxi, Vint maxj)
- set element topology - Specify the topology of a 2D solid element. If maxi is set to 3 then a quadratic element form is assumed. If maxi is set to 4 then a cubic element form is assumed. The default topology is - SYS_SHAPEQUADwith maxi = maxj = 0. Serendipity elements with missing mid-side nodes are also supported. See Computational Cells for a description on how to define element topologies with missing mid-side nodes.- Errors
- SYS_ERROR_VALUEis generated if an improper maxi or maxj is specified.
- SYS_ERROR_ENUMis generated if an improper shape is input.
 
 - Parameters
- p – Pointer to Solid2d object. 
- shape – The topological shape for the element - x=SYS_SHAPETRI Triangle =SYS_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 <= 4 and 2 <= maxi <= 4, then a Lagrange finite element is assumed. If maxi != 0 then maxj = 0 or maxj = maxi. 
 
 
- 
void vfe_Solid2DSetPropPtr(vfe_Solid2D *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 reference temperature is 0. and the depth is 1. - Errors
- SYS_ERROR_ENUMis generated if an improper type is specified.
 - Parameters
- p – Pointer to Solid2d object. 
- type – Type of element property - x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_TEMPREF Reference temperatures =VFE_PROP_DEPTH Depth 
- propptr – Pointer to start of element nodal properties 
 
 
- 
void vfe_Solid2DStiff(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- kl – [out] Degree of freedom stiffness matrix 
 
 
- 
void vfe_Solid2DStrsStrn(vfe_Solid2D *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 global Cartesian coordinate system. The stress and strain values are ordered first by the 6 tensor components followed by the number of element nodes. For example, a total of 24 values would be returned in strs and 24 values in strn for a 4 node quadrilateral element. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DReact(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of displacements 
- r – [out] Degree of freedom reaction vector 
 
 
- 
void vfe_Solid2DReactStiff(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DSetMatlSystem(vfe_Solid2D *p, Vint type, Vdouble vec[], Vdouble angle)
- set material axes direction - Specify the convention to be used to construct the orientation of the material axes. By default the material axes orientation is - SYS_ELEMSYS_GLOBALwith angle set to 0. If the material axes are not aligned with the element axes (the global Cartesian system), this function may be used to specify the orientation of the material axes.- For a description of element coordinate systems, type, and associated orientation vector data, please see Element Coordinate Systems - Parameters
- p – Pointer to Solid2d object. 
- type – Material system convention 
- vec – Orientation vector data 
- angle – Angle to rotate material x’,y’ axes about the z’ axis in degrees. 
 
 
- 
void vfe_Solid2DInitHist(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
 - Parameters
- p – Pointer to Solid2d object. 
 
- 
void vfe_Solid2DNumDof(vfe_Solid2D *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_Solid2DDofMap()to return the location and type of each degree of freedom.- Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid2d 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_Solid2DNumIntPnt(vfe_Solid2D *p, Vint analysistype, Vint *nepnts)
- query number of element integration points - Query for number of element integration points nepnts. - Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid2d 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_Solid2DDofMap(vfe_Solid2D *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_Solid2DNumDof()to return the number of element degrees of freedom.- Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid2d 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_Solid2DSetHistPtr(vfe_Solid2D *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_Solid2DNumIntPnt(). 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 Solid2d object. 
- oldhist – Pointer to start of material history at previous step 
- newhist – Pointer to start of material history at current step 
 
 
- 
void vfe_Solid2DMass(vfe_Solid2D *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_Solid2DMassDiag()to compute a diagonal mass matrix.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- m – [out] Degree of freedom consistent mass matrix 
 
 
- 
void vfe_Solid2DGeomStiff(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of displacements 
- kg – [out] Degree of freedom geometric stiffness matrix 
 
 
- 
void vfe_Solid2DDistLoad(vfe_Solid2D *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 along a 2D solid element edge. The vector q contains values for the load type for each node in the element. If the loadtype is - VFE_DISTLOAD_PRESthen q contains a scalar pressure at each element node. If the loadtype is- VFE_DISTLOAD_TRACthen q contains a vector traction at each element node. If the loadtype is- VFE_DISTLOAD_TANGFORCEthen q contains a scalar force/area at each element node. Note that the input array of node locations, x, contains the coordinate locations for all nodes in the element - not just on the edge of interest. Correspondingly the output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element.- Errors
- SYS_ERROR_ENUMis generated if an improper enttype or loadtype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero edge Jacobian is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_PRES Load directed along inward normal to element edge =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_Solid2DConcLoad(vfe_Solid2D *p, Vdouble x[][3], Vint enttype, Vint no, Vint loadtype, Vdouble xc[3], Vdouble v[], Vdouble f[])
- concentrated load vector - Compute the consistent degree of freedom loads given a concentrated force q on the edge of the shell. The load location, xc, is projected to the surface defined by the element node locations. The projected point location is used as the point of application of the concentrated force even if it falls outside of the boundary of the element. - The vector q contains three values if loadtype is - VFE_CONCLOAD_TRAC, or a single value if loadtype is- VFE_CONCLOAD_PRES. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The load is in units of force for both- VFE_CONCLOAD_TRACand- VFE_CONCLOAD_PRES.- Errors
- SYS_ERROR_ENUMis generated if an improper enttype or loadtype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero Jacobian is computed.
- SYS_ERROR_COMPUTEis generated if the location xc cannot be projected onto the surface defined by the element node locations.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- enttype – Entity type on which load is applied - =SYS_EDGE Element edge 
- no – Element edge number 
- loadtype – Concentrated load type - x=VFE_CONCLOAD_TRAC Concentrated force vector =VFE_CONCLOAD_PRES Concentrated force along normal direction
- xc – Location of concentrated load 
- v – Concentrated load vector 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid2DElemLoad(vfe_Solid2D *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 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_COMPUTEis generated if a non-positive Jacobian is computed.
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- q – Array of node accelerations 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid2DDistHeat(vfe_Solid2D *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 2D solid 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 area. - Errors
- SYS_ERROR_ENUMis generated if an improper enttype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero edge Jacobian is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DElemHeat(vfe_Solid2D *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_COMPUTEis generated if a non-positive Jacobian is computed.
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- q – Array of node heat fluxes 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid2DPower(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- r – [out] Degree of freedom power vector 
 
 
- 
void vfe_Solid2DPowerCond(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DCond(vfe_Solid2D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- kl – [out] Degree of freedom conductance matrix 
 
 
- 
void vfe_Solid2DHFlxTGrd(vfe_Solid2D *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 fluxes and gradients are computed in the global Cartesian coordinate system. The flux and gradient values are ordered first by the 3 vectoral components followed by the number of element nodes. For example, a total of 12 values would be returned in hflx and 12 values in tgrd for a 4 node quadrilateral element. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DCap(vfe_Solid2D *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_Solid2DCapDiag()to compute a diagonal capacitance matrix. This calculation requires material density and specific heat.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- c – [out] Degree of freedom consistent capacitance matrix 
 
 
- 
void vfe_Solid2DCapDiag(vfe_Solid2D *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_Solid2DCap()to compute a consistent capacitance matrix. This calculation requires material density and specific heat.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- cd – [out] Degree of freedom diagonal capacitance matrix 
 
 
- 
void vfe_Solid2DStrsAdapt(vfe_Solid2D *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 2. These quantities are useful for computing new characteristic element length for mesh adaptation. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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_Solid2DHFlxAdapt(vfe_Solid2D *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 2. These quantities are useful for computing new characteristic element length for mesh adaptation. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid2d 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 
 
 
- 
void vfe_Solid2DShapeGrad(vfe_Solid2D *p, Vdouble x[][3], Vint centflag, Vdouble h[], Vdouble phr[][2], Vdouble dj[], Vdouble phx[][3])
- shape functions and gradients - Return shape functions and gradients at element nodes or centroid only. The shape functions, h, contain shape functions for each node for each node or centroid. The gradients, phr and phx contains gradients for each node for each node or centroid. The determinant of the Jacobian, dj, is a single value for each node or centroid. If the determinant is zero then the corresponding partials wrt. coordinates are undefined. - Parameters
- p – Pointer to Solid2d object. 
- x – Array of node locations. 
- centflag – Centroid flag - x=SYS_ON Centroid location only =SYS_OFF Element node locations 
- h – [out] Shape functions 
- phr – [out] Gradient of shape function wrt. natural coordinates 
- dj – [out] Determinant of Jacobian 
- phx – [out] Gradient of shape function wrt. coordinates 
 
 
2.3. 3D Solid Elements - Solid3D
The Solid3D module is used to simulate 3D continuum.
The methods associated with a Solid3D object are the following.
- Begin and end an instance of an object, generic object functions - vfe_Solid3DBegin()- create an instance of a Solid3D object
- vfe_Solid3DEnd()- destroy an instance of a Solid3D object
- vfe_Solid3DError()- return Solid3D object error flag
 
- Attributes and Parameters - vfe_Solid3DSetHistPtr()- set pointers to material history
- vfe_Solid3DSetMatlSystem()- set material axes direction
- vfe_Solid3DSetObject()- set attribute object
- vfe_Solid3DSetParami()- set element formulation parameters
- vfe_Solid3DSetPropPtr()- set pointer to element nodal properties
- vfe_Solid3DSetTopology()- set input element topology
 
- Degree of Freedom and Integration Information - vfe_Solid3DDofMap()- query element degree of freedom map
- vfe_Solid3DNumDof()- query number of element degrees of freedom
- vfe_Solid3DNumIntPnt()- query number of element integration points
 
- General - vfe_Solid3DShapeGrad()- shape functions and gradients
 
- Structural analysis computations - vfe_Solid3DBMatrix()- B(beta) matrix
- vfe_Solid3DConcLoad()- concentrated load vector
- vfe_Solid3DDistLoad()- distributed load vector.
- vfe_Solid3DElemLoad()- body force vector
- vfe_Solid3DGeomStiff()- geometric stiffness matrix
- vfe_Solid3DInitHist()- initialize material history
- vfe_Solid3DMass()- consistent mass matrix
- vfe_Solid3DMassDiag()- diagonal mass matrix
- vfe_Solid3DReact()- reaction vector
- vfe_Solid3DReactStiff()- reaction vector, stiffness matrix
- vfe_Solid3DStiff()- linear stiffness matrix
- vfe_Solid3DStrsAdapt()- stress based error analysis
- vfe_Solid3DStrsStrn()- stress and strain
 
- Thermal analysis computations - vfe_Solid3DDistHeat()- distributed heat loads.
- vfe_Solid3DElemHeat()- body heat generation
- vfe_Solid3DCap()- consistent capacitance matrix
- vfe_Solid3DCapDiag()- diagonal capacitance matrix
- vfe_Solid3DPower()- thermal power
- vfe_Solid3DPowerCond()- thermal power, conductance matrix
- vfe_Solid3DCond()- conductance matrix
- vfe_Solid3DHFlxAdapt()- heat flux based error analysis
- vfe_Solid3DHFlxTGrd()- heat flux and temperature gradient
 
Instance a Solid3D object using vfe_Solid3DBegin().
Use the function
vfe_Solid3DStrsAdapt() 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_Solid3DHFlxAdapt() performs a similar
computation for heat transfer analysis.
The element technology to be used is set using the
function vfe_Solid3DSetParami() with type VFE_TECH.
Standard fully integrated, VFE_TECH_ISOP,
and underintegrated, VFE_TECH_URED exist for all elements.
The default technology is VFE_TECH_ISOP which will work in all cases.
However for some element types, alternative technologies should be
considered for improved performance.
The recommended technology varies with element topology and order and
is listed below.  The number of integration points used for fully
integrated technologies appears in paraenthesis.  Underintegrated
technologies use a single integration point for the linear element
forms and for the higher order element forms use the fully integrated
rule for the one-less order.
- Linear tetrahedron, VFE_TECH_ISOP (1) 
- Parabolic tetrahedron, VFE_TECH_ISOP (4) 
- Cubic tetrahedron, VFE_TECH_ISOP (11) 
 
- Linear pyramid, VFE_TECH_ISOP (5) 
- Parabolic pyramid, VFE_TECH_ISOP (13) 
- Cubic pyramid, VFE_TECH_ISOP (27) 
 
- Linear pentahedron, VFE_TECH_ISOP (6) Use VFE_TECH_ENHANCED if bending deformations are particularly important with a coarse mesh. 
- Parabolic pentahedron, VFE_TECH_URED if adjacent elements exist across 4 or more faces else VFE_TECH_ISOP (9). 
- Cubic pentahedron, VFE_TECH_ISOP (18) 
 
- Linear hexahedron, VFE_TECH_MIXED or VFE_TECH_ENHANCED. Use VFE_TECH_MIXED unless bending deformations are particularly important with a coarse mesh. The fully integrated technology, VFE_TECH_ISOP (8), is overly stiff. 
- Parabolic hexahedron - if Lagrangian VFE_TECH_ISOP (27), if Serendipity VFE_TECH_URED if adjacent elements exist across 5 or more faces else VFE_TECH_ISOP. An enhanced formulation, VFE_TECH_ENHANCED, also is supported. 
- Cubic hexahedron - VFE_TECH_ISOP (64) 
2.4. Function Descriptions
The currently available Solid3D functions are described in detail in this section.
- 
vfe_Solid3D *vfe_Solid3DBegin(void)
- create an instance of a Solid3D object - Create an instance of a Solid3D object. Memory is allocated for the object private data and the pointer to the object is returned. Default topology is the 8-noded Serendipity hexahedron with isoparametric technology. - Destroy an instance of a Solid3D object using - void vfe_Solid3DEnd (vfe_Solid3D *solid3d) - Return the current value of a Solid3D object error flag using - Vint vfe_Solid3DError (vfe_Solid3D *solid3d) - Returns
- The function returns a pointer to the newly created Solid3D object. If the object creation fails, NULL is returned. 
 
- 
void vfe_Solid3DEnd(vfe_Solid3D *p)
- create an instance of a Solid3D object 
- 
Vint vfe_Solid3DError(vfe_Solid3D *p)
- return the current value of a Solid3D object error flag 
- 
void vfe_Solid3DSetObject(vfe_Solid3D *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_OBJECTTYPEis generated if an improper objecttype is specified.
 - Parameters
- p – Pointer to Solid3d object. 
- objecttype – The object type identifier - x=VFE_MATLFUN MatlFun object 
- object – Pointer to the object to be set. 
 
 
- 
void vfe_Solid3DSetParami(vfe_Solid3D *p, Vint type, Vint iparam)
- set element formulation parameters - Set element formulation technology parameters. The basic element technology is set with the - VFE_TECHparameter. By default- VFE_TECHis set to- VFE_TECH_ISOP.- The parameter - VFE_CHECKTOTALVOLtoggles of Jacobian checking performed during element volume integrations. If this parameter is not enabled, then an error is generated if the Jacobian contribution of any integration point is non-positive. If this parameter is enabled then an error is generated only if the total Jacobian contribution over all integration points is non-positive. By default- VFE_CHECKTOTALVOLis set to- SYS_OFF.- Set element strain type using - VFE_STRAINTYPEwith a value of either- VFE_LARGESTRAINto enable large strain or- VFE_SMALLSTRAINto enable small strains. By default- VFE_STRAINTYPEis set to- VFE_SMALLSTRAIN.- The parameter - VFE_TEMPMATLAVGtoggles 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_TEMPMATLAVGis set to- SYS_ON.- Errors
- SYS_ERROR_ENUMis generated if an improper type is specified.
- SYS_ERROR_VALUEis generated if an improper iparam is specified.
 
 - Parameters
- p – Pointer to Solid3d object. 
- type – Type of formulation parameter to set - x=VFE_TECH Element technology =VFE_CHECKTOTALVOL Total Jacobian checking =VFE_TEMPMATLAVE Average material temperature flag =VFE_STRAINTYPE Element strain type 
- iparam – Integer parameter value. - x=SYS_ON Enable =SYS_OFF Disable =VFE_TECH_ISOP Standard isoparametric technology =VFE_TECH_URED Uniform reduced technology =VFE_TECH_MIXED Mixed technology =VFE_TECH_ENHANCED Enhanced technology =VFE_SMALLSTRAIN Small strain =VFE_LARGESTRAIN Large strain 
 
 
- 
void vfe_Solid3DSetTopology(vfe_Solid3D *p, Vint shape, Vint maxi, Vint maxj, Vint maxk)
- set element topology - Specify the topology of a 3D solid element. If maxi is set to 3 then a quadratic element form is assumed. If maxi is set to 4 then a cubic element form is assumed. The default topology is - SYS_SHAPEHEXwith maxi = maxj = maxk = 0. Serendipity elements with missing mid-side nodes are also supported. See Computational Cells for a description on how to define element topologies with missing mid-side nodes.- Errors
- SYS_ERROR_VALUEis generated if an improper maxi, maxj or maxk is specified.
- SYS_ERROR_ENUMis generated if an improper shape is input.
 
 - Parameters
- p – Pointer to Solid3d object. 
- shape – The topological shape for the element - x=VIS_SHAPETET Tetrahedron =VIS_SHAPEPYR Pyramid =VIS_SHAPEWED Wedge =VIS_SHAPEHEX Hexahedron 
- 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 <= 4 and 2 <= maxi <= 4, then a Lagrange finite element is assumed. If maxi != 0 then maxj = 0 or maxj = maxi. 
- maxk – The number of points along the k direction. Currently maxk must be equal to maxj. 
 
 
- 
void vfe_Solid3DSetPropPtr(vfe_Solid3D *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. and reference temperature is 0. In thermal analysis computations, - VFE_PROP_VOLFACT_MASSwill impact the capacitance,- VFE_PROP_VOLFACT_STIFFthe conductance, and- VFE_PROP_VOLFACT_LOADthe body heat loads.- Errors
- SYS_ERROR_ENUMis generated if an improper type is specified.
 - Parameters
- p – Pointer to Solid3d object. 
- type – Type of element property - x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_TEMPREF Reference temperatures =VFE_PROP_VOLFACT All volume scale factors =VFE_PROP_VOLFACT_MASS Volume scale factor for the mass only =VFE_PROP_VOLFACT_STIFF Volume scale factor for the stiffness only =VFE_PROP_VOLFACT_LOAD Volume scale factor for the body loads only 
- propptr – Pointer to start of element nodal properties 
 
 
- 
void vfe_Solid3DStiff(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- kl – [out] Degree of freedom stiffness matrix 
 
 
- 
void vfe_Solid3DStrsStrn(vfe_Solid3D *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 global Cartesian coordinate system. The stress and strain values are ordered first by the 6 tensor components followed by the the number of element nodes. For example, a total of 48 values would be returned in strs and 48 values in strn for a 8 node hexahedral element. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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_Solid3DReact(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of displacements 
- r – [out] Degree of freedom reaction vector 
 
 
- 
void vfe_Solid3DReactStiff(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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_Solid3DNumDof(vfe_Solid3D *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_Solid3DDofMap()to return the location and type of each degree of freedom.- Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid3d 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_Solid3DNumIntPnt(vfe_Solid3D *p, Vint analysistype, Vint *nepnts)
- query number of element integration points - Query for number of element integration points nepnts. - Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid3d 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_Solid3DDofMap(vfe_Solid3D *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_Solid3DNumDof()to return the number of element degrees of freedom.- Errors
- SYS_ERROR_ENUMis generated if an improper analysistype is specified.
 - Parameters
- p – Pointer to Solid3d 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_Solid3DSetHistPtr(vfe_Solid3D *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_Solid3DNumIntPnt(). 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 Solid3d object. 
- oldhist – Pointer to start of material history at previous step 
- newhist – Pointer to start of material history at current step 
 
 
- 
void vfe_Solid3DMass(vfe_Solid3D *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_Solid3DMassDiag()to compute a diagonal mass matrix.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- m – [out] Degree of freedom consistent mass matrix 
 
 
- 
void vfe_Solid3DMassDiag(vfe_Solid3D *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_Solid3DMass()to compute a consistent mass matrix.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- md – [out] Degree of freedom diagonal mass vector 
 
 
- 
void vfe_Solid3DGeomStiff(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of displacements 
- kg – [out] Degree of freedom geometric stiffness matrix 
 
 
- 
void vfe_Solid3DDistLoad(vfe_Solid3D *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 solid element edge or face. The vector q contains values for the load type for each node in the element. If the loadtype is - VFE_DISTLOAD_PRESthen q contains a scalar pressure at each element node. If the loadtype is- VFE_DISTLOAD_TRACthen q contains a vector traction at each element node. If the traction is applied to an edge the units are force/length, if applied to a face the units are force/area. If the loadtype is- VFE_DISTLOAD_TANGFORCEthen q contains a scalar force/length at each element node.- Note that the input array of node locations, x, contains the coordinate locations for all nodes in the element - not just on the face of interest. Correspondingly the output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. - Errors
- SYS_ERROR_ENUMis generated if an improper enttype or loadtype is specified.
- SYS_ERROR_OPERATIONis generated if an invalid combination of enttype and loadtype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero edge or face Jacobian is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- enttype – Entity type on which load is applied - =SYS_EDGE Element edge =SYS_FACE Element face 
- no – Element face number 
- loadtype – Distributed load type - x=VFE_DISTLOAD_TRAC Load directed along vector =VFE_DISTLOAD_PRES Load directed along inward normal to element face =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_Solid3DConcLoad(vfe_Solid3D *p, Vdouble x[][3], Vint enttype, Vint no, Vint loadtype, Vdouble xc[3], Vdouble v[], Vdouble f[])
- concentrated load vector - Compute the consistent degree of freedom loads given a concentrated force q on the edge of the shell. The load location, xc, is projected to the surface defined by the element node locations. The projected point location is used as the point of application of the concentrated force even if it falls outside of the boundary of the element. - The vector q contains three values if loadtype is - VFE_CONCLOAD_TRAC, or a single value if loadtype is- VFE_CONCLOAD_PRES. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The load is in units of force for both- VFE_CONCLOAD_TRACand- VFE_CONCLOAD_PRES.- Errors
- SYS_ERROR_ENUMis generated if an improper enttype or loadtype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero Jacobian is computed.
- SYS_ERROR_COMPUTEis generated if the location xc cannot be projected onto the surface defined by the element node locations.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- enttype – Entity type on which load is applied - =SYS_FACE Element face 
- no – Element face number 
- loadtype – Concentrated load type - x=VFE_CONCLOAD_TRAC Concentrated force vector =VFE_CONCLOAD_PRES Concentrated force along normal direction
- xc – Location of concentrated load 
- v – Concentrated load vector 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid3DElemLoad(vfe_Solid3D *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_COMPUTEis generated if a non-positive Jacobian is computed.
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- q – Array of node accelerations 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid3DDistHeat(vfe_Solid3D *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 solid element face. 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 area. - Errors
- SYS_ERROR_ENUMis generated if an improper enttype is specified.
- SYS_ERROR_VALUEis generated if an improper no is specified.
- SYS_ERROR_COMPUTEis generated if a zero face Jacobian is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- enttype – Entity type on which load is applied - =SYS_FACE Element face 
- no – Element face number 
- q – Vector of distributed load values 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid3DElemHeat(vfe_Solid3D *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_COMPUTEis generated if a non-positive Jacobian is computed.
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- q – Array of node heat fluxes 
- f – [out] Degree of freedom vector of consistent loads. 
 
 
- 
void vfe_Solid3DPowerCond(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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_Solid3DPower(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- r – [out] Degree of freedom power vector 
 
 
- 
void vfe_Solid3DCond(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- kl – [out] Degree of freedom conductance matrix 
 
 
- 
void vfe_Solid3DHFlxTGrd(vfe_Solid3D *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 fluxes and gradients are computed in the global Cartesian coordinate system. The flux and gradient values are ordered first by the 3 vectoral components followed by the number of element nodes. For example, a total of 24 values would be returned in hflx and 24 values in tgrd for a 8 node hexahedral element. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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_Solid3DCap(vfe_Solid3D *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_Solid3DCapDiag()to compute a diagonal capacitance matrix. This calculation requires material density and specific heat.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- c – [out] Degree of freedom consistent capacitance matrix 
 
 
- 
void vfe_Solid3DCapDiag(vfe_Solid3D *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_Solid3DCap()to compute a consistent capacitance matrix. This calculation requires material density and specific heat.- Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- u – Degree of freedom vector of temperatures 
- cd – [out] Degree of freedom diagonal capacitance matrix 
 
 
- 
void vfe_Solid3DInitHist(vfe_Solid3D *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_NULLOBJECTis generated if a MatlFun attribute object has not been set.
 - Parameters
- p – Pointer to Solid3d object. 
 
- 
void vfe_Solid3DSetMatlSystem(vfe_Solid3D *p, Vint type, Vdouble vec[], Vdouble angle)
- set material axes direction - Specify the convention to be used to construct the orientation of the material axes. By default the material axes orientation is - SYS_ELEMSYS_GLOBALwith angle set to 0. If the material axes are not aligned with the element axes (the global Cartesian system), this function may be used to specify the orientation of the material axes.- For a description of element coordinate systems, type, and associated orientation vector data, please see Element Coordinate Systems - Parameters
- p – Pointer to Solid3d object. 
- type – Material system convention 
- vec – Orientation vector data 
- angle – Angle to rotate material x’,y’ axes about the z’ axis in degrees. 
 
 
- 
void vfe_Solid3DStrsAdapt(vfe_Solid3D *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 3. These quantities are useful for computing new characteristic element length for mesh adaptation. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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_Solid3DHFlxAdapt(vfe_Solid3D *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 3. These quantities are useful for computing new characteristic element length for mesh adaptation. - Errors
- SYS_ERROR_NULLOBJECTis generated if a MatlFun attribute object has not been set.
- SYS_ERROR_COMPUTEis generated if a negative Jacobian transformation is computed.
 
 - Parameters
- p – Pointer to Solid3d 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 
 
 
- 
void vfe_Solid3DShapeGrad(vfe_Solid3D *p, Vdouble x[][3], Vint centflag, Vdouble h[], Vdouble phr[][3], Vdouble dj[], Vdouble phx[][3])
- shape functions and gradients - Return shape functions and gradients at element nodes or centroid only. The shape functions, h, contain shape functions for each node for each node or centroid. The gradients, phr and phx contains gradients for each node for each node or centroid. The determinant of the Jacobian, dj, is a single value for each node or centroid. If the determinant is zero then the corresponding partials wrt. coordinates are undefined. - The strain-displacement matrix for 3D isoparametric solid elements (B-Matrix) relates the 3D strain measures to the nodal displacements. For a 3D element defined by N nodes, the displacements are given by shape functions h_i(r,s,t), with 1=1,M, as follows: where ux_i, uy_i, and uz_i are the components of the nodal displacement vector.ui.- | ux | | ux_i | | uy | = SUM_(i=1,N) h_i(r,s,t) * | uy_i | | uz | | uz_i | - The 3D strains E are given by We rewrite the above equations as- | e_xx | | ux,x | | h_i,x 0 0 | | e_yy | | uy,y | | 0 h_i,y 0 | | ux_i | E = | e_zz | = | uz,z | = SUM_(i=1,N) | 0 0 h_i,z | * | uy_i | | g_xy | | ux,y + uy,x | | h_i,y h_i,x 0 | | uz_i | | g_yz | | uy,z + uz,y | | 0 h_i,z h_i,y | | g_zx | | uz,x + ux,z | | h_i,z 0 h_i,x |where- E = SUM_(i=1,N) B_i * U_i The resulting B-Matrix is given by- | h_i,x 0 0 | | 0 h_i,y 0 | B_i = | 0 0 h_i,z | | h_i,y h_i,x 0 | | 0 h_i,z h_i,y | | h_i,z 0 h_i,x |- B = | B_1 B_2 ... B_N | - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- centflag – Centroid flag - x=SYS_ON Centroid location only =SYS_OFF Element node locations 
- h – [out] Shape functions 
- phr – [out] Gradient of shape function wrt. natural coordinates 
- dj – [out] Determinant of Jacobian 
- phx – [out] Gradient of shape function wrt. coordinates 
 
 
- 
void vfe_Solid3DBMatrix(vfe_Solid3D *p, Vdouble x[][3], Vint *nipt, Vdouble rst[][3], Vdouble h[], Vdouble phr[][3], Vdouble phx[][3], Vdouble w[], Vdouble dj[], Vdouble tm[][3][3], Vdouble dm[][21], Vdouble b[][6])
- B(beta) matrix. - Compute the B(beta) matrix relating the strain components (exx,eyy,ezz,gxy,gyz,gzx) to the element degrees of freedom. Use - vfe_Solid3DNumDof()to determine the number of element degrees of freedom, and- vfe_Solid3DDofMap()to determine the meaning of each degree of freedom.- This function also returns the number of integration points needed to form the element stiffness matrix, nipt; the (r,s,t) natural coordinates of the element at the integration points; the nodal shape functions at each integration point in the order h_(i1,n1), h_(i1,n2), … , h_(i1,nen), h_(i2,n1), h_(i2,n2), … h_(nipt,nen) where i1, i2, …, nipt are integration point numbers, and n1, n2, … , nen are node numbers; the gradient phr of the shape functions h with respect to the natural coordinates (r,s,t); the gradient phx of the shape functions h with respect to the global coordinates (x,y,z); the integration point weights w at each integration point; the determinant of the Jacobian transformation d(x,y,z)/d(r,s,t); the orthogonal material matrices tm at each integration point; the material matrix dm at each integration point, where the 21 entries refer to the lower triangle of the 6x6 matrix relating the strains to the stresses; and the B(beta) matrix at each integration point relating the strains to the element degrees of freedom. The rows of B are defined for all element degrees of freedom for the first integration point, followed by the values at the second integration point, etc. - Errors
- SYS_ERROR_OPERATIONis generated if element is improperly configured.
- SYS_ERROR_COMPUTEis generated if- vfe_Solid3DSetPropPtr()has been called with- VFE_PROP_VOLFACTor- VFE_PROP_VOLFACT_STIFFand the nodal volume factors interpolate to a negative value at an integration point, or if a negative determinant of the Jacobian is detected.
 
 - Parameters
- p – Pointer to Solid3d object. 
- x – Array of node locations. 
- nipt – [out] Number of integration points 
- rst – [out] Natural coordinates (r,s,t) at integration points 
- h – [out] Nodal shape functions at integration points 
- phr – [out] Gradient of h w.r.t. natural coordinates at integration points 
- phx – [out] Gradient of h w.r.t. global coordinates at integration points 
- w – [out] Weights at integration points 
- dj – [out] Determinant of Jacobian matrix at integration points 
- tm – [out] Material coordinate system at integration points 
- dm – [out] Material matrix at integration points 
- b – [out] B(beta) matrix relating element dofs to strains 
 
 
