Field#

A field container with pre-defined fields is created with:

FieldsMixed(region[, n, values, ...])

A mixed field based on a region and returns a FieldContainer instance.

A Field-container holds one or more fields.

FieldContainer(fields)

A container for fields based on a list or tuple of Field instances.

Available kinds of fields:

Field(region[, dim, values])

A Field on points of a region with dimension dim and initial point values.

FieldAxisymmetric(region[, dim, values])

A axisymmetric Field on points of a two-dimensional region with dimension dim (default is 2) and initial point values (default is 0).

FieldPlaneStrain(region[, dim, values])

A plane strain Field on points of a two-dimensional region with dimension dim (default is 2) and initial point values (default is 0).

Detailed API Reference

class felupe.FieldsMixed(region, n=3, values=(0, 0, 1, 0), axisymmetric=False, planestrain=False, offset=0, npoints=None, mesh=None, **kwargs)[source]#

A mixed field based on a region and returns a FieldContainer instance.

copy()#

Return a copy of the field.

extract(grad=True, sym=False, add_identity=True)#

Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.

Parameters:
  • grad (bool, optional (default is True)) – Flag for gradient evaluation.

  • sym (bool, optional (default is False)) – Flag for symmetric part if the gradient is evaluated.

  • add_identity (bool, optional (default is True)) – Flag for the addition of the identity matrix if the gradient is evaluated.

Returns:

(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.

Return type:

array

Link value array of other field.

plot(*args, **kwargs)#

Plot the first field of the container.

See also

felupe.Scene.plot

Plot method of a scene.

screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)#

Take a screenshot of the first field of the container.

See also

pyvista.Plotter.screenshot

Take a screenshot of a PyVista plotter.

values()#

Return the field values.

view(point_data=None, cell_data=None, cell_type=None)#

View the field with optional given dicts of point- and cell-data items.

Parameters:
  • point_data (dict or None, optional) – Additional point-data dict (default is None).

  • cell_data (dict or None, optional) – Additional cell-data dict (default is None).

  • cell_type (pyvista.CellType or None, optional) – Cell-type of PyVista (default is None).

Returns:

A object which provides visualization methods for felupe.FieldContainer.

Return type:

felupe.ViewField

See also

felupe.ViewField

Visualization methods for felupe.FieldContainer.

class felupe.FieldContainer(fields)[source]#

A container for fields based on a list or tuple of Field instances.

copy()[source]#

Return a copy of the field.

extract(grad=True, sym=False, add_identity=True)[source]#

Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.

Parameters:
  • grad (bool, optional (default is True)) – Flag for gradient evaluation.

  • sym (bool, optional (default is False)) – Flag for symmetric part if the gradient is evaluated.

  • add_identity (bool, optional (default is True)) – Flag for the addition of the identity matrix if the gradient is evaluated.

Returns:

(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.

Return type:

array

Link value array of other field.

plot(*args, **kwargs)[source]#

Plot the first field of the container.

See also

felupe.Scene.plot

Plot method of a scene.

screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)[source]#

Take a screenshot of the first field of the container.

See also

pyvista.Plotter.screenshot

Take a screenshot of a PyVista plotter.

values()[source]#

Return the field values.

view(point_data=None, cell_data=None, cell_type=None)[source]#

View the field with optional given dicts of point- and cell-data items.

Parameters:
  • point_data (dict or None, optional) – Additional point-data dict (default is None).

  • cell_data (dict or None, optional) – Additional cell-data dict (default is None).

  • cell_type (pyvista.CellType or None, optional) – Cell-type of PyVista (default is None).

Returns:

A object which provides visualization methods for felupe.FieldContainer.

Return type:

felupe.ViewField

See also

felupe.ViewField

Visualization methods for felupe.FieldContainer.

class felupe.Field(region, dim=1, values=0, **kwargs)[source]#

A Field on points of a region with dimension dim and initial point values. A slice of this field directly accesses the field values as 1d-array.

The interpolation method returns the field values evaluated at the numeric integration points q of all cells c in the region.

\[u^i_{(pc)} = \hat{u}_a^i h_{a(pc)}\]

The gradient method returns the gradient of the field values w.r.t. the undeformed mesh point coordinates, evaluated at the integration points of all cells in the region.

