3. Shell and Membrane Elements - Shell2D, Shell3D, Mem3D
Shell elements are a specialization of a volume element in which a given direction, the shell thickness, may be assumed to be small with respect to the other directions - the shell surface. Shell elements generally assume that the stress normal to the shell surface is zero. There are two modules for the simulation of shells: Shell2D for plane stress, plane strain, and axisymmetric 2D analysis and Shell3D for full 3D analysis.
Membrane elements are used to represent thin surfaces which have no bending stiffness, only stiffness in the plane of the surface is retained. Use the Mem3D module to model membranes for 3D analysis.
3.1. 2D Shell Elements - Shell2D
The Shell2D module is used to simulate 2D shells. 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 Shell2D functions as a three component position vector, in all cases the third component is ignored.
All symmetric shell resultants are given in the order Nxx, Nzz, Nxz, Mxx, Mzz. Mxz. Qxy. Qzy. Note that in plane strain, the ezz, gxz and kzz, kxz and gzy components of strain are identically zero, while the Nxz and Mxx, Mxz and Qzy components of stress are zero. In plane stress, the gxz and kzz, kxz and gzy components of strain are identically zero, while the Nzz, Nxz and Mxx, Mxz and Qzy components of stress are zero.
2D Shell elements are one-dimensional elements in 2D space which support in-plane deformation, bending and tranverse shear through the thickness of the shell. The stress normal to the shell surface is assumed to be zero. The methods associated with a Shell2D object are the following.
Begin and end an instance of an object, generic object functions
vfe_Shell2DBegin()
- create an instance of a Shell2D objectvfe_Shell2DEnd()
- destroy an instance of a Shell2D objectvfe_Shell2DError()
- return Shell2D object error flag
Attributes and Parameters
vfe_Shell2DSetHistPtr()
- set pointers to material historyvfe_Shell2DSetObject()
- set attribute objectvfe_Shell2DSetParamd()
- set element formulation parametersvfe_Shell2DSetParami()
- set element formulation parametersvfe_Shell2DSetPropPtr()
- set pointer to element nodal propertiesvfe_Shell2DSetTopology()
- set input element topology
Degree of Freedom and Integration Information
vfe_Shell2DDofMap()
- query element degree of freedom mapvfe_Shell2DNumDof()
- query number of element degrees of freedomvfe_Shell2DNumIntPnt()
- query number of element integration points
Structural analysis computations
vfe_Shell2DConcLoad()
- concentrated load vectorvfe_Shell2DDirCos()
- compute shell local direction cosinesvfe_Shell2DDistLoad()
- distributed load vectorvfe_Shell2DElemLoad()
- body force vectorvfe_Shell2DGeomStiff()
- geometric stiffness matrixvfe_Shell2DInitHist()
- initialize material historyvfe_Shell2DMass()
- consistent mass matrixvfe_Shell2DMassDiag()
- diagonal mass matrixvfe_Shell2DReact()
- reaction vectorvfe_Shell2DReactStiff()
- reaction vector, stiffness matrixvfe_Shell2DStiff()
- linear stiffness matrixvfe_Shell2DStrsAdapt()
- stress based error analysisvfe_Shell2DStrsStrn()
- resultant stress and strain/curvaturevfe_Shell2DStressStrain()
- three-dimensional stress and strain
3.2. Element Geometry
The basic shell element geometry is defined by the node coordinates of the shell
reference surface.
The shell wall is defined by the node offsets of the shell midsurface
from the reference surface, shell thicknesses and shell normals.
These properties which are set using vfe_Shell2DSetPropPtr()
are illustrated below in Figure 3-1.
Other element quantities which may be set using
this function are
shell temperature and effective shear factor.
If the node offsets are zero, the shell midsurface is coincident with the shell reference surface. In general, large midsurface offset values, which are much greater than the shell thickness, should be avoided. By default the shell node normals are computed to be normal to the shell reference surface. The user may optionally enter node normals to achieve a “miter” joint between adjacent shell elements. Shell thickness and offset are always measured in the direction of the normal to the shell reference surface regardless of the direction of the shell node normals. The element natural coordinate, r, lies on the shell midsurface with the positive direction in the direction of the element connectivity. The element local coordinate system for stress and material orientation is the local x’ axis along the element local r direction and the local z’. direction along the global z axis. The local y’ axis is normal to the shell reference surface. The local axes form a right handed orthogonal system.
The element technology to be used is set using the
function vfe_Shell2DSetParami()
with type VFE_TECH.
Standard fully integrated, VFE_TECH_ISOP, formulations exist for all orders.
The default technology is VFE_TECH_ISOP.
Linear VFE_TECH_ISOP (1). Use a Kirckhoff formulation, VFE_TECH_KIRCHHOFF, for cubic bending.
Parabolic VFE_TECH_ISOP (2).
The exact nature of the shell wall construction which couples element
and primitive material properties,
is managed using the ShellProp module.
The ShellProp module provides for simple monocoque,
vfe_ShellPropSetMonocoque()
,
and general integrated shell walls, vfe_ShellPropSetIntegrate()
, as well as
preintegrated shell wall constructions input as a general matrix,
vfe_ShellPropSetMatrix()
, relating
midsurface strains and curvatures to stress and moment resultants.
For monocoque and general integrated shell walls the ShellProp module
allows primitive materials, such as defined by the LinMat
module, to be specified at the integration points.
Note that when using the general integrated shell wall defined by
vfe_ShellPropSetIntegrate()
,
the actual distance of the shell wall integration point from the
midsurface is computed internally by multiplying the integration point
coordinates input in vfe_ShellPropSetIntegrate()
by one-half of the
thickness specified using vfe_Shell2DSetPropPtr()
.
Therefore, if the actual shell thickness is specified using
vfe_Shell2DSetPropPtr()
, the shell wall integration point
coordinates input in vfe_ShellPropSetIntegrate()
should be in the range
[-1.,1.] from the bottom to the top shell surface, with the coordinate
0. lying on the midsurface.
Use the function
vfe_Shell2DStrsAdapt()
to aid in computing element strain energy, strain
energy error and other useful quantities to aid in solution error estimation
and mesh adaptation.
3.3. Function Descriptions
The currently available Shell2D functions are described in detail in this section.
-
vfe_Shell2D *vfe_Shell2DBegin(void)
create an instance of a Shell2D object
Create an instance of a Shell2D object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a Shell2D object using
void vfe_Shell2DEnd (vfe_Shell2D *shell2d)
Return the current value of a Shell2D object error flag using
Vint vfe_Shell2DError (vfe_Shell2D *shell2d)
- Returns:
The function returns a pointer to the newly created Shell2D object. If the object creation fails, NULL is returned.
-
void vfe_Shell2DEnd(vfe_Shell2D *p)
destroy an instance of a Shell2D object
-
Vint vfe_Shell2DError(vfe_Shell2D *p)
return the current value of a Shell2D object error flag
-
void vfe_Shell2DSetObject(vfe_Shell2D *p, Vint objecttype, Vobject *object)
set attribute object
Set a pointer to an attribute object. Users must supply a MatlFun object prior to computing any quantity that requires a material model definition.
- Errors
SYS_ERROR_OBJECTTYPE
is generated if an improper objecttype is specified.
- Parameters:
p – Pointer to Shell2D object.
objecttype – The object type identifier
x=VFE_MATLFUN MatlFun object
object – Pointer to the object to be set.
-
void vfe_Shell2DSetParami(vfe_Shell2D *p, Vint type, Vint iparam)
set element formulation parameters
Set element technology and formulation parameters.
The basic element technology is set with the
VFE_TECH
parameter. By defaultVFE_TECH
is set toVFE_TECH_ISOP
.Toggle large rotations using
VFE_LARGEROTATION
. By defaultVFE_LARGEROTATION
is set toSYS_OFF
.The addition of rotary inertias to the element mass matrix is set with the
VFE_ROTINERTIA
parameter. By defaultVFE_ROTINERTIA
is set toSYS_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 defaultVFE_2D
is set toVFE_PLANESTRESS
.The parameter
VFE_TEMPMATLAVG
toggles the method for computing the temperature used for evaluating temperature dependent material properties. If enabled, the temperature used for temperature dependent material properties is the average of the element node point temperatures. If disabled, the temperature is isoparametrically interpolated from the node point temperatures at each element integration point. By defaultVFE_TEMPMATLAVG
is set toSYS_ON
.The parameter
VFE_CHECKCORNNORM
toggles a comparison between computed normals at integration points in shells and a computed normal interpolated from corner node values. If this parameter is enabled then the dot product of these two normals is computed and, if negative, an error is flagged. By defaultVFE_CHECKCORNNORM
is set toSYS_ON
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper iparam is specified.
- Parameters:
p – Pointer to Shell2D object.
type – Type of formulation parameter to set
x=VFE_TECH Element technology =VFE_LARGEROTATION Toggle for large rotations =VFE_ROTINERTIA Enable rotary inertias =VFE_TEMPMATLAVE Average material temperature flag =VFE_2D 2D Approximation =VFE_CHECKCORNNORM Corner normal checking
iparam – Integer parameter value.
x=VFE_TECH_ISOP Standard isoparametric technology =VFE_TECH_KIRCHHOFF Kirchhoff technology =VFE_AXISYMMETRIC 2D axisymmetric =VFE_PLANESTRAIN 2D plane strain =VFE_PLANESTRESS 2D plane stress =SYS_OFF Disable =SYS_ON Enable
-
void vfe_Shell2DSetParamd(vfe_Shell2D *p, Vint type, Vdouble dparam)
set element formulation parameters
Set element formulation parameters.
Use
VFE_MAXPROJANG
to set the maximum angle in degrees between the shell normal at an element integration point and the shell normal at a node. By defaultVFE_MAXPROJANG
is set to 90. degrees.Use
VFE_ADDEDMASS
to add non-structural mass/area to the shell element. By defaultVFE_ADDEDMASS
is set to 0.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper dparam is specified.
- Parameters:
p – Pointer to Shell2D object.
type – Type of formulation parameter to set
x=VFE_MAXPROJANG Maximum nodal projection angle =VFE_ADDEDMASS Additional mass/area
dparam – Double precision parameter value.
-
void vfe_Shell2DSetTopology(vfe_Shell2D *p, Vint maxi)
set element topology
- 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 Shell2D object.
maxi – The number of points along the i direction. Specify the topology of a 2D shell element. If maxi is set to 3 then a quadratic element form is assumed.
-
void vfe_Shell2DSetPropPtr(vfe_Shell2D *p, Vint type, Vdouble *propptr)
set pointer to element nodal properties
Set a pointer to the start of a specified type of element nodal properties. Note that the properties are not copied by this function, only the pointer itself is copied. If a property pointer is not set the element assumes a default value for the associated property. By default the temperature is 0., the thickness is 2., the offset is 0. and the effective shear factor is 5/6. If the element node normals are not specified, the normals are computed to be normal to the shell reference surface. By default the reference temperature is 0.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to Shell2D object.
type – Type of element property
x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_TEMPREF Reference temperatures =VFE_PROP_THICKNESS Thicknesses =VFE_PROP_NORMAL Normals =VFE_PROP_OFFSET Offsets =VFE_PROP_KS Effective shear factor
propptr – Pointer to start of element nodal properties
-
void vfe_Shell2DSetHistPtr(vfe_Shell2D *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 total number of integration points in the element. The number of history variables depends on the underlying material model and may be queried using
vfe_MatlFunNumHist()
. The total number of integration points is the product of the number of element integration points returned usingvfe_Shell2DNumIntPnt()
and the number of section integration points returned usingvfe_MatlFunNumStressStrain()
. 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 Shell2D object.
oldhist – Pointer to start of material history at previous step
newhist – Pointer to start of material history at current step
-
void vfe_Shell2DStiff(vfe_Shell2D *p, Vdouble x[][3], Vdouble kl[])
linear stiffness matrix
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
-
void vfe_Shell2DGeomStiff(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vdouble kg[])
geometric stiffness matrix
Compute the geometric stiffness matrix, kg, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the geometric stiffness is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kg – [out] Degree of freedom geometric stiffness matrix
-
void vfe_Shell2DReact(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
reaction vector
Compute the reaction vector, r, given the node coordinates, x, and the degree of freedom displacement vector, u.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
r – [out] Degree of freedom reaction vector
-
void vfe_Shell2DReactStiff(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
reaction vector, stiffness matrix
Compute the reaction vector, r, and optionally the stiffness matrix, k, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D 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_Shell2DStrsStrn(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[])
resultant stress and strain/curvature
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 shell local coordinate system. The actual direction cosine matrices of the shell local systems at the element output locations may be returned using
vfe_Shell2DDirCos()
.The strs is composed of 8 stress resultants and the strn is composed of 8 associated midsurface strains and curvatures at each output location. The stress and strain values are ordered first by the 8 stress resultant or strain components followed by the the number of element nodes. For example, for a 3 node shell element a total of 8*3 = 24 stress values will be returned in strs and 24 strain values returned in strn.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D 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_Shell2DStressStrain(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[], Vdouble stress[], Vdouble strain[])
three-dimensional stress and strain
Compute nodal stresses and strains, strs and strn, and the three-dimensional stresses and strains, stress and strain, given the node coordinates, x, and the degree of freedom displacement vector, u.
The stress resultants and strain/curvatures are the same as those obtained with
vfe_Shell2DStrsStrn()
, and the same remarks regarding local coordinate systems apply.The three-dimensional stresses and strains are computed at a series of points in the beam’s cross section. The number of cross sectional points is obtained using
vfe_MatlFunNumStressStrain()
and depend on the type of cross section used. See ShellProp for more details on each cross section type.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
strs – [out] Array of nodal stress resultants
strn – [out] Array of nodal strains/curvatures
stress – [out] Array of nodal three dimensional stresses at recovery points
strain – [out] Array of nodal three dimensional strains at recovery points
-
void vfe_Shell2DInitHist(vfe_Shell2D *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 Shell2D object.
-
void vfe_Shell2DNumDof(vfe_Shell2D *p, Vint analysistype, Vint *nepnts)
query number of element integration points
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_Shell2DDofMap()
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 Shell2D object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
nepnts – [out] Number of element integration points
-
void vfe_Shell2DDofMap(vfe_Shell2D *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, rotation, temperature, etc.
The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_Shell2DNumDof()
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 Shell2D object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
loc – [out] Vector of degree of freedom locations
tag – [out] Vector of degree of freedom types
-
void vfe_Shell2DMass(vfe_Shell2D *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_Shell2DMassDiag()
to compute a diagonal mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
m – [out] Degree of freedom consistent mass matrix
-
void vfe_Shell2DMassDiag(vfe_Shell2D *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_Shell2DMass()
to compute a consistent mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
md – [out] Degree of freedom diagonal mass vector
-
void vfe_Shell2DDistLoad(vfe_Shell2D *p, Vdouble x[][3], Vint enttype, Vint loadtype, Vdouble q[], Vdouble f[])
distributed load vector
Compute the consistent degree of freedom loads given a distributed load, q along a 2D shell element edge. The vector q contains values for the load type for each node in the element. If the loadtype is
VFE_DISTLOAD_PRES
then q contains a scalar pressure at each element node. If the loadtype isVFE_DISTLOAD_TRAC
then q contains a vector traction at each element node. The pressure and traction units are force/length.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_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_COMPUTE
is generated if a zero edge Jacobian is computed.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
enttype – Entity type on which load is applied
=SYS_EDGE Element edge
loadtype – Distributed load type
x=VFE_DISTLOAD_TRAC Load directed along vector =VFE_DISTLOAD_PRES Load directed along inward normal to element face
q – Vector of distributed load values
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Shell2DConcLoad(vfe_Shell2D *p, Vdouble x[][3], Vint enttype, 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 isVFE_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 bothVFE_CONCLOAD_TRAC
andVFE_CONCLOAD_PRES
.- Errors
SYS_ERROR_ENUM
is generated if an improper enttype or loadtype is specified.SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.SYS_ERROR_COMPUTE
is generated if the location xc cannot be projected onto the surface defined by the element node locations.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
enttype – Entity type on which load is applied
=SYS_EDGE Element edge
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_Shell2DElemLoad(vfe_Shell2D *p, Vdouble x[][3], Vdouble q[][3], Vdouble f[])
body force vector
Compute the consistent degree of freedom body loads given acceleration load vector, q, on an element. The vector q contains an acceleration vector for for each node in the element. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The input element loads are in the units of force per unit mass. Note that the computation of consistent loads uses the material density.
- Errors
SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.
- Parameters:
p – Pointer to Shell2D object.
x – Array of node locations.
q – Array of node accelerations
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Shell2DNumIntPnt(vfe_Shell2D *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 Shell2D object.
analysistype – The type of analysis
x=VFE_ANALYSIS_STRUCTURAL Structural analysis
nepnts – [out] Number of element integration points
-
void vfe_Shell2DDirCos(vfe_Shell2D *p, Vdouble x[][3], Vdouble u[], Vdouble tm[][3][3])
compute shell local direction cosines
Compute the direction cosine matrices of the shell 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_Shell2DStrsStrn()
orvfe_Shell2DStressStrain()
are expressed. Given that X’,Y’ and Z’ are three orthonormal vectors indicating the direction of the local coordinate axes in the global coordinate system (x,y,z), then the direction cosine matrix, tm for this local coordinate system is defined as:The element degrees of freedom, u, are only required if large rotations have been enabled. Otherwise this argument is ignored.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
- Parameters:
p – Pointer to Shell2D object.
x – Array of point locations defining shell surface.
u – Degree of freedom vector of displacements
tm – [out] Array of direction cosine matrices at the shell nodes.
-
void vfe_Shell2DStrsAdapt(vfe_Shell2D *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 stress resultants, strss. In addition useful quantities such as the characteristic length, effective polynomial order and dimension of the element are returned. The element dimension, d, is 1. These quantities are useful for computing new characteristic element length for mesh adaptation.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell2D 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
3.4. 3D Shell Elements - Shell3D
The Shell3D module is used to simulate 3D shells. 3D Shell elements are two-dimensional elements in 3D space which support in-plane deformation, bending and tranverse shear through the thickness of the shell. The stress normal to the shell surface is assumed to be zero.
The methods associated with a Shell3D object are the following.
Begin and end an instance of an object, generic object functions
vfe_Shell3DBegin()
- create an instance of a Shell3D objectvfe_Shell3DEnd()
- destroy an instance of a Shell3D objectvfe_Shell3DError()
- return Shell3D object error flag
Attributes and Parameters
vfe_Shell3DSetHistPtr()
- set pointers to material historyvfe_Shell3DSetLocalSystem()
- set stress axes directionvfe_Shell3DSetMatlSystem()
- set material axes directionvfe_Shell3DSetObject()
- set attribute objectvfe_Shell3DSetParamd()
- set element formulation parametersvfe_Shell3DSetParami()
- set element formulation parametersvfe_Shell3DSetPropPtr()
- set pointer to element nodal propertiesvfe_Shell3DSetTopology()
- set input element topology
Degree of Freedom and Integration Information
vfe_Shell3DDofMap()
- query element degree of freedom mapvfe_Shell3DNumDof()
- query number of element degrees of freedomvfe_Shell3DNumIntPnt()
- query number of element integration points
Structural analysis computations
vfe_Shell3DBMatrix()
- B(beta) matrixvfe_Shell3DConcLoad()
- concentrated load vectorvfe_Shell3DDirCos()
- compute shell local direction cosinesvfe_Shell3DDistLoad()
- distributed load vector.vfe_Shell3DElemLoad()
- body force vectorvfe_Shell3DGeomStiff()
- geometric stiffness matrixvfe_Shell3DInitHist()
- initialize material historyvfe_Shell3DMass()
- consistent mass matrixvfe_Shell3DMassDiag()
- diagonal mass matrixvfe_Shell3DReact()
- reaction vectorvfe_Shell3DReactStiff()
- reaction vector, stiffness matrixvfe_Shell3DStiff()
- linear stiffness matrixvfe_Shell3DStrsAdapt()
- stress based error analysisvfe_Shell3DStrsStrn()
- resultant stress and strain/curvaturevfe_Shell3DStressStrain()
- three-dimensional stress and strain
3.5. Element Geometry
The basic shell element geometry is defined by the node coordinates of the shell
reference surface.
The shell wall is defined by the node offsets of the shell midsurface
from the reference surface, shell thicknesses and shell normals.
These properties which are set using vfe_Shell3DSetPropPtr()
are illustrated below in Figure 3-2.
Other element quantities which may be set using this function are
shell temperature and effective shear factor.
If the node offsets are zero, the shell midsurface is coincident with the shell reference surface. In general, large midsurface offset values, which are much greater than the shell thickness, should be avoided. By default the shell node normals are computed to be normal to the shell reference surface. The user may optionally enter node normals to achieve a “miter” joint between adjacent shell elements. Shell thickness and offset are always measured in the direction of the normal to the shell reference surface regardless of the direction of the shell node normals. The element natural coordinates, r, s, lie in the shell midsurface.
The element technology to be used is set using the
function vfe_Shell2DSetParami()
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 VFE_TECH_ISOP (1),
Parabolic triangle VFE_TECH_ANS. The fully integrated technology VFE_TECH_ISOP (3), is overly stiff.
Linear quadrilateral, VFE_TECH_MIXED or VFE_TECH_ANS. 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.
The exact nature of the shell wall construction which couples element
and primitive material properties,
is managed using the ShellProp module.
The ShellProp module provides for simple monocoque,
vfe_ShellPropSetMonocoque()
,
and general integrated shell walls, vfe_ShellPropSetIntegrate()
, as well as
preintegrated shell wall constructions input as a general matrix,
vfe_ShellPropSetMatrix()
, relating
midsurface strains and curvatures to stress and moment resultants.
For monocoque and general integrated shell walls the ShellProp module
allows primitive materials, such as defined by the LinMat
module, to be specified at the integration points.
Note that when using the general integrated shell wall defined by
vfe_ShellPropSetIntegrate()
,
the actual distance of the shell wall integration point from the
midsurface is computed internally by multiplying the integration point
coordinates input in vfe_ShellPropSetIntegrate()
by one-half of the
thickness specified using vfe_Shell3DSetPropPtr()
.
Therefore, if the actual shell thickness is specified using
vfe_Shell3DSetPropPtr()
, the shell wall integration point
coordinates input in vfe_ShellPropSetIntegrate()
should be in the range
[-1.,1.] from the bottom to the top shell surface, with the coordinate
0. lying on the midsurface.
Use the function
vfe_Shell3DStrsAdapt()
to aid in computing element strain energy, strain
energy error and other useful quantities to aid in solution error estimation
and mesh adaptation.
3.6. Function Descriptions
The currently available Shell3D functions are described in detail in this section.
-
vfe_Shell3D *vfe_Shell3DBegin(void)
create an instance of a Shell3D object
Create an instance of a Shell3D object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a Shell3D object using
void vfe_Shell3DEnd (vfe_Shell3D *shell3d)
Return the current value of a Shell3D object error flag using
Vint vfe_Shell3DError (vfe_Shell3D *shell3d)
- Returns:
The function returns a pointer to the newly created Shell3D object. If the object creation fails, NULL is returned.
-
void vfe_Shell3DEnd(vfe_Shell3D *p)
destroy an instance of a Shell3D object
-
Vint vfe_Shell3DError(vfe_Shell3D *p)
return the current value of a Shell3D object error flag
-
void vfe_Shell3DSetObject(vfe_Shell3D *p, Vint objecttype, Vobject *object)
set attribute object
Set a pointer to an attribute object. Users must supply a MatlFun object prior to computing any quantity that requires a material model definition.
- Errors
SYS_ERROR_OBJECTTYPE
is generated if an improper objecttype is specified.
- Parameters:
p – Pointer to Shell3D object.
objecttype – The object type identifier
x=VFE_MATLFUN MatlFun object
object – Pointer to the object to be set.
-
void vfe_Shell3DSetParami(vfe_Shell3D *p, Vint type, Vint iparam)
set element formulation parameters
Set element technology and formulation parameters.
The basic element technology is set with the
VFE_TECH
parameter. By defaultVFE_TECH
is set toVFE_TECH_ISOP
.Toggle large rotations using
VFE_LARGEROTATION
. By defaultVFE_LARGEROTATION
is set toSYS_OFF
.The addition of drilling stiffness to the shell element stiffness matrix is set with the
VFE_DRILLSTIFF
parameter. By defaultVFE_DRILLSTIFF
is set toSYS_ON
.The addition of rotary inertias to the element mass matrix is set with the
VFE_ROTINERTIA
parameter. By defaultVFE_ROTINERTIA
is set toSYS_OFF
.The parameter
VFE_TEMPMATLAVG
toggles the method for computing the temperature used for evaluating temperature dependent material properties. If enabled, the temperature used for temperature dependent material properties is the average of the element node point temperatures. If disabled, the temperature is isoparametrically interpolated from the node point temperatures at each element integration point. By defaultVFE_TEMPMATLAVG
is set toSYS_ON
.The parameter
VFE_CHECKCORNNORM
toggles a comparison between computed normals at integration points in shells and a computed normal interpolated from corner node values. If this parameter is enabled then the dot product of these two normals is computed and, if negative, an error is flagged. By defaultVFE_CHECKCORNNORM
is set toSYS_ON
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper iparam is specified.
- Parameters:
p – Pointer to Shell3D object.
type – Type of formulation parameter to set
x=VFE_TECH Element technology =VFE_LARGEROTATION Toggle for large rotations =VFE_DRILLSTIFF Enable drill stiffness =VFE_ROTINERTIA Enable rotary inertias =VFE_TEMPMATLAVE Average material temperature flag =VFE_CHECKCORNNORM Corner normal checking
iparam – Integer parameter value.
x=VFE_TECH_ISOP Standard isoparametric technology =VFE_TECH_URED Uniform reduced technology =VFE_TECH_MIXED Mixed technology =VFE_TECH_ANS Assumed natural strains technology =VFE_TECH_KIRCHHOFF Kirchhoff technology =SYS_OFF Disable =SYS_ON Enable
-
void vfe_Shell3DSetParamd(vfe_Shell3D *p, Vint type, Vdouble dparam)
set element formulation parameters
Set element formulation parameters.
The factor applied to the drilling stiffness in the shell element stiffness matrix is set with the
VFE_DRILLFACTOR
parameter. By defaultVFE_DRILLFACTOR
is set to 1. If unique shell normal directions are supplied then the drill factor can be reduced to a value of .01 . The addition of drill stiffness is enabled usingvfe_Shell3DSetParami()
.Use
VFE_MAXPROJANG
to set the maximum angle in degrees between the shell normal at an element integration point and the shell normal at a node. By defaultVFE_MAXPROJANG
is set to 90. degrees.Use
VFE_ADDEDMASS
to add non-structural mass/area to the shell element. By defaultVFE_ADDEDMASS
is set to 0..- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper dparam is specified.
- Parameters:
p – Pointer to Shell3D object.
type – Type of formulation parameter to set
x=VFE_DRILLFACTOR Drill stiffness factor =VFE_MAXPROJANG Maximum nodal projection angle =VFE_ADDEDMASS Additional mass/area
dparam – Double precision parameter value.
-
void vfe_Shell3DSetTopology(vfe_Shell3D *p, Vint shape, Vint maxi, Vint maxj)
set element topology
Specify the topology of a 3D shell element. If maxi is set to 3 then a quadratic element form is assumed. The default topology is shape =
SYS_SHAPEQUAD
, 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_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 Shell3D 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 <= 3 and 2 <= maxi <= 3, then a Lagrange finite element is assumed. If maxi != 0 then maxj = 0 or maxj = maxi.
-
void vfe_Shell3DSetPropPtr(vfe_Shell3D *p, Vint type, Vdouble *propptr)
set pointer to element nodal properties
Set a pointer to the start of a specified type of element nodal properties. Note that the properties are not copied by this function, only the pointer itself is copied. If a property pointer is not set the element assumes a default value for the associated property. By default the temperature is 0., the thickness is 2., the offset is 0. and the effective shear factor is 5/6. If the element node normals are not specified, the normals are computed to be normal to the shell reference surface. By default the reference temperature is 0.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to Shell3D object.
type – Type of element property
x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_TEMPREF Reference temperatures =VFE_PROP_THICKNESS Thicknesses =VFE_PROP_NORMAL Normals =VFE_PROP_OFFSET Offsets =VFE_PROP_KS Effective shear factor
propptr – Pointer to start of element nodal properties
-
void vfe_Shell3DSetHistPtr(vfe_Shell3D *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 total number of integration points in the element. The number of history variables depends on the underlying material model and may be queried using
vfe_MatlFunNumHist()
. The total number of integration points is the product of the number of element integration points returned usingvfe_Shell3DNumIntPnt()
and the number of section integration points returned usingvfe_MatlFunNumStressStrain()
. 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 Shell3D object.
oldhist – Pointer to start of material history at previous step
newhist – Pointer to start of material history at current step
-
void vfe_Shell3DStiff(vfe_Shell3D *p, Vdouble x[][3], Vdouble kl[])
linear stiffness matrix
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
-
void vfe_Shell3DGeomStiff(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vdouble kg[])
geometric stiffness matrix
Compute the geometric stiffness matrix, kg, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the geometric stiffness is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kg – [out] Degree of freedom geometric stiffness matrix
-
void vfe_Shell3DReact(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
reaction vector
Compute the reaction vector, r, given the node coordinates, x, and the degree of freedom displacement vector, u.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
r – [out] Degree of freedom reaction vector
-
void vfe_Shell3DReactStiff(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
reaction vector, stiffness matrix
Compute the reaction vector, r, and optionally the stiffness matrix, k, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D 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_Shell3DStrsStrn(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[])
resultant stress and strain/curvature
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 shell local coordinate system. The convention used to generate local coordinate systems is specified using
vfe_Shell3DSetLocalSystem()
. The actual direction cosine matrices of the shell local systems at the element output locations may be returned usingvfe_Shell3DDirCos()
.The strs is composed of 8 stress resultants and the strn is composed of 8 associated midsurface strains and curvatures at each output location. The stress and strain values are ordered first by the 8 stress resultant or strain components followed by the the number of element nodes. For example, for a 4 node shell element a total of 8*4 = 32 stress values will be returned in strs and 32 strain values returned in strn.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D 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_Shell3DStressStrain(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vdouble strs[], Vdouble strn[], Vdouble stress[], Vdouble strain[])
three-dimensional stress and strain
Compute nodal stresses and strains, strs and strn, and the three-dimensional stresses and strains, stress and strain, given the node coordinates, x, and the degree of freedom displacement vector, u.
The stress resultants and strain/curvatures are the same as those obtained with
vfe_Shell3DStrsStrn()
, and the same remarks regarding local coordinate systems apply.The three-dimensional stresses and strains are computed at a series of points in the beam’s cross section. The number of cross sectional points is obtained using
vfe_MatlFunNumStressStrain()
and depend on the type of cross section used. See ShellProp for more details on each cross section type.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
strs – [out] Array of nodal stress resultants
strn – [out] Array of nodal strains/curvatures
stress – [out] Array of nodal three dimensional stresses at recovery points
strain – [out] Array of nodal three dimensional strains at recovery points
-
void vfe_Shell3DInitHist(vfe_Shell3D *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 Shell3D object.
-
void vfe_Shell3DNumDof(vfe_Shell3D *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_Shell3DDofMap()
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 Shell3D 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_Shell3DDofMap(vfe_Shell3D *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, rotation, temperature, etc.
The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_Shell3DNumDof()
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 Shell3D 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_Shell3DMass(vfe_Shell3D *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_Shell3DMassDiag()
to compute a diagonal mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
m – [out] Degree of freedom consistent mass matrix
-
void vfe_Shell3DMassDiag(vfe_Shell3D *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_Shell3DMass()
to compute a consistent mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
md – [out] Degree of freedom diagonal mass vector
-
void vfe_Shell3DDistLoad(vfe_Shell3D *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 3D shell 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_PRES
then q contains a scalar pressure at each element node. The inward normal to an edge is formed by the cross product of the tangent vector along the edge with the inward normal to the shell. If the loadtype isVFE_DISTLOAD_TRAC
then q contains a vector traction at each element node. If the loadtype isVFE_DISTLOAD_TANGFORCE
then q contains a scalar force/length at each element node. If the loadtype isVFE_DISTLOAD_TANGMOMENT
then q contains a scalar moment/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 edge or face of interest. Correspondingly the output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The distributed loads for edges are in units of force per unit length. The distributed loads for faces are in units of force per unit area.- Errors
SYS_ERROR_ENUM
is generated if an improper enttype or loadtype is specified.SYS_ERROR_OPERATION
is generated if an invalid combination of enttype and loadtype is specified.SYS_ERROR_VALUE
is generated if an improper no is specified.SYS_ERROR_COMPUTE
is generated if a zero edge or face Jacobian is computed.
- Parameters:
p – Pointer to Shell3D 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 number
loadtype – Distributed load type
x=VFE_DISTLOAD_TRAC Load directed along vector to element face or element edge =VFE_DISTLOAD_PRES Load directed along inward normal to element face or element edge =VFE_DISTLOAD_TANGFORCE Load directed along element edge =VFE_DISTLOAD_TANGMOMENT Moment directed about element edge
q – Vector of distributed load values
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Shell3DConcLoad(vfe_Shell3D *p, Vdouble x[][3], Vint enttype, 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 face 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 isVFE_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 bothVFE_CONCLOAD_TRAC
andVFE_CONCLOAD_PRES
.- Errors
SYS_ERROR_ENUM
is generated if an improper enttype or loadtype is specified.SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.SYS_ERROR_COMPUTE
is generated if the location xc cannot be projected onto the surface defined by the element node locations.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
enttype – Entity type on which load is applied
=SYS_FACE Element face
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_Shell3DElemLoad(vfe_Shell3D *p, Vdouble x[][3], Vdouble q[][3], Vdouble f[])
body force vector
Compute the consistent degree of freedom body loads given acceleration load vector, q, on an element. The vector q contains an acceleration vector for for each node in the element. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The input element loads are in the units of force per unit mass. Note that the computation of consistent loads uses the material density.
- Errors
SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.
- Parameters:
p – Pointer to Shell3D object.
x – Array of node locations.
q – Array of node accelerations
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Shell3DNumIntPnt(vfe_Shell3D *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 Shell3D 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_Shell3DSetLocalSystem(vfe_Shell3D *p, Vint type, Vdouble vec[], Vdouble angle)
set local coordinate system convention
Specify the convention to be used to construct the orientation of the shell element local x’,y’,z’ coordinate system with respect to the shell reference surface. For stress and strain computation for output using
vfe_Shell3DStrsStrn()
orvfe_Shell3DStressStrain()
, 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 z’ axis is always constructed to be normal to the shell reference surface. The orientation of the x’ and y’ axes in the plane of the shell reference surface is determined by type. The vec array is only used if the specified type requires position or direction vectors. An additional rotation of the x’,y’ axes about the z’ 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 Shell3D object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate shell x’,y’ axes about the shell z’ axis in degrees.
-
void vfe_Shell3DDirCos(vfe_Shell3D *p, Vdouble x[][3], Vdouble u[], Vdouble tm[][3][3])
compute shell local direction cosines
Compute the direction cosine matrices of the shell 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_Shell3DStrsStrn()
orvfe_Shell3DStressStrain()
are expressed. Given that X’,Y’ and Z’ are three orthonormal vectors indicating the direction of the local coordinate axes in the global coordinate system (x,y,z), then the direction cosine matrix, tm for this local coordinate system is defined as:The element degrees of freedom, u, are only required if large rotations have been enabled. Otherwise this argument is ignored.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_Shell3DSetLocalSystem()
.- Parameters:
p – Pointer to Shell3D object.
x – Array of point locations defining shell surface.
u – Degree of freedom vector of displacements
tm – [out] Array of direction cosine matrices at the shell nodes.
-
void vfe_Shell3DSetMatlSystem(vfe_Shell3D *p, Vint type, Vdouble vec[], Vdouble angle)
set material system convention
Specify the convention to be used to construct the orientation of the shell element local material x’,y’,z’ coordinate system with respect to the shell reference surface. This local system is computed at each integration point location on the shell surface and is assumed to be the coordinate system in which the material properties of the element at the integration point are expressed.
The z’ axis is always constructed to be normal to the shell reference surface. The orientation of the x’ and y’ axes in the plane of the shell reference surface is determined by type. The vec array is only used if the specified type requires position or direction vectors. An additional rotation of the x’,y’ axes about the z’ 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 Shell3D object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate shell x’,y’ axes about the shell z’ axis in degrees.
-
void vfe_Shell3DStrsAdapt(vfe_Shell3D *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 stress resultants, 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_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Shell3D 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_Shell3DBMatrix(vfe_Shell3D *p, Vdouble x[][3], Vint *nipt, Vdouble rst[][2], Vdouble h[], Vdouble phr[][2], Vdouble w[], Vdouble dj[], Vdouble tm[][3][3], Vdouble dm[][36], Vdouble b[][8])
B(beta) matrix.
Compute the B(beta) matrix relating the strain components (exx,eyy,gxy,kxx,kyy,kxy,g_xz,g_yz) to the element degrees of freedom. Use
vfe_Shell3DNumDof()
to determine the number of element degrees of freedom, andvfe_Shell3DDofMap()
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) wheree 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 integration point weights w at each integration point; the determinant of the Jacobian transformation d(x,y,z)/d(r,s); 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 8x8 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_OPERATION
is generated if element is improperly configured.SYS_ERROR_COMPUTE
is generated if a negative determinant of the Jacobian is detected.
- Parameters:
p – Pointer to Shell3D 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
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
3.7. 3D Membrane Elements - Mem3D
The Mem3D module is used to simulate 3D membranes. 3D membrane elements are two-dimensional elements in 3D space which support in-plane deformation only. The stress normal to the membrane surface is assumed to be zero. The transverse shear stresses are also zero.
The methods associated with a Mem3D object are the following.
Begin and end an instance of an object, generic object functions
vfe_Mem3DBegin()
- create an instance of a Mem3D objectvfe_Mem3DEnd()
- destroy an instance of a Mem3D objectvfe_Mem3DError()
- return Mem3D object error flag
Attributes and Parameters
vfe_Mem3DSetHistPtr()
- set pointers to material historyvfe_Mem3DSetLocalSystem()
- set stress axes directionvfe_Mem3DSetMatlSystem()
- set material axes directionvfe_Mem3DSetObject()
- set attribute objectvfe_Mem3DSetParamd()
- set element formulation parametersvfe_Mem3DSetParami()
- set element formulation parametersvfe_Mem3DSetPropPtr()
- set pointer to element nodal propertiesvfe_Mem3DSetTopology()
- set input element topology
Degree of Freedom and Integration Information
vfe_Mem3DDofMap()
- query element degree of freedom mapvfe_Mem3DNumDof()
- query number of element degrees of freedomvfe_Mem3DNumIntPnt()
- query number of element integration points
Structural analysis computations
vfe_Mem3DConcLoad()
- concentrated load vectorvfe_Mem3DDirCos()
- compute membrane local direction cosinesvfe_Mem3DDistLoad()
- distributed load vector.vfe_Mem3DElemLoad()
- body force vectorvfe_Mem3DGeomStiff()
- geometric stiffness matrixvfe_Mem3DInitHist()
- initialize material historyvfe_Mem3DMass()
- consistent mass matrixvfe_Mem3DMassDiag()
- diagonal mass matrixvfe_Mem3DReact()
- reaction vectorvfe_Mem3DReactStiff()
- reaction vector, stiffness matrixvfe_Mem3DStiff()
- linear stiffness matrixvfe_Mem3DStrsAdapt()
- stress based error analysisvfe_Mem3DStrsStrn()
- stress and strain
Thermal analysis computations
vfe_Mem3DDistHeat()
- distributed heat loads.vfe_Mem3DElemHeat()
- body heat generationvfe_Mem3DCap()
- consistent capacitance matrixvfe_Mem3DCapDiag()
- diagonal capacitance matrixvfe_Mem3DPower()
- thermal powervfe_Mem3DPowerCond()
- thermal power, conductance matrixvfe_Mem3DCond()
- conductance matrixvfe_Mem3DHFlxAdapt()
- heat flux based error analysisvfe_Mem3DHFlxTGrd()
- heat flux and temperature gradient
Instance a Mem3D object using vfe_Mem3DBegin()
.
Once a Mem3D is instanced,
set the material function attribute object MatlFun using
vfe_Mem3DSetObject()
.
The current topology of the element is specified using vfe_Mem3DSetTopology()
.
The thickness of the membrane element is specified using vfe_Mem3DSetPropPtr()
.
Query the element degree of freedom map using vfe_Mem3DDofMap()
and
vfe_Mem3DNumDof()
.
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_Mem3DSetHistPtr()
and vfe_Mem3DInitHist()
.
Since all element geometry is available directly to the element in terms
of node coordinates and thicknesses,
the material function should point to a primitive material.
3.8. Element Geometry
The basic membrane element geometry is defined
by the node coordinates of the membrane reference surface and the node thicknesses
There are no offset properties available for the membrane element, therefore
the membrane midsurface is coincident with the reference surface.
The node thickness properties are set using vfe_Mem3DSetPropPtr()
and are illustrated below in Figure 3-3.
Other element quantities which may be set using
this function are node temperature and reference temperature.
The z’ component of stress is zero as well as the transverse shear
stresses in the y’z’ and z’x’ planes.
The element natural coordinates, r, s, lie in the membrane midsurface.
Use the function
vfe_Mem3DStrsAdapt()
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_Mem3DHFlxAdapt()
performs a similar
computation for heat transfer analysis.
3.9. Function Descriptions
The currently available Mem3D functions are described in detail in this section.
-
vfe_Mem3D *vfe_Mem3DBegin(void)
create an instance of a Mem3D object
Create an instance of a Mem3D object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a Mem3D object using
void vfe_Mem3DEnd (vfe_Mem3D *mem3d)
Return the current value of a Mem3D object error flag using
Vint vfe_Mem3DError (vfe_Mem3D *mem3d)
- Returns:
The function returns a pointer to the newly created Mem3D object. If the object creation fails, NULL is returned.
-
void vfe_Mem3DEnd(vfe_Mem3D *p)
destroy an instance of a Mem3D object
See
vfe_Mem3DBegin()
-
Vint vfe_Mem3DError(vfe_Mem3D *p)
return the current value of a Mem3D object error flag
See
vfe_Mem3DBegin()
-
void vfe_Mem3DSetObject(vfe_Mem3D *p, Vint objecttype, Vobject *object)
set attribute object
Set a pointer to an attribute object. Users must supply a MatlFun object prior to computing any quantity that requires a material model definition.
- Errors
SYS_ERROR_OBJECTTYPE
is generated if an improper objecttype is specified.
- Parameters:
p – Pointer to Mem3D object.
objecttype – The object type identifier
x=VFE_MATLFUN MatlFun object
object – Pointer to the object to be set.
-
void vfe_Mem3DSetParami(vfe_Mem3D *p, Vint type, Vint iparam)
set element formulation parameters
Set element technology and formulation parameters.
The basic element technology is set with the
VFE_TECH
parameter. By defaultVFE_TECH
is set toVFE_TECH_ISOP
.Set element strain type using
VFE_STRAINTYPE
with a value of eitherVFE_LARGESTRAIN
to enable large strain orVFE_SMALLSTRAIN
to enable small strains. By defaultVFE_STRAINTYPE
is set toVFE_SMALLSTRAIN
. If a small strain material is used underVFE_LARGESTRAIN
then it is assumed that the material law relates the Green-Lagrange strain tensor with the second second Piola-Kirchhoff stress tensor. Otherwise, the material law relates the deformation gradient with the Cauchy stress.The parameter
VFE_TEMPMATLAVG
toggles the method for computing the temperature used for evaluating temperature dependent material properties. If enabled, the temperature used for temperature dependent material properties is the average of the element node point temperatures. If disabled, the temperature is isoparametrically interpolated from the node point temperatures at each element integration point. By defaultVFE_TEMPMATLAVG
is set toSYS_ON
.The parameter
VFE_CHECKCORNNORM
toggles a comparison between computed normals at integration points in shells and a computed normal interpolated from corner node values. If this parameter is enabled then the dot product of these two normals is computed and, if negative, an error is flagged. By defaultVFE_CHECKCORNNORM
is set toSYS_ON
.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper iparam is specified.
- Parameters:
p – Pointer to Mem3D object.
type – Type of formulation parameter to set
x=VFE_TECH Element technology =VFE_STRAINTYPE Element strain type =VFE_TEMPMATLAVE Average material temperature flag =VFE_CHECKCORNNORM Corner normal checking
iparam – Integer parameter value.
x=VFE_TECH_ISOP Standard isoparametric technology =VFE_TECH_URED Uniform reduced technology =VFE_SMALLSTRAIN Small strain =VFE_LARGESTRAIN Large strain =SYS_OFF Disable =SYS_ON Enable
-
void vfe_Mem3DSetParamd(vfe_Mem3D *p, Vint type, Vdouble dparam)
set element formulation parameters
Set element formulation parameters.
Use
VFE_MAXPROJANG
to set the maximum angle in degrees between the membrane surface normal at an element integration point and the membrane normal at a node. By defaultVFE_MAXPROJANG
is set to 90. degrees.Use
VFE_ADDEDMASS
to add non-structural mass/area to the membrane element. By defaultVFE_ADDEDMASS
is set to 0..- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_VALUE
is generated if an improper dparam is specified.
- Parameters:
p – Pointer to Mem3D object.
type – Type of formulation parameter to set
x=VFE_MAXPROJANG Maximum nodal projection angle =VFE_ADDEDMASS Additional mass/area
dparam – Double precision parameter value.
-
void vfe_Mem3DSetTopology(vfe_Mem3D *p, Vint shape, Vint maxi, Vint maxj)
set element topology
Specify the topology of a 3D membrane element. If maxi is set to 3 then a quadratic element form is assumed. The default topology is shape =
SYS_SHAPEQUAD
, maxi = maxj = 0.- 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 Mem3D 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 <= 3 and 2 <= maxi <= 3, then a Lagrange finite element is assumed. If maxi != 0 then maxj = 0 or maxj = maxi.
-
void vfe_Mem3DSetPropPtr(vfe_Mem3D *p, Vint type, Vdouble *propptr)
set pointer to element nodal properties
Set a pointer to the start of a specified type of element nodal properties. Note that the properties are not copied by this function, only the pointer itself is copied. If a property pointer is not set the element assumes a default value for the associated property. By default the thickness is 1., the temperature is 0. and the reference temperature is 0. If a pre-stress is set it is added to the computed stress resulting from the material model.
- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.
- Parameters:
p – Pointer to Mem3D object.
type – Type of element property
x=VFE_PROP_TEMPERATURE Temperatures =VFE_PROP_TEMPREF Reference temperatures =VFE_PROP_THICKNESS Thicknesses =VFE_PROP_PRESTRESS Pre-stress
propptr – Pointer to start of element nodal properties
-
void vfe_Mem3DSetHistPtr(vfe_Mem3D *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 total number of integration points in the element. The number of history variables depends on the underlying material model and may be queried using
vfe_MatlFunNumHist()
. The total number of integration points is the product of the number of element integration points returned usingvfe_Mem3DNumIntPnt()
and the number of section integration points returned usingvfe_MatlFunNumStressStrain()
. 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 Mem3D object.
oldhist – Pointer to start of material history at previous step
newhist – Pointer to start of material history at current step
-
void vfe_Mem3DStiff(vfe_Mem3D *p, Vdouble x[][3], Vdouble kl[])
linear stiffness matrix
Compute the linear stiffness matrix, kl, given the node coordinates, x. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
kl – [out] Degree of freedom stiffness matrix
-
void vfe_Mem3DCond(vfe_Mem3D *p, Vdouble x[][3], Vdouble kl[])
thermal conductance matrix
Compute the linear conductance matrix, kl, given the node coordinates, x. The lower triangle of the conductance matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
kl – [out] Degree of freedom conductance matrix
-
void vfe_Mem3DGeomStiff(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vdouble kg[])
geometric stiffness matrix
Compute the geometric stiffness matrix, kg, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the geometric stiffness is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
kg – [out] Degree of freedom geometric stiffness matrix
-
void vfe_Mem3DReact(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
reaction vector
Compute the reaction vector, r, given the node coordinates, x, and the degree of freedom displacement vector, u.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
u – Degree of freedom vector of displacements
r – [out] Degree of freedom reaction vector
-
void vfe_Mem3DPower(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vdouble r[])
thermal power vector
Compute the power vector, r, given the node coordinates, x, and the degree of freedom temperature vector, u.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
u – Degree of freedom vector of temperatures
r – [out] Degree of freedom power vector
-
void vfe_Mem3DReactStiff(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
reaction vector, stiffness matrix
Compute the reaction vector, r, and optionally the stiffness matrix, k, given the node coordinates, x, and the degree of freedom displacement vector, u. The lower triangle of the stiffness matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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_Mem3DPowerCond(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vint kflag, Vdouble r[], Vdouble k[])
thermal power, conductance matrix
Compute the power vector, r, and optionally the conductance matrix, k, given the node coordinates, x, and the degree of freedom temperature vector, u. The lower triangle of the conductance matrix is returned.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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_Mem3DStrsStrn(vfe_Mem3D *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 membrane local coordinate system. The convention used to generate local coordinate systems is specified using
vfe_Mem3DSetLocalSystem()
. The actual direction cosine matrices of the membrane local systems at the element output locations may be returned usingvfe_Mem3DDirCos()
.The strs is composed of 6 stresses and the strn is composed of 6 associated strains at each output location. The stress and strain values are ordered first by the 6 tensor components followed by the the number of element nodes. For example, for a 4 node membrane element a total of 6*4 = 24 stress values will be returned in strs and 24 strain values returned in strn.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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_Mem3DHFlxTGrd(vfe_Mem3D *p, Vdouble x[][3], Vdouble u[], Vdouble hflx[], Vdouble tgrd[])
heat flux and thermal gradient
Compute nodal heat fluxes and temperature gradients, hflx and tgrd, given the node coordinates, x, and the degree of freedom temperature vector, u. The flux and gradients are computed in the membrane local system. The convention used to generate local coordinate systems is specified using
vfe_Mem3DSetLocalSystem()
.The flux and gradient values are ordered first by the 3 vectoral components followed by the number of element nodes.
- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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_Mem3DInitHist(vfe_Mem3D *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 Mem3D object.
-
void vfe_Mem3DNumDof(vfe_Mem3D *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_Mem3DDofMap()
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 Mem3D 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_Mem3DDofMap(vfe_Mem3D *p, Vint analysistype, Vint loc[], Vint tag[])
query element degree of freedom map
Query for element degree of freedom map. The degree of freedom map consists of a location index, loc and type, tag for each degree of freedom used by the element.
The location index is either a positive node index into the element connectivity indicating a nodal freedom or a zero value indicating an elemental degree of freedom. The tag indicates the type of the degree of freedom. Tag values are one of a set of enumerated types which indicate whether the degree of freedom is a translation or temperature.
The length of the loc and tag vectors is equal to the number of element degrees of freedom. Use
vfe_Mem3DNumDof()
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 Mem3D 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_Mem3DMass(vfe_Mem3D *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_Mem3DMassDiag()
to compute a diagonal mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
m – [out] Degree of freedom consistent mass matrix
-
void vfe_Mem3DMassDiag(vfe_Mem3D *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_Mem3DMass()
to compute a consistent mass matrix.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
md – [out] Degree of freedom diagonal mass vector
-
void vfe_Mem3DCapDiag(vfe_Mem3D *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_Mem3DCap()
to compute a consistent capacitance matrix. This calculation requires material density and specific heat.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
u – Degree of freedom vector of temperatures
cd – [out] Degree of freedom diagonal capacitance matrix
-
void vfe_Mem3DCap(vfe_Mem3D *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_Mem3DCapDiag()
to compute a diagonal capacitance matrix. This calculation requires material density and specific heat.- Errors
SYS_ERROR_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
u – Degree of freedom vector of temperatures
c – [out] Degree of freedom consistent capacitance matrix
-
void vfe_Mem3DDistLoad(vfe_Mem3D *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 3D membrane 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_PRES
then q contains a scalar pressure at each element node. If the loadtype is The inward normal to an edge is formed by the cross product of the tangent vector along the edge with the inward normal to the shell. If the loadtype isVFE_DISTLOAD_TRAC
then q contains a vector traction at each element node. If the loadtype isVFE_DISTLOAD_TANGFORCE
then q contains a scalar force/length at each element node. If the loadtype is Note that the input array of node locations, x, contains the coordinate locations for all nodes in the element - not just on the edge or face of interest. Correspondingly the output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The distributed loads for edges are in units of force per unit length. The distributed loads for faces are in units of force per unit area.- Errors
SYS_ERROR_ENUM
is generated if an improper enttype or loadtype is specified.SYS_ERROR_OPERATION
is generated if an invalid combination of enttype and loadtype is specified.SYS_ERROR_VALUE
is generated if an improper no is specified.SYS_ERROR_COMPUTE
is generated if a zero edge or face Jacobian is computed.
- Parameters:
p – Pointer to Mem3D 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 number
loadtype – Distributed load type
x=VFE_DISTLOAD_TRAC Load directed along vector to element face or element edge =VFE_DISTLOAD_PRES Load directed along inward normal to element face or 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_Mem3DConcLoad(vfe_Mem3D *p, Vdouble x[][3], Vint enttype, 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 face 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 isVFE_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 bothVFE_CONCLOAD_TRAC
andVFE_CONCLOAD_PRES
.- Errors
SYS_ERROR_ENUM
is generated if an improper enttype or loadtype is specified.SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.SYS_ERROR_COMPUTE
is generated if the location xc cannot be projected onto the surface defined by the element node locations.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
enttype – Entity type on which load is applied
=SYS_FACE Element face
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_Mem3DDistHeat(vfe_Mem3D *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 membrane 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_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 Mem3D 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_Mem3DElemLoad(vfe_Mem3D *p, Vdouble x[][3], Vdouble q[][3], Vdouble f[])
body force vector
Compute the consistent degree of freedom body loads given acceleration load vector, q, on an element. The vector q contains an acceleration vector for for each node in the element. The output array of consistent degree of freedom loads, f, contains loads for all degrees of freedom in the element. The input element loads are in the units of force per unit mass. Note that the computation of consistent loads uses the material density.
- Errors
SYS_ERROR_COMPUTE
is generated if a zero Jacobian is computed.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
q – Array of node accelerations
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Mem3DElemHeat(vfe_Mem3D *p, Vdouble x[][3], Vdouble q[], Vdouble f[])
body heat generation
Compute the consistent degree of freedom body heat generation given nodal heat generation vector, q on an element. The vector q contains heat generation per volume for each node in the element. The output array of consistent degree of freedom loads, f, contains the heat generation in the element. The input element loads are in the units of power per unit volume.
- Errors
SYS_ERROR_COMPUTE
is generated if a non-positive Jacobian is computed.
- Parameters:
p – Pointer to Mem3D object.
x – Array of node locations.
q – Array of node heat fluxes
f – [out] Degree of freedom vector of consistent loads.
-
void vfe_Mem3DNumIntPnt(vfe_Mem3D *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 Mem3D 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_Mem3DSetLocalSystem(vfe_Mem3D *p, Vint type, Vdouble vec[], Vdouble angle)
set local coordinate system convention
Specify the convention to be used to construct the orientation of the membrane element local x’,y’,z’ coordinate system with respect to the membrane reference surface. For stress and strain computation for output using
vfe_Mem3DStrsStrn()
, 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 z’ axis is always constructed to be normal to the membrane reference surface. The orientation of the x’ and y’ axes in the plane of the membrane reference surface is determined by type. The vec array is only used if the specified type requires position or direction vectors. An additional rotation of the x’,y’ axes about the z’ 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 Mem3D object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate membrane x’,y’ axes about the membrane z’ axis in degrees.
-
void vfe_Mem3DDirCos(vfe_Mem3D *p, Vdouble x[][3], Vdouble tm[][3][3])
compute membrane local direction cosines
Compute the direction cosine matrices of the membrane 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_Mem3DStrsStrn()
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_Mem3DSetLocalSystem()
.- Parameters:
p – Pointer to Mem3D object.
x – Array of point locations defining membrane surface.
tm – [out] Array of direction cosine matrices at the membrane nodes.
-
void vfe_Mem3DSetMatlSystem(vfe_Mem3D *p, Vint type, Vdouble vec[], Vdouble angle)
set material system convention
Specify the convention to be used to construct the orientation of the membrane element local material x’,y’,z’ coordinate system with respect to the membrane reference surface. This local system is computed at each integration point location on the membrane surface and is assumed to be the coordinate system in which the material properties of the element at the integration point are expressed.
The z’ axis is always constructed to be normal to the membrane reference surface. The orientation of the x’ and y’ axes in the plane of the membrane reference surface is determined by type. The vec array is only used if the specified type requires position or direction vectors. An additional rotation of the x’,y’ axes about the z’ 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 Mem3D object.
type – Local system convention
vec – Orientation vector data
angle – Angle to rotate membrane x’,y’ axes about the membrane z’ axis in degrees.
-
void vfe_Mem3DStrsAdapt(vfe_Mem3D *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 stress resultants, 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_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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_Mem3DHFlxAdapt(vfe_Mem3D *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_NULLOBJECT
is generated if a MatlFun attribute object has not been set.SYS_ERROR_COMPUTE
is generated if a negative Jacobian transformation is computed or the maximum nodal projection angle is exceeded.
- Parameters:
p – Pointer to Mem3D 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