1. Introduction

1.1. Module Summary

Each solution results data file or library is viewed as a collection of named datasets which logically group associated information on the library. Each dataset in turn is composed of three parts; parameters, named attributes, and actual data. These components together provide data description and data storage mechanisms for finite element analysis data. There are two broad classes of libraries, native and external.

  • Native libraries, data files which are created by and maintained directly by CEETRON Access (using the Natlib module). These libraries are accessible for read and write operations.

  • External libraries, data files which have been created by an external software system such as a commercial FEA program. These libraries are accessible for read operations only.

External libraries are accessed by a mechanism called direct binary access. This mechanism allows data to be accessed from external libraries as quickly as the FEA programs which created them. Direct binary access, as the term implies, accesses data randomly on external libraries in binary mode with little or no translation processing involved. The modules currently available in CEETRON Access may be divided into 4 categories: 1) library device components, 2) library interfaces 3) data functions and 4) high level library management.

  • Library Device Components

    Library     Generic library information
    Dataset     Generic dataset information
    Attribute   Generic attribute information
    
  • Library Devices

    ABAFil       Read/Write ABAQUS input file
    ABALib       Read/Write ABAQUS .odb output data base
    AdamsLib     Read MSC/Adams
    ANSLib       Read ANSYS .rst, .rth and .rfl results files
    ANSFil       Read/Write ANSYS input file
    AUTODYNLib   Read AUTODYN results file
    CFXLib       Read ANSYS/CFX result file
    CGNSVLib     Read/Write CGNS data base
    COMSOLLib    Read COMSOL mph and section files
    D3DFil       Read/Write LSTC LS-DYNA input file
    D3DLib       Read LSTC State and time history databases
            DPFLib       Ansys DPF framework for accessing CAE files
    EnSightLib   Read/Write ANSYS/EnSight Casefile format
    FEMAPLib     Read FEMAP Neutral file
    FDILib       Read FIDAP Neutral file
    FLUENTLib    Read/Write FLUENT case and data file format
    GMVLib       Read GMV file
    HMAFil       Read Altair HyperMesh ASCII file
    H3DLib       Read Altair H3D file
    MarcLib      Read MSC/MARC .t16 and .t19 post files
    MemLib       Read and write native memory objects
    NASFil       Read/Write MSC/NASTRAN bulk data input file
    NASLib       Read/Write MSC/NASTRAN OUTPUT2, XDB and HDF5 output file
    NatLib       Read and write native libraries
    OBJFil       Read/Write OBJ file
    OpenFOAMLib  Read OpenFOAM databases, Write mesh entities
    PAMLib       Read ESI/PAM-CRASH DAISY and ERF files
    PatLib       Read/Write Patran neutral files
    PERMASLib    Read INTES/PERMAS result file
    PLOT3DLib    Read NASA/PLOT3D grid, solution and function files
    PluginLib    Read plugin-compatible files
    POLYFLOWLib  Read POLYFLOW mesh and solution files
    RASFil       Write PTC/Mechanica FEM Neutral Format
    RASLib       Read PTC/Mechanica design study
    SAMCEFLib    Read Siemens SAMCEF des files
    SDRCLib      Read/Write SDRC/IDEAS universal file
    STARCCMLib   Read STAR-CCM results file
    STLFil       Read/Write 3D Systems/STL ASCII text and binary files
    TecplotLib   Read/Write Tecplot file format
    VTKLib       Read Vtk legacy file format
    
  • Data Function

    DataFun     Define direct binary access function pointers
    DataIPC     Implement data functions with interprocesss communication
    
  • Library Management

    LMan        Implement overall library management functions.
    

Each instance of a module is termed an object. In order to access the finite element model and results data from an external library a series of steps must be followed to open the data file, declare its format, query for the amount and type of information on the file and finally request the information to be read into memory. Using the SDRC Universal file format as an example, the methodology for reading the base model and results information can be described by the following steps:

  • For each external library type (SDRC Universal file, NASTRAN Ouput 2 file, etc.) a library device interface object is instantiated. In this case a SDRCLib object is created. The SDRCLib module is just one of many library device interfaces available in CEETRON Access. In order to provide a single programming interface for opening, querying and reading from an external library regardless of its specific type, the data function module, DataFun, is used. This module contains a set of function pointers which define the abstract interface between the user and any CEETRON Access library device. Each library device interface module, such as SDRCLib, has a function which fills a DataFun object with the function pointers particular to that library device interface module. Once a DataFun object has been configured in this way all further data access of the library device is handled through the data function object. This procedure is illustrated below:

    vdm_SDRCLib *sdrclib;
    vdm_DataFun *df;
                /* create SDRC Universal File device object */
    sdrclib = vdm_SDRCLibBegin();
                /* create data function */
    df = vdm_DataFunBegin();
                /* load data functions for SDRCLib device */
    vdm_SDRCLibDataFun (sdrclib,df);
    
  • To open a library device use vdm_DataFunOpen(). Note that the path name passed to vdm_DataFunOpen() is assumed to be UTF-8 encoded. There are exceptions for some interfaces on Windows platforms. Opening a library device is often the most time consuming operation performed on a library. In general the library file must be completely scanned to determine the number and type of datasets, number of nodes and elements, etc. Any data read into memory in order to scan the file is not retained at this time. The opened file must contain the basic node and element model information.

    Vint ierr;
                /* open library on file "example.unv" */
    vdm_DataFunOpen (df,0,"example.unv",VDM_SDRC_UNIVERSAL);
    ierr = vdm_DataFunError(df);
                /* error opening file */
    if(ierr) {
    }
    
  • It is possible, for some library devices, to append the results of one or more additional files. The appended files do not need to be of the same type as the initial opened file. Any node and element model related datasets on the appended file are not generated. Only results datasets are generated. Not all opened library devices support appended files and only certain file types are appendable. The return error code should be tested.

                /* append library on file "append.op2" */
    vdm_DataFunAppend (df,"append.op2",VDM_NASTRAN_OUTPUT2);
    ierr = vdm_DataFunError(df);
                /* error appending file */
    if(ierr) {
    }
    
  • The library device is now ready for access. CEETRON Access allows the user to query for certain global library device parameters. The number of entities such as finite elements and nodes may be retrieved using vdm_DataFunGetNumEntities(). In order to begin to query for the existence of datasets, the internal Library object used by the library device must be retrieved. Use vdm_DataFunGetLibrary() to get the Library object. The finite element model information consists of many different datasets for node coordinates, element topology, element connectivity, etc. Each model dataset is given a name on the library. Each dataset on a library is assigned a dataset index ranging from zero to the number of datasets on the library minus one. The procedure for querying for the dataset index, idst, of the node coordinates, dataset name X.N, is illustrated below. The number of coordinate datasets found, ndsts, should equal 1.

    Vint numnp, numel;
    Vint idst, ndsts;
    
                /* get number of nodes and element */
    vdm_DataFunGetNumEntities (datafun,SYS_NODE,&numnp);
    vdm_DataFunGetNumEntities (datafun,SYS_ELEM,&numel);
                /* get library object */
    vdm_DataFunGetLibrary (datafun,&library);
                /* search for node coordinate dataset */
    vdm_LibrarySearchDataset (library,"X.N",1,&idst,&ndsts);
    
  • The model datasets may be read individually to reconstruct the node and element model information. However, it is strongly recommended to use the LMan module and the vdm_LManLoadModel() function to read all model associated datasets and populate a VisTools Model object. Similarly, it is recommended that all results datasets be read using the vdm_LManLoadState() to read the datasets and populate a VisTools State object.

  • In order to ultimately read the contents of the node coordinate dataset into memory, the user must query for the dataset parameters. The dataset parameters include dataset name, length, number of rows, number of columns and data type. This process illustrates the use of the component modules Library and Dataset. The dataset parameters allow the user to determine the amount of memory required to read the contents of the dataset.

    vdm_Dataset *dataset;
    Vchar name[33];
    Vlong lrec;
    Vint nrow, ncol, ntyp;
    
                /* get dataset object of idst dataset */
    vdm_LibraryGetDataset (library,idst,&dataset);
                /* inquire dataset parameters */
    vdm_DatasetInq (dataset,name,&lrec,&nrow,&ncol,&ntyp);
    
  • Now the data may be read into memory. Memory must be allocated to hold the data. Use the dataset length, lrec, and data type, ntyp, to determine how may bytes to allocate. The returned value of ntyp should equal 2 which indicates single precision floating point data. The data is read using vdm_DataFunReadDataset(). The node coordinates are now in memory in a simple rectangular array format.

    Vfloat *x;
    
                /* allocate memory */
    x = (Vfloat*)malloc(lrec*sizeof(Vfloat));
                /* read data */
    vdm_DataFunReadDataset (df,idst,(void*)x);
    
  • To close a library device and terminate CEETRON Access, use vdm_DataFunClose() to close the library. Free any user memory allocated and delete all objects.

                /* free user memory */
    free(x);
                /* close library device */
    vdm_DataFunClose (df);
                /* delete objects */
    vdm_SDRCLibEnd (sdrclib);
    vdm_DataFunEnd (df);
    
  • Other finite element model datasets may be accessed in a similar manner. The procedure for reading finite element results data is also similar except that there are many more dataset names for results datasets and they are parameterized in general by load case, time step, etc. There are two basic techniques for finding results datasets. One is to search for specific results datasets using the function vdm_LibrarySearchDataset() with wild card characters in the dataset names. The second is to traverse each dataset in the library one by one using vdm_LibraryGetDataset() and examining dataset names for results specific parameters using vdm_DatasetResult() and related functions.

1.2. Interface Support Information

The following section is meant to summarize the following information for each interface module: 1) the vendor of the commercial analysis system producing the result files and 2) the release levels of the analysis systems supported by the interface module.

1.3. Compiling and Linking a CEETRON Access Application

