GeMA
The GeMA main application
Accessor object methods

Accessor objects are the abstraction used by GeMA for querying and updating values from data sets associated with nodes, cells, properties, Gauss points, boundary conditions, contact boundary conditions and discontinuity sets. The following advantages arise from using accessors:

  • Access to the data is done in a simple and unique way, independent of the underlying storage data structure and also independent of the data format (scalar, vector or matrix data);
  • Provides a way for evaluating user functions transparently (the user of an accessor does not need to know if the returned data came from a constant stored value or from the automatic evaluation of a user function);
  • Provides unit conversion also in a transparent way (if the requested unit is different from the data unit, values are converted automatically - see comments on the Working with units page).

There are six types of accessors, differing basically on the type of information needed for indexing the accessor. They are:

  • ValueAccessor: An accessor object used for retrieving data associated with mesh nodes: state variables, node coordinates and node attributes. Indexed by a node number. It is also used for retrieving property data directly from a property set (and not properties associated with a cell). In that case, the accessor is indexed by property set row indices.
  • CellAccessor: An accessor object used for retrieving data associated with mesh cells: cell attributes and cell properties. Indexed by a cell object.
  • GaussAccessor: An accessor object used for retrieving data associated with Gauss points. Indexed by a cell object and an integration point number.
  • BcAccessor: An accessor object used for retrieving data associated with boundary conditions. Indexed by a boundary condition line number and the list index of the node / cell inside the row application point (needed to identify the node / cell for function evaluation if the boundary condition is given by a user function) for value querying methods (value() and valueStr()) or just by the boundary condition line number for the remaining methods.
  • ContactBcAccessor: An accessor object used for retrieving property values associated with a contact boundary condition. Indexed by a property value line number (usually returned by a call to cbc:contactIndex()), the boundary condition line number and the list index of the node / cell inside the row application point (needed to identify the node / cell for function evaluation if the boundary condition is given by a user function) for value querying methods (value() and valueStr()) or just by the property value line number for the remaining methods.
  • DiscontinuityAccessor: An accessor object used for retrieving data associated with discontinuities: attributes and properties. Indexed by a discontinuity object. Important: discontinuity sets do not support function evaluation, neither for attributes nor for properties.

IMPORTANT: Following the Lua spirit, all indices used when calling accessor methods are one based, and not zero based as in the C++ API.

When querying values associated with cell objects, including Gauss attributes and (contact) boundary conditions attached to cells, an evaluation coordinate might be needed if the evaluated data is provided by an user function and one of the function parameters needs interpolation or is a coordinate value. If the coordinate is surely not needed (for node based values, for example), a nil value can be passed as coordinate, but otherwise, the required coordinate should be provided, usually an integration point coordinate.

When working with value accessors referencing a value set containing ghost nodes (affectedNodes equal to 'ghost' or 'both' - see the Data options page), some care is needed while indexing the "ghost part" of the accessor. See the comments on the Working with ghost nodes section.

Accessor objects are retrieved through the following methods from mesh, property set, boundary condition, contact boundary condition and discontinuity set objects: mesh:nodeCoordAccessor(), mesh:nodeValueAccessor(), mesh:cellAttributeAccessor(), mesh:cellPropertyAccessor(), mesh:gaussAttributeAccessor(), bc:propertyAccessor(), cbc:propertyAccessor(), cbc:materialAccessor(), pset:propertyAccessor(), dset:attributeAccessor() and dset:propertyAccessor().

Example:

local m = modelData:mesh('myMeshName')
local coordAc = m:nodeCoordAccessor()
local svAc = m:nodeValueAccessor('stateVarId', 'cm')
local naAc = m:nodeValueAccessor('nodeAttrId')
local caAc = m:cellAttributeAccessor('cellAttrId')
local pAc = m:cellPropertyAccessor('propertyId')
local gAc = m:gaussAttributeAccessor('gaussAttrId', 'N/m2')
local bcAc = modelData:boundaryCondition('bcName'):propertyAccessor('bcPropertyId')
local cbcAc = modelData:contactBoundaryCondition('cbcName'):propertyAccessor('cbcPropertyId')
local mAc = modelData:contactBoundaryCondition('cbcName'):materialAccessor()
local psAc = modelData:propertySet('psName'):propertyAccessor('propertyId')
local dsaAc = modelData:discontinuitySet('dsName'):attributeAccessor('dsAttributeId')
local dspAc = modelData:discontinuitySet('dsName'):propertyAccessor('dsPropertyId')

