OOF2: The Manual

Name

Flux — Base class for Fluxes

Synopses

Only those methods useful when extending OOF2 are listed here.

C++ Synopsis

#include "engine/flux.h"
class Flux {
  const std::string& name() const;
  int ndof() const;
  int divergence_dim() const;
  void negate();
  Components* components(Planarity planarity) const;
  Components* divergenceComponents() const;
  Components* outOfPlaneComponents() const;
  virtual FieldIndex* getIndex(const std::string& name) const;
  ArithmeticOutputValue output(const FEMesh* mesh,
                               const Element* element,
                               const MasterPosition& pos) const;

}

Python Synopsis

from ooflib.SWIG.engine import Flux
class Flux:
  def name(self)
  def ndof(self)
  def divergence_dim(self)
  def negate(self)
  def components(self, planarity)
  def divergenceComponents(self)
  def outOfPlaneComponents(self)
  def getIndex(self, name)

Source Files

  • SRC/engine/flux.h: C++ headers
  • SRC/engine/flux.C: C++ source code
  • SRC/engine/flux.swg: SWIG source code
  • SRC/engine/flux.spy: python code included in flux.swg

Description

Flux is the base class for all Flux objects. Like Fields and Equations, Fluxes are global objects. There is only one stress flux object, even though stress may be computed on many different meshes. Flux objects store information about the physical flux, but do not store its values.

Fluxes should only be created by calling the Python constructors of the derived classes. As with Fields and Equations references to Flux objects are kept in the main OOF2 namespace in a variable whose name is the name passed to the derived class constructor. See Section 8.3. It is also possible to retrieve a Flux by name using the getFlux function.

Methods

std::string &name() const

This returns the name that was assigned to this Flux when it was created. See Section 8.3.

int ndof() const

ndof() returns the number of floating point numbers required to represent a value of the Flux.

int divergence_dim() const

divergence_dim() returns the number of in-plane components in the Flux's divergence. This is equal to the number of Equations required to solve the Flux's divergence equation.

void negate()

OOF2 attempts to be generic where possible, making it easier to maintain and extend. However, sometimes tradition gets in the way. For example, stress is given by an elastic modulus times a derivative (strain):

σij=Cijklεkl
(179)

but heat flux 𝐉 is the negative of the thermal conductivity modulus κ times a derivative (temperature gradient):

𝐉i=κijTxj.
(180)

and the electric displacement (aka electric flux density) 𝐃 also contains a minus sign:

𝐃i=ϵijφxj.
(181)

OOF2 needs to know when to insert the minus signs into the divergence equations (2.8) and (2.9). This is done by calling Flux::negate() after creating a Flux. If negate() is called, then the Flux will be treated like Stress. If it is not called, it will be treated like the Heat FLux.

const Components* components(Planarity planarity) const

components() returns a Components pointer that can be used to loop over and refer to the components of the Flux. The given Planarity determines whether the iteration should include the in-plane components, out-of-plane components, or both.

Components* divergenceComponents() const

divergenceComponents() returns a Components pointer that can be used to loop over the in-plane components of the divergence of the Flux.

Components* outOfPlaneComponents() const

outOfPlaneComponents() returns a Components pointer that can be used to loop over the out-of-plane components of the Flux.

[Caution] Possible Confusion

The difference between Flux::components(OUT_OF_PLANE) and Flux::outOfPlaneComponents() needs some explanation. They both refer to the same components of the Flux, in the same order. But the first is used to locate the out-of-plane components in a list that contains all of the components and the second locates them in a list that contains only the out-of-plane components. The difference is illustrated by using the FieldIndex::integer() method to examine the list indices of the components.

For example, here are all the components of the Stress. (SymTensorIndex is the FieldIndex subclass for symmetric tensor indices):

>>> list(Stress.components())
[SymTensorIndex(0,0), SymTensorIndex(1,1), SymTensorIndex(2,2), SymTensorIndex(1,2), SymTensorIndex(0,2), SymTensorIndex(0,1)] 

FieldIndex::integer converts a FieldIndex into the equivalent integer position in an array or list:

>>> list(ij.integer() for ij in Stress.components())
[0, 1, 2, 3, 4, 5] 

components(OUT_OF_PLANE) simply returns the out-of-plane parts of this set of indices, with the same integer positions:

>>> list(Stress.components(OUT_OF_PLANE))
[SymTensorIndex(2,2), SymTensorIndex(1,2), SymTensorIndex(0,2)]
>>> list(ij.integer() for ij in Stress.components(OUT_OF_PLANE))
[2, 3, 4] 

outOfPlaneComponents return the same indices, using a different FieldIndex subclass, and the integer indices are appropriate for a list or array containing only the out-of-plane components:

>>> list(Stress.outOfPlaneComponents())
[OutOfPlaneSymTensorIndex(2, 2), OutOfPlaneSymTensorIndex(1, 2), OutOfPlaneSymTensorIndex(0, 2)]
>>> list(ij.integer() for ij in Stress.outOfPlaneComponents())
[0, 1, 2] 

A similar issue arises with Fields.

FieldIndex *getIndex(const std::string& name) const

getIndex() returns a FieldIndex object, given the name of the desired index. Different Flux subclasses expect different strings, and return different subclasses of FieldIndex. For example, a SymmetricTensorFlux expects the name to be "xx", "xy", "yz", etc, and returns a SymTensorIndex.

In C++ the returned value is a pointer to a newly allocated object. The caller is responsible for deleting it. An easy way to do that is to wrap it in an IndexP:

Flux *flux = ...;
IndexP index = IndexP(flux->getIndex("xx")); 

ArithmeticOutputValue output(const FEMesh*, const Element*, const MasterPosition&) const

output() is similar to Field::output() — it returns the value of the Flux at a given point. It differs in that it is computed at an arbitrary point within an Element, specified by a MasterPosition, instead of at a Node.[69]

The return value is an OutputValue, which wraps the appropriate OutputVal derived class.

The arguments are:

const FEMesh* mesh

The finite element mesh on which the Flux is to be computed.

const Element* element

The element containing the point at which the Flux is to be computed.

const MasterPosition& pos

The position of the output point, specified as a point in the element's master coordinate space. Master space coordinates can be converted to physical coordinates by Element::from_master().



[69] This is because Fluxes depend on Field gradients, which are computed within Elements, whereas Field values are stored at Nodes.