In addition to the overall C preprocessor definitions described in section DevTools, Compiling and Linking a DevTools Application, another C preprocessor definition specifies external library interface APIs. The current external library APIs used in CEETRON Access are to support the ABAQUS .odb, PAM-CRASH .DSY, PAM-CRASH .erfh5, MSC/NASTRAN .h5, STAR-CCM results file, CGNS data base, AUTODYN results file, Altair H3D file, CFX results file, SAMCEF results file interfaces and any HDF5 based interface.

VKI_LIBAPI_ABAODB      ABAQUS .odb interface
VKI_LIBAPI_AUTODYN     AUTODYN interface
VKI_LIBAPI_CFX         CFX results file
VKI_LIBAPI_CGNS        CGNS data base
VKI_LIBAPI_FEMZIP      FEMZIP LS-DYNA d3plot interface
VKI_LIBAPI_FEMZIP_     FEMZIP LS-DYNA d3plot interface appended
VKI_LIBAPI_HDF5        PAM-CRASH ERF, MSC/NASTRAN .h5, Native HDF5
VKI_LIBAPI_H3D         Altair H3D interface
VKI_LIBAPI_PAMDSY      PAM-CRASH .DSY interface
VKI_LIBAPI_STARCCM     STAR-CCM .ccm results file
VKI_LIBAPI_SAMCEF      SAMCEF .des results file
VKI_LIBAPI_ZLIB        ANSYS 17 compressed results files

1.3.1. ABAQUS .odb interface

Since the .odb API is implemented in C++, the .odb interface requires the use of a C++ compiler and an installation of ABAQUS. The compiler required on a specific system must match the compiler used by ABAQUS. See the following URL.

http://www.simulia.com

The file vdm/vdmodb.cxx must be compiled with the required C++ compiler. All other CEETRON Access source code may be compiled with a compatible C compiler. Compiler options and linker options and required libraries are dependent upon the specific version of ABAQUS and computer platform.

1.3.2. AUTODYN interface

The AUTODYN interface is implemented using the AUTODYN toolkit provided by ANSYS AUTODYN. Calls to the AUTODYN functions are conditionally compiled by defining VKI_LIBAPI_AUTODYN, then the application must be linked with the AUTODYN libraries. See the following URL.

http://www.ansys.com

1.3.3. FEMZIP LS-DYNA d3plot interface

FEMZIP is used to compress a family of LS-DYNA d3plot files. The FEMZIP API is implemented using a toolkit provided by Fraunhofer SCAI. Apparently FEMZIP is written in Fortran so that the function names may be appended by an underscore on some platforms. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_FEMZIP_ or VKI_LIBAPI_FEMZIP to access functions with or without an appended underscore respectively. See the following URL.

http://www.sidact.com/femzip.html

1.3.4. CFX results file

The CFX results file API is implemented using a toolkit provided by ANSYS CFX. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_CFX, then the application must be linked with the CFX libraries. See the following URL.

http://www.ansys.com

1.3.5. CGNS data base

The CFD General Notation System (CGNS) data base is a collection of conventions and software for the storage and retrieval of CFD data. Calls to the CGNS functions are conditionally compiled by defining VKI_LIBAPI_CGNS. If CGNS is used by compiled with VKI_LIBAPI_CGNS then the application must be linked with the CGNS libraries. The source code for the underlying CGNS library are publicly available. See the following URL.

http://cgns.sourceforge.net

1.3.6. Altair H3D interface

The H3D interface is implemented using the H3d toolkit provided by Altair Engineering. Calls to the H3D functions are conditionally compiled by defining VKI_LIBAPI_H3D then the application must be linked with the H3D libraries. See the following URL.

http://www.altair.com

1.3.7. PAM-CRASH .DSY interface

The .DSY API is implemented using a toolkit provided by ESI-GROUP. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_PAMDSY. See the following URL,

http://www.esi-group.com

1.3.8. PAM-CRASH .erfh5 interface, MSC/NASTRAN .h5, Native .vh5

The above file formats are implemented on top of the HDF5 interface. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_HDF5. See the following URL,

http://www.hdfgroup.org

1.3.9. STAR-CCM .ccm results file interface

The STAR-CCM API is implemented using a toolkit provided by CD-Adapco. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_STARCCM. See the following URL,

http://www.cd-adapco.com

1.3.10. SAMCEF .des results file interface

The SAMCEF API is implemented using a toolkit provided by Siemens. Calls to the toolkit functions are conditionally compiled by defining VKI_LIBAPI_SAMCEF. See the following URL,

http://www.siemens.com

The SAMCEF interface requires the definition of two environment variables: SAM_ZONE refers to the memory usage in the SAMCEF library. SAM_EXE identifies the SAMCEF library installation path.

1.4. Dataset Specification Syntax

Datasets are identified by name. A dataset specification has the following form,

dataset_name:id1:id2:id3

1.4.1. Dataset Name, dataset_name

The dataset name is a character string which may contain any letter or digit including the special characters period, dollar sign and underscore. The period is most often used in dataset names to separate the name into two or more fields. By convention the last field, delimited by a period, indicates the data structure. The length of a dataset specification name may not exceed 256 characters. Below are examples of some legal dataset specifications.

V.N:1
S.E:10
D.N:1:5
TEMP.[Total].E
SET.ELEM.T

1.4.2. Numeric Identifiers, :id1:id2:id3

The numeric identifiers are positive integers. They are delimited from the dataset name and other numeric identifiers by the colon, :. The numeric identifiers are optional in the dataset specification, however one or more of the numeric identifiers is usually present. Their exact interpretation depends upon the context in which they are used. They are used in many places to signify the load step or time step number, increment number, vibration, or buckling mode number, etc. For example, on an external ABAQUS .fil database the node displacements at load step 3, increment 2 are resident in dataset D.N:3:2

1.4.3. Dataset Specification Expansion

It is useful in many instances to be able to refer to a group of related datasets with a single dataset specification or a single dataset with an abbreviated but unique character string. These tasks are accomplished with dataset name expansion. The dataset name and numerical identifiers may be expanded with a syntax which resembles UNIX file name expansion. The following wild card substitutions are used.

  • The asterisk, *, matches any string of characters including a null character string.

  • The question mark, ?, matches any single character.

  • Parentheses, ( and )., enclose a set of matching characters. If the first character in the set is a caret, it matches a character NOT in the set. The shorthand S-E is used to specify a set of characters S upto E, inclusive.

For example,

D*

matches all datasets starting with the letter D. The dataset specification,

\*:(2-4)

matches all datasets associated with solution step 2, 3 and 4.

In a departure from standard UNIX file expansion syntax, dataset specifications may also contain numeric identifier ranges. Ranges are specified in the form, FiTjBk, representing From the ith To the jth By increments of k. For example,

D.N:F1T49B2

matches all datasets named D.N:1, D.N:3, D.N:5, …, D.N:49. The highest and lowest numeric identifier may be referenced using the letters H and L respectively. This is useful in quickly selecting the first or last solution on a library. For example,

S.EL:H

matches the dataset S.EL.n, where n is the highest numbered solution step. The asterisk wild card may be combined with identifier ranges. For example,

D.N:H:*

matches displacement datasets at the last iteration of each load step.

1.5. Dataset Parameters and Attributes

Associated with each dataset on a library are a set of dataset parameters and attributes. The dataset parameters describe the basic physical characteristics, data structure and data type of the dataset. Dataset attributes are given a name and a value. They are used to help describe the contents of a dataset or refine its data structure.

Dataset parameters include record length, number of rows and columns, data type and number of dataset attributes. The total size of the dataset in (bytes) is the product of record length and the size of the data type. The data type of a dataset should be homogeneous. The basic data structure of a dataset is that of a dimensioned array. The dimensionality of the array is specified by the number of rows and columns. Variations on the basic dimensioned array data type are accomplished using special reserved dataset attributes. These data structures are described in more detail later in this chapter. The description of each dataset parameter is as follows.

  • lrec, record length in units of size of data type

  • nrow, number of rows

  • ncol, number of columns

  • ntyp, data type

    =1  integer (Vint) (SYS_INTEGER)
    =2  real, single precision (Vfloat) (SYS_FLOAT)
    =3  Hollerith (Vuint) (SYS_HOLLERITH)
    =4  real, double precision (Vdouble) (SYS_DOUBLE)
    =6  complex, single precision (Vfloat[2]) (SYS_COMPLEX)
    =7  complex, double precision (Vdouble[2]) (SYS_DOUBLECOMPLEX)
    
  • natt, number of dataset attributes

Dataset attributes are meant to contain information about a dataset. An attribute has an associated name, length, data type and value. An individual attribute is identified by specifying a dataset name and attribute name. Attribute names are character strings containing 16 characters or less. The period is often used to delimit an attribute name into two or more fields. For example, attributes which have other dataset names as values are named as follows,

Link.Size
Link.Index
Link.Complex
Link.Cid
Link.Section
Link.Layers

The data type and length of an attribute describe the data type and length of an attribute value. An attribute value cannot contain more than 16 integer or floating point numbers, 8 double precision numbers or 256 characters. An attribute value must be of a single data type.

1.6. Data Structures