TIP: Whenever a Lua accessor object is created, a C++ accessor is also created. If an accessor will be used several times, it should be cached in a local variable to avoid needless creation and destruction of accessors. On the example bellow, the second option is much more resource friendly than the first.

local m = modelData:mesh('myMeshName')
-- Bad code ahead (although correct)
for i = 1, m:numNodes() do
m:nodeValueAccessor('Attr'):setValue(i, 0.0)
end
-- Better code (creates the accessor only once)
local ac = m:nodeValueAccessor('Attr')
for i = 1, m:numNodes() do
ac:setValue(i, 0.0)
end

Index:

Value metadata methods

ac:info()
Description: Returns an object with metadata about the values that can be queried/updated by this accessor.
Parameters: None.
Returns: Returns the accessor value info object.

Example:

local info = ac:info()


ac:size()
Description: Returns the number of values (lines) accessible through this accessor. Its interpretation depends on the actual data type.
- ValueAccessors: The returned size is either the number of geometry nodes in the mesh, the number of ghost nodes or their sum, depending of the actual affectedNodes property for this attribute or state variable (see Data options). It can also be the number of lines in a property set or the number of contacts in a contact boundary condition (for an accessor returned by the cbc:materialAccessor() method).
- CellAccessors: The returned size will be the number of cells in the mesh.
- GaussAccessors: The returned size will be the number of integration points in the mesh, taking into consideration that meshes can have heterogeneous element types, each one with its own integration rule.
- BcAccessors: The returned size will be the number of lines in the boundary condition definition.
- ContactBcAccessors: The returned size will be the number of property values (material contact properties) in the contact boundary condition definition.
- DiscontinuityAccessors: The returned size will be the number of discontinuity objects in the set.
Parameters: None.
Returns: Returns the number of entries in the associated data set.

Example:

local nlines = ac:size()


ac:valueSize()
Description: Returns the dimension of each individual value stored in the associated data set. The returned value will be equal to 1 for scalar values and to the data matrix dimension (lines * columns) for vector or matrix values.
Parameters: None.
Returns: Returns the dimension for each data entry in the associated data set.

Example:

local vsize = ac:valueSize()


ac:isScalar()
Description: Returns true if the value returned by the accessor is a scalar value, false otherwise. Can be false even if ac:valueSize() is equal to 1 (an 1x1 matrix, for example).
Parameters: None.
Returns: Returns true if the value returned by the accessor is a scalar value.

Example:

if ac:isScalar() then
...
end


ac:unit()
Description: Returns the unit in which data values are stored. See also the Working with units page.
Parameters: None.
Returns: Returns a string with the value unit.

Example:

local acunit = coord1Ac:unit()


ac:defValue()
Description: Returns the default data value converted to the acessor unit. Meaningfull only when the default value is NOT a function.
Parameters: None.
Returns: Returns a number for scalar data and a table for non scalar data. Matrix values are linearized using column major format (FORTRAN style).

Example:

local defVal = ac:defValue()


Value querying methods