\[\left( \frac{\partial u^i}{\partial X^J} \right)_{(pc)} = \hat{u}^i_{a(pc)} \left( \frac{\partial h_a}{\partial X^J} \right)_{(pc)}\]
Parameters:
  • region (Region) – The region on which the field will be created.

  • dim (int (default is 1)) – The dimension of the field.

  • values (float (default is 0.0) or array) – A single value for all components of the field or an array of shape (region.mesh.npoints, dim)`.

  • kwargs (dict, optional) – Optional keyword arguments of the field.

copy()[source]#

Return a copy of the field.

extract(grad=True, sym=False, add_identity=True)[source]#

Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.

Parameters:
  • grad (bool, optional (default is True)) – Flag for gradient evaluation.

  • sym (bool, optional (default is False)) – Flag for symmetric part if the gradient is evaluated.

  • add_identity (bool, optional (default is True)) – Flag for the addition of the identity matrix if the gradient is evaluated.

Returns:

(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.

Return type:

array

fill(a)[source]#

Fill all field values with a scalar value.

grad(sym=False)[source]#

Gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated.

Parameters:

sym (bool, optional (default is False)) – Calculate the symmetric part of the gradient.

Returns:

Gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.

Return type:

array

interpolate()[source]#

Interpolate field values at points and evaluate them at the integration points of all cells in the region.

class felupe.FieldAxisymmetric(region, dim=2, values=0)[source]#

A axisymmetric Field on points of a two-dimensional region with dimension dim (default is 2) and initial point values (default is 0).

  • component 1 = axial component

  • component 2 = radial component

 x_2 (radial direction)

  ^
  |        _
  |       / \
--|-----------------> x_1 (axial rotation axis)
          \_^

This is a modified Field class in which the radial coordinates are evaluated at the numeric integration points. The grad()-method is modified in such a way that it does not only contain the in-plane 2d-gradient but also the circumferential stretch as shown in Eq.(1).

            |  dudX(2d) :   0   |
dudX(axi) = | ..................|                  (1)
            |     0     : u_r/R |
region#

The region on which the field will be created.

Type:

Region

dim#

The dimension of the field.

Type:

int (default is 2)

values#

A single value for all components of the field or an array of shape (region.mesh.npoints, dim)`.

Type:

float (default is 0.0) or array

copy()#

Return a copy of the field.

extract(grad=True, sym=False, add_identity=True)#

Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.

Parameters:
  • grad (bool, optional (default is True)) – Flag for gradient evaluation.

  • sym (bool, optional (default is False)) – Flag for symmetric part if the gradient is evaluated.

  • add_identity (bool, optional (default is True)) – Flag for the addition of the identity matrix if the gradient is evaluated.

Returns:

(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.

Return type:

array

fill(a)#

Fill all field values with a scalar value.

grad(sym=False)[source]#

3D-gradient as partial derivative of field values at points w.r.t. the undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is returned.

            |  dudX(2d) :   0   |
dudX(axi) = | ..................|
            |     0     : u_r/R |
Parameters:

sym (bool, optional (default is False)) – Calculate the symmetric part of the gradient.

Returns:

Full 3D-gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.

Return type:

array

interpolate()[source]#

Interpolate field values at points and evaluate them at the integration points of all cells in the region.

class felupe.FieldPlaneStrain(region, dim=2, values=0)[source]#

A plane strain Field on points of a two-dimensional region with dimension dim (default is 2) and initial point values (default is 0).

This is a modified Field class in which the radial coordinates are evaluated at the numeric integration points. The grad()-method is modified in such a way that it does not only contain the in-plane 2d-gradient but also the circumferential stretch as shown in Eq.(1).

                    |  dudX(2d) :   0   |
dudX(planestrain) = | ..................|                  (1)
                    |     0     :   0   |
region#

The region on which the field will be created.

Type:

Region

dim#

The dimension of the field.

Type:

int (default is 2)

values#

A single value for all components of the field or an array of shape (region.mesh.npoints, dim)`.

Type:

float (default is 0.0) or array

copy()#

Return a copy of the field.

extract(grad=True, sym=False, add_identity=True)#

Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.

Parameters:
  • grad (bool, optional (default is True)) – Flag for gradient evaluation.

  • sym (bool, optional (default is False)) – Flag for symmetric part if the gradient is evaluated.

  • add_identity (bool, optional (default is True)) – Flag for the addition of the identity matrix if the gradient is evaluated.

Returns:

(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.

Return type:

array

fill(a)#

Fill all field values with a scalar value.

grad(sym=False)[source]#

3D-gradient as partial derivative of field values at points w.r.t. the undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is returned.

                    |  dudX(2d) :   0   |
dudX(planestrain) = | ..................|
                    |     0     :   0   |
Parameters:

sym (bool, optional (default is False)) – Calculate the symmetric part of the gradient.

Returns:

Full 3D-gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.

Return type:

array

interpolate()[source]#

Interpolate field values at points and evaluate them at the integration points of all cells in the region.