*vfs_SysVectorBegin - create an instance of a SysVector object vfs_SysVectorPre - specify real or complex vector vfs_SysVectorEnd - destroy an instance of a SysVector object vfs_SysVectorError - return SysVector object error flag vfs_SysVectorCopy - copy a SysVector object
vfs_SysVectorAdd - add two vectors vfs_SysVectorAssem - assemble an element vector vfs_SysVectorThreadAssem - assemble an element vector in thread vfs_SysVectorDef - define number of global dof Inq - inquire number of global dof vfs_SysVectorDot - compute dot product of two vectors vfs_SysVectorExpand - expand a system vector vfs_SysVectorGather - gather values from vector vfs_SysVectorThreadGather - gather values from vector in thread vfs_SysVectorInit - initialize vector to a given value vfs_SysVectorInverse - compute the inverse of a vector vfs_SysVectorMaskAdd - add two vectors with mask vfs_SysVectorMaskGather - gather values from vector with mask vfs_SysVectorMaskScatter - scatter values to vector with mask vfs_SysVectorMove - move vector values to another vector vfs_SysVectorNorm - compute the norm of a vector vfs_SysVectorLInfNorm - compute the L-Infinity norm of a vector vfs_SysVectorProd - compute the product of two vectors vfs_SysVectorAddProd - compute and add the product of two vectors vfs_SysVectorRandom - generate a random vector vfs_SysVectorReduce - reduce a system vector vfs_SysVectorScale - scale a vector by a value vfs_SysVectorScatter - scatter values to vector vfs_SysVectorThreadScatter - scatter values to vector in thread vfs_SysVectorScatterVal - scatter a constant value to vector vfs_SysVectorSetComplexMode - set complex mode vfs_SysVectorSetThreadComplexMode - set complex mode in thread vfs_SysVectorSetNumThreads - set number of threads vfs_SysVectorSquareRoot - computes square root of vector
All values in a SysVector object may be initialized to some user defined real value using vfs_SysVectorInit. Values may be set into a SysVector object at selected degrees of freedom using vfs_SysVectorScatter. Values may be added into a SysVector object at selected degrees of freedom using vfs_SysVectorAssem. Note that the values should be initialized to zero before using vfs_SysVectorAssem.
Random values may be generated using vfs_SysVectorRandom.
Various operations can be performed upon system vectors such as dot product addition, scaling, etc using vfs_SysVectorDot, vfs_SysVectorAdd, vfs_SysVectorScale, etc.
Values may be selectively retrieved from a system vector at selected global dof numbers using vfs_SysVectorGather.
SysVector objects are entered as arguments to a SysMatrix object for performing matrix-vector products, back substitution on factored system matrices, etc. Destroy an instance of a SysVector object using vfs_SysVectorEnd.
If Gather/Scatter/Assem operations are performed in parallel over several threads then the complex mode can be specified on each thread using vfs_SysVectorSetThreadComplexMode.
*vfs_SysVectorBegin - create an instance of a SysVector object
vfs_SysVector *vfs_SysVectorBegin ()
None
Destroy an instance of a SysVector object using
void vfs_SysVectorEnd (vfs_SysVector *sysvector)
Return the current value of a SysVector object error flag using
Vint vfs_SysVectorError (vfs_SysVector *sysvector)
Make a copy of a SysVector object. The private data from the fromsysvector object is copied to the sysvector object. Any previous private data in sysvector is lost.
void vfs_SysVectorCopy (vfs_SysVector *sysvector, vfs_SysVector *fromsysvector)
vfs_SysVectorPre - specify real or complex vector
void vfs_SysVectorPre (vfs_SysVector *sysvector, Vint complextype)
sysvector Pointer to SysVector object. complextype Enumerated type for real or complex matrix =SYS_DOUBLE Real matrix =SYS_DOUBLECOMPLEX Complex matrix
None
SYS_ERROR_ENUM is issued if complextype is invalid.
This call must be made before vfs_SysVectorDef is called. Defaults to real.
vfs_SysVectorAdd - add two vectors
void vfs_SysVectorAdd (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vdouble s, vfs_SysVector *sysvectory)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x. s Scale factor sysvectory Pointer to SysVector object y.
None
vfs_SysVectorAssem - assemble an element vector
void vfs_SysVectorAssem (vfs_SysVector *sysvector, Vint nedofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. nedofs Number of element degrees of freedom lm Vector of element global dof numbers v Element vector
None
vfs_SysVectorThreadAssem - assemble an element vector in a thread
void vfs_SysVectorThreadAssem (vfs_SysVector *sysvector, Vint ithread, Vint nedofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. ithread Thread number starting from 1. nedofs Number of element degrees of freedom lm Vector of element global dof numbers v Element vector
None
If a complex mode is specified with vfs_SysVectorSetThreadComplexMode then the array v refers to the mode specified.
vfs_SysVectorDef - define number of global dof
void vfs_SysVectorDef (vfs_SysVector *sysvector, Vint neq)
sysvector Pointer to SysVector object. neq Number of global dof
None
Inquire of a defined neq as an output argument using
void vfe_SysVectorInq (vfe_SysVector *sysvector, Vint *neq)
vfs_SysVectorDot - compute dot product of two vectors
void vfs_SysVectorDot (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vdouble *s)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x.
s Dot product
vfs_SysVectorExpand - expand a system vector
void vfs_SysVectorExpand (vfs_SysVector *sysvector, Vint nvec, vfs_SysVector *sysvectort[], vfs_SysVector *sysvectorq)
sysvector Pointer to SysVector object. nvec Number of system vectors 0 < nvec sysvectort Pointer to array of nvec SysVector objects sysvectorq Pointer to expanded SysVector object.
None
sysvector = [ a_1 a_2 ... a_N ]
then
sysvectorq = [ a_1 * sysvectort[0] + a_2 * sysvectort[1] + ... + a_N * sysvectort[N-1]
vfs_SysVectorGather - gather values from vector
void vfs_SysVectorGather (vfs_SysVector *sysvector, Vint ndofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. ndofs Number of degrees of freedom lm Vector of global dof numbers
v Gathered values
vfs_SysVectorThreadGather - gather values from vector in thread
void vfs_SysVectorThreadGather (vfs_SysVector *sysvector, Vint ithread, Vint ndofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. ithread Thread number starting from 1. ndofs Number of degrees of freedom lm Vector of global dof numbers
v Gathered values
If a complex mode is specified with vfs_SysVectorSetThreadComplexMode then the array v refers to the mode specified.
vfs_SysVectorInit - initialize vector to a given value
void vfs_SysVectorInit (vfs_SysVector *sysvector, Vdouble v)
sysvector Pointer to SysVector object. v Initialization value
None
vfs_SysVectorInverse - compute the inverse of a vector
void vfs_SysVectorInverse (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x.
None
vfs_SysVectorMaskAdd - add two vectors with mask
void vfs_SysVectorMaskAdd (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vdouble s, vfs_SysVector *sysvectory, Vint n, Vint mask[])
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object containing first vector in sum. s Scaling parameter. sysvectory Pointer to SysVector object containing second vector in sum, whose components given by the 0-based mask array will be scaled by s during the sum. n Number of components to copy mask Vector of 0-based components to use when referencing system vector sysvectory.
None
sysvector = [ a_1 a_2 ... a_N ] sysvectorx = [ x_1 x_2 ... x_N ] sysvectory = [ y_1 y_2 ... y_n ]
then
a_1 = x_1 + s * y_mask[0] a_2 = x_2 + s * y_mask[1] ... a_[n-1] = x_[n-1] + s * y_mask[n-1]
vfs_SysVectorMaskGather - gather values from vector with mask
void vfs_SysVectorMaskGather (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vint n, Vint mask[])
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object where data will be copied from. n Number of components to copy mask Vector of 0-based components to use when referencing system vector sysvectorx.
None
sysvector = [ a_1 a_2 ... a_N ] sysvectorx = [ x_1 x_2 ... x_N ]
then
a_1 = x_mask[0] a_2 = x_mask[1] ... a_[n-1] = x_mask[n-1]
vfs_SysVectorMaskScatter - scatter values to vector with mask
void vfs_SysVectorMaskScatter (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vint n, Vint mask[])
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object where data will be copied from. n Number of components to copy mask Vector of 0-based components to use when referencing system vector sysvectorx.
None
sysvector = [ a_1 a_2 ... a_N ] sysvectorx = [ x_1 x_2 ... x_N ]
then
a_mask[0] = x_1 a_mask[1] = x_2 ... a_mask[n-1] = x_[n-1]
vfs_SysVectorMove - move vector values to another vector
void vfs_SysVectorMove (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x.
None
vfs_SysVectorNorm - compute the norm of a vector
void vfs_SysVectorNorm (vfs_SysVector *sysvector, Vint type, Vdouble *value)
sysvector Pointer to SysVector object. type Type of vector norm to compute =SYSVECTOR_NORM_L1 Sum of absolute values =SYSVECTOR_NORM_L2 Square root of sum of squares =SYSVECTOR_NORM_LINF Maximum component absolute value
value Norm value
vfs_SysVectorLInfNorm - compute the L-Infinity norm of a vector
void vfs_SysVectorLInfNorm (vfs_SysVector *sysvector, Vdouble *value, Vint *dof)
sysvector Pointer to SysVector object.
value L-Infinity norm value dof Global dof number
vfs_SysVectorProd - compute the product of two vectors
void vfs_SysVectorProd (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, vfs_SysVector *sysvectory)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x. sysvectory Pointer to SysVector object y.
None
vfs_SysVectorAddProd - compute and add the product of two vectors
void vfs_SysVectorAddProd (vfs_SysVector *sysvector, Vdouble s, vfs_SysVector *sysvectorx, vfs_SysVector *sysvectory)
sysvector Pointer to SysVector object. s Scale factor sysvectorx Pointer to SysVector object x. sysvectory Pointer to SysVector object y.
None
vfs_SysVectorRandom - generate a random vector
void vfs_SysVectorRandom (vfs_SysVector *sysvector)
sysvector Pointer to SysVector object.
None
vfs_SysVectorReduce - reduce a system vector
void vfs_SysVectorReduce (vfs_SysVector *sysvector, Vint nvec, vfs_SysVector *sysvectort[], vfs_SysVector *sysvectorq)
sysvector Pointer to SysVector object. nvec Number of system vectors 0 < nvec sysvectort Pointer to array of nvec SysVector objects sysvectorq Pointer to reduced SysVector object.
None
sysvectorq = [ a_1 a_2 ... a_N ]
then
a_1 = sysvector * sysvectort[0] a_2 = sysvector * sysvectort[1] ... a_N = sysvector * sysvectort[N-1]
vfs_SysVectorScale - scale a vector by a value
void vfs_SysVectorScale (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx, Vdouble s)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object x. s Scale factor
None
vfs_SysVectorScatter - scatter values to vector
void vfs_SysVectorScatter (vfs_SysVector *sysvector, Vint ndofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. ndofs Number of degrees of freedom lm Vector of global dof numbers v Scattered values
None
sysvector = [ a_1 a_2 ... a_N ]
then
a_lm[0] = v[0] a_lm[1] = v[1] ... a_lm[ndofs-1] = v[ndofs-1]
vfs_SysVectorThreadScatter - scatter values to vector in thread
void vfs_SysVectorThreadScatter (vfs_SysVector *sysvector, Vint ithread, Vint ndofs, Vint lm[], Vdouble v[])
sysvector Pointer to SysVector object. ithread Thread number starting from 1. ndofs Number of degrees of freedom lm Vector of global dof numbers v Scattered values
None
If a complex mode is specified with vfs_SysVectorSetThreadComplexMode then the array v refers to the mode specified.
vfs_SysVectorScatterVal - scatter a constant value to vector
void vfs_SysVectorScatter (vfs_SysVector *sysvector, Vint ndofs, Vint lm[], Vdouble v)
sysvector Pointer to SysVector object. ndofs Number of degrees of freedom lm Vector of global dof numbers v Constant scattered value
None
vfs_SysVectorSetComplexMode - set complex mode
void vfs_SysVectorSetComplexMode (vfs_SysVector *sysvector, Vint complexmode)
sysvector Pointer to SysVector object. complexmode The type of data being passed through the Vdouble[] argument list. =SYS_COMPLEX_REAL Data is real =SYS_COMPLEX_IMAGINARY Data is imaginary =SYS_COMPLEX_REALIMAGINARY Data is complex
None
The methods vfs_SysVectorAssem, vfs_SysVectorGather, and vfs_SysVectorScatter all take arrays of floating point numbers. The values entered or returned can be declared as SYS_COMPLEX_REAL, SYS_COMPLEX_IMAGINARY, or SYS_COMPLEX_REALIMAGINARY. The data in the argument then takes the specified mode.
This is useful if only real or only imaginary data is being processed. The mode SYS_COMPLEX_REAL indicates that only one value, the real part, is being entered for each a_i in the array. In this case, the imaginary entries in the matrix remain unchanged. Likewise, the mode SYS_COMPLEX_IMAGINARY only affects the imaginary part. If mode is set to SYS_COMPLEX_REALIMAGINARY then two consecutive numbers in the array stand for the real and imaginary parts of the data. Defaults to SYS_COMPLEX_REAL.
vfs_SysVectorSetThreadComplexMode - set complex mode in thread
void vfs_SysVectorSetComplexMode (vfs_SysVector *sysvector, Vint ithread, Vint complexmode)
sysvector Pointer to SysVector object. ithread Thread number starting from 1. complexmode The type of data being passed through the Vdouble[] argument list. =SYS_COMPLEX_REAL Data is real =SYS_COMPLEX_IMAGINARY Data is imaginary =SYS_COMPLEX_REALIMAGINARY Data is complex
None
The methods vfs_SysVectorAssem, vfs_SysVectorGather, and vfs_SysVectorScatter all take arrays of floating point numbers. The values entered or returned can be declared as SYS_COMPLEX_REAL, SYS_COMPLEX_IMAGINARY, or SYS_COMPLEX_REALIMAGINARY. The data in the argument then takes the specified mode. Each thread in parallel operations can take different values for the complexmode value.
This is useful if only real or only imaginary data is being processed. The mode SYS_COMPLEX_REAL indicates that only one value, the real part, is being entered for each a_i in the array. In this case, the imaginary entries in the matrix remain unchanged. Likewise, the mode SYS_COMPLEX_IMAGINARY only affects the imaginary part. If mode is set to SYS_COMPLEX_REALIMAGINARY then two consecutive numbers in the array stand for the real and imaginary parts of the data. Defaults to SYS_COMPLEX_REAL.
vfs_SysVectorSquareRoot - computes square root of vector
void vfs_SysVectorSquareRoot (vfs_SysVector *sysvector, vfs_SysVector *sysvectorx)
sysvector Pointer to SysVector object. sysvectorx Pointer to SysVector object to take square root of.
none
vfs_SysVectorSetNumThreads - set number of threads
void vfs_SysVectorSetNumThreads (vfs_SysVector *sysvector, Vint num)
sysvector Pointer to SysVector object. num Number of threads =1 Default, serial execution. >1 Number of threads for parallel execution
None
*vfs_SysMatrixBegin - create an instance of a SysMatrix object vfs_SysMatrixEnd - destroy an instance of a SysMatrix object vfs_SysMatrixPre - specify real or complex matrix vfs_SysMatrixError - return SysMatrix object error flag
vfs_SysMatrixAbort - abort a sparse direct factorization vfs_SysMatrixAssem - assemble a symmetric element matrix vfs_SysMatrixAssemDiag - assemble a diagonal element matrix vfs_SysMatrixAssemFull - assemble a full element matrix vfs_SysMatrixColumn - extract a column vfs_SysMatrixDef - define number of global dof and structure Inq - inquire number of global dof and structure vfs_SysMatrixDiag - extract the diagonal vfs_SysMatrixEigenAll - extract all eigenvalues and eigenvectors vfs_SysMatrixEigen - extract set of eigenvalues and eigenvectors vfs_SysMatrixEigenMode - retrieve a computed eigenvector vfs_SysMatrixEigenFree - free all computed eigenvector storage vfs_SysMatrixFactor - factorize vfs_SysMatrixGather - gather non-zero values of a column vfs_SysMatrixGetDouble - get double precision matrix information vfs_SysMatrixGetInteger - get integer matrix information vfs_SysMatrixGetLong - get long matrix information vfs_SysMatrixGetSysVector - retrieve internal SysVector objects vfs_SysMatrixInit - initialize matrix to a given value vfs_SysMatrixInner - form matrix as inner product of vectors vfs_SysMatrixPreProcess - establish matrix structure vfs_SysMatrixProcess - process for factorization vfs_SysMatrixProd - multiply by a system vector vfs_SysMatrixReact - multiply by a system vector vfs_SysMatrixReduce - reduce a matrix vfs_SysMatrixSetComplexMode - set complex mode vfs_SysMatrixSetThreadComplexMode - set complex mode in thread vfs_SysMatrixSetFunction - set user function vfs_SysMatrixSetInitSol - set initial solution vector vfs_SysMatrixSetNumThreads - set number of threads vfs_SysMatrixSetObject - set attribute objects vfs_SysMatrixSetOocFile - set out-of-core file name vfs_SysMatrixSetParamd - set solution parameters vfs_SysMatrixSetParami - set solution parameters vfs_SysMatrixSetTemp - set temporary restraints vfs_SysMatrixSetNumThreadsAssem - set number of threads for assembly vfs_SysMatrixSetValue - set matrix value vfs_SysMatrixSolve - solve given a system vector vfs_SysMatrixSolveMult - solve given an array of system vectors vfs_SysMatrixSum - sum a system matrix vfs_SysMatrixThreadAssem - assemble a symmetric element matrix vfs_SysMatrixThreadAssemDiag - assemble a diagonal element matrix vfs_SysMatrixThreadAssemFull - assemble a full element matrix
The SysMatrix module utilizes the BLAS routines for performance reasons. Use VKI_LIBAPI_BLASMKL_SEQUENTIAL or VKI_LIBAPI_BLASMKL_THREAD when compiling with Intel's MKL BLAS routines for optimal performance. Use VKI_LIBAPI_BLASMKL_ILP64 if your application uses the 64-bit integer API for the MKL BLAS library. Use VKI_LIBAPI_OPENBLAS if the OpenBLAS are to be used.
See section VfsTools, Compiling and Linking a VfsTools Application. The complete, external BLAS functions are required for the Block Lanczos and AMLS eigensolvers.
The direct sparse solver used in SysMatrix may be changed from the default native solver. Two external solver packages are supported, Intel Pardiso and MUMPS. See section VfsTools, Compiling and Linking a VfsTools Application.
The assembled matrix is established in a SysMatrix object by assembly of individual finite element matrices and optional multipoint constraints. Multipoint constraints in general may be enforced using a variety of methods such as penalty functions, Lagrange multipliers or direct matrix transformation. If multipoint constraints are enforced using penalty functions or Lagrange multipliers then they are assembled as symmetric "element" matrices in the conventional way. Note that the direct sparse solvers currently do not support full pivoting which is necessary for robust treatment of Lagrange multipliers. If they are to be enforced by direct matrix transformation, then the multipoint constraints (MPCs) must have been specified using vfs_DofTabSetMPC for each MPC equation before any element matrix information is entered into SysMatrix. SysMatrix will transform the matrix internally and eliminate all dependent degrees of freedom from the equation system. Note that transformations can not be applied to a diagonal matrix type.
To determine the structure of a sparse matrix, the SysMatrix object uses information which has been previously entered in a DofTab object. The DofTab object must be set using vfs_SysMatrixSetObject. Once this is done use vfs_SysMatrixPreProcess to establish the structure of the matrix. The matrix must then be initialized using vfs_SysMatrixInit, and then assembled by calling vfs_SysMatrixAssem, vfs_SysMatrixAssemDiag, or vfs_SysMatrixAssemFull for each element with the actual element matrix values. If threaded assembly is used then the assembly functions vfs_SysMatrixThreadAssem, vfs_SysMatrixThreadAssemDiag and vfs_SysMatrixThreadAssemFull are used.
The system may be solved for a given right hand side vector using either a direct or iterative solution procedure. Use vfs_SysMatrixSetParami to specify the solution procedure and preconditioner type if necessary. vfs_SysMatrixProcess performs any pre-processing operations required before factorization occurs. vfs_SysMatrixFactor either factorizes the matrix or computes the preconditioner depending upon the solution procedure. vfs_SysMatrixSolve solves the system and produces a solution vector given a right hand side vector. vfs_SysMatrixSolveMult solves the system for multiple right hand side vectors producing multiple solution vectors.
Other than boundary conditions, a mechanism is available to eliminate singularities. Use vfs_SysMatrixSetParami with parameter SYSMATRIX_AUTOSPC set to SYS_ON. If enabled, whenever a pivot falls below the singularity threshold, the value SYSMATRIX_AUTOSPC_FACTOR multiplied by the matrix largest entry is set to the pivot to eliminate the singularity. SYSMATRIX_AUTOSPC default to SYS_OFF. SYSMATRIX_AUTOSPC_FACTOR defaults to 1.e+5.
Use vfs_SysMatrixColumn and vfs_SysMatrixDiag to extract a column or the diagonal of the assembled matrix into a SysVector object. Perform a matrix vector product using vfs_SysMatrixProd. Project a matrix to a subspace defined by a set of system vectors using vfs_SysMatrixReduce.
A function can be set with vfs_SysMatrixSetFunction to monitor the progress during a direct, sparse factorization procedure. Progress is measured by retrieving information from the SysMatrix module with vfs_SysMatrixGetInteger, vfs_SysMatrixGetLong, or vfs_SysMatrixGetDouble.
The factorization can be aborted with a call to vfs_SysMatrixAbort.
Destroy an instance of a SysMatrix object using vfs_SysMatrixEnd.
The type of solver to use depends upon the following:
Iterative PCG/GMRES solvers may provide substantial improvements in CPU time and memory. This is particularly true for 3D solid dominated models. The use of higher order elements generally improves the performance of iterative solvers. Use vfs_SysMatrixSetInitSol to set an optional initial guess solution vector. The iterative solver preconditioner must be factorized and during iteration, back substitutions on the preconditioner are performed. Any of the linear equation solvers may be used for this purpose.
There are three convergence tolerances used in the iterative solver, relative solution error, residual error and energy error. Tolerances are specified using vfs_SysMatrixSetParamd with parameters SYSMATRIX_ITER_UTOL, SYSMATRIX_ITER_FTOL and SYSMATRIX_ITER_ETOL respectively.
The following code fragment illustrates the basic framework of using the SysMatrix module to assemble and solve of a set of symmetric finite element matrices. The framework allows for the optional use of the iterative solver.
Vint i, j; /* declare objects */ vfs_SysMatrix *sysmatrix; vfs_SysVector *sysvectorb, sysvectorx; vfs_DofTab *doftab; /* number of nodes, elements */ Vint numnp, numel; /* number of equations, global degrees of freedom */ Vint neq; /* number of degrees of freedom in element, element degree of freedom numbers */ Vint nedofs, lm[MAX_DOF_ELEMENT]; /* degree of freedom types within an element */ Vint loc[MAX_DOF_ELEMENT], tag[MAX_DOF_ELEMENT]; /* lower triangle of element stiffness by rows */ Vdouble ke[MAX_DOF_ELEMENT*(MAX_DOF_ELEMENT+1)/2]; /* element topology and order */ Vint shape, maxi, maxj, maxk; /* number of nodes in element and connectivity */ Vint nix, ix[MAX_NODE_ELEMENT]; /* node coordinates */ Vdouble x[3]; /* create sysmatrix object */ sysmatrix = vfs_SysMatrixBegin(); doftab = vfs_DofTabBegin(); vfs_DofTabDef (doftab,numnp,numel); for(eid = 1; eid <= numel; eid++ { vfs_DofTabSetMap (doftab,nedofs,loc,tag); /* localize element connectivity ix[] from host application */ . . vfs_DofTabElemUse (doftab,eid,nix,ix); } /* modification for iterative solver */ if(iterative_solver) { vfs_DofTabSetParami (doftab,DOFTAB_REDUCE,DOFTAB_REDUCE_ITER); /* gather and set node coordinates */ for(j = 1; j <= numnp; j++) ... vfs_DofTabSetCoords (doftab,j,x); } /* gather and set element topology, for(i = 1; i <= numel; i++) { ... vfs_DofTabSetTopology (doftab,shape,maxi,maxj,maxk); } vfs_SysMatrixDef (sysmatrix,SYSMATRIX_SYMM_ITER,neq); } else { vfs_SysMatrixDef (sysmatrix,SYSMATRIX_SYMM_SPARSE,neq); } /* register DofTab with SysMatrix object */ vfs_SysMatrixSetObject (sysmatrix,VFS_DOFTAB,doftab); /* establish matrix structure */ vfs_SysMatrixPreProcess (sysmatrix); /* Initialize matrix values to zero */ vfs_SysMatrixInit (sysmatrix,0.); /* assembly for each element */ for(i = 1; i <= numel; i++) { /* form element matrix and gather degrees of freedom */ ... vfs_SysMatrixAssem (sysmatrix,nedofs,lm,ke); } /* equation ordering and/or symbolic factorization */ vfs_SysMatrixProcess (sysmatrix); /* factorize for direct solution */ /* form preconditioner for iterative solution */ vfs_SysMatrixFactor (sysmatrix); /* create load and displacement sysvector objects */ sysvectorb = vfs_SysVectorBegin(); vfs_SysVectorDef (sysvectorb,neq); sysvectorx = vfs_SysVectorBegin(); vfs_SysVectorDef (sysvectorx,neq); /* fill load vector using vfs_SysVectorAssem or vfs_SysVectorScatter, etc. */ /* solve */ vfs_SysMatrixSolve (sysmatrix,sysvectorb,sysvectorx); /* access displacement vector using vfs_SysVectorGather */ /* delete objects */ vfs_DofTabEnd (doftab); vfs_SysMatrixEnd (sysmatrix); vfs_SysVectorEnd (sysvectorb); vfs_SysVectorEnd (sysvectorx);Note that multiple load cases may be solved by calling vfs_SysMatrixSolve for each load case.
The blocked Lanczos method is specifically designed for extracting large numbers of vibration modes through a frequency range. In this case, the Lanczos method will use an automatic shifting strategy, if required, to move through the requested number of frequencies and/or frequency range from low to high frequencies. Several shifts may be required to complete the extraction of the requested vibration modes. The automatic shifting strategy is activated with any eigenvalue request which asks for the lowest modes. The initial shift in these cases is always to some point at or below the lowest frequency bound. Additional shifts are applied if all the requested eigenvalues are not found with the initial shift. Any request asking for eigenvalues nearest to a frequency will be extracted using a single shift. The user must be aware of this fact when designing the eigenvalue request.
The block size used by the Lanczos method is set using vfs_SysMatrixSetParami with type SYSMATRIX_LANCZOS_BLOCKSIZE. In general the block size should be set to be greater than the expected eigenvalue multiplicity. For structural vibration problems with rigid body modes a value of 7 is suggested.
The AMLS eigensolution method is applicable to large scale vibration problems in which large numbers of eigenvalues are desired. The eigensolution is approximate in the sense that at successive substructures in which intermediate eigensolutions are performed, the eigensolution is truncated to some degree. The method can be several times faster than an equivalent block Lanczos solution. The AMLS method requires a factorization of the mass matrix in addition to a factorization of the stiffness matrix. Out-of-core storage of both matrices may be enabled usng vfs_SysMatrixSetParami with type SYSMATRIX_FACTOR_OOC. The names of the out-of-core files may be specified using vfs_SysMatrixSetOocFile. The degree of eigenvalue truncation affects the accuracy and runtime of the AMLS method. Use SYSMATRIX_AMLS_CUTFACT1 to set the frequency cutoff factor. This factor is multiplied by the high frequency limit of the requested frequency spectrum. Use SYSMATRIX_AMLS_CUTFACT2 to control a secondary cutoff frequency used in the reduced problem. This parameter is particularly effective in increasing the accuracy of the highest frequencies.
The subspace method should be used only for buckling problems. There is no automatic shifting strategy available for this method due to the normally low number of modes requested.
The function vfs_SysMatrixEigen is designed to compute a subset of eigenvalues and eigenvectors for large systems using the subspace, Lanczos or AMLS method. The function vfs_SysMatrixEigenAll is designed to compute all eigenvalues and eigenvectors for the generalized eigenproblem for small systems.
The computation of eigenvectors is performed during the eigenvalue extraction process for the subspace and block Lanczos eigensolvers. It is performed, on demand, as a post processing step for the AMLS eigensolver. Therefore, an option exists for out-of-core storage of the eigenvectors to reduce memory usage for the subspace and Lanczos eigensolvers, use vfs_SysMatrixSetParami with type SYSMATRIX_EIGEN_OOC. The file name used for the out-of-core vectors may be specified using vfs_SysMatrixSetOocFile. Use the function vfs_SysMatrixEigenMode to compute/retrieve eigenvectors. For the subspace and Lanczos methods, the pre-computed eigenvectors are retrieved from memory or read from the specified out-of-core file. For the AMLS method the eigenvectors are computed when requested in memory blocks for computational efficiency. The requested eigenvector is selected from the currently computed block. If there are memory constraints, use the function vfs_SysMatrixSetParami with type SYSMATRIX_AMLS_RECOVERSIZE to set the number of vectors in the computational block. Since eigenvectors are computed on-demand for the AMLS method, the vfs_SysMatrixEigenMode function call can take considerable time if the computation of a new block of eigenvectors is triggered by the call. After all desired eigenvectors are retrieved, free memory and delete any out-of-core file using vfs_SysMatrixEigenFree.
*vfs_SysMatrixBegin - create an instance of a SysMatrix object
vfs_SysMatrix *vfs_SysMatrixBegin ()
None
Destroy an instance of a SysMatrix object using
void vfs_SysMatrixEnd (vfs_SysMatrix *sysmatrix)
Return the current value of a SysMatrix object error flag using
Vint vfs_SysMatrixError (vfs_SysMatrix *sysmatrix)
vfs_SysMatrixPre - specify real or complex matrix
void vfs_SysMatrixPre (vfs_SysMatrix *sysmatrix, Vint complextype)
sysmatrix Pointer to SysMatrix object. complextype Enumerated type for real or complex matrix =SYS_DOUBLE Real matrix =SYSMATRIX_DOUBLECOMPLEX Complex matrix
None
SYS_ERROR_ENUM is issued if complextype is invalid.
This call must be made before vfs_SysMatrixDef is called. Defaults to real.
vfs_SysMatrixAbort - abort a sparse direct factorization
void vfs_SysMatrixAbort (vfs_SysMatrix *sysmatrix)
sysmatrix Pointer to SysMatrix object.
None
None
vfs_SysMatrixAssem,vfs_SysMatrixAssemDiag,vfs_SysMatrixAssemFull - assemble an element matrix
void vfs_SysMatrixAssem (vfs_SysMatrix *sysmatrix, Vint nedofs, Vint lm[], Vdouble s[]) void vfs_SysMatrixAssemDiag (vfs_SysMatrix *sysmatrix, Vint nedofs, Vint lm[], Vdouble sd[]) void vfs_SysMatrixAssemFull (vfs_SysMatrix *sysmatrix, Vint nedofs, Vint lm[], Vdouble sf[])
sysmatrix Pointer to SysMatrix object. nedofs Number of element degrees of freedom lm Vector of element global dof numbers s Element matrix sd Element diagonal matrix sf Element full matrix
None
vfs_SysMatrixDef - define number of global dof and structure
void vfs_SysMatrixDef (vfs_SysMatrix *sysmatrix, Vint type, Vint neq)
sysmatrix Pointer to SysMatrix object. type Type of matrix =SYSMATRIX_SYMM_SPARSE Symmetric sparse =SYSMATRIX_NSYMM_SPARSE Non-symmetric sparse =SYSMATRIX_SYMM_FULL Symmetric full =SYSMATRIX_NSYMM_FULL Non-symmetric full =SYSMATRIX_SYMM_ITER Symmetric iterative =SYSMATRIX_NSYMM_ITER Non-symmetric iterative =SYSMATRIX_DIAG Diagonal neq Number of global dof
None
Inquire of a defined type and neq as output arguments using
void vfe_SysMatrixInq (vfe_SysMatrix *sysmatrix, Vint *type, Vint *neq)
vfs_SysMatrixDiag - extract the diagonal
void vfs_SysMatrixDiag (vfs_SysMatrix *sysmatrix, vfs_SysVector *sysvector)
sysmatrix Pointer to SysMatrix object. sysvector Pointer to SysVector object
None
vfs_SysMatrixColumn - extract a column
void vfs_SysMatrixColumn (vfs_SysMatrix *sysmatrix, Vint lmcol, vfs_SysVector *sysvector)
sysmatrix Pointer to SysMatrix object. lmcol Column number to extract. sysvector Pointer to SysVector object
None
vfs_SysMatrixEigenAll - extract all eigenvalues and eigenvectors
void vfs_SysMatrixEigenAll (vfs_SysMatrix *sysmatrix, vfs_SysMatrix *sysmatrixb, vfs_SysVector *sysvectore, vfs_SysVector *sysvectorx[])
sysmatrix Pointer to SysMatrix object for matrix A. sysmatrixb Pointer to SysMatrix object for matrix B. sysvectore Pointer to SysVector object of eigenvalues sysvectorx Pointer to array of SysVector objects of eigenvectors X.
None
The type of eigenproblem solved is dependent upon the SYSMATRIX_EIGEN parameter set using vfs_SysMatrixSetParami. If the eigenproblem type is SYSMATRIX_EIGEN_STAN then the matrix B is assumed to be the identity and sysmatrixb is ignored. If the eigenproblem type is SYSMATRIX_EIGEN_GENERAL then the matrix B must be input and the eigenvalues as defined above are returned. If the eigenproblem type is SYSMATRIX_EIGEN_BUCK then the matrix B is assumed to be the geometric stiffness matrix. The eigenvalues returned in sysvectore will be set to -lambda. If the eigenproblem type is SYSMATRIX_EIGEN_VIBE then the matrix B is assumed to be a mass matrix. The eigenvalues returned in sysvectore will be set to sqrt(lambda)/(2*pi). If the eigenproblem type is SYSMATRIX_EIGEN_RIGID then the matrix B is assumed to be the identity and sysmatrixb is ignored.
vfs_SysMatrixEigen - extract set of eigenvalues and eigenvectors
void vfs_SysMatrixEigen (vfs_SysMatrix *sysmatrix, vfs_SysMatrix *sysmatrixb, Vint *neigen, Vdouble eigenval[])
sysmatrix Pointer to SysMatrix object for matrix A. sysmatrixb Pointer to SysMatrix object for matrix B.
neigen Number of eigenvalues found eigenval Vector of neigen computed eigenvalues
The eigenvalue/eigenvector pairs computed by this function depend on the values of integer and double precision parameters set, respectively, with vfs_SysMatrixSetParami and vfs_SysMatrixSetParamd. The supported combinations are outlined in the table below.
Supported eigenvalue parameter combinations
Case |
SYSMATRIX_EIGEN_TYPE | SYSMATRIX_EIGEN_LOWERTEST | SYSMATRIX_EIGEN_UPPERTEST | SYSMATRIX_EIGEN_NUM |
1 |
SYS_EIGEN_ALL | SYS_ON | SYS_ON | N/A |
2 |
SYS_EIGEN_ALL | SYS_OFF | SYS_ON | N/A |
3 |
SYS_EIGEN_LOWEST | SYS_ON | SYS_ON | 0 |
4 |
SYS_EIGEN_LOWEST | SYS_ON | SYS_ON | n |
5 |
SYS_EIGEN_LOWEST | SYS_OFF | SYS_ON | 0 |
6 |
SYS_EIGEN_LOWEST | SYS_OFF | SYS_ON | n |
7 |
SYS_EIGEN_LOWEST | SYS_ON | SYS_OFF | 0 |
8 |
SYS_EIGEN_LOWEST | SYS_ON | SYS_OFF | n |
9 |
SYS_EIGEN_NEAREST | N/A | N/A | 0 |
10 |
SYS_EIGEN_NEAREST | N/A | N/A | n |
Case 1 will compute all eigenvalues found in the interval defined by [SYSMATRIX_EIGEN_LOWER,SYSMATRIX_EIGEN_UPPER]. The number of requested eigenvalues is ignored, and the user should ensure that the eigenval argument is sufficiently large to hold all eigenvalues encountered in the interval.
Case 2 is similar to case 1, except that the interval is replaced by [0,SYSMATRIX_EIGEN_UPPER].
Case 3 computes the single lowest eigenvalue in the interval [SYSMATRIX_EIGEN_LOWER,SYSMATRIX_EIGEN_UPPER].
Case 4 is similar to Case 3, except that the n lowest eigenvalues in the interval are computed.
Cases 5 and 6 are similar to Cases 3 and 4, respectively, except that the interval is replaced by [0,SYSMATRIX_EIGEN_UPPER].
Cases 7 and 8 are similar to cases 3 and 4, respectively, except that the interval is replaced by [SYSMATRIX_EIGEN_LOWER,infinity].
Case 9 computes the single eigenvalue closest to the shift SYSMATRIX_EIGEN_SHIFT, while Case 10 computes the n closest eigenvalues to the shift SYSMATRIX_EIGEN_SHIFT.
When using the AMLS eigensolver, only Case 4 and Case 6 are supported. The number of eigenvalues and upper bound must be specified. The lower bound is optional.
vfs_SysMatrixEigenMode - retrieve a computed eigenvector
void vfs_SysMatrixEigenMode (vfs_SysMatrix *sysmatrix, Vint ieigen, vfs_SysVector *sysvector,
sysmatrix Pointer to SysMatrix object for matrix ieigen Eigen mode to retrieve sysvector Pointer to SysVector for eigenvector
None
vfs_SysMatrixEigenFree - free all computed eigenvector storage
void vfs_SysMatrixEigenFree (vfs_SysMatrix *sysmatrix)
sysmatrix Pointer to SysMatrix object for matrix
None
vfs_SysMatrixFactor - factorize
void vfs_SysMatrixFactor (vfs_SysMatrix *sysmatrix)
sysmatrix Pointer to SysMatrix object.
None
vfs_SysMatrixGather - gather non-zero values of a column
void vfs_SysMatrixGather (vfs_SysMatrix *sysmatrix, Vint lmcol, Vint *ndof, Vint lm[], Vdouble s[])
sysmatrix Pointer to SysMatrix object. lmcol Column number to extract.
ndofs Number of non-zero degrees of freedom lm Vector of global dof numbers of non-zero values s Gathered values
vfs_SysMatrixGetDouble - get double precision matrix information
void vfs_SysMatrixGetDouble (vfs_SysMatrix *sysmatrix, Vint type, Vdouble param[])
sysmatrix Pointer to SysMatrix object. type Type of integer information to query =SYSMATRIX_ITER_FCONV Iterative solver residual error =SYSMATRIX_ITER_UCONV Iterative solver solution error =SYSMATRIX_ITER_ECONV Iterative solver energy error =SYSMATRIX_EIGENSHIFT The currently used eigenvalue shift =SYSMATRIX_EIGENVAL Eigenvalues in subspace =SYSMATRIX_EIGENERR Error in each eigenvalue in subspace =SYSMATRIX_FACTOR_OPS Number of FLOPS for factor =SYSMATRIX_MANTISSA Mantissa of matrix determinant =SYSMATRIX_MINPIVOTVAL The smallest pivot reduction =SYSMATRIX_MOD_DIAG The pivot's algorithmically updated value =SYSMATRIX_NEW_DIAG The factorized value of the diagonal pivot =SYSMATRIX_OLD_DIAG The original value of the diagonal pivot =SYSMATRIX_SOLVE_OPS Number of FLOPS for solve
param Returned double precision information
The double precision queries SYSMATRIX_FACTOR_OPS and SYSMATRIX_SOLVE_OPS require the system matrix to have been processed using vfs_SysMatrixProcess. The double precision query SYSMATRIX_MANTISSA requires the system matrix to have been factorized using vfs_SysMatrixFactor. The matrix determinant is represented by a double precision mantissa times ten raised to an exponent. The mantissa is returned using the type SYSMATRIX_MANTISSA. Query for the exponent using vfs_SysMatrixGetInteger with type SYSMATRIX_EXPONENT.
The double precision queries SYSMATRIX_ITER_FCONV, SYSMATRIX_ITER_UCONV and SYSMATRIX_ITER_ECONV return the normalized residual, solution and energy convergence numbers. These numbers are updated at every iteration, and can be queried in a monitor function. They can be directly compared with the tolerances specified with vfs_SysMatrixSetParamd using the types SYSMATRIX_ITER_FTOL, SYSMATRIX_ITER_UTOL and SYSMATRIX_ITER_ETOL respectively.
The queries SYSMATRIX_OLD_DIAG and SYSMATRIX_NEW_DIAG return the pivot values before and after factorization. If SYSMATRIX_AUTOSPC is set to SYS_ON, then SYSMATRIX_MOD_DIAG gives the value of the algorithmically modified pivot. All three values are best used in a monitor function along with the integer parameters SYSMATRIX_FACTOR_EQN and SYSMATRIX_FACTOR_UEQN.
Use SYSMATRIX_EIGENSHIFT to retrieve the currently used eigenvalue shift.
Use SYSMATRIX_EIGENVAL and SYSMATRIX_EIGENERR to return an array dimensioned the subspace size with the value of each eigenvalue in the subspace, converged or not, within the required limits or not, and its error, respectively. The subspace size can be obtained by using vfs_SysMatrixGetInteger and using the parameter SYSMATRIX_SUBSPACE_SIZE.
After factorization, use SYSMATRIX_MINPIVOTVAL to retrieve the minimum relative pivot value drop.
vfs_SysMatrixGetInteger - get integer matrix information
void vfs_SysMatrixGetInteger (vfs_SysMatrix *sysmatrix, Vint type, Vint param[])
sysmatrix Pointer to SysMatrix object. type Type of integer information to query =SYSMATRIX_EIGENCONV Number of converged eigenvalues =SYSMATRIX_EIGENPHASE Phase in an eigen solution =SYSMATRIX_EIGENSTEP The subspace iteration number =SYSMATRIX_EXPONENT Exponent base ten of matrix determinant =SYSMATRIX_FACTORIZED Factorized flag =SYSMATRIX_FACTOR_EQN Current factor column from 1 to neq =SYSMATRIX_FACTOR_UEQN Current user-numbered factor equation =SYSMATRIX_ITERSTATUS Status of the iterative solution =SYSMATRIX_MAXSUPERSIZE Size of largest supernode =SYSMATRIX_MINPIVOTID User-numbered minimum pivot equation =SYSMATRIX_MOD_DIAG Toggle for modified factored pivot =SYSMATRIX_NUMITER Number of iterative solver iterations =SYSMATRIX_NUMNEG Number negative roots =SYSMATRIX_PHASE Phase object is in =SYSMATRIX_REQEIGEN Number of requested eigenvalues =SYSMATRIX_SINGULAR Matrix singular flag =SYSMATRIX_STACK_SIZE Stack size for out-of-core solver =SYSMATRIX_SUBSPACE_SIZE Computed eigenvalue subspace size =SYSMATRIX_ZEROPIVOT Zero pivot degree of freedom
param Returned integer information
The integer queries SYSMATRIX_NUMNEG and SYSMATRIX_EXPONENT require the system matrix to have been factorized using vfs_SysMatrixFactor. SYSMATRIX_MAXSUPERSIZE is the size of the largest supernode used during sparse factorization. Use SYSMATRIX_MINPIVOTID to retrieve the user equation number that experienced the largest relative drop in its value during factorization.
The matrix determinant is represented by a double precision mantissa times ten raised to an exponent. The exponent is returned using the type SYSMATRIX_EXPONENT. Query for the mantissa using vfs_SysMatrixGetDouble with type SYSMATRIX_MANTISSA.
Query for the status of the factorization using SYSMATRIX_FACTORIZED. If the function vfs_SysMatrixFactor fails due to a singular matrix, query for this condition using SYSMATRIX_SINGULAR. If the matrix is singular, the degree of freedom of the singular pivot may be returned using SYSMATRIX_ZEROPIVOT.
Query for the number of iterations with type SYSMATRIX_NUMITER. During an iterative solution, the value of this parameter is updated at every iteration, making this an ideal parameter to track the iteration procedure with a monitor function.
Use a monitor function to retrieve updated parameters during computations.
During vfs_SysMatrixFactor, if a sparse solver is being used, SYSMATRIX_FACTOR_EQN reflects the monotonically increasing column number being factorized. Likewise, SYSMATRIX_FACTOR_UEQN returns the equation number, in user space, for the current column. These two numbers differ because of reordering that is done for efficiency purposes. If SYSMATRIX_AUTOSPC has been set, then SYSMATRIX_MOD_DIAG returns SYS_ON if the singularity removal algorithm was invoked for this equation.
During vfs_SysMatrixEigen use SYSMATRIX_EIGENPHASE to determine the phase in process. You can also determine that vfs_SysMatrixEigen has been called by querying for SYSMATRIX_PHASE and verifying that its value is SYSMATRIX_PHASE_EIGENSOLUTION. If that's the case, all eigensolver-related queries are available. Depending on the phase, certain parameters may be available. SYSMATRIX_SUBSPACE_SIZE is the size of the subspace used by the eigensolver; SYSMATRIX_EIGENCONV is the number of eigenvalues converged up to the current iteration; SYSMATRIX_EIGEN_STEP is the current subspace iteration number.
The out-of-core solver requires auxiliary memory in the form of a stack. Use SYSMATRIX_STACK_SIZE to query for the number of double precision works used in the stack.
Use SYSMATRIX_PHASE to determine the phase SysMatrix is currently executing. Only two phases are meaningful: SYSMATRIX_PHASE_EIGENSOLUTION and SYSMATRIX_PHASE_SOLVE. That's because these two phases provide the user with additional information.
During vfs_SysMatrixSolve - which can be determined by verifying that SYSMATRIX_PHASE takes the value SYSMATRIX_PHASE_SOLVE - you can query for additional iteration parameters depending on the iteration iteration status - SYSMATRIX_ITERSTATUS. If the status is set to SYSMATRIX_ITERSTATUS_INIT then the iteration is about to start; if set to SYSMATRIX_ITERSTATUS_ITER then SysMatrix is iterating to find a solution. If set to SYSMATRIX_ITERSTATUS_TERM then the iteration has terminated. During the iterative solution you also query for SYSMATRIX_NUMITER, the number of iterations. You can also use vfs_SysMatrixGetDouble to query for specific convergence norms: SYSMATRIX_ITER_FCONV, for residual convergence; SYSMATRIX_ITER_UCONV, for solution convergence; and SYSMATRIX_ITER_ECONV, for energy convergence. If the GMRES solver is used only SYSMATRIX_ITER_FCONV is available.
vfs_SysMatrixGetLong - get long matrix information
void vfs_SysMatrixGetLong (vfs_SysMatrix *sysmatrix, Vint type, Vlong lparam[])
sysmatrix Pointer to SysMatrix object. type Type of integer information to query =SYSMATRIX_ASSEM_NONZ Number non zeros in assembled matrix =SYSMATRIX_FACTOR_NONZ Number Non zeros in factored matrix
param Returned integer information
The number of non-zero entries in a sparse matrix can be retrieved with SYSMATRIX_ASSEM_NONZ.
The query SYSMATRIX_FACTOR_NONZ requires the system matrix to have been processed using vfs_SysMatrixProcess.
vfs_SysMatrixGetSysVector - retrieve internal SysVector objects
void vfs_SysMatrixGetSysVector (vfs_SysMatrix *sysmatrix, Vint type, vfs_SysVector **sysvector)
sysmatrix Pointer to SysMatrix object. type Type of SysVector to retrieve =SYSMATRIX_VEC_PRODA Vector a =SYSMATRIX_VEC_PRODB Vector b
None
vfs_SysMatrixInit - initialize matrix to a given value
void vfs_SysMatrixInit (vfs_SysMatrix *sysmatrix, Vdouble v)
sysmatrix Pointer to SysMatrix object. v Initialization value
None
vfs_SysMatrixInner - form matrix as inner product of vectors
void vfs_SysMatrixInner (vfs_SysMatrix *sysmatrix, Vint nvec, vfs_SysVector *sysvectorx[], vfs_SysVector *sysvectory[])
sysmatrix Pointer to SysMatrix object A. sysvectorx Pointer to vector of SysVector objects X sysvectory Pointer to vector of SysVector objects Y v Initialization value
None
vfs_SysMatrixPreProcess - establish matrix structure
void vfs_SysMatrixPreProcess (vfs_SysMatrix *sysmatrix)
sysmatrix Pointer to SysMatrix object.
None
vfs_SysMatrixProcess - process before factorization
void vfs_SysMatrixProcess (vfs_SysMatrix *sysmatrix)
sysmatrix Pointer to SysMatrix object.
None
vfs_SysMatrixProd - multiply by a system vector
void vfs_SysMatrixProd (vfs_SysMatrix *sysmatrix, vfs_SysVector *sysvectorx, vfs_SysVector *sysvectory)
sysmatrix Pointer to SysMatrix object. sysvectorx Pointer to SysVector object x sysvectorx Pointer to SysVector object y
None
vfs_SysMatrixReact - compute reaction forces from applied displacements
void vfs_SysMatrixReact (vfs_SysMatrix *sysmatrix, vfs_SysVector *sysvectorx, vfs_SysVector *sysvectory)
sysmatrix Pointer to SysMatrix object. sysvectorx Pointer to SysVector object x
sysvectorx Pointer to SysVector object y
vfs_SysMatrixReduce - reduce a matrix
void vfs_SysMatrixReduce (vfs_SysMatrix *sysmatrix, Vint nvec, Vint nattach, Vint nvibe, Vint lmattach[], vfs_SysVector *sysvector[], vfs_SysMatrix *sysmatrixq)
sysmatrix Pointer to SysMatrix object. nvec Number of system vectors 0 < nvec nattach Number of attachment modes. nvibe Number of vibration modes. lmattach Vector of attachment mode global dof numbers sysvector Pointer to array of nvec SysVector objects sysmatrixq Pointer to SysMatrix object of reduced matrix
None
If the matrix to be reduced is the stiffness matrix then there is a considerable computational efficiency in specifying those modes which are attachment modes. Attachment modes are those modes due to applied unit displacement at a specific degree of freedom. These are the "retained" degrees of freedom. These nattach modes must appear first in the set of sysvector modes and nattach related global degrees of freedom must be specified in lmattach. If vibration modes are present they must come after the attachment points. If the matrix to be reduced is not a stiffness matrix (eg. mass matrix), set nattach to zero and lmattach may be input NULL.
In matrix form, given
sysvector[0] = [ a_00 a_01 ... a_0N ] sysvector[1] = [ a_10 a_11 ... a_1N ] ... sysvector[n] = [ a_n0 a_n1 ... a_nN ]
then
-- -- -- -- | a_00 a_01 ... a_0N | | a_00 a_10 ... a_n0 | sysmatrixq = | a_10 a_11 ... a_1N | sysmatrix | a_01 a_11 ... a_n1 | | ... | | ... | | a_n0 a_n1 ... a_nN | | a_0N a_1N ... a_nN | -- -- __ __
vfs_SysMatrixSetComplexMode - set complex mode
void vfs_SysMatrixSetComplexMode (vfs_SysMatrix *sysmatrix, Vint complexmode)
sysmatrix Pointer to SysMatrix object. complexmode The type of data being passed through the Vdouble[] argument list. =SYS_COMPLEX_REAL Data is real =SYS_COMPLEX_IMAGINARY Data is imaginary =SYS_COMPLEX_REALIMAGINARY Data is complex
None
The methods vfs_SysMatrixAssem, vfs_SysMatrixAssemDiag, vfs_SysMatrixAssemFull, and all take arrays of floating point numbers as lower triangular, full, or diagonal element matrices. In addition, vfs_SysMatrixColumn returns an array of floating point numbers for the column values. The values entered or returned can be declared as SYS_COMPLEX_REAL, SYS_COMPLEX_IMAGINARY, or SYS_COMPLEX_REALIMAGINARY. The data in the argument then takes the specified mode.
This is useful if a real element stiffness matrix is being assembled. The mode SYS_COMPLEX_REAL indicates that only one value is being entered for each a_ij in the array, and this value is assembled into the matrix. In this case, the imaginary entries in the matrix remain unchanged. Likewise, if an imaginary damping matrix is being assembled, the mode SYS_COMPLEX_IMAGINARY can be specified so that only the imaginary parts of SysMatrix are updated. If mode is set to SYS_COMPLEX_REALIMAGINARY then two consecutive numbers in the array stand for the real and imaginary parts of the data. Defaults to SYS_COMPLEX_REAL.
vfs_SysMatrixSetThreadComplexMode - set complex mode in thread
void vfs_SysMatrixSetThreadComplexMode (vfs_SysMatrix *sysmatrix, Vint ithread, Vint complexmode)
sysmatrix Pointer to SysMatrix object. ithread Thread number starting from 1. complexmode The type of data being passed through the Vdouble[] argument list. =SYS_COMPLEX_REAL Data is real =SYS_COMPLEX_IMAGINARY Data is imaginary =SYS_COMPLEX_REALIMAGINARY Data is complex
None
The methods vfs_SysMatrixThreadAssem, vfs_SysMatrixThreadAssemDiag, and vfs_SysMatrixThreadAssemFull all take arrays of floating point numbers as lower triangular, full, or diagonal element matrices. In addition, vfs_SysMatrixColumn returns an array of floating point numbers for the column values. The values entered or returned can be declared as SYS_COMPLEX_REAL, SYS_COMPLEX_IMAGINARY, or SYS_COMPLEX_REALIMAGINARY. The data in the argument then takes the specified mode.
This is useful if a real element stiffness matrix is being assembled. The mode SYS_COMPLEX_REAL indicates that only one value is being entered for each a_ij in the array, and this value is assembled into the matrix. In this case, the imaginary entries in the matrix remain unchanged. Likewise, if an imaginary damping matrix is being assembled, the mode SYS_COMPLEX_IMAGINARY can be specified so that only the imaginary parts of SysMatrix are updated. If mode is set to SYS_COMPLEX_REALIMAGINARY then two consecutive numbers in the array stand for the real and imaginary parts of the data. Defaults to SYS_COMPLEX_REAL.
vfs_SysMatrixSetFunction - set a function pointer
void vfs_SysMatrixSetFunction (vfs_SysMatrix *sysmatrix, Vint functype, void (*function)(vfs_SysMatrix*,Vobject*) Vobject *object)
sysmatrix Pointer to SysMatrix object. functype The function type identifier =SYSMATRIX_FUN_MONITOR Monitoring function =SYSMATRIX_FUN_PROD External function that multiplies the system matrix by a vector function Pointer to the function to be called object Pointer to a user object to be passed to the monitoring function.
None
The method vfs_SysMatrixAbort can be called from within the monitoring function if necessary.
The available data SysMatrix can be queried for depends on the function called. Functions for which a monitor function is called are vfs_SysMatrixFactor, vfs_SysMatrixSolve, and vfs_SysMatrixEigen.
A monitor function can be used with vfs_SysMatrixFactor if vfs_SysMatrixDef is defined with SYSMATRIX_SYMM_SPARSE, SYSMATRIX_NSYMM_SPARSE, SYSMATRIX_SYMM_FULL, or SYSMATRIX_SYMM_FULL, SYSMATRIX_TECH is set to SYSMATRIX_TECH_SPARSE (i.e., a direct, sparse factorization is being performed). The function vfs_SysMatrixGetInteger can be called to query for the parameter SYSMATRIX_FACTOR_EQN, the column number currently being factored. This parameter varies from 1 to the number of equations and can be used to measure progress in the factorization. The user equation can also be retrieved via the parameter SYSMATRIX_FACTOR_UEQN. You can also query for three different double precision parameters with vfs_SysMatrixGetDouble: SYSMATRIX_OLD_DIAG, SYSMATRIX_NEW_DIAG, and SYSMATRIX_UPD_DIAG. The first parameter gives you the original value of the diagonal for the equation being processed; the second gives its value after factorization. If SYSMATRIX_AUTOSPC is set to SYS_ON, the parameter SYSMATRIX_UPD_DIAG gives the value of the updated pivot. In addition, the integer parameter SYSMATRIX_MOD_DIAG flags whether the pivot has been modified by the presence of SYSMATRIX_AUTOSPC.
During factorization, the frequency with which the monitor function is called can be controlled with the integer parameter SYSMATRIX_MONITOR_FREQ, which defaults to 1000.
When used along with vfs_SysMatrixEigen, use vfs_SysMatrixGetInteger with the parameter SYSMATRIX_EIGENPHASE to determine the current phase in the eigensolver. The following phases are available:
Phase |
Description |
SYSMATRIX_EIGENPHASE_FACTOR |
Matrix being factored |
SYSMATRIX_EIGENPHASE_SETUP |
Factorization complete |
SYSMATRIX_EIGENPHASE_ITER |
Subspace iteration |
During factorization, SYSMATRIX_FACTOR_EQN can be queried for. During the setup phase, vfs_SysMatrixGetInteger can be used to query for the subspace size with SYSMATRIX_SUBSPACE_SIZE. During the subspace iteration process, vfs_SysMatrixGetInteger can be used to query for the iteration number with SYSMATRIX_EIGENSTEP, and for the number of converged eigenvalues at the step, with SYSMATRIX_EIGENCONV. At any time during vfs_SysMatrixEigen you can retrieve the number of eigenvalues requested using vfs_SysMatrixGetInteger with the parameter SYSMATRIX_NUM_REQEIGEN.
If SYSMATRIX_FUN_PROD is used, the method vfs_SysMatrixGetSysVector must be called from within the product function; the product is constructed as sysvectorb = A sysvectora. If there are MPCs set, these will be taken care of internally by SysMatrix.
vfs_SysMatrixSetInitSol - set initial solution vector
void vfs_SysMatrixSetInitSol (vfs_SysMatrix *sysmatrix, vfs_SysVector *initsol)
sysmatrix Pointer to SysMatrix object. initsol Pointer to SysVector object initial guess solution
None
vfs_SysMatrixSetNumThreads - set number of threads
void vfs_SysMatrixSetNumThreads (vfs_SysMatrix *sysmatrix, Vint num)
sysmatrix Pointer to SysMatrix object. num Number of threads =0 Default, serial execution. >=1 Number of threads for parallel execution
None
vfs_SysMatrixSetObject - set attribute object
void vfs_SysMatrixSetObject (vfs_SysMatrix *sysmatrix, Vint objecttype, Vobject *object)
sysmatrix Pointer to SysMatrix object. objecttype The object type identifier =VFS_DOFTAB DofTab object object Pointer to the object to be set.
None
vfs_SysMatrixSetOocFile - set out-of-core filename
void vfs_SysMatrixSetOocFile (vfs_SysMatrix *sysmatrix, Vint type, Vchar *file[])
sysmatrix Pointer to SysMatrix object. type Parameter type to set. =SYSMATRIX_KFACTOR_FILE Out-of-core stiffness factor =SYSMATRIX_MFACTOR_FILE Out-of-core mass factor =SYSMATRIX_EIGEN_FILE Out-of-core eigenvectors file Specifies the name of the file to be used
none
none
Use SYSMATRIX_MFACTOR_FILE to specify the name of the file to be used to store the mass matrix factor during an out-of-core matrix factorization. Defaults to vki.moc. A mass matrix factor is generated by the AMLS eigensolution procedure.
Use SYSMATRIX_EIGEN_FILE to specify the name of the file to be used to store the eigenvectors during eigen extraction. Defaults to vki.eoc.
vfs_SysMatrixSetParamd - set solution parameters
void vfs_SysMatrixSetParamd (vfs_SysMatrix *sysmatrix, Vint type, Vdouble param)
sysmatrix Pointer to SysMatrix object. type Parameter type to set. =SYSMATRIX_ITER_UTOL Convergence tolerance displacement =SYSMATRIX_ITER_FTOL Convergence tolerance force =SYSMATRIX_ITER_ETOL Convergence tolerance energy =SYSMATRIX_ITER_ATOL Convergence tolerance absolute =SYSMATRIX_EIGEN_LOWER Lower limit of eigenvalue interval =SYSMATRIX_EIGEN_UPPER Upper limit of eigenvalue interval =SYSMATRIX_EIGEN_SHIFT Given eigenvalue shift =SYSMATRIX_AMLS_CUTFACT1 Overall frequency cutoff factor =SYSMATRIX_AMLS_CUTFACT2 Reduced frequency cutoff factor =SYSMATRIX_LANCZOS_UTOL Lanczos eigenvalue tolerance =SYSMATRIX_SUBSPACE_UTOL Subspace eigenvalue tolerance =SYSMATRIX_AUTOSPC_FACTOR Factor singularity penalty multiplier =SYSMATRIX_FACTOR_SINGTOL Relative singularity tolerance =SYSMATRIX_FACTOR_ABSZERO Absolute singularity tolerance param Specifies the double value that type will be set to.
none
Specify the lower limit of the eigenvalue interval using SYSMATRIX_EIGEN_LOWER and the upper limit using SYSMATRIX_EIGEN_UPPER. The defaults are [0.,1.].
Specify a given eigenvalue shift using SYSMATRIX_EIGEN_SHIFT. The default is 0.
Specify an overall cutoff factor at substructure levels using SYSMATRIX_AMLS_CUTFACT1. The frequencies at each substructure level are truncated at this factor times the upper frequency limit of requested spectrum. This value will affect the accuracy of all modes. Suggest specified values be no larger than 10. By default SYSMATRIX_AMLS_CUTFACT1 is 5.
Specify a reduced cutoff factor using SYSMATRIX_AMLS_CUTFACT2. Higher values will have a tendency to increase the accuracy of higher modes. By default SYSMATRIX_AMLS_CUTFACT2 is 1.7 Suggest specified values be no larger than 2.5. Note that 1. < SYSMATRIX_AMLS_CUTFACT2 < SYSMATRIX_AMLS_CUTFACT1.
Specify a subspace eigenvalue convergence tolerance using SYSMATRIX_SUBSPACE_UTOL. The default is .0001 . Specify a Lanczos eigenvalue convergence tolerance using SYSMATRIX_LANCZOS_UTOL. The default is 3.7e-11 It is not recommended to alter these values.
Specify a factorization singularity tolerance using SYSMATRIX_FACTOR_SINGTOL. This tolerance is used to compare the original diagonal pivot to its reduced value. If the integer parameter SYSMATRIX_AUTOSPC is set to SYS_ON, if a singularity is encountered then the reduced pivot is replaced by the value given by the parameter SYSMATRIX_AUTOSPC_FACTOR multiplied by the matrix largest entry. This procedure is equivalent to introducing a penalty parameter that multiplies the matrix diagonal in order to enforce a constraint. Users should exert care in using this procedure, especially if there is a load applied in the direction of the singularity, or if the system is indefinite. By default SYSMATRIX_FACTOR_SINGTOL is 1.E-13 and SYSMATRIX_AUTOSPC_FACTOR is 1.e+5.
Specify an absolute singularity tolerance using SYSMATRIX_FACTOR_ABSZERO. The system is considered singular if a pivot's absolute value falls below this threshold.
vfs_SysMatrixSetParami - set solution parameters
void vfs_SysMatrixSetParami (vfs_SysMatrix *sysmatrix, Vint type, Vint param)
sysmatrix Pointer to SysMatrix object. type Integer parameter type to set. =SYSMATRIX_EIGEN Eigenproblem procedure =SYSMATRIX_EIGEN_LOWERTEST Test for lower bound on eigenvalue =SYSMATRIX_EIGEN_UPPERTEST Test for upper bound on eigenvalue =SYSMATRIX_EIGEN_MAXITER Maximum iterations per shift =SYSMATRIX_EIGEN_METHOD Eigensolver method for vibration =SYSMATRIX_EIGEN_NORM Eigenvector normalization =SYSMATRIX_EIGEN_NUM Number of eigenvalues requested =SYSMATRIX_EIGEN_TYPE Eigensolution type =SYSMATRIX_EIGEN_ZEROMODES Rigid body modes flag =SYSMATRIX_EIGEN_INITVECTYPE Vibration eigenvector initialization type =SYSMATRIX_EIGEN_OOC Out-of-core eigenvector flag =SYSMATRIX_AMLS_RECOVERSIZE AMLS vector recovery block size =SYSMATRIX_LANCZOS_BLOCKSIZE Lanczos block size =SYSMATRIX_TECH Solution procedure =SYSMATRIX_CACHESIZE Cache size (KiloBytes) =SYSMATRIX_WORKSIZE_FACTOR Out of core factor memory (MB) =SYSMATRIX_SOLVETEMP Toggle applied displacements during solve =SYSMATRIX_MAXCGITER Max iterative solver iterations =SYSMATRIX_SKIPZERODIAG Toggle to ignore zero entries in diagonal matrix =SYSMATRIX_FACTOR_OOC Out-of-core matrix factorization flag =SYSMATRIX_AUTOSPC Toggle to penalize singularities =SYSMATRIX_MONITOR_FREQ Monitoring frequency for factorization =SYSMATRIX_SOLVERTYPE Linear equation solver param Specifies the integer value that type will be set to. =SYSMATRIX_EIGEN_BUCK Buckling eigenproblem =SYSMATRIX_EIGEN_GENERAL Generalized eigenproblem =SYSMATRIX_EIGEN_RIGID Rigid body eigenproblem =SYSMATRIX_EIGEN_STAN Standard eigenproblem =SYSMATRIX_EIGEN_VIBE Vibration eigenproblem =SYS_EIGEN_AMLS AMLS method for vibration =SYS_EIGEN_LANCZOS Lanczos method for vibration =SYS_EIGEN_SUBSPACE Subspace method for vibration =SYS_EIGEN_ALL Get all eigenvalues in range =SYS_EIGEN_LOWEST Get lowest eigenvalues in range =SYS_EIGEN_NEAREST Get nearest eigenvalues =SYS_EIGEN_NORMMASS Mass normalization =SYS_EIGEN_NORMMAX Max component normalization =SYS_ON Enable =SYS_OFF Disable =SYS_SOLVERTYPE_LL Left looking solver =SYS_SOLVERTYPE_MF Multifrontal solver =SYS_SOLVERTYPE_MFP Multifrontal parallel solver =SYS_SOLVERTYPE_PARDISO Intel Pardiso solver =SYS_SOLVERTYPE_MUMPS MUMPS solver
None
Specify the eigenproblem procedure using SYSMATRIX_EIGEN. By default SYSMATRIX_EIGEN is SYSMATRIX_EIGEN_STAN.
Specify the number of eigenvalues sought using SYSMATRIX_EIGEN_NUM. The number of requested eigenvalues must be less than or equal to the number of global dof set using vfs_SysMatrixDef. The default is one.
Set the eigensolution type using SYSMATRIX_EIGEN_TYPE. This parameter determines whether the nearest eigenvalues about a given shift are sought, SYS_EIGEN_NEAREST; the lowest eigenvalues in an interval are sought, SYS_EIGEN_LOWEST; or all the eigenvalues in an interval are sought, SYS_EIGEN_ALL. In any case the eigensolution will terminate if the specified number of eigenvalues set using the SYSMATRIX_EIGEN_NUM parameter are computed. The interval limits are set using vfs_SysMatrixSetParamd. By default SYSMATRIX_EIGEN_TYPE is SYS_EIGEN_LOWEST.
Specify the maximum number of eigenvalue iterations taken at each shift using SYSMATRIX_EIGEN_MAXITER. The default is 50.
Specify the normalization of eigenvectors using SYSMATRIX_EIGEN_NORM. If the parameter value is SYS_EIGEN_NORMMAX, then the eigenvector is normalized to unit value of the largest component. If the parameter value is SYS_EIGEN_NORMMASS, then the eigenvector is normalized to unit value of the generalized mass. By default SYSMATRIX_EIGEN_NORM is SYS_EIGEN_NORMMAX.
Specify the existence of rigid body modes (zero eigenvalues) using SYSMATRIX_EIGEN_ZEROMODES. If it is known that no rigid body modes exist in the eigensolution then it is computationally advantageous to disable this flag. By default SYSMATRIX_EIGEN_ZEROMODES is SYS_ON.
Use SYSMATRIX_EIGEN_OOC to indicate that eigenvectors are to be stored out of core.
The toggle SYSMATRIX_SKIPZERODIAG allows for the factorization and solution of a diagonal matrix where some entries are zero. If set to SYS_ON then zero entries will be skipped and the matrix is assumed non-singular. Defaults to SYS_OFF.
The parameter SYSMATRIX_CACHESIZE sets the central processor cache size in kilobytes. This value should be set to the hardware cache size as closely as possible for optimum performance of factor and solve operations. By default SYSMATRIX_CACHESIZE is 256.
The parameter SYSMATRIX_WORKSIZE_FACTOR sets the block size in megabytes to be used with the out of core solver option. By default SYSMATRIX_WORKSIZE_FACTOR is 256.
Use SYSMATRIX_FACTOR_OOC to indicate that an out-of-core solver is to be used.
The toggle SYSMATRIX_AUTOSPC enables an algorithm that penalizes a detected singularity so that factorization can proceed. This parameter should be used in conjunction with the floating point parameters SYSMATRIX_FACTOR_SINGTOL and SYSMATRIX_AUTOSPC_FACTOR. By default SYSMATRIX_AUTOSPC is set to SYS_OFF.
Use SYSMATRIX_MONITOR_FREQ to specify the frequency to be used to call a monitoring function during a sparse factorization procedure. Defaults to 1000.
Use SYSMATRIX_EIGEN_INITVECTYPE to set the algorithm that will be employed during the initialization phase of a vibration problem. If set to 0, the first eigenvector is set equal to the scaled diagonal of the mass matrix, while all subsequent vectors are randomly initialized with values between -1.0 and 1.0. If set to 1, then the first vector a random value between -0.1 and 0.1 is added to all degrees of freedom for which the diagonal of the mass matrix is larger than a tolerance specified with vfs_SysMatrixSetParamd with the parameter SYSMATRIX_EIGEN_MASSTOL. An array is then constructed with the ratio of the diagonal of the mass and stiffness matrices is computed. If the stiffness matrix entry is smaller than the parameter SYSMATRIX_EIGEN_STIFFTOL then the ratio for the degree of freedom is set to zero. This array is sorted. For the second vector, the degree of freedom associated with the largest ratio is set to 1.0, while all others for which the mass diagonal is larger than SYSMATRIX_EIGEN_MASSTOL are set to a random value between -0.1 and 0.1. For the third vector, the degree of freedom with the second largest ratio is chosen and the algorithm is repeated for all other vectors. Defaults to 0.
Use SYSMATRIX_EIGEN_METHOD to set the eigensolver type used for vibration analysis in vfs_SysMatrixEigen. Method can be SYS_EIGEN_SUBSPACE, SYS_EIGEN_LANCZOS or SYS_EIGEN_AMLS. Defaults to SYS_EIGEN_LANCZOS. The SYS_EIGEN_AMLS method can only be used if the solver type is SYS_SOLVERTYPE_MFP. Note the block Lanczos and AMLS eigensolvers require the MKL BLAS or OpenBLAS.
Use SYSMATRIX_AMLS_RECOVERSIZE to set the eigenvector block size used during the eigenvector recovery computation in the AMLS eigensolver. This should only be changed to a smaller number if there are memory constraints. The default is set to 128.
Use SYSMATRIX_LANCZOS_BLOCKSIZE to select the block size used by the Lanczos eigensolver. If set to 0, a simplified unblocked Lanczos algorithm is used which avoids the need for an external BLAS library. Note that performance is significantly degraded if an external BLAS library is not used. By default, if VKI_LIBAPI_BLASMKL_SEQUENTIAL, VKI_LIBAPI_BLASMKL_THREAD or VKI_LIBAPI_OPENBLAS is defined, then SYSMATRIX_LANCZOS_BLOCKSIZE defaults to 7, else it defaults to 0. An error will be issued if the block size is set to a value greater than zero and some version of the MKL BLAS has not been defined. Generally the block size should be greater than or equal to the largest expected eigenvalue multiplicity. For structural vibration problems with rigid body modes a value of 7 is appropriate.
Use SYSMATRIX_TECH to set the linear equation technology to be used. The technology defaults to SYSMATRIX_TECH_SPARSE for any of the direct methods, to SYSMATRIX_TECH_CONJGRAD for SYSMATRIX_SYMM_ITER, and to SYSMATRIX_TECH_GMRES for SYSMATRIX_NSYMM_ITER. Use SYSMATRIX_TECH_GMRES for matrices of type SYSMATRIX_SYMM_ITER if the matrix is expected to be indefinite.
Use SYSMATRIX_SOLVERTYPE to set the specific direct sparse equation solver package to be used. Set SYS_SOLVERTYPE_LL to use the serial left-looking, in-core only solver, SYS_SOLVERTYPE_MF to use the legacy serial multi-frontal solver and SYS_SOLVERTYPE_MFP to use the parallel multi-frontal solver. Note the parallel multi-frontal solver requires the MKL BLAS. Set SYS_SOLVERTYPE_PARDISO to use Intel's MKL Pardiso solver. Pardiso selection requires that DevTools be compiled with VKI_LIBAPI_PARDISO. Set SYS_SOLVERTYPE_MUMPS to use the NUMPS solver. MUMPS selection requires that DevTools be compiled with VKI_LIBAPI_MUMPS. Defaults to SYS_SOLVERTYPE_MFP.
vfs_SysMatrixSetTemp - set temporary restraints
void vfs_SysMatrixSetTemp (vfs_SysMatrix *sysmatrix, Vint ndofs, Vint lm[])
sysmatrix Pointer to SysMatrix object. ndofs Number of degrees of freedom lm Vector of global dof numbers to restrain
None
vfs_SysMatrixSetNumThreadsAssem - set number of threads for assembly
void vfs_SysMatrixSetNumThreadsAssem (vfs_SysMatrix *sysmatrix, Vint nthreads)
sysmatrix Pointer to SysMatrix object. nthreads Maximum number of threads to be used.
None
vfs_SysMatrixSetValue - set matrix value
void vfs_SysMatrixSetValue (vfs_SysMatrix *sysmatrix, Vint i, Vint j, Vdouble value[])
sysmatrix Pointer to SysMatrix object. i 0-based column number. j 0-based row number. value Value to set matrix entry to.
None
The column and row values, i and j, must come from the matrix structure specified by the user with vfs_DofTabPointersIndicesPtr. SysMatrix will ensure that, if any reordering is performed, the specified value will be assigned to the properly permuted matrix.
vfs_SysMatrixSolve - solve given a system vector
void vfs_SysMatrixSolve (vfs_SysMatrix *sysmatrix, vfs_SysVector *sysvectorx, vfs_SysVector *sysvectory)
sysmatrix Pointer to SysMatrix object. sysvectorx Pointer to SysVector object x sysvectory Pointer to SysVector object y
None
vfs_SysMatrixSolveMult - solve given an array of system vectors
void vfs_SysMatrixSolveMult (vfs_SysMatrix *sysmatrix, Vint nvec, vfs_SysVector *sysvectorx[], vfs_SysVector *sysvectory[])
sysmatrix Pointer to SysMatrix object. nvec Number of system vectors 0 < nvec sysvector Pointer to array of nvec SysVector objects x sysvectory Pointer to array of nvec SysVector objects y
None
vfs_SysMatrixSum - sum a system matrix
void vfs_SysMatrixSum (vfs_SysMatrix *sysmatrix, Vdouble s, vfs_SysMatrix *sysmatrixb)
sysmatrix Pointer to SysMatrix object. sysmatrixb Pointer to SysMatrix object b
None
vfs_SysMatrixThreadAssem,vfs_SysMatrixThreadAssemDiag,vfs_SysMatrixThreadAssemFull - assemble an element matrix
void vfs_SysMatrixThreadAssem (vfs_SysMatrix *sysmatrix, Vint ithread, Vint nedofs, Vint lm[], Vdouble s[]) void vfs_SysMatrixThreadAssemDiag (vfs_SysMatrix *sysmatrix, Vint ithread, Vint nedofs, Vint lm[], Vdouble sd[]) void vfs_SysMatrixThreadAssemFull (vfs_SysMatrix *sysmatrix, Vint ithread, Vint nedofs, Vint lm[], Vdouble sf[])
sysmatrix Pointer to SysMatrix object. ithread Thread number starting from 1. nedofs Number of element degrees of freedom lm Vector of element global dof numbers s Element matrix sd Element diagonal matrix sf Element full matrix
None