The basic data structures supported by the data management system are the single data type rectangular array and single data type variable row array. Dataset attributes are used to logically link datasets to other datasets to augment these basic data structures.

  • Rectangular array. This basic data structure is described by the dataset parameters. The dimensions of the array are nrow by ncol. Data is physically stored in the array in row major order. The data type is given by ntyp. The dataset may be optionally stored in sparse format in which the actual column indices appearing in the dataset are contained in an associated dataset. Columns are numbered starting from 1. If a column is not specified it means that the associated data is undefined. The associated dataset is specified by the value of the following dataset attribute. If this dataset attribute is not present, the dataset is not stored in sparse format.

    Link.Index     integer vector of column indices.
    

    For example, the sparse array, nrow=2, ncol=4, would be represented as follows:

    Non sparse format     Sparse format    Indices
    --         --         --         --    -- --
    | 0.  | 0.  |         | 2.1 | 1.5 |    | 2 |
    | 2.1 | 1.5 |         | 4.0 | 6.0 |    | 4 |
    | 0.  | 0.  |         --         --    -- --
    | 4.0 | 6.0 |
    --         --
    

    A variation of the rectangular array is used to specifically store element face and element edge data. In this case the dataset is always stored in sparse format and the Link.Index points to an integer 2-row vector of element indices and associated element face or edge number.

    Link.Index     integer 2-row vector of element indices and entity number.
    
  • Variable row array. The variable row array is an array with variable number of rows per column. If the rows contain a data quantity such as stress tensors, the number of components in the data quantity is given by nrow. If the rows contain element connectivity indices, the maximum number of rows in a column is given by nrow. The number of columns is given by ncol. The dataset attribute, Structure must be set to either VariableRow, ElementNode, or ElementIP. If the structure is VariableRow then the length of each column in the dataset is contained in an associated dataset.

    Link.Size      integer vector of column lengths, length ncol
    

    If the structure is ElementNode then the dataset specifically contains data at element node locations. The length of each column in the dataset is computed using the number of dataset rows, nrow times the number of nodes in each element in dataset ELEM.SIZE.E times the number of element sections in the optional contents of the Link.Section dataset attribute. If the structure is ElementIP then the dataset specifically contains data at element integration point locations. The length of each column in the dataset is computed using the number of dataset rows, nrow times the number of integration points in each element in dataset ELEM.EIPS.E times the number of element sections in the optional contents of the Link.Section dataset attribute.

    Link.Section   integer vector of section lengths.
    

    If the structure is ElementFaceNode or ElementEdgeNode then the dataset specifically contains data at element face node or element edge node locations. The length of each column in the dataset is computed algoritmically using the number of dataset rows, nrow times the number of element face or edge nodes. In this case the dataset is always stored in sparse format and the Link.Index points to an integer 2-row vector of element indices and associated element face or edge number.

    Link.Index     integer 2-row vector of element indices and entity number.
    

    If the structure is set to LowerTriangle or Diagonal, then the data is associated with a set of generalized degrees of freedom. The number of columns is the number of degrees of freedom in the set. A lower triangular matrix is stored by rows.

CEETRON Access uses a dataset naming convention to indicate the data structure and data class of each dataset. These conventions are referred to as formats. The following formats are used.

  • Node array format, .N, node array datasets are rectangular arrays of any data type. The data type is given by ntyp. The row dimension is the number of items or components associated with each node. A row dimension of 1 implies a scalar value is stored in the dataset for each node. A row dimension of 3 implies a vector is stored at each node, etc. The column dimension is equal to the number of nodes in the finite element model. Node array datasets contain node coordinates, displacement vectors, node normals, etc.

  • Element array format, .E, element array datasets are analogous to node array datasets except that element data is stored rather than node data. Element array datasets are used to store element temperatures, element stresses, element shape information, etc.

  • Element Node array format, .EL, element node arrays contain variable column length data. Element node arrays are used to store element connectivity, stress results, etc.

  • Element Integration Point array format, .EIP, element integration point arrays contain variable column length data. Element integration point arrays contain the same type of information as element node arrays except that the data locations are at the element integration points rather than the element nodes.

  • Element Face, Edge, Face Node and Edge Node format, .EF, .EE, .EFN and .EEN contain data and element faces, element edges, element face nodes and element edge nodes respectively.

  • Degree of freedom format, .D, datasets contain data associated with system wide degrees of freedom. The dataset attribute, Structure may be set to LowerTriangle or Diagonal. A lower triangular matrix is stored by rows. This structure is used to store reduced symmetric stiffness matrices, etc.

  • Mode format, .MOD, mode format datasets are rectangular arrays of any data type with data associated with a mode. A mode is normally the result of an eigensolution. The types of results at modes are typically frequency, generalized mass, etc.

  • Particle format, .PCL, particle format datasets are rectangular arrays of any data type with data located at so-called particles. Particles are point objects and may be either physical objects or arbitrary locations in space used for visualization, etc. The column dimension is equal to the number of particles. Particle datasets are never indexed.

  • Contact Pair format, .CPR, contact pair format datasets are rectangular arrays of any data type with data associated with a contact pair.

  • Panel format, .PNL, panel format datasets are rectangular arrays of any data type with data associated with a panel.

  • Table format, .T, table format datasets are rectangular arrays of miscellaneous objects such as coordinate system definitions, material properties, etc. The row dimension is the number of items stored per object. The column dimension is the number of objects.

  • Table variable row array format, .TL, variable row arrays contain variable column length data for generic entities.

  • Library dataset, .L, contains a sublibrary which in turn contains other datasets and sublibraries. These datasets are analogous to a LINUX directory. Use vdm_DataFunLibDataset() to make this library current. Library datasets are used, for example, to contain models and associated results from adaptive meshing or superelement recovery. The following naming conventions are used.

    REMESH.L:ith   adaptive mesh model, ith is sequential remeshing number
    SEMESH.L:seid  superelement recovery mesh model, seid is superelement id
    GEOMESH.L      tesselated geometry model
    

Complex results quantities may be stored as real-imaginary or magnitude-phase pairs of datasets. The real or magnitude part of the complex result will be stored in a dataset which is named as if the result were purely real. The imaginary or phase part of the complex result is stored in a dataset which adds a .I immediately after the dataset name (including optional dataset qualifier) for the real or magnitude part. A Complex attribute is added with the character value Real, Magnitude, Imaginary or Phase. Note that phase is always in degrees. The real or magnitude dataset will have an additional Link.Complex attribute which names the associated dataset containing the imaginary or phase part of the complex result. If only the real or magnitude part exists, then the Link.Complex Attribute will be omitted. For example, a complex displacement result containing both magnitude and phase parts would result in the following datasets and attributes.

D.N:1:1
 Complex = Magnitude
 Link.Complex = D.I.N:1:1

 D.I.N:1:1
 Complex = Phase

1.7. Finite Element Model Dataset Names and Descriptions

Datasets appearing on native and external libraries follow a naming convention and adhere to one of the basic data structures detailed in the preceding section. Each model dataset will contain a “Model” attribute. The character value of the attribute will give the class of model data that the dataset contains.

  • Possible Model attributes

    Mesh
    Element Data Arrays
    Node Sets
    Element Sets
    Element Entity Sets
    Load Cases
    Restraint Cases
    Multipoint Constraint Cases
    Initial Condition Cases
    Global Properties
    Material Properties
    Element Properties
    Solution Properties
    Function Properties
    Units
    Coordinate Systems
    Analytical Surfaces
    Contact Pairs

Each model dataset also contains a “Description” attribute whose character value is a basic description of the dataset contents. A summary of the “Mesh” datasets is given below.

1.7.1. Model Datasets

COLORID.E COLORID.N

Integer vectors of object colors for each object type. ncol=maximum possible number of objects of the given type, nrow=1.

COLORMAP.T

Real array of RGB floating point triples for each color index. nrow=3, ncol=maximum number of color map entries defined. Column j contains the red, green and blue intensities in the range [0.,1.] of color map identifier j.

CSYS.T

Real array of the origin and orientation of defined coordinate system objects. nrow=12, ncol=maximum number of coordinate systems defined. Column j contains the definition of coordinate system internal identifier j.

row_1-3

Global coordinates of origin of coordinate system

row_4-12

Direction cosines of coordinate system axes

CSYS.TYPE.T

Integer vector of types for each coordinate system. ncol=maximum number of coordinate systems, nrow=1.

row_1        Coordinate system type, type
            =1          Rectangular        (SYS_CARTESIAN)
            =2          Cylindrical        (SYS_CYLINDRICAL)
            =3          Spherical, type a  (SYS_SPHERICAL)
            =4          Spherical, type b  (SYS_SPHERICAL_ALT)
            =5          Toroidal           (SYS_TORODIAL)

CSYS.ID.T

Integer vector of user assigned identifiers for each coordinate system. ncol=maximum number of coordinate systems, nrow=1.

DOF.CID.N

Integer vector of coordinate system identifiers in which nodal degrees of freedom are expressed for each node. ncol=maximum number of nodes, nrow=1.

EID.E

Integer vector of user assigned identifiers for each element. ncol=maximum number of elements, nrow=1.

ELEM.NODE.EL

Integer variable row array of element node connectivity for each element. ncol=maximum number of elements, nrow=maximum length of element node connectivity. Connectivity references internal node number ordering. The datasets .ELEM.NODE.PNTR.E and .ELEM.NODE.SIZE.E are the integer pointer and size of each column in dataset ELEM.NODE.EL.

ELEM.SHAP.E

Integer vector of element shape of each element. ncol=maximum number of elements, nrow=1. See section VisTools, Computational Cells for a description of element topologies.

shape        Element shape
            =SYS_SHAPEPOINT
            =SYS_SHAPELINE
            =SYS_SHAPETRI
            =SYS_SHAPEQUAD
            =SYS_SHAPETET
            =SYS_SHAPEPYR
            =SYS_SHAPEWED
            =SYS_SHAPEHEX
            =SYS_SHAPEPOLYGON
            =SYS_SHAPEPOLYHED

ELEM.MIJK.E

Integer array of element order for each element. ncol=maximum number of elements, nrow=3. Column j contains 3 form parameters for element j.

maxi         Number of nodes in I direction
maxj         Number of nodes in J direction
maxk         Number of nodes in K direction

ELEM.SPEC.E

Integer vector of element specific type of each element. The specific type is designed to distinguish variations within the general element types contained in dataset ELEM.TYPE.E. See section VisTools, Element Types for a description of element specific types. ncol=maximum number of elements, nrow=1.

