OOF2: The Manual

Name

FluxData — Store data during flux matrix computation

Synopses

C++ Synopsis

#include "engine/flux.h"
class FluxData {
  double& matrix_element(const FieldIndex& fluxcomponent,
                         const Field* field,
                         const FieldIndex& fieldcomponent,
                         const ElementFuncNodeIterator& node);

  double& matrix_element(const FieldIndex& fluxcomponent,
                         const Field* field,
                         const ElementFuncNodeIterator& node);

  double offset_element(const FEMesh* mesh,
                        const FieldIndex& fluxcomponent) const;

  double& offset_element(const FEMesh* mesh,
                         const FieldIndex& fluxcomponent);

  double& rhs_element(int );
}

Python Synopsis

from oof2.SWIG.engine import flux
class FluxData:
  def add_matrix_element(self, fluxIndex, field, fieldIndex, node, increment)
  def add_offset_element(self, fluxIndex, increment)
  def add_rhs_element(self, index, increment)

Source Files

  • SRC/engine/flux.h: C++ headers
  • SRC/engine/flux.C: C++ source code
  • SRC/engine/flux.swg: SWIG source code

Description

FluxData objects are temporary containers used to accumulate data when building the finite element stiffness matrix and right-hand-side vector. They're relevant to the authors of OOF2 extensions because they're passed in to Property::fluxmatrix() and Property::fluxrhs(). Those functions must use FluxData methods to make their contributions to the flux matrix.

Methods

Only methods useful to the authors of OOF2 extensions are listed here. Note that the C++ and Python functions are distinct — the Python equivalents of the C++ functions have different names and argument lists, so they're documented separately (unlike elsewhere in the manual).

double& matrix_element(const FieldIndex&, const Field*, const FieldIndex&, const ElementNodeIterator&)

matrix_element returns a reference to an element of the flux matrix. Recall that the fluxmatrix is the matrix that gives the components of a Flux when multiplied by a vector of degrees of freedom of an Element. Each call to Property::fluxmatrix (which must be written by the extension author) adds the contribution from one Node's shapefunction.

The arguments to matrix_element are:

const FieldIndex& fluxIndex

This is the component of the Flux that's being computed. The calling program should loop over all values of fluxIndex and call FluxData::matrix_element for each one. It should use Flux::iterator() to do the iteration. (The IteratorP returned by Flux::iterator will automatically be converted to a FieldIndex&.)

const Field*

The Field whose contribution to the Flux is being computed. This is half of the information needed to specify the degree of freedom.

const FieldIndex& fieldIndex

This specifies the component of the Field whose contribution to the Flux is being computed. It's the other half of the information required for specifying the degree of freedom. The FieldIndex should be obtained from Field::iterator().

const ElementFuncNodeIterator& node

This indicates the Node whose shapefunction is being used. It's given in the form of a ElementFuncNodeIterator because that class is more useful than Node in the fluxmatrix context.

An alternative version of matrix_element omits the fieldIndex argument and can be used for scalar Fields.

See Property::fluxmatrix for an example of the use of matrix_element.

double& offset_element(const FieldIndex&)

Properties that make Field-independent contributions to a Flux must make those contributions by calling offset_element from Property::fluxrhs. offset_element returns a reference to a double, to which fluxrhs must add a Property-dependent value.

The argument to offset_element is a FieldIndex specifying which Flux component is being computed. The FieldIndex can be obtained from Flux::iterator().

See Property::fluxrhs for an example of the use of offset_element.

double& rhs_element(int)

rhs_element is similar to offset_element, but it's used when computing a direct contribution to the right hand side of a divergence equation, rather than a direct contribution to a flux. The argument is the component of the equation. For an example, see Property:fluxrhs.

add_matrix_element(fluxIndex, field, fieldIndex, node, increment)

This is the Python equivalent of matrix_element. add_matrix_element adds the given increment to the flux matrix element specified by the first five arguments, which are interpreted as they are in matrix_element.

add_offset_element(fluxIndex, increment)

add_offset_element is the Python equivalent of offset_element. It adds the given increment to the Flux component specified by fluxIndex.

add_rhs_element(index, increment)

add_rhs_element is the Python equivalent of rhs_element. It adds the given increment to the divergence equation component specified by index.