2. Degree of Freedom and Connection Management - DofTab
The DofTab module is designed to manage the assignment of global degrees of freedom (dofs) arising from the assembly of individual finite elements. It is also responsible for storing model information such as nodal coordinates and element connectivity that are necessary for the iterative solver in SysMatrix. The DofTab module will manage both nodal and elemental (internal element) degrees of freedom. Once the global dof have been assigned then the global dof at a node or an element may be queried. This information is used to assemble element matrices and vectors into the global equation system for solution and to gather the resultant solution quantities for output. The methods associated with a DofTab object are the following.
Begin and end an instance of an object, generic object functions
vfs_DofTabBegin()
- create an instance of a DofTab objectvfs_DofTabEnd()
- destroy an instance of a DofTab objectvfs_DofTabError()
- return DofTab object error flag
Operations
vfs_DofTabDef()
- define number of nodes and elementsvfs_DofTabInq()
- inquire number of nodes and elementsvfs_DofTabElemDof()
- query element global dofvfs_DofTabElemSup()
- query element global suppressed dofvfs_DofTabElemUse()
- specify element degree of freedom usevfs_DofTabGetUseTag()
- query used degree of freedom tagsvfs_DofTabInvDof()
- query inverse degree of freedom informationvfs_DofTabMaxElemDof()
- query for maximum size of element global dofvfs_DofTabNodeDirCos()
- set direction cosine matrix for preconditionervfs_DofTabNodeDof()
- query node global dofvfs_DofTabNodeDofSpec()
- flag degree of freedom for preconditionervfs_DofTabNodeSpec()
- flag node for preconditionervfs_DofTabNodeSup()
- query node global suppressed dofvfs_DofTabNodePerm()
- set nodal degree of freedom constraintsvfs_DofTabPointersIndicesPtr()
- set pointers and indicesvfs_DofTabProcess()
- process and assign global dofvfs_DofTabSetCoords()
- set the nodal coordinatesvfs_DofTabSetFeaType()
- set the element typevfs_DofTabSetAssemBlock()
- set running assembling block idvfs_DofTabSetMode()
- set operating modesvfs_DofTabSetMap()
- set current element degree of freedom mapvfs_DofTabSetMPC()
- set multipoint constraintsvfs_DofTabSetNumThreads()
- set number of threadsvfs_DofTabSetParami()
- set integer parametersvfs_DofTabSetTopology()
- set element topologyvfs_DofTabThreadAssemNum()
- query number of assembly groupsvfs_DofTabThreadAssemGroups()
- query assembly groups
Instance a DofTab object using vfs_DofTabBegin()
.
Once a DofTab object is instanced,
set the number of nodes and elements for which degrees of freedom are
to be managed using vfs_DofTabDef()
.
Initially all node and element freedoms are assumed to be constrained.
Degrees of freedom are activated by declaring their use by a finite element.
By looping through all elements and setting the freedoms used by each, all
degrees of freedom required by the assembly of all elements will be activated.
This is accomplished by setting the degree of freedom map for each element
using vfs_DofTabSetMap()
and then declaring the degrees of freedom used
by calling vfs_DofTabElemUse()
.
After the degree of freedom use for all elements has been established,
permanent degree of freedom suppressions at a node are applied
using vfs_DofTabNodePerm()
. These suppressions reflect restraints
at nodes. Note that there is no mechanism in DofTab to suppress
internal elemental freedoms.
if multipoint constraints (MPCs) are present and they are to be
enforced by direct matrix transformation, then the
MPCs must be specified using vfs_DofTabSetMPC()
for each MPC equation.
At this point the DofTab object will have a representation
of all active degrees of freedom required by the finite element model.
The global active dof can now be numbered by calling
vfs_DofTabProcess()
. The global dof are numbered starting from
1. The total number of global active dof and the total number of global
suppressed dof are returned by this function.
The function vfs_DofTabProcess()
also separately numbers the global suppressed dof starting from 1.
The degree of freedom numbering is done to minimize the fill for matrix
factoring. The current fill reducing algorithm uses Metis.
Once the global dof have been assigned to the active freedoms,
the user can query for the global dof numbers at a single node using
vfs_DofTabNodeDof()
or for an entire
element using vfs_DofTabElemDof()
.
Calls to vfs_DofTabElemDof()
must be preceded by vfs_DofTabSetMap()
.
The complete global dof for an element consist of the
global dof assigned
to the nodal freedoms at each element node plus the global dof assigned
to the internal elemental freedoms for the element.
A global dof number of 0 returned by these functions represents a
suppressed degree of freedom.
The global suppressed dof at a node are returned using vfs_DofTabNodeSup()
.
The global suppressed dof for an element are returned using vfs_DofTabElemSup()
.
A global dof number of 0 returned by this function represents an active degree
of freedom.
The suppressed dof numbers are useful for assembling reactions.
The function, vfs_DofTabGetUseTag()
, returns the number and types, or “tags”,
of all the degrees of freedom used at nodes. For example in a finite element
model using only 3D solid elements which use translational degrees of freedom
at nodes, the number of degree of freedom types is 3 and the degree of
freedom tags are SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TZ.
Use the function vfs_DofTabInvDof()
to find the node or element ids and
degree of freedom tags associated with a set of global degrees of freedom.
Destroy an instance of a DofTab object using vfs_DofTabEnd()
.
2.1. Parallel Assembly
The DofTab and SysMatrix modules have been extensively
parallelized internally for performance. One important parallelization
opportunity requires parallelization at the application level and
that is assembly of element matrices in a SysMatrix object.
The currently supported parallelization methodology involves
computing a series of groups of elements which are disjoint
using an element coloring algorithm. These groups are computed
by DofTab if enabled using vfs_DofTabSetParami()
with
type DOFTAB_THREADASSEM.
All of the elements in each group can be safely assembled
in parallel. The groups are computed by the function vfs_DofTabProcess()
.
Once this function is complete,
use vfs_DofTabThreadAssemNum()
to query the number
of element groups and vfs_DofTabThreadAssemGroups()
to get
the number of elements and list of elements in each group.
It suggested to use OpenMP to parallelize the loop over elements
in each group.
The basic procedure is as follows where numel is the number of elements.
Vint num, hasserial;
Vint *numelem, *blockid, *elems;
vfs_DofTabThreadAssemNum (doftab,&num,&hasserial);
/* allocate arrays */
numelem = (Vint*)malloc(num*sizeof(Vint);
blockid = (Vint*)malloc(num*sizeof(Vint);
elems = (Vint*)malloc(numel*sizeof(Vint);
vfs_DofTabThreadAssemGroups (doftab,&num,numelem,blockid,elems);
The number of groups is num. The variable hasserial is non-zero if the last group is meant to be processed serially.
The array numelems contains the number of elements in each group, and the array elems contains the elements in each group. The element numbers in elems are one based. The first numelem[0] elements in elems are in the first group, the next numelem[1] elements are in the second group, etc. Note that if DOFTAB_THREADASSEM has not be enabled, num is returned as 1 and hasserial is enabled.
Assembly of elements can then proceed as follows. In this example 8 threads are used.
Vint n, m;
Vint in, istart;
Vint numthreads;
maxthreads = 8;
/* setup SysMatrix for parallel assembly */
vfs_SysMatrixSetNumThreadsAssem (sysmatrix,maxthreads);
/* loop over groups */
istart = 0;
for(n = 0; n < num; n++) {
/* special case of serial last group */
if(blockid[n] <= 0) {
numthreads = 1;
} else {
numthreads = maxthreads;
}
#pragma omp parallel for num_threads(numthreads)
/* loop over elements in group */
for(m = 0; m < numelem[n]; m++) {
Vint ithread;
Vint in, nedofs, lm[maximum_element_dofs];
Vint nix, ix[maximum_element_connectivity];
Vdouble s[maximum_element_matrix_size];
in = elems[istart+m];
/* get element connectivity, ix[nix] */
...
/* get number of element freedoms and freedom list */
vfs_DofTabElemDof (doftab,in,nix,ix,&nedofs,lm);
/* gather/compute element matrix s element in */
...
/* get thread number, increment by one */
ithread = omp_get_thread_num()+1;
/* assemble */
vfs_SysMatrixThreadAssem (sysmatrix,ithread,nedofs,lm,s);
}
istart += numelem[n];
}
/* assembly complete, free memory */
free(numelem);
free(elems);
2.2. Special Considerations for the Iterative Solver
If the iterative solver is to be used in
SysMatrix, all elements must have their topology specified with
vfs_DofTabSetTopology()
, and all nodes must have their coordinates specified
with vfs_DofTabSetCoords()
.
If there are local coordinate systems for unknowns at nodes, their
direction cosine matrices must be entered using vfs_DofTabNodeDirCos()
.
If there are temporary restraints to be applied in SysMatrix using
vfs_SysMatrixSetTemp()
then this information must be supplied
to DofTab using the function vfs_DofTabNodeDofSpec()
.
This function is a subset of the function vfs_DofTabNodeSpec()
.
The effect is that the degree of freedom(s) is/are exactly represented
in the preconditioner matrix. Any freedom which is temporarily constained
for any reason, eg. to apply a non-zero displacement, can be treated
exactly.
In addition vfs_DofTabSetParami()
must be used to set DOFTAB_REDUCE to DOFTAB_REDUCE_ITER.
A feature of the iterative solver is to flag particular nodes
to be given special treatment - the “stiffness” associated
with the node shall be represented exactly in the preconditioner at
the expense of making the preconditioner somewhat more heavy weight.
Use the function vfs_DofTabNodeSpec()
to specify such nodes.
Examples of candidate nodes would be nodes attached to distorted
elements, nodes in areas of very high solution gradients, nodes
on discrete contacting surfaces, etc.
All of this information must be supplied before vfs_DofTabProcess()
is called.
The iterative solver preconditioner technology has two “aggressiveness”
settings, DOFTAB_ITERMEMLEVEL and DOFTAB_ITEREXPLEVEL, which
can be changed from the default values if a model is relatively
large and contains highly refined regions. These parameters
are set using vfs_DofTabSetParami()
.
2.3. Special Considerations for Rigid Mode Analysis
If a rigid body analysis is to be performed by SysMatrix,
all elements must have their topology specified with
vfs_DofTabSetTopology()
, and all nodes must have their coordinates specified
with vfs_DofTabSetCoords()
.
If there are local coordinate systems for unknowns at nodes, their
direction cosine matrices must be entered using vfs_DofTabNodeDirCos()
.
If there are temporary restraints to be applied in SysMatrix using
vfs_SysMatrixSetTemp()
then this information should be supplied
to DofTab using the function vfs_DofTabNodeDofSpec()
.
This information is identical to the requirements of the iterative solver.
Also element types should be set using vfs_DofTabSetFeaType()
.
The efficiency of the rigid mode analysis is dependent upon the
higher dimensionality element types being properly set.
For example, 3D elements of type SYS_ELEM_SOLID and
SYS_ELEM_SHELL and 2D elements of type SYS_ELEM_SOLID are
particularly important.
In addition vfs_DofTabSetParami()
must be used to set DOFTAB_REDUCE to DOFTAB_REDUCE_RIGID.
2.4. Function Descriptions
The currently available DofTab functions are described in detail in this section.
-
vfs_DofTab *vfs_DofTabBegin(void)
create an instance of a DofTab object
Create an instance of a DofTab object. Memory is allocated for the object private data and the pointer to the object is returned.
Destroy an instance of a DofTab object using
void vfs_DofTabEnd (vfs_DofTab *doftab )
Return the current value of a DofTab object error flag using
Vint vfs_DofTabError (vfs_DofTab *doftab )
- Returns:
The function returns a pointer to the newly created DofTab object. If the object creation fails, NULL is returned.
-
void vfs_DofTabEnd(vfs_DofTab *p)
destroy an instance of a DofTab object
-
Vint vfs_DofTabError(vfs_DofTab *p)
return the current value of a DofTab object error flag
-
void vfs_DofTabDef(vfs_DofTab *p, Vint numnp, Vint numel)
define number of nodes and elements
Specify the number of nodes and elements to be managed in the DofTab object. All degrees of freedom are initialized as constrained.
Inquire of defined numnp and numel as output arguments using
void vfs_DofTabInq (vfs_DofTab *doftab, Vint *numnp, Vint *numel)
- Errors
SYS_ERROR_VALUE
is generated if either numnp or numel is less than or equal to zero.
- Parameters:
p – Pointer to DofTab object.
numnp – Number of nodes
numel – Number of elements
-
void vfs_DofTabInq(vfs_DofTab *p, Vint *numnp, Vint *numel)
inquire of defined numnp and numel as output arguments
-
void vfs_DofTabSetParami(vfs_DofTab *p, Vint type, Vint param)
set integer parametersset integer parameters
Set processing parameters.
The parameter
DOFTAB_REDUCE
enables the processing required by the iterative solver or the rigid mode solver. The iterative solver is selected usingDOFTAB_REDUCE_ITER
, the rigid mode solver is selected usingDOFTAB_REDUCE_RIGID
. By defaultDOFTAB_REDUCE
is set toDOFTAB_REDUCE_NONE
.The parameter
DOFTAB_ITERMEMLEVEL
sets the aggressiveness level for the iterative solver used for higher order element contributions to the system matrix. By defaultDOFTAB_ITERMEMLEVEL
is set to zero. A level of 1 would be used for models containing large, highly refined higher order element regions.The parameter
DOFTAB_ITEREXPLEVEL
sets the aggressiveness level for the iterative solver used for linear order element contributions to the system matrix. By defaultDOFTAB_ITEREXPLEVEL
is set to 3. A level of 4 or 5 would be used for models containing large, highly refined linear order element regions.The parameter
DOFTAB_THREADASSEM
is used to enable the calculation of disjoint element groups using a coloring algorithm. The number of groups and group contents are queried usingvfs_DofTabThreadAssemNum()
andvfs_DofTabThreadAssemGroups()
. These disjoint groups are to be used with the SysMatrix object using a threaded assembly procedure. By defaultDOFTAB_THREADASSEM
is disabled.- Errors
SYS_ERROR_ENUM
is generated if an improper type is specified.SYS_ERROR_ENUM
is generated if an improper param is specified.
- Parameters:
p – Pointer to DofTab object.
type – Integer parameter type to set.
x=DOFTAB_REDUCE Select reduced solver =DOFTAB_ITERMEMLEVEL Iterative low memory level =DOFTAB_ITEREXPLEVEL Iterative expansion level =DOFTAB_THREADASSEM Threaded assembly flag
param – Specifies the integer value that type will be set to.
=DOFTAB_REDUCE_NONE No reduction =DOFTAB_REDUCE_ITER Iterative solver reduction =DOFTAB_REDUCE_RIGID Rigid mode reduction =SYS_OFF Disable =SYS_ON Enable
-
void vfs_DofTabSetMap(vfs_DofTab *p, Vint nedofs, Vint loc[], Vint tag[])
set current element degree of freedom map
Set the current number of element degrees of freedom and the degree of freedom map. The loc vector consists of nedofs one-based node indices into the element connectivity. A zero entry in loc denotes an internal elemental degree of freedom. The tag vector consists of nedofs defined constants which specify the degree of freedom type, eg. x-translation, y-rotation, etc. See Degrees of Freedom for a definition of the degree of freedom types. This function must be called before element usage can be declared with
vfs_DofTabElemUse()
.For example, the degree of freedom map for a 4 node tetrahedron, with translational degrees of freedom at each node, is set as follows.
Vint loc[12] = {1,1,1, 2,2,2, 3,3,3, 4,4,4}; Vint tag[12] = {SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TZ, SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TZ, SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TZ, SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TZ}; vfs_DofTabSetMap (doftab,12,loc,tag);
As another example, the degree of freedom map for a homogeneous multipoint constraint equation of the general form a*tx(node 20) + b*ty(node 20) + c*ty(node 1) = 0. to be enforced using a Lagrange multiplier is set as follows.
Vint loc[4] = {1,1,2,0}; Vint tag[12] = {SYS_DOF_TX, SYS_DOF_TY, SYS_DOF_TY, SYS_DOF_LAGM}; vfs_DofTabSetMap (doftab,4,loc,tag);
- Errors
SYS_ERROR_VALUE
is generated if an improper nedofs or loc is specified.SYS_ERROR_ENUM
is generated if an improper tag is specified.
- Parameters:
p – Pointer to DofTab object.
nedofs – Number of element degrees of freedom
loc – Vector of degree of freedom locations
tag – Vector of degree of freedom types
-
void vfs_DofTabElemUse(vfs_DofTab *p, Vint eid, Vint nix, Vint ix[])
specify element degree of freedom use
Declare all degrees of freedom used by an element. The degrees of freedom for an element consist of the nodal degrees of freedom used by nix element nodes in vector ix, and the internal elemental freedoms used by element eid. The degree of freedom map for the element must have been previously specified using
vfs_DofTabSetMap()
.- Errors
SYS_ERROR_VALUE
is generated if an improper eid or nix is specified.SYS_ERROR_OPERATION
is generated if the global dof have already been processed byvfs_DofTabProcess()
or if no element degree of freedom map has been specified.
- Parameters:
p – Pointer to DofTab object.
eid – Element number
nix – Number of nodes per element
ix – Vector of element node connectivity
-
void vfs_DofTabSetCoords(vfs_DofTab *p, Vint nid, Vdouble x[3])
set the nodal coordinates
Set the nodal coordinate so that a preconditioner can be built for SysMatrix.
- Errors
SYS_ERROR_OPERATION
is generated ifvfs_DofTabDef()
has not been called.SYS_ERROR_VALUE
is generated if nid is less than 1 or nid is greater than the number of nodes.
- Parameters:
p – Pointer to DofTab object.
nid – Node id
x – Vector of nodal (x,y,z) coordinates
-
void vfs_DofTabNodePerm(vfs_DofTab *p, Vint nid, Vint ntags, Vint tag[])
set nodal degree of freedom constraints
Specify permanent constraints at node nid for a set of ntags degree of freedom types tag. The specified degrees of freedom are marked as inactive and will not be included in the global active dof when
vfs_DofTabProcess()
is called. This function should be called after all element usage has been specified usingvfs_DofTabElemUse()
.- Errors
SYS_ERROR_VALUE
is generated if an improper nid or ntags is specified.SYS_ERROR_ENUM
is generated if an improper tag is specified.SYS_ERROR_OPERATION
is generated if the global dof have already been processed.
- Parameters:
p – Pointer to DofTab object.
nid – Node number
ntags – Number of degrees of freedom types at node
tag – Vector of degree of freedom types
-
void vfs_DofTabElemDof(vfs_DofTab *p, Vint eid, Vint nix, Vint ix[], Vint *nedofs, Vint lm[])
query element global dof
Query for the global dof assigned to element eid. The complete global dof for an element consist of the global dof assigned to the nodal freedoms at each element node plus the global dof assigned to the internal elemental freedoms for the element. The element number, eid is required to assign any elemental freedoms. The element node connectivity, ix, is required to assign nodal freedoms. The global dof must have been previously processed using
vfs_DofTabProcess()
. Usevfs_DofTabElemSup()
to return the global suppressed dof for an element.- Errors
SYS_ERROR_VALUE
is generated if an improper eid entry, nix or ix entry is specified.SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
eid – Element number
nix – Number of nodes per element
ix – Vector of element node connectivity
nedofs – [out] Number of elements dofs
lm – [out] Vector of element global dof numbers
-
void vfs_DofTabElemSup(vfs_DofTab *p, Vint eid, Vint nix, Vint ix[], Vint *nedofs, Vint lm[])
query element global suppressed dof
Query for the global suppressed dof assigned to element eid. The global dof must have been previously processed using
vfs_DofTabProcess()
. Usevfs_DofTabElemDof()
to return the global active dof for an element.- Errors
SYS_ERROR_VALUE
is generated if an improper eid entry, nix or ix entry is specified.SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
eid – Element number
nix – Number of nodes per element
ix – Vector of element node connectivity
nedofs – [out] Number of element dofs
lm – [out] Vector of element global dof numbers
-
void vfs_DofTabNodeDof(vfs_DofTab *p, Vint nid, Vint ntags, Vint tag[], Vint lm[])
query node global dof
Query for the global dof assigned to a node nid for a set of ntags degree of freedom types in vector tag. Use
vfs_DofTabGetUseTag()
to determine the number and type of degrees of freedom used at nodes. The global dof must have been previously processed usingvfs_DofTabProcess()
. Usevfs_DofTabNodeSup()
to return the global suppressed dof for a node.- Errors
SYS_ERROR_VALUE
is generated if an improper nid or ntags is specified.SYS_ERROR_ENUM
is generated if an improper tag is specified.SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
nid – Node number
ntags – Number of degrees of freedom types used at nodes
tag – Vector of degree of freedom types at nodes
lm – [out] Vector of global dof numbers at node
-
void vfs_DofTabNodeSup(vfs_DofTab *p, Vint nid, Vint ntags, Vint tag[], Vint lm[])
query node global suppressed dof
Query for the global suppressed dof assigned to a node nid for a set of ntags degree of freedom types in vector tag. Use
vfs_DofTabGetUseTag()
to determine the number and type of degrees of freedom used at nodes. The global dof must have been previously processed usingvfs_DofTabProcess()
. Usevfs_DofTabNodeDof()
to return the global active dof for a node.- Errors
SYS_ERROR_VALUE
is generated if an improper nid or ntags is specified.SYS_ERROR_ENUM
is generated if an improper tag is specified.SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
nid – Node number
ntags – Number of degrees of freedom types used at nodes
tag – Vector of degree of freedom types at nodes
lm – [out] Vector of global suppressed dof numbers at node
-
void vfs_DofTabProcess(vfs_DofTab *p, Vint *neq, Vint *nre)
process and assign global dof
Number the global dof and global suppressed dof. The global dof consist of all degrees of freedom used by elements excluding any explicit suppressions at nodes. The global suppressed dof consist of explicitly suppressed degrees of freedom at nodes. This function must be called before the global dof can be queried using
vfs_DofTabElemDof()
,vfs_DofTabNodeDof()
, etc. If the degrees of freedom have already been processed this function returns with the values of neq and nre without reprocessing the degrees of freedom.- Parameters:
p – Pointer to DofTab object.
neq – [out] Number of global (active) dof
nre – [out] Number of global suppressed dof
-
void vfs_DofTabGetUseTag(vfs_DofTab *p, Vint *ntags, Vint tag[])
query used degree of freedom tags
Return the number, ntags, and types, tag, of all the degrees of freedom used at nodes. The global dof must have been previously processed using
vfs_DofTabProcess()
. This information is useful for gathering computed displacements, etc. at nodes from the global system solution vector. Usevfs_DofTabNodeDof()
to query for the global dof numbers at a node given ntags and tag returned by this function.- Errors
SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
ntags – [out] Number of degrees of freedom types used at nodes
tag – [out] Vector of degree of freedom types
-
void vfs_DofTabMaxElemDof(vfs_DofTab *p, Vint *maxelemdof)
query for maximum size of element global dof
Return the maximum number of degrees of freedom occurring in any element.
- Parameters:
p – Pointer to DofTab object.
maxelemdof – [out] Maximum number of element degrees of freedom
-
void vfs_DofTabInvDof(vfs_DofTab *p, Vint ndofs, Vint lm[], Vint ids[], Vint tag[])
query inverse degree of freedom information
Return the node or element ids and degree of freedom types, tag, for a set of ndofs global dof numbers lm. The global dof must have been previously processed using
vfs_DofTabProcess()
.- Errors
SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
ndofs – Number of degrees of freedom
lm – Vector of global dof numbers
ids – [out] Node or element ids associated with global dof
tag – [out] Vector of degree of freedom types
-
void vfs_DofTabSetTopology(vfs_DofTab *p, Vint id, Vint shape, Vint maxi, Vint maxj, Vint maxk)
set the element topology
Set the element topology so that a preconditioner can be subsequently built. See section Computational Cells for a description of element topology conventions. For example, the topology of element number 350, a 10 node tetrahedron, with Serendipity node connectivity (corner nodes first followed by midsides) is set as follows.
vfs_DofTabSetTopology (doftab,350,SYS_SHAPETET,3,0,0)
- Errors
SYS_ERROR_OPERATION
is generated ifvfs_DofTabDef()
has not been called.SYS_ERROR_ENUM
is generated if an improper shape is specified. AndSYS_ERROR_VALUE
is generated if id is greater than the number of elements, or if maxi, maxj, or maxk are invalid.
- Parameters:
p – Pointer to DofTab object.
id – Element id
shape – Element shape parameter
x=SYS_SHAPEPOINT Point =SYS_SHAPELINE Line =SYS_SHAPETRI Triangle =SYS_SHAPEQUAD Quadrilateral =SYS_SHAPETET Tetrahedron =SYS_SHAPEWED Wedge =SYS_SHAPEHEX Hexahedron
maxi – The number of points along the i direction. If maxi = 0 then the linear element form of the specified shape is assumed.
maxj – The number of points along the j direction. If maxj = 0 and 2 <= maxi <= 4, then a Serendipity finite element is assumed. If 2 <= maxj and 2 <= maxi, then a Lagrange finite element is assumed. For triangle, tetrahedron and wedge shapes, set either maxj = 0 or maxj = maxi.
maxk – The number of points along the k direction. For tetrahedron shapes, set either maxk = 0 or maxk = maxi For hexahedron and wedge shapes, set either maxk = 2 or maxk = maxi
-
void vfs_DofTabNodeSpec(vfs_DofTab *p, Vint level, Vint nix, Vint ix[])
flag all of a node’s degrees of freedom for preconditioner
Flag all active degrees of freedom on a list of nodes to be treated with special care by the preconditioner. Typical candidates are degrees of freedom where temporary restraints have been specified with
vfs_SysMatrixSetTemp()
, or where the user suspects potential sources of ill-conditioning may arise (e.g., very distorted elements). If level is set to 1, then all degrees of freedom at the node are flagged; if level is set to 2, then all degrees of freedom of the node and its neighbors are flagged. This function may be called as many times as necessary.- Errors
SYS_ERROR_VALUE
is generated if level is less than 1 or greater than 2; or if nix is less than 1; or if any entry in ix is less than 1 or greater than the number of nodes specified invfs_DofTabDef()
; or if any entry in tags is not a valid nodal degree of freedom.SYS_ERROR_OPERATION
is generated ifvfs_DofTabProcess()
has already been called, or ifvfs_DofTabDef()
has not yet been called.
- Parameters:
p – Pointer to DofTab object.
level – Propagation level for flagging node
nix – Number of nodes being flagged
ix – List of nodes being flagged
-
void vfs_DofTabNodeDirCos(vfs_DofTab *p, Vint nid, Vdouble tm[3][3])
set direction cosine matrix at nodes for iterative solver
Set the direction cosine matrix necessary for the iterative solver
- Parameters:
p – Pointer to DofTab object.
nid – Node id
tm – Direction cosine matrix at node
-
void vfs_DofTabNodeDofSpec(vfs_DofTab *p, Vint nix, Vint ix[], Vint ntags, Vint tags[])
flag degree of freedom for preconditioner
Flag a list of degrees of freedom on a list of nodes to be treated with special care by the preconditioner. Typical candidates are degrees of freedom where temporary restraints have been specified with
vfs_SysMatrixSetTemp()
, or where the user suspects potential sources of ill-conditioning may arise (e.g., very distorted elements).- Errors
SYS_ERROR_VALUE
is generated if nix is less than 1; or if any entry in ix is less than 1 or greater than the number of nodes specified invfs_DofTabDef()
; or if any entry in tags is not a valid nodal degree of freedom.SYS_ERROR_OPERATION
is generated ifvfs_DofTabProcess()
has already been called, or ifvfs_DofTabDef()
has not yet been called.
- Parameters:
p – Pointer to DofTab object.
nix – Number of nodes being flagged
ix – List of nodes being flagged
ntags – Number of nodal degree of freedom tags being flagged
tags – List of degrees of freedom being flagged
-
void vfs_DofTabSetFeaType(vfs_DofTab *p, Vint eid, Vint featype)
set the element type
Set the element type so a rigid body analysis can be efficiently performed by SysMatrix. In particular, 2D and 3D elements of type
SYS_ELEM_SOLID
and 3D elements of typeSYS_ELEM_SHELL
are relevant to a rigid body analysis.- Errors
SYS_ERROR_OPERATION
is generated ifvfs_DofTabDef()
has not been called.SYS_ERROR_VALUE
is generated if eid is less than 1 or eid is greater than the number of elements.
- Parameters:
p – Pointer to DofTab object.
eid – Element id
featype – The element type
-
void vfs_DofTabSetMode(vfs_DofTab *p, Vint mode, Vint param)
set operating mode
Set operating modes.
Specify the shape point iterative solver preconditioner treatment for shape point elements using
DOFTAB_SHAPEPOINTMODE
. If enabled shape point element connected nodes will be treated internally as specified nodes using the functionvfs_DofTabNodeSpec()
. If disabled they will not be given this treatment. If a shape point element nodes represent very discrete topologies then enable this mode. If the shape point nodes represent a distributed connection along an edge or particularly a face disable this mode. By defaultDOFTAB_SHAPEPOINTMODE
mode is toSYS_ON
.- Errors
SYS_ERROR_ENUM
is generated if an improper mode is input.
- Parameters:
p – Pointer to DofTab object.
mode – Type of mode to set
x=DOFTAB_SHAPEPOINTMODE shape point preconditioner treatment
param – Mode parameter
x=SYS_ON Enable mode =SYS_OFF Disable mode
-
void vfs_DofTabSetMPC(vfs_DofTab *p, Vint ndofs, Vint ix[], Vint tag[], Vdouble c[], Vdouble rhs)
set multipoint constraints
Set the multipoint node connectivity, degree of freedom type, coefficients and right hand side for a multipoint constraint. The number of multipoint degrees of freedom is ndofs. The first degree of freedom is assumed to be the dependent degree of freedom. The multipoint constraints must be set before
vfs_DofTabProcess()
is called.If all coefficients in c are identically zero, the value rhs will still be associated with the first degree of freedom, and forces will be updated with this value during a solution procedure. However, in this case, and in this case only, the degree of freedom will not be considered dependent.
- Errors
SYS_ERROR_VALUE
is generated if an improper nedofs is specified.
- Parameters:
p – Pointer to DofTab object.
ndofs – Number of multipoint degrees of freedom
ix – Vector of node connectivity
tag – Vector of degree of freedom types
c – Vector of multipoint coefficients
rhs – Right hand side
-
void vfs_DofTabSetNumThreads(vfs_DofTab *p, Vint num)
set number of threads
Set the number of threads used for parallel execution of DofTab functions. If num is 0, then DofTab functions are run serially. If num is greater than or equal to 1 then num threads are used for parallel execution of various algorithms in
vfs_DofTabProcess()
including, for example, equation numbering using Metis for factorization.- Errors
SYS_ERROR_VALUE
is generated if num is less than 0.
- Parameters:
p – Pointer to DofTab object.
num – Number of threads
x=0 Default, serial execution. >=1 Number of threads for parallel execution
-
void vfs_DofTabThreadAssemNum(vfs_DofTab *p, Vint *num)
query number of assembly groups
Return the number of assembly groups. Use
vfs_DofTabThreadAssemGroups()
to get the number of elements and elements in each group.- Errors
SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
num – [out] Number of assembly groups
-
void vfs_DofTabThreadAssemGroups(vfs_DofTab *p, Vint *num, Vint numelem[], Vint blockid[], Vint elems[])
query assembly groups
Return the number of assembly groups, the number of elements in each group, the blockid in each group, and the list of elements in each group. The arrays numelems and blockid have a length of num and contains the number of elements and block id in each group. The array elems is length total number of elements and contains the elements in each group. The element numbers in elems are one based. The first numelem[0] elements are in the first group, the next numelem[1] elements are in the second group, etc. If a negative block id exists it will immediately follow its corresponding positive value, if it exists.
- Errors
SYS_ERROR_OPERATION
is generated if the global dof have not been processed.
- Parameters:
p – Pointer to DofTab object.
num – [out] Number of assembly groups
numelem – [out] Vector of number of elements in each group
blockid – [out] Vector of block ids in each group
elems – [out] List of elements in each group
-
void vfs_DofTabPointersIndicesPtr(vfs_DofTab *p, Vlong pointers[], Vint indices[])
set pointers and indices
Both pointers and indices are 0-based arrays. The pointers[i] array entry points to a location in the indices array where the indices for the i-th column start. The number of nonzeros in the i-th column is given by pointers[i+1]-pointers[i] .
This function is used in conjunction with
vfs_SysMatrixSetValue()
. In this function, entry a_ij refers to the matrix value for column i, row j.- Errors
SYS_ERROR_OPERATION
is generated if DofTabProcess has already been called.
- Parameters:
p – Pointer to DofTab object.
pointers – Location in indices for first row in each column
indices – List of row indices for each column
-
void vfs_DofTabSetAssemBlock(vfs_DofTab *p, Vint blockid)
set running assembling block id
Set the blockid for all subsequent elements declared with
vfs_DofTabElemUse()
.This function is optional, and is only relevant for when users want to assemble elements in
vfs_SysMatrix()
in parallel. When this is the case, DofTab will group together elements that can be assembled in parallel.Certain applications group elements in a particular way. As a result, the assembly parallelization must be done only for the group of elements available at the time. By assigning a blockid to each element users are guaranteed that a group of elements that can be assembled in parallel belongs to the same blockid.
In general, blockid values are positive. However, users can also specify negative blockid values. When this is the case, it is expected that all elements associated with this block will be assembled in serial mode and need not be grouped.
DofTab guarantees that if both blockid and -blockid are present, then the positive values will be listed immediately before the negative values. Hence, for each positive/negative pair, the elements that can be assembled in parallel will be listed first, immediately followed by the corresponding negative blockid values.
If this function is called for any element, then it must be called for all elements. By default, all elements have blockid set to zero. If any element has an assigned non-zero blockid, then the presence of a zero-valued blockid will generate an error during
vfs_DofTabProcess()
.- Errors
SYS_ERROR_VALUE
is generated if blockid is zero.
- Parameters:
p – Pointer to DofTab object.
blockid – Running block id