row_1        Element specific type code

ELEM.TYPE.E

Integer vector of internal element type of each element. Use dataset ELEM.SPEC.E to determine element variations within a particular element type. See section VisTools, Element Types for a description of element types. ncol=maximum number of elements, nrow=1.

row_1        Element type code
            =SYS_ELEM_SOLID
            =SYS_ELEM_SHELL
            =SYS_ELEM_MEMBRANE
            =SYS_ELEM_BEAM
            =SYS_ELEM_TRUSS
            =SYS_ELEM_INFINITE
            =SYS_ELEM_GAP
            =SYS_ELEM_JOINT
            =SYS_ELEM_SPRINGDASHPOT
            =SYS_ELEM_RIGID
            =SYS_ELEM_CONSTRAINT
            =SYS_ELEM_PLOT
            =SYS_ELEM_MASS
            =SYS_ELEM_INTER
            =SYS_ELEM_SUPER
            =SYS_ELEM_REINFORCEMENT

ELEM.TYPE.EXT.E

Integer vector of external element type of each element. The external element type number is the identifying element type number which is assigned by the originating, external finite element system. For example, for ANSYS, this number is the ANSYS element routine number. ncol=maximum number of elements, nrow=1.

ELEM.TYPE.HOL.E

Character valued external element name for each element. Each element name consists of 8 characters. Character strings are stored in Hollerith form, 4 characters per integer word. ncol=maximum number of elements, nrow=2.

row_1        First four characters
row_2        Second four characters

ERS.CID.E

Integer vector of element local coordinate system identifiers. This dataset may be referenced by element result datasets to indicate the coordinate system in which element vector and tensor results are expressed. ncol=maximum number of elements, nrow=1. User defined coordinate systems are located in the CSYS.T model dataset. Vector or position orientations are located in the ERS.VEC.E model dataset. Direction cosine matrices are located in the DIR_COS result dataset. See section VisTools, Element Coordinate Systems for a description of element coordinate systems.

row_1       Element local coordinate system identifier
            >0                           User defined coordinate system
            =SYS_ELEMSYS_GLOBAL          Global coordinate system
            =SYS_ELEMSYS_STANDARD        Standard
            =SYS_ELEMSYS_POSITION        Position
            =SYS_ELEMSYS_GLOBALPROJECT   Global project
            =SYS_ELEMSYS_VECTOR          Vector
            =SYS_ELEMSYS_VECTORELEMNODE  Element node vector
            =SYS_ELEMSYS_BISECTOR        Bisector
            =SYS_ELEMSYS_ROTANG          Rotation angles
            =SYS_ELEMSYS_ROTANGELEMNODE  Element node rotation angles
            =SYS_ELEMSYS_FIRSTEDGE       First edge
            =SYS_ELEMSYS_CYLINDRICAL_ALT Alternate Cylindrical

ERS.VEC.E

Real array of global coordinates of element orientation positions or vectors. ncol=maximum number of elements, nrow=3. Each column j contains the x,y,z global coordinates of a orientation position or vector for the jth element. These orientations are required for certain types of element local systems. Element local coordinate system identifiers are located in the ERS.CID.E model dataset.

INDX.ELEM.T:id
INDX.NODE.T:id
INDX.MODE.T:id
INDX.PANEL.T:id

Integer vector of element, node, mode or panel indices. ncol=number of element, node, mode or panel indices, nrow=1.

INDX.ELEM.FACE.T:id
INDX.ELEM.EDGE.T:id

Integer 2 row vector of element indices and associated element face or edge number. ncol=number of element indices, nrow=2.

MID.E

Integer vector of user assigned material property identifiers for each element. ncol=maximum number of elements, nrow=1.

NID.N

Integer vector of user assigned identifiers for each node. ncol=maximum number of nodes, nrow=1.

NORMAL.N

Real array of the components of normal vectors at nodes. ncol=maximum number of nodes, nrow=3. Each column j contains the components nx,ny,nz of the normal vector at the jth node expressed in the global coordinate system.

PARAMETER.HOL.T

Character valued active model parameters. Character strings are stored in Hollerith form, 4 characters per word.

col_1-20

Title

col_21-24

Originating code name

col_25-28

Originating code version

col_29-32

Origination time

col_33-36

Origination date

PARAMETER.IJK.T

Integer array of grid dimensions for a model containing one or more structured grids. ncol=number of grids in structured grid model, nrow=3. Column j contains the number of vertices in each of the I,J,K directions.

row_1        Number of vertices in I direction
row_2        Number of vertices in J direction
row_3        Number of vertices in K direction

PARAMETER.INT.T

Integer valued active model parameters.

col_1        Number of nodes
col_2        Number of elements
col_3        Dimensionality of model
col_4        Number of degrees of freedom per node
col_5        Number of material properties
col_6        Number of physical properties
col_7        Total number of active degrees of freedom
col_8        Total number of suppressed degrees of freedom
col_9        Analysis type
            =0          Unknown
            =1          Structural   (SYS_ANALYSIS_STRUCTURAL)
            =2          Thermal      (SYS_ANALYSIS_THERMAL)
            =3          Electric     (SYS_ANALYSIS_ELECTRIC)
            =4          Magnetic     (SYS_ANALYSIS_MAGNETIC)
            =5          Fluid        (SYS_ANALYSIS_FLUID)
            =6          Acoustic     (SYS_ANALYSIS_ACOUSTIC)
            =7          Diffusion    (SYS_ANALYSIS_DIFFUSION)
col_10       Solution type
            =0          Unknown
            =1          Static              (SYS_SOL_STATIC)
            =2          Vibration           (SYS_SOL_VIBRATION)
            =3          Buckling            (SYS_SOL_BUCKLING)
            =4          Transient           (SYS_SOL_TRANSIENT)
            =5          Superelement        (SYS_SOL_SUPERELEMENT)
            =6          Frequency response  (SYS_SOL_FREQRESPONSE)
            =7          Complex eigenvalue  (SYS_SOL_COMPLEXEIGEN)
col_11       2D type
            =0          Unknown
            =1          Plane strain          (SYS_PLANESTRAIN)
            =2          Plane stress          (SYS_PLANESTRESS)
            =3          Axisymmetric          (SYS_AXISYMMETRIC)
            =4          Axisymmetric Fourier  (SYS_AXISYMFOURIER)
col_12       Nonlinear flag
            =0          Linear
            =1          Nonlinear
col_13       Number of grids in structured grid model
col_14       Library type

PARTID.E

Integer vector of user assigned part identifiers for each element. Part identifiers are determined to be the most useful integer partitioning quantity for a finite element model. ncol=maximum number of elements, nrow=1.

PARTID.N

Integer vector of user assigned part identifiers for each node. ncol=maximum number of nodes, nrow=1.

PID.E

Integer vector of user assigned physical property identifiers for each element. ncol=maximum number of elements, nrow=1.

SET.ELEM.T:id

Integer vector of element indices assigned to element set identifier id. ncol=number of elements, nrow=1.

SET.ELEM.EDGE.T:id

Integer vector of element indices and edge numbers assigned to element edge set identifier id. ncol=number of elements, nrow=2.

row_1        Element index
row_2        Edge number

SET.ELEM.FACE.T:id

Integer vector of element indices and face numbers assigned to element face set identifier id. ncol=number of elements, nrow=2.

row_1        Element index
row_2        Face number

SET.NODE.T:id

Integer vector of node indices assigned to node set identifier id. ncol=number of nodes, nrow=1.

X.N

Real array of global coordinates of nodes. ncol=maximum number of nodes, nrow=3. Each column j contains the x,y,z global coordinates of the jth node.

1.8. Finite Element Result Dataset Names and Descriptions

Results data may be resident on native and/or external libraries. If results data is resident on a native library, the native library interface uses binary direct access I/O to read data directly to memory resident data structures without transformation. If the results data resides on an external library, the data is meant to be able to be accessed logically in the same manner as native libraries. The external library interface contains several functions which have been designed for this purpose. The basic concept is to establish a dataset naming convention which maps dataset names to finite element results data residing on the external library. The physical layer of the external library interface must be written to access data from the external library in response to dataset access requests and perform any transformations necessary to reformat the data into the specified CEETRON Access data structure.

CEETRON Access supports finite element results data located at nodes, element centroids and element nodes. The results may be represented as conventional scalar, vector or symmetric tensor quantities. The special cases of stress resultants and strains in shell and beam elements are handled as well.

The data management system is designed to access, with any single request for results data, a result at all nodes or elements in the finite element model. For example, a single request for nodal temperatures expects a column of scalars representing the temperature at each node in the finite element model to be returned. Dataset names are mapped by convention to specific finite element results quantities which reside on the library. Dataset attributes are mapped to appropriate global variables associated with each dataset. An example of a global variable is the frequency of a vibration mode. As discussed earlier, CEETRON Access dataset names have the following syntax,

dataset_name:id1:id2:id3

The dataset name is an alphanumeric string consisting of two or more fields delimited by periods. The first field, termed the dataset root name, is chosen to reflect the contents of the datasets, such as TEMP for temperature and, SE_DENSITY for strain energy density. The last field must be one of several reserved strings which indicates the format of the dataset. Formats applicable to results quantities are .N (node array), .E (element array), and .EL (element node array). The dataset numeric identifiers may be used to indicate any appropriate identifiers such as load case, time step, load step, iteration, vibration or buckling mode, etc.

Very often there may be several variations of a specific result type such as temperature. For example there may be a quantity described as “Temperature Total”. In this case the specific result type is temperature and “Total” is considered to be a qualifier. Qualifier strings are appended to the dataset root name enclosed in square brackets. The “Temperature Total” at elements for step one would then result in the following dataset name.

TEMP.[TOTAL].E:1