valueAc:value(index)
cellAc:value(cell, coord, ip)
gaussAc:value(cell, ip, coord)
bcAc:value(bcIndex, listIndex, coord)
cbcAc:value(contactIndex, bcIndex, listIndex, coord)
dsAc:value(discontinuity)
Description: Returns the data set value at the supplied index, converted to the accessor unit if needed.
Parameters: index For ValueAccessors only, provides the node number for the queried node attribute or state variable. Should be a value from 1 to ac:size(). When querying a property table directly (rare since cell properties should be queried through a cell accessor), the index is the row number in the property table.
cell For CellAccessors and GaussAccessors only, provides the cell object whose cell, property or Gauss attribute is being queried.
ip For CellAccessors and GaussAccessors only, provides the index of the integration point being queried fro Gauss accessors. Should be a value from 1 to the number of integration points for the cell. For cell accessors, the ip number is needed only when querying a value from a cell attribute associated with a function that takes as parameter a Gauss attribute. In those special cases, the ip index is necessary for evaluating the function parameter and, if absent, the Gauss attribute parameter will be interpolated at the evaluation coordinate instead.
bcIndex For BcAccessors and ContactBcAccessors only, provides the boundary condition line index (a value from 1 to ac:size(). For a contact boundary condition accessor, used only for function evaluations.
listIndex For BcAccessors and ContactBcAccessors only, provides the list index of the node / cell inside the row application point (used for function evaluations only).
contactIndex For ContactBcAccessors only, provides the contact boundary condition property values line index (a value from 1 to ac:size().
discontinuity For DiscontinuityAccessors only, provides the discontinuity object whose attribute or prperty is being queried.
coord The set of coordinates needed when evaluating cell user functions expecting coordinate parameters or interpolated node parameters. Also used for cell user functions with Gauss attributes as parameters to interpolate the attribute value when the evaluation point is not an integration point or when the integration point index was not given by the user.
If needed, coordinates should be given by either a table or a matrix object. Otherwise, this parameter can be omitted or set to nil. If the mesh hosting the cell is an element mesh, coordinates should be given as natural element coordinates. Otherwise (cell mesh), coordinates should be cartesian coordinates.
Returns: Returns a number for scalar values or a matrix object for vectors / matrices.

Example:

-- Node coordinate accessor, returning a 3d point coordinate for the first
-- mesh node, with values converted to 'cm'
local ac = mesh:nodeCoordAccessor('cm')
local coord = ac:value(1)
local y = coord(2)
-- Traverses all integration points for the first mesh cell, querying the
-- value of the 'rho' scalar Gauss attribute
local ac = mesh:gaussAttributeAccessor('rho')
local rs = ac:info():ruleSet() -- Get the rule set that this accessor is bound to
local c = mesh:cell(1)
local ir = mesh:elementIntegrationRule(c:type(), rs)
for j=1, ir:numPoints() do
-- Get integration point coordinates
local ip, w = ir:integrationPoint(j)
-- Get current value (the ip parameter can be skipped if the rho attribute does not supports user function)
local rho = ac:value(c, j, ip)
print(rho)
end


valueAc:valueStr(index, options)
cellAc:valueStr(cell, coord, options, ip)
gaussAc:valueStr(cell, ip, coord, options)
bcAc:valueStr(bcIndex, listIndex, coord, options)
cbcAc:valueStr(contactIndex, bcIndex, listIndex, coord, options)
dsAc:valueStr(discontinuity, options)
Description: Returns the same value returned by ac:value(...), converted to a string.
Parameters: index, cell, ip, bcIndex, listIndex, contactIndex, discontinuity The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation.
coord The (optional) set of evaluation coordinates. See details on the ac:value() method documentation.
options An optional table with formatting options. Can contain the following fields (missing fields will assume their defaults):
- evalFunctions: Boolean flag defining whether functions should be evaluated or have its name printed. Default = false.
- defNil: Boolean flag defining whether default values should be printed as nil. Default = false.
- format: A string following the printf format specifying the numeric format that should be used. Default is equal to the value returned by ac:info():format().
Returns: Returns the data value converted to a string.

Example:

-- Prints the accessor value at the given node index using four decimal places
print(ac1:valueStr(nodeIndex, {format = '.4f'}))
-- Prints the accessor value at the given cell, evaluating cell functions at the element centroid
-- (assumes a quad element whose centroid natural coordinates are 0, 0).
print(ac2:valueStr(cell, {0, 0}, {evalFunctions = true})


valueAc:isDefValue(index)
cellAc:isDefValue(cell)
gaussAc:isDefValue(cell, ip)
bcAc:isDefValue(bcIndex)
cbcAc:isDefValue(contactIndex)
dsAc:isDefValue(discontinuity)
Description: Returns true if the value at the specified index is equal to the default value. Works even if the deafault value is a user function (in that case, will check if the data set value at the index is a refrence to the default user function).
Parameters: index, cell, ip, bcIndex, contactIndex, discontinuity - The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation. For Bc and ContactBc accessors, only the first index parameter is required, since no function evaluation is needed.
Returns: Returns true if the value at the specified index is equal to the default value, false otherwise.

Example:

-- Checks if the cell attribute value at the given element is equal to the default value or not
if ac:isDefValue(elem) then
...
end


valueAc:isFunction(index)
cellAc:isFunction(cell)
gaussAc:isFunction(cell, ip)
bcAc:isFunction(bcIndex)
cbcAc:isFunction(contactIndex)
dsAc:isFunction(discontinuity)
Description: Returns true if the value at the specified index is equal to a user function.
Parameters: index, cell, ip, bcIndex, contactIndex, discontinuity - The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation. For Bc and ContactBc accessors, only the first index parameter is required, since no function evaluation is needed.
Returns: Returns true if the value at the specified index is a user function, false otherwise.

Example:

-- Checks if the Gauss point attribute value at the first integration point of the given
-- element is equal to the default value or not.
if ac:isFunction(elem, 1) then
...
end


valueAc:functionId(index)
cellAc:functionId(cell)
gaussAc:functionId(cell, ip)
bcAc:functionId(bcIndex)
cbcAc:functionId(contactIndex)
dsAc:functionId(discontinuity)
Description: Returns the function name (id) of the value at the specified index.
Parameters: index, cell, ip, bcIndex, contactIndex, discontinuity - The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation. For Bc and ContactBc accessors, only the first index parameter is required, since no function evaluation is needed.
Returns: Returns the function name (id) or an empty string if the value is not a function.

Example:

local fname = ac:functionId(nodeIndex)


gaussAc:numPoints(cell)
Description: Returns the number of integration points of the specified cell. Available only for Gauss accessors.
Parameters: cell - The cell object.
Returns: Returns the number of integration points for the cell.

Example:

local npoints = ac:numPoints(cell)


valueAc:adjustLinearIndex(index, mesh)
valueAc:adjustLinearIndex(index, firstGhostIndex)
Description: Converts a linear index (from 1 to mesh:totalNumNodes()) into an index that can be used with every data-set affectedNodes type, including ghost only data-sets. See the comments on the Working with ghost nodes section for more details. Can be called with a mesh object as second parameter or with the linear index of the first mesh ghost node.
Parameters: index The linear index (from 1 to mesh:totalNumNodes()).
mesh The mesh object to which the data set is attached.
firstGhostIndex The linear index of the first mesh ghost node. This alternative form can be used if the mesh object is not available.
Returns: Returns the adjusted index or nil if the given index cannot be used with this accessor (geometry index for a ghost only accessor or ghost index for a geometry only accessor).

Example:

-- The second ghost node
local index = mesh:numNodes() + 2
-- Adjust index so that if ac supports geometry and ghost nodes, index remains unchanged, but if ac
-- only supports ghost nodes, index will revert to 2 and so can be used for indexing the accessor.
index = ac:adjustLinearIndex(index, mesh)
local val = ac:value(index)


Value updating methods

valueAc:setValue(index, value) / valueAc:setValue(index, dim, value)
cellAc:setValue(cell, value) / cellAc:setValue(cell, dim, value)
gaussAc:setValue(cell, ip, value) / gaussAc:setValue(cell, ip, dim, value)
bcAc:setValue(bcIndex, value) / bcAc:setValue(bcIndex, dim, value)
cbcAc:setValue(contactIndex, value) / cbcAc:setValue(contactIndex, dim, value)
dsAc:setValue(discontinuity, value) / dsAc:setValue(discontinuity, dim, value)
Description: Updates the data set value at the supplied index with the given value. Can also be used to update only one dimension of a vector / matrix value.
Parameters: index, cell, ip, bcIndex, contactIndex, discontinuity The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation. For Bc and ContactBc accessors, only the first index parameter is required, since no function evaluation is needed.
dim An optional dimension parameter. If given, should be a value between 1 and ac:valueSize() and implies that the given value will be a scalar used to update only the referenced dimension of the underlying vector/matrix data.
value The value to be set. Should be a number for scalar accessors or if dim has been given. Should be a table or a matrix object for vector / matrix accessors. If the value is given as a Lua table, this table can be a list with subtables, each containing the values for one matrix line, or a flat table linearized in column major order (FORTRAN style). The value can also be a string with a function name if the data set allows functions, or a constant name if the data set has an associated constant map (if the data set supports both, the given string will first be compared to constant names and then to function names).
IMPORTANT: Values should be given in the accessor unit (the unit used when creating the accessor or the data unit if no unit was given on the creation).
Returns: Nothing.

Example:

-- Updates the 2D coordinate of the fourth mesh node with the new {10.0, 10.0} coordinate
local ac = mesh:nodeCoordAccessor()
ac:setValue(4, {10.0, 10.0})
ac:setValue(4, Vector{10.0, 10.0})
-- Updates only the Y coordinate of the fourth mesh node to 20.0
ac:setValue(4, 2, 20.0)
-- Updates the value of the 'dummy' cell attribute, of the second mesh cell, with the 'dummyF' user function
local ac = mesh:cellAttributeAccessor('dummy')
local c = mesh:cell(2)
ac:setValue(c, 'dummyF')
-- Traverses all integration points for the first mesh cell, updating the
-- value of the 'rho' scalar Gauss attribute with a value equal to the
-- current value + 10%
local ac = mesh:gaussAttributeAccessor('rho')
local rs = ac:info():ruleSet() -- Get the rule set that this accessor is bound to
local c = mesh:cell(1)
local ir = mesh:elementIntegrationRule(c:type(), rs)
for j=1, ir:numPoints() do
-- Get integration point coordinates
local ip, w = ir:integrationPoint(j)
-- Get current value (the ip parameter can be skipped if the rho attribute does not supports user function)
local rho = ac:value(c, j, ip)
-- Update value
ac:setValue(c, j, rho*1.1)
end


valueAc:setValueAsDef(index)
cellAc:setValueAsDef(cell)
gaussAc:setValueAsDef(cell, ip)
bcAc:setValueAsDef(bcIndex)
cbcAc:setValueAsDef(contactIndex)
dsAc:setValueAsDef(discontinuity)
Description: Updates the data set value at the supplied index to its default value.
Parameters: index, cell, ip, bcIndex, contactIndex, discontinuity - The values needed as an accessor index for the specific accessor type. See details on the ac:value() method documentation. For Bc and ContactBc accessors, only the first index parameter is required, since no function evaluation is needed.
Returns: Nothing.

Example:

ac:setValueAsDef(nodeIndex)


Working with ghost nodes

When a value accessor references a data set containing ghost nodes (affectedNodes equal to "ghost" or "both"), some care must be taken when providing indices for methods to query or update values, especially if the data set stores values for ghost nodes only (affectedNodes equal to "ghost").

Indices passed to the accessor methods can be either:

  • For affectedNodes equal to "both": a 'ghost index', i.e., a value with the highest bit set, (as described by isMeshGhostNode()) or a 'linear' index, a value between 1 and mesh:totalNumNodes(), encompassing both geometry and ghost nodes.
  • For affectedNodes equal to "ghost": a 'ghost index' or a value between 1 and mesh:numGhostNodes().

In that way, querying all node values can be easily done in several use cases, as seen below for a coordinate accessor (that will always have affectedNodes equal to "both"):

local ac = mesh:nodeCoordAccessor()
-- Continuous linear indexing. Will access first regular geometric nodes (i <= mesh:numNodes())
-- and then ghost nodes (i > mesh:numNodes())
for i = 1, mesh:totalNumNodes() do
local coord = ac:value(i)
end
-- Indexing nodes in a cell. Regular nodes have an index between 1 and mesh:numNodes()
-- while ghost nodes have a value between 1 and mesh:numGhostNodes() with the highest bit set
local c = mesh:cell(cellId);
for i = 1, c:totalNumNodes() do
local coord = ac:value(c:nodeIndex(i))
end

The catchy part comes when dealing with ghost only accessors mixed with other accessors in a general setting. A loop traversing all available linear node indices (from 1 to mesh:totalNumNodes()) can cause unexpected errors since the linear index can not be used to access the ghost only accessor directly, as it first traverses the geometry part, that has no data for this accessor.

A solution for this case is given by using the ac:adjustLinearIndex() method as seen below.

for i = 1, mesh:totalNumNodes() do
-- Check and adjust i for usage with possibly ghost only accessor ac
-- If i is not a valid index for ac, returns nil. Otherwise returns an adjusted index
local index = ac:adjustLinearIndex(i, mesh)
if index then
local val = ac:value(index)
end
end

Keep in mind that in general, when working with a single known accessor, the problem above is a no issue since one knows the kind of indexing that the accessor needs. Generic code that can possibly handle with mixed affectedNodes types, on the other hand, should follow the guides above.