If the dataset contains the Complex attribute value of Imaginary or Phase, then the field, I, is appended to the dataset root name or the optional dataset qualifier. This is done in order to distinguish the dataset name from the dataset name containing the Complex attribute value of Real or Magnitude.

In some cases, a result quantity will appear on an external library for which there is no current predefined result type. If the external library contains a descriptive auxiliary string for the result quantity then the first two fields in the dataset name will be constructed consisting of the dataset root string UNKNOWN. followed by the descriptive auxiliary string treated as a qualifier. For example, an unknown element quantity with descriptive auxiliary string ViscosityMolecular would have the following dataset name.

UNKNOWN.[ViscosityMolecular].E:1

The Contents attribute contains an alphanumeric string which is a short description of the contents of the dataset.

The Title, Subtitle, Label and Sublabel attributes contain alphanumeric strings which are short descriptions of the dataset. These attributes are meant to range from the more general to the more specific descriptions of the dataset. More specifically, the Title and Subtitle attributes are associated with the overall model and the Label and Sublabel attributes are associated with solution steps. Some interfaces will also define up to four additional subtitles, ie. Subtitle1, Subtitle2, Subtitle3 or Subtitle4.

The DataSource attribute contains an alphanumeric string which indicates the primary source of the data contained in the dataset. This string may be a file name or a data block or record number, etc. depending on the technical nature of the interface.

The DataType attribute contains an alphanumeric string which specifies the mathematical data type of the contents of the dataset.

  • Possible DataType attributes

    Scalar

    Vector

    x, y, z

    Tensor

    xx, yy, zz, xy, yz, zx

    GeneralTensor

    xx, xy, xz, yx, yy, yz, zx, zy, zz

    ElementResult

    Beams - Nxx, Myy, Mzz, Torque, Qxy, Qzx, TBimom
    Beams - Exx, Kyy, Kzz, Twist, Txy, Tzx, TBicur
    Shells - Nxx, Nyy, Nxy, Mxx, Myy, Mxy, Qxz, Qyz
    Shells - Exx, Eyy, Exy, Kxx, Kyy, Kxy, Txz, Tyz
    Spring - Fx, Fy, Fz, Mx, My, Mz
    Spring - Ex, Ey, Ez, Kx, Ky, Kz

    Scalars

    SixDof

    Tx, Ty, Tz, Rx, Ry, Rz

The Category attribute contains an alphanumeric string which qualifies the contents of the dataset.

  • Possible Category attributes

    Buckling

    buckling mode

    Vibration

    vibration mode

    Static

    static solution

    Transient

    transient solution

    Constraint

    single degree of freedom applied unit value mode

    ConcentratedLoad

    single degree of freedom concentrated load mode

    DistributedLoad

    distributed load mode

The Complex attribute contains an alphanumeric string which specifies the interpretation of the components of complex data types.

  • Possible Complex attributes

    Real

    Real part of real-imaginary results

    Magnitude

    Magnitude part of magnitude-phase results

    Imaginary

    Imaginary part of real-imaginary results

    Phase

    Phase part of magnitude-phase results

The Link.Complex attribute which names the associated dataset containing the imaginary or phase part of the complex result. If only the real or magnitude part exists, then the Link.Complex Attribute will be omitted. Note that phase is always in degrees.

The StrainType attribute contains an alphanumeric string which specifies the specific interpretation of strain results. This attribute will only appear on results datasets containing strain. The shear strain components of engineering strain are twice the value of the corresponding components of tensor strain.

  • Possible StrainType attributes

    Tensor

    Tensor strain

    Engineering

    Engineering strain

If the contents of a results dataset are expressed in a coordinate system other than the global coordinate system, then a Link.Cid global variable will be present and will point to a dataset of associated coordinate system identifiers for each node or element.

A complete list of supported result types are listed below. A defined constant is associated with each result type. The convention for defined constants is to prepend “SYS_RES_” to the dataset root. For example the defined constant SYS_RES_A refers to acceleration results. A data type is associated with each result type. It is possible for some result types with an associated vector data type in the table below to actually be either a vector (SYS_VECTOR) or sixdof (SYS_SIXDOF). data type. This is particularly true of Displacement type results. It is also possible for any result type to be either a scalar (SYS_SCALAR) or scalars (SYS_SCALARS) if the number of components representing the result does not match the default number of components for the data type.

The results dimensions use the following symbols for basic physical dimensions.

  • Physical dimension symbols and conventions

    None, Dimensionless
    Unknown, Unknown dimensions
    Variable, Variable dimensions
    M, Mass
    L, Length
    T, Time
    Q, Charge
    K, Temperature
    %, Percent
    Rad, Radians
    Deg, Degrees
    Cyc, Cycles

Specific Results Types

Dataset Root

Data Type

Contents

Dimensions

A

Vector

Acceleration

L/T2,Rad/T2

AREA

Scalar

Area

L2

BODY_FORCE

Vector

Body Force

M/L2T2

CAP_MAT

Matrix

Capacitance Matrix

ML2/T2K

CHEM_SHRINKAGE

Scalar

Chemical Shrinkage

None

CLOSURE

Scalar

Closure

L

CODE_CHECK

Scalar

Code Check Value

deprecate

CONC

Scalar

Concentration

Unknown

COND

Scalar

Thermal Conductivity

ML/T3K

COND_MAT

Matrix

Thermal Conductivity Matrix

ML/T3K

CONV_COEF

Scalar

Convection Coefficient

M/T3K

CP

Scalar

Specific Heat

M2/KT2

CRACK_DENSITY

Scalar

Crack Density

1/L

D

Vector

Displacement

L,Rad

DAMAGE

Scalars

Damage

DAMP

Scalar

Damping

M/T

DECIBEL

Scalar

Decibel

None

DELETED

Scalar

Deleted Entities

None

DENS

Scalar

Density

M/L3

DENS_GRAD

Vector

Density Gradient

M/L4

DIR

Vector

Direction Vector

None

DIR_COS

GeneralTensor

Direction Cosine Matrix

None

DIST

Scalar

Distance

L

DOM_FLUID_PHASE

Scalar

Dominant Fluid Phase

D_MAT

Matrix

Damping Matrix

M/T

E

Tensor

Strain

None

ELEC_FIELD

Vector

Electric Field

ML/T2Q

ELEC_FLUX

Vector

Electric Flux

ML3/T2Q

ELEC_POT

Scalar

Electric Potential

ML2/T2Q

ENERGY

Scalar

Energy

ML2/T2

ENERGY_DENSITY

Scalar

Energy Density

M/LT2

ENTROPY

Scalar

Entropy

ML2/T2K

EN_FLUX

Scalar

Element Nodal Heat Flux

ML2/T3

EN_FORC

Vector

Element Nodal Force

ML/T2,ML2/T2

E_RATE

Tensor

Strain Rate

1/T

FACTOR

Scalar

Factor

None

FAIL_INDEX

Scalars

Failure Index

None

FAT_DAMAGE

Scalar

Fatigue Damage

FAT_DAMAGE_DIR

Vector

Fatigue Damage Direction

None

FAT_LIFE

Scalar

Fatigue Life

FAT_SAFE_FACT

Scalar

Fatigue Safety Factor

None

FILM_COEF

Scalar

Film Coefficient

M/T3K

FLUENCE

Scalar

Fluence

M/T2

FORCE

Vector

Force

ML/T2

FORCE_MOMENT

SixDof

Force Moment

ML/T2,ML2/T2

FRACTION

Scalar

Fraction

None

FREQ

Scalar

Frequency

Cyc/T

GAP

Scalar

Gap

L

H

Scalar

Enthalpy

ML2/T2

HEAT

Scalar

Heat Generated

ML2/T2

HEAT_FLOW

Scalar

Heat flow

ML2/T3

HEAT_FLUX

Vector

Heat Flux

M/T3

HEAT_GRAD

Vector

Heat Gradient

ML/T2

H_DOT

Scalar

Enthalpy Rate

ML2/T3

ID

Scalar

Identifier

None

INERTIA

Tensor

Rotary Inertia

ML/Rad

INTENSITY

Vector

Intensity

M/T3

J

Scalar

Current

Q/T

J_DENSITY

Vector

Current Density

Q/TL2

KE

Scalar

Kinetic Energy

ML2/T2

KE_DENSITY

Scalar

Kinetic Energy Density

M/LT2

KE_PERCENT

Scalar

Kinetic Energy Percent

%

K_MAT

Matrix

Stiffness Matrix

M/T2

LENGTH

Scalar

Length

L

LOAD_FACT

Scalar

Load Factor

None

L_VEC

Vector

Load Vector

ML/T2

MACH

Scalar

Mach Number

None

MAG_FIELD

Vector

Magnetic Field

M/TQ

MAG_FLUX

Vector

Magnetic Flux

M/TQ

MAG_POT

Scalar

Magnetic Potential

ML/TQ

MARG_SAFE

Scalar

Margin of Safety

None

MASS

Scalar

Mass

M

MASS_FLOW

Scalar

Mass Flow

M/T

MASS_FLUX

Vector

Mass Flux

M/L2T1

MOMENT

Vector

Moment

ML2/T2

MU_LAMB

Scalar

Mu Lambda

MU_TURB

Scalar

Mu Turbulent

M_MAT

Matrix

Mass Matrix

M

NUMBER

Scalar

Number

None

ORDER

Scalar

Order

None

P

Vector

Momentum

ML/T

PENE_CONTACT

Scalar

Contact Penetration

L

PHASE_DIAMETER

Scalar

Phase Diameter

POISSONS_RATIO

Scalar

Poisson’s Ratio

None

POROSITY

Scalar

Porosity

None

POWER

Scalar

Power

ML2/T3

POWER_DENSITY

Scalar

Power Density

M/LT3

PRANDTL

Scalar

Prandtl Number

None

PRES

Scalar

Pressure

M/LT2

PRES_COEF

Scalar

Pressure Coefficient

None

PRES_DOT

Scalar

Pressure 1st Time Derivative

M/LT3

PRES_DOTDOT

Scalar

Pressure 2nd Time Derivative

M/LT4

PRES_GRAD

Scalar

Pressure Gradient

M/L2T2

PROBABILITY

Scalar

Probability

None

Q

Scalar

Electric

Q

QUAL_INDEX

Scalar

Quality Index

deprecate

R

Vector

Reaction Force

ML/T2,ML2/T2

RADIANCE

Scalar

Radiance

M/T3

RADIUS

Scalar

Radius

L

RC_PROD

Scalar

RC Product

T

REACTION_PROGRESS

Scalar

Reaction Progress

REYNOLDS

Scalar

Reynolds Number

None

ROTATION

Vector

Rotation

Rad

ROT_ANG

Vector

Rotation Angle Vector

Deg

ROUGHNESS

Scalar

Roughness

L

R_HEAT_FLOW

Scalar

Reaction Heat Flow

ML2/T3

R_J

Scalar

Reaction Current

Q/T

R_MASS_FLOW

Scalar

Reaction Mass Flow

M/T

R_Q

Scalar

Reaction Charge

Q

S

Tensor

Stress

M/LT2

SAFE_FACT

Scalar

Safety Factor

None

SCALARS

Scalars

Scalars

Unknown

SD

Vector

Element Displacement

L

SDV

Scalars

Solution Dependent Variable

Unknown

SE

Scalar

Strain Energy

ML2/T2

SEK

ElementResult

Strain and Curvature

Variable

SEP

Scalar

Separation

L

SE_DENSITY

Scalar

Strain Energy Density

M/LT2

SE_ERROR

Scalar

Strain Energy Error

ML2/T2

SE_PERCENT

Scalar

Strain Energy Percent

%

SF

Vector

Element Force

ML/T2

SFM

ElementResult

Stress and Moment Resultant

Variable

SHEAR_MODULUS

Scalar

Shear Modulus

M/LT2

SOUND_LEVEL

Scalar

Sound Level

SDRC

SOUND_MODEL

Scalar

Sound Wave Model Source

STAT

Scalar

Status

None

STATE

Scalar

State

None

STIFF

Scalar

Stiffness

M/T2

STREAM

Scalar

Stream Function

L2/T

STRENGTH_RATIO

Scalars

Strength Ratio

None

STRENGTH_SAFE_FACT

Scalar

Strength Safety Factor

None

S_ERROR_DENSITY

Scalar

Stress Error Density

deprecate

TE

Scalar

Thermal Energy

ML2/T2

TEMP

Scalar

Temperature

K

TEMP_DOT

Scalar

Temperature 1st Time Derivative

K/T

TEMP_GRAD

Vector

Temperature Gradient

T/L

TEXP_COEF

Scalar

Thermal Expansion Coefficient

1/K

TE_DENSITY

Scalar

Thermal Energy Density

M/LT2

TE_ERROR

Scalar

Thermal Energy Error

ML2/T2

TF

Vector

Total Residual Force

ML/T2

THICKNESS

Scalar

Element Thickness

L

TIME

Scalar

Time

T

TRAC

Vector

Traction

M/LT2

TRANSLATION

Vector

Translation

L

TURB_DIST

Scalar

Turbulent Distance

L

TURB_ED

Scalar

Turbulent Energy Dissipation

TURB_KE

Scalar

Turbulent Kinetic Energy

L2/T2

TURB_SD

Scalar

Turbulent Specific Dissipation

UNKNOWN

Scalars

Unknown

Unknown

USER

Scalar

User Selected Variable

Unknown

UTAU

Scalar

U Tau

V

Vector

Velocity

L/T,Rad/T

VIEW_FACT

Scalar

View Factor

None

VISC

Scalar

Viscosity

M/LT

VISC_EDDY

Scalar

Eddy Viscosity

M/LT

VOF

Scalar

Volume of Fluid

VOID_RATIO

Scalar

Void Ratio

None

VOLT

Scalar

Voltage

ML2/T2Q

VOLUME

Scalar

Volume

L3

VOLUME_FLUX

Scalar

Volume Flux

L/T

VORTICITY

Vector

Vorticity

1/T

V_DIV

Scalar

Velocity Divergence

1/T

V_GRAD

GeneralTensor

Velocity Gradient

1/T

WALL_SHEAR

Vector

Wall Shear

M/LT2

WATER_ACCUM

Scalar

Water Accumulation

Maya

WEIGHT

Scalar

Weight

ML/T2

X

Vector

Position

L

XF

Vector

External Applied Force

ML/T2,ML2/T2

XF_HEAT_FLOW

Scalar

External Applied Heat Flow

ML2/T3

XF_J

Scalar

External Applied Current

Q/T

XF_MASS_FLOW

Scalar

External Applied Mass Flow

M/T

XF_Q

Scalar

External Applied Charge

Q

YOUNGS_MODULUS

Scalar

Young’s Modulus

M/LT2

YPLUS

Scalar

Y Plus

None

Note that DELETED is defined as follows:

0. active, not deleted
1. deleted

Note that STATE is defined as follows:

0. open and not near contact or simply open
1. open and near contact
2. closed and sliding or simply closed
3. closed and sticking
4. sealed
5. leaked
6. crushed

A complete list of supported specific qualifier types are listed below. Like result types, a defined constant is associated with each qualifier type. The convention for qualifier defined constants is to prepend “SYS_QUA_” to the dataset qualifier string. For example the defined constant SYS_QUA_MAG refers to the MAG qualifier. A dataset name may contain several qualifiers.

Specific Qualfier Types

Dataset Qualifier

Contents

X

X Component

Y

Y Component

Z

Z Component

MAG

Magnitude

XX

XX Component

YY

YY Component

ZZ

ZZ Component

XY

XY Component

YZ

YZ Component

ZX

ZX Component

VONVISES

Von Mises

NXX

NXX Component

NYY

NYY Component

NXY

NXY Component

MXX

MXX Component

MYY

MYY Component

MXY

MXY Component

QXY

QXY Component

QZX

QZX Component

TX

X Translation

TY

Y Translation

TZ

Z Translation

RX

X Rotation

RY

Y Rotation

RZ

Z Rotation

MINIMUM

Minimum

MAXIMUM

Maximum

INTERMEDIATE

Intermediate

EFF

Effective

TOT

Total

INC

Incremental

REL

Relative

ABS

Absolute

EQUIV

Equivalent

PLAST

Plastic

CREEP

Creep

THERMAL

Thermal

ELAST

Elastic

INELAST

Inelastic

REYNOLDS

Reynolds

STAG

Stagnation

CAUCHY

Cauchy

PK

Piola-Kirchhoff

LOG

Logarithmic

TRESCA

Tresca

YIELD

Yield

MEAN

Mean

NONLIN

Nonlinear

TORSION

Torsion

SWELLING

Swelling

CRACKING

Cracking

NORMAL

Normal

SHEAR

Shear

PSD

Power Spectral Density

RMS

Root Mean Square

NOM

Nominal

COMPONENT

Component

INTERLAMINAR

Interlaminar

PLY

Ply

BOND

Bond

CONSTRAINT

Constraint

PORE

Pore

GASKET

Gasket

CONTACT

Contact

INITIAL

Initial

FINAL

Final

TENS

Tension

COMP

Compression

INFRARED

Infrared

DIFFUSE

Diffuse

COLLIMATED

Collimated

SOLAR

Solar

SOUND

Sound

RADIATIVE

Radiative

CONDUCTIVE

Conductive

CONVECTIVE

Convective

RESIDUAL

Residual

ADJUSTED

Adjusted

PSIDE

Plus Side

MSIDE

Minus Side

LOCAL

Local

INTERNAL

Internal

BULK

Bulk

STATIC

Static

DYNAMIC

Dynamic

SLIP

Slip

FREE

Free

FORCED

Forced

FLUID

Fluid

HARMONIC

Harmonic

LORENTZ

Lorentz

BEARING

Bearing

RADIAL

Radial

TANG

Tangential

AXIAL

Axial

GREEN

Green

MEANPRES

Mean Pressure

DEFORM

Deformation

GPF

Element Grid Point

VISC

Viscous

GLUE

Glue

FRICTION

Friction

SCALAR

Scalar

PRINCIPAL

Principal

NOZ

Number of Zero Crossings

ATO

Auto-correlation

CRM

CRMS

SOURCE

Source

SUM

Sum

NET

Net

ACOUSTIC

Acoustic

RADIOSITY

Radiosity

IRRADIANCE

Irradiance

TRANSMITTED

Transmitted

REFLECTED

Reflected

INCIDENT

Incident

ABSORBED

Absorbed

SOLID

Solid

DISSIPATED

Dissipated

OPT

Optimization

CONC

Concentrated

DIST

Distributed

ATTACHMENT

Attachment

INERTIA

Inertia

REDUCED

Reduced

MODAL

Modal

PARTICLE

Particle

NEIGHBORS

Neighbors

SHELL

Shell

ERROR

Error

NORM

Norm

EIP

Element Integration Point

PANEL

Panel

LCR

Level Crossing

PEAK

Peak

DIFFERENCE

Difference

CYCLIC

Cyclic

SIN

Sine

COS

Cosine

MPC

Multi-Point Constraint

TURB

Turbulent

SYM

Symmetric

ASYM

ASymmetric

DAMAGE

Damage

COHESIVE

Cohesive

BOLT

Bolt

CHOCKING

Chocking

MECH

Mechanical

ENTHALPY

Enthalpy

FLAME

Flame

GLASSTRANS

Glass Transition

EMITTED

Emitted

COUPLED

Coupled

UNCOUPLED

Uncoupled

STRUCT

Structural

TRANSPORT

Transport

WELD

Weld

KIRCH

Kirchhoff

NORMALIZED

Normalized

APPLIED

Applied

GYROSCOPIC

Gyroscopic

DAMPING

Damping

CIRCULATORY_FORCES

Circulatory Forces

TSAI_HILL

Tsai Hill

TSAI_WU

Tsai Wu

HASHIN

Hashin

STRESS

Stress

STRAIN

Strain

RICE_TRACEY

Rice Tracey

HOFFMAN

Hoffman

PUCK

Puck

LARC04

Larc04

TRACTION

Traction

PLANE_STRESS

Plane Stress

FIBER

Fiber

MATRIX

Matrix

RATIO

Ratio

VARIANT

Variant

CRITICAL

Critical

FULL

Full

MATERIAL

Material

AWEIGHT

A Weighted

INACTIVE

Inactive

OVERHEATING

Overheating

UNCONVERGED

Unconverged

VOLUME

Volume

VIBE

Vibration

PENETRATION

Penetration

BIOT

Biot

ENG

Engineering

CG

Center of Gravity

PARTICIPATION

Participation

AERO

Aeroelastic

MASS

Mass

RIGID

Rigid

UNIT

Unit

RESTRAINED

Restrained

UNRESTRAINED

Unrestrained

HOMOGENIZED

Homogenized

A complete list of supported valued qualifier types are listed below. The convention for valued qualifier defined constants is to prepend “SYS_QUAVAL_” to the dataset qualifier string before the = sign in capital letters. For example the defined constant SYS_QUAVAL_PHASE refers to the Phase= qualifier. A dataset name may contain several qualifiers.

Specific Valued Qualfier Types

SYS_QUAVAL_xxxxx

Dataset Qualifier/Contents

Value

CONTACTPAIR

ContactPair=integer

Contact Pair

DESIGNCYCLE

DesignCycle=integer

Design cycle

HARMONIC

Harmonic=integer

Harmonic number

INTPNT

IntPnt=integer

Element integration point number

MATERIAL

Material=integer

Material number

NODE

Node=integer

Node number

PHASE

Phase=integer

Chemical phase

ROTORSPEED

RotorSpeed=real

Rotor speed

SET

Set=integer

Set number

SECTION

Section=integer

Section number or BOT, MID, TOP

SPECIES

Species=integer

Species number

BOLTSEQUENCE

BoltSequence=integer

Bolt sequence number

STEP

Step=integer

Step identifier

1.9. Finite Element Results at Section Points

The term “section point” refers to a set of results located at a point within an element which is not spanned by the connecting nodes. This type of point is used primarily for shell and beam type elements. For example stress data may be positioned at layers in a laminated shell element or points within the cross section of a beam element. In this case the dataset will have a Link.Section attribute which points to an integer element dataset containing the layer position type and number of section points per element. The dataset may also have a Link.Layers attribute which points to an integer dataset containing the layer position type and layer number of each section. The Link.Layers attribute is only present if the layer numbers in the element are not numbered sequentially from one. For example, a stress dataset with section data would appear as follows:

S.EL:1:1
Link.Section = ELEM.SECT.E:1

If the layer numbers are not numbered sequentially from one, then a Link.Layers attribute will appear.

S.EL:1:1
Link.Section = ELEM.SECT.E:1
Link.Layers = ELEM.LAYS.T:1

An entry in the Link.Section dataset for an element consists of the layer position type and the number of sections. These two items are represented by 8 bit and 24 bit numbers and are concatenated together into a single integer. The layer position type is in the upper 8 bits and the number of sections is in the lower 24 bits. Possible layer position types are:

  • SYS_LAYERPOSITION_NONE, 0, layer position unknown, generic section

  • SYS_LAYERPOSITION_MID, mid layer

  • SYS_LAYERPOSITION_BOTTOP, bottom and top of layer

  • SYS_LAYERPOSITION_BOTMIDTOP, bottom, middle and top of layer

  • SYS_LAYERPOSITION_INTPNT, integration point

  • SYS_LAYERPOSITION_BOTMID, bottom and middle of layer

  • SYS_LAYERPOSITION_MIDTOP, middle and top of layer

  • SYS_LAYERPOSITION_BOT, bottom of layer

  • SYS_LAYERPOSITION_TOP, top of layer

  • SYS_LAYERPOSITION_B1M, between bottom and middle of layer

  • SYS_LAYERPOSITION_M1T, between middle and top of layer

  • SYS_LAYERPOSITION_B5T, bottom, between bottom and middle, middle, between middle and top and top of layer

For example if the layer position type is SYS_LAYERPOSITION_BOTTOP and the number of sections is 10, this implies that there are 5 layers. The first section is the bottom of layer 1, the second section is the top of layer 1, the third section is the bottom of layer 2 and so on until the 10th section is the top of layer 5.

Elements which generally do not have multiple section points, such as solid elements, are given a layer position type of SYS_LAYERPOSITION_NONE and the number of sections is 1.

If there is a Link.Layers attribute, then there is an entry in the Link.Layers dataset for each section point at each element in ascending element index order. This implies that the length of this dataset is the sum of the number of section points indicated in the Link.Section dataset. An entry in the Link.Layers dataset consists of a layer position and a layer number. The possible layer position types are a subset of the layer position types listed above.

  • SYS_LAYERPOSITION_NONE, 0, layer position unknown, generic section

  • SYS_LAYERPOSITION_MID, mid layer

  • SYS_LAYERPOSITION_INTPNT, integration point

  • SYS_LAYERPOSITION_BOT, bottom of layer

  • SYS_LAYERPOSITION_TOP, top of layer

  • SYS_LAYERPOSITION_B1M, between bottom and middle of layer

  • SYS_LAYERPOSITION_M1T, between middle and top of layer

These two items are represented in a manner similar to the Link.Section dataset. The layer position type is in the upper 8 bits and the layer number is in the lower 24 bits. This limits the layer number to a non-negative value less than approximately 16 million. This dataset must be traversed serially using the number of sections for each element obtained from the Link.Section dataset.

If the layer position for an element in the Link.Section dataset is SYS_LAYERPOSITION_NONE, then the layer position is assumed to be in the Link.Layers dataset entry for the element.

1.10. Finite Element Results at System Degrees of Freedom

The term system degree of freedom refers to a result associated with a specific degree of freedom (unknown) at a node or element and/or a generalized degree of freedom associated with a mode of some type, ie. vibration mode or residual load mode. Examples of system degrees of freedom are the x translation at a node or the pressure in an element. A dataset with results associated with system degrees of freedom is given a name with a .D suffix and requires a list of associated entities and degree of freedom types. In this case the dataset will have a Link.EntDof attribute which points to an integer dataset with two rows containing the entity index and degree of freedom type for each system degree of freedom. These datasets are useful for storing reduced mass and stiffness matrices resulting from some form of superelement generation.

For example, a superelement stiffness dataset with coefficients associated with specific degrees of freedom would appear as follows:

K_RED.D:1
Link.EntDof = ENT.DOF.T:1

1.11. Finite Element Results History Datasets

A results history dataset (or simply history dataset) contains a result for a small subset of a finite element model for many steps of a solution. This is in contrast to a full model results dataset which generally contains a result for an entire finite element model at a single step. For example, a single history dataset may contain the velocity for nodes 100,200 and 300 for every 10th time step from time steps 10 to 1000.

The rules for constructing a history dataset name and attributes for a particular result quantity are built upon the methodology used for a full model results dataset. All history dataset names start with the string HIST. The remaining fields are constructed in a manner identical to the syntax for a full model dataset for the same quantity (TEMP, SE_DENSITY, etc.) and format (.N, .E, .EL, etc.). The dataset numeric identifiers for a history dataset generally contain one less identifier than the corresponding full model results dataset. This is due to the fact that the last identifier in a full model results dataset is the single “step” value for which the corresponding history dataset will contain values for many steps.

A number of additional history datasets are required to hold the node or element indices of the history subset, the step numbers and any dataset attributes such as Time or Frequency. It is also possible that some components of a multi-dimensional data type, such as a vector or tensor data type, may be defined and others not defined. The results history dataset will contain link attributes which point to these additional history datasets. The following links are used.

Link.Index       node or element indices for subset
Link.Component   component defined/undefined flags
Link.Step        step numbers
Link.Time        step times
Link.Frequency   step frequencies
Link.LoadFactor  step load factors

For example, if a series of 100 results datasets for displacement at nodes for a transient response analysis were to exist as follows:

D.N:2:1
D.N:2:2
...
D.N:2:100

A corresponding time history dataset containing displacements for a subset of nodes for all 100 steps would be named as follows:

HIST.D.N:2
 Link.Index = INDX.NODE.T
 Link.Time = HIST.TIME.T:2
 Link.Step = HIST.STEP.T:2

For node array format, .N, or element array format, .E, ncol=number of nodes or elements and nrow=number of data type components. In general, the number of steps is must be computed as the length of the Link.Step dataset. However, in the case of rectangular array format the number of steps can be computed by dividing the record length, lrec , by ncol*nrow. The order of the data in history datasets is the data components of the first entity at all steps followed by the data components of the second entity at all steps, etc.

The node indices of the subset are placed in the following dataset.

INDX.NODE.T

The associated dataset attributes such as Time, Frequency which may exist for the history of steps are placed in separate history datasets.

HIST.TIME.T:2
 Link.Step = HIST.STEP.T:2
HIST.FREQUENCY.T:2
 Link.Step = HIST.STEP.T:2

Certain quantities such as strain energy for the whole model are placed in separate history datasets.

HIST.SE.T:2
 Link.Time = HIST.TIME.T:2
 Link.Step = HIST.STEP.T:2

The step numbers themselves, in this case 1,2,…,100, are placed in a separate integer valued history dataset.

HIST.STEP.T:2

If only a subset of the vector components are defined for certain nodes, then the defined component flags for the velocity components of each node are placed in a dataset and a link attribute is created as follows.

Link.Component = COMP.N:1

The dataset, COMP.N:1, contains an integer flag for each component of the vector data for each node. The value will be 0 if the component is undefined. The component flags dataset itself is, in general, a sparse dataset and will therefore have a link to the defined indices.

COMP.N:1
  Link.Index = INDX.NODE.T

1.12. Monitor Functions

The user may define a monitor function to track the progress of an interface as it opens a file or reads data from the file. The file opening operation, vdm_DataFunOpen(), can be particularly time consuming. The opening operation generally reads and develops the basic simulation model first and then continues to build a list of the simulation results.

Using a monitor function, the progress may be monitored by querying for the progress phase and source. The phase may be one of three integer values, VDM_PHASE_OPENMODEL, VDM_PHASE_OPENMODELCOMPLETE and VDM_PHASE_OPENRESULT. The source is an interface dependent string which describes the current position in the file which is being read. The phase is queried using the function vdm_DataFunGetInteger() with type VDM_PHASE and the source is queried using the function vdm_DataFunGetString() with type VDM_SOURCE. When the phase value is VDM_PHASE_OPENMODELCOMPLETE the model construction is complete and the user may load the model and return to allow the opening of the simulation results.

Because of the appearance of the interface object as the first argument in the monitor function, a separate callback function must be defined for each interface that requires monitoring. For example, the NASLib interface monitor function is set in NASLib using vdm_NASLibSetFunction(). A similar methodology is used for all interface objects.

1.13. Multiple Domains

Some file types support multiple domains, usually associated with a domain decomposition for a parallel solution. It is possible to query for the number of domains and then open each domain separately, instancing an interface object and an associated DataFun object for each one. This allows the individual domains to be opened and read in parallel.

The vdm_DataFunNumDomains() function may be called to perform a scan to determine the number of domains. The number of domains is returned. If a library device does not support multiple domains, then the number of domains will be zero. If a non-zero number of domains is returned then vdm_DataFunOpen() may then be called with mode set to zero to open all domains or with mode > 0 indicating the specific domain number to open. If the number of domains is returned as zero, then vdm_DataFunOpen() must be called with mode set to zero to open the library device.

  • For example, the OpenFOAMLib interface supports multiple domains. vdm_DataFunNumDomains() is called to determine the number of domains, numdomains.

    vdm_OpenFOAMLib *openfoamlib;
    vdm_DataFun *df;
    Vint numdomains;
                /* create OpenFOAM device object */
    openfoamlib = vdm_OpenFOAMLibBegin();
                /* create data function */
    df = vdm_DataFunBegin();
                /* load data functions for OpenFOAMLib device */
    vdm_OpenFOAMLibDataFun (openfoamlib,df);
                /* query library on file "example/system/controlDict" */
    vdm_DataFunNumDomains (df,"example/system/controlDict",SYS_OPENFOAM,
                        &numdomains);
    
  • It is possible that the value of numdomains is returned as zero. In this case the interface does not support or the specific file does not contain multiple domains. Note that even if the number of domains is non-zero, vdm_DataFunOpen() may be called with a zero. In this case all domains will be opened.

    Vint i;
    Vint ierr;
                /* open each individual domain */
    if(numdomains) {
    dfdom =(vdm_DataFun**)malloc(numdomains*sizeof(vdm_DataFun*));
    ofdom =(vdm_OpenFOAMLib**)malloc(numdomains*sizeof(vdm_OpenFOAMLib*));
    for(i = 0; i < numdomains; i++) {
        ofdom[i] = vdm_OpenFOAMLibBegin();
        dfdom[i] = vdm_DataFunBegin();
        vdm_OpenFOAMLibDataFun (ofdom[i],dfdom[i]);
        vdm_DataFunOpen (dfdom[i],i+1,"example/system/controlDict",
                        SYS_OPENFOAM);
        if(ierr) {
        }
    }
                /* open all */
    } else {
    vdm_DataFunOpen (df,0,"example/system/controlDict",SYS_OPENFOAM);
    ierr = vdm_DataFunError(df);
    if(ierr) {
    }
    }
    
                /* close and delete */
    if(numdomains) {
    for(i = 0; i < numdomains; i++) {
        vdm_DataFunClose (dfdom[i]);
        vdm_OpenFOAMLibEnd (ofdom[i]);
        vdm_DataFunEnd (dfdom[i]);
    }
    free(dfdom);
    free(ofdom);
    } else {
    vdm_DataFunClose (df);
    }
                /* delete common objects */
    vdm_OpenFOAMLibEnd (openfoamlib);
    vdm_DataFunEnd (df);
    

1.14. A First Program

As an example of a simple CEETRON Access application the following program illustrates reading the node coordinates from an SDRC Universal file. Initially a SDRCLib library device object and an associated DataFun object are created. The DataFun object is loaded with function pointers by the SDRCLib object. At this point all further data access of the library device is handled through the data function object.

The library device is opened and attached to a SDRC Universal file residing on file bumper.unv. The library name, “Bumper”, is completely arbitrary. Opening a library device is often the most time consuming operation performed on a library. In general the library file must be completely scanned to determine the number and type of datasets, number of nodes and elements, etc. At this point the library device may be queried for dataset information. First the number of nodes are queried. Next, in order to begin to query for the existence of datasets, the internal Library object used by the library device must be retrieved. The Library object is searched for the reserved dataset name for node coordinates, X.N. The dataset indices and number of datasets which satisfy the search are returned. In this case only a single dataset is found. Every library device should contain the node coordinate dataset, X.N, as part of the basic finite element model description. The node coordinates are read into memory and the first and last node coordinates are printed.

Finally the library device is closed and the DataFun and SDRCLib objects are destroyed. Note that the Library object does not “belong” to the user (was not instanced by the user) and therefore must not explicitly destroyed by the user.

#include <stdio.h>
#include <stdlib.h>
#include "base/base.h"
#include "vdm/vdm.h"

/*----------------------------------------------------------------------
                     Read and Print Finite Element Node Coordinates
----------------------------------------------------------------------*/
int
main()
{
    vdm_SDRCLib* sdrclib;
    vdm_DataFun* datafun;
    vdm_Library* library;
    Vint filetype;
    Vint numnp;
    Vint ndsts;
    Vint idst;
    Vfloat* px;
    Vint ind;

    /* create data function object */
    datafun = vdm_DataFunBegin();

    /* create SDRC Universal File library device */
    filetype = VDM_SDRC_UNIVERSAL;
    sdrclib = vdm_SDRCLibBegin();
    vdm_SDRCLibDataFun(sdrclib, datafun);

    /* open library device */
    vdm_DataFunOpen(datafun, 0, "bumper.unv", filetype);

    /* get number of nodes */
    vdm_DataFunGetNumEntities(datafun, SYS_NODE, &numnp);

    /* get library object */
    vdm_DataFunGetLibrary(datafun, &library);

    /* search for node coordinate dataset */
    vdm_LibrarySearchDataset(library, "X.N", 1, &idst, &ndsts);

    /* allocate memory for coordinates */
    px = (Vfloat*)malloc(3 * numnp * sizeof(Vfloat));

    /* read node coordinates */
    vdm_DataFunReadDataset(datafun, idst, px);

    /* print first and last node coordinates */
    printf("Node Coordinates\n");
    ind = 0;
    printf("%10d %12f %12f %12f\n", 1, px[ind], px[ind + 1], px[ind + 2]);
    ind = 3 * (numnp - 1);
    printf("%10d %12f %12f %12f\n", numnp, px[ind], px[ind + 1], px[ind + 2]);
    /* free memory */
    free(px);

    /* close library device */
    vdm_DataFunClose(datafun);

    /* free objects */
    vdm_DataFunEnd(datafun);
    vdm_SDRCLibEnd(sdrclib);
    return 0;
}

The output of this example program appears below.

Node Coordinates
        1    30.826799     3.968360    -1.003940
    652   -47.007900    15.984300    -2.007870

1.15. A First Program - C++ Version

The following program is a listing of the C++ version of the same “A First Program” listed above which used C language bindings.

#include <stdio.h>
#include <stdlib.h>
#include "base/base.h"
#include "vdm/vdm.h"

/*----------------------------------------------------------------------
                     Read and Print Finite Element Node Coordinates
----------------------------------------------------------------------*/
int
main()
{
    vdm_SDRCLib* sdrclib;
    vdm_DataFun* datafun;
    vdm_Library* library;
    Vint filetype;
    Vint numnp;
    Vint ndsts;
    Vint idst;
    Vfloat* px;
    Vint ind;

    /* create data function object */
    datafun = new vdm_DataFun();

    /* create SDRC Universal File library device */
    filetype = VDM_SDRC_UNIVERSAL;
    sdrclib = new vdm_SDRCLib();
    sdrclib->DataFun(datafun);

    /* open library device */
    datafun->Open(0, "bumper.unv", filetype);

    /* get number of nodes */
    datafun->GetNumEntities(SYS_NODE, &numnp);

    /* get library object */
    datafun->GetLibrary(&library);

    /* search for node coordinate dataset */
    library->SearchDataset("X.N", 1, &idst, &ndsts);

    /* allocate memory for coordinates */
    px = (Vfloat*)malloc(3 * numnp * sizeof(Vfloat));

    /* read node coordinates */
    datafun->ReadDataset(idst, px);

    /* print first and last node coordinates */
    printf("Node Coordinates\n");
    ind = 0;
    printf("%10d %12f %12f %12f\n", 1, px[ind], px[ind + 1], px[ind + 2]);
    ind = 3 * (numnp - 1);
    printf("%10d %12f %12f %12f\n", numnp, px[ind], px[ind + 1], px[ind + 2]);
    /* free memory */
    free(px);

    /* close library device */
    datafun->Close();

    /* free objects */
    delete datafun;
    delete sdrclib;
    return 0;
}