FELiCS.Fields.Field#
Module Contents#
Classes#
Class representing a finite element field. |
Attributes#
- FELiCS.Fields.Field.logger#
- class FELiCS.Fields.Field.Field(FEMSpace, mesh, name=None, isStateVector=False, m=None)#
Class representing a finite element field.
This class provides methods for handling finite element fields, including coefficient manipulation, boundary condition application, and expression evaluation. It supports complex-valued fields and can optionally represent a state vector in a mixed formulation. The class also supports an optional spectral dimension via a wave number.
The class interacts with
dolfinx.fem.Functionfor finite element operations and includes utilities for working with PETSc vectors and UFL expressions.Initialize the Field object
- Parameters:
FEMSpace (
dolfinx.fem.FunctionSpace) – The finite element function space.mesh (
FELiCS.SpaceDisc.FELiCSMesh) – FELiCS mesh associated with the function space.name (
str, optional) – Name of the field. If not provided, a default name is chosen based on the field type (scalar, vector, or mixed).isStateVector (
bool, optional) – If True, the field is interpreted as a state vector in a mixed space, and subfield names may be taken fromspace.stateVectorNameswhere available. Default is False.m (
int, optional) – Wave number associated with a spectral spatial dimension. If this is set (even to zero), the field is assumed to have one spectral spatial dimension, andhasSpectralDimensionis set to True.
- property name#
- getNamesOfSubFields()#
Return the list of subfield names for vector or mixed spaces.
For scalar fields, an empty list is returned and a warning is logged. For vector fields, names are generated from the base field name and the mesh axis names (e.g.
u_x,u_y). For mixed fields, names are taken from previously set values, fromstateVectorNames(whenisStateVectoris True), or reasonable defaults such asscalar1,vector1, etc.- Returns:
List of subfield names. Empty for pure scalar fields.
- Return type:
listofstr
- setNamesOfSubFields(nameList)#
Assign explicit names to the subfields of a mixed or vector field.
- Parameters:
nameList (
listofstr) – List of subfield names. The length must match the number of subspaces in the underlying function space; otherwise, a warning is logged and the names are not changed.
- getTensor()#
Wrap the field as a tensor-aware object.
- Returns:
Tensor wrapper around the underlying finite element function, constructed with the mesh coordinate system and the field’s spectral settings (
mandhasSpectralDimension).- Return type:
Notes
Mixed functions are not yet fully supported and are treated as a single tensor object.
- isReal()#
Check whether the field is (numerically) real-valued.
- Returns:
True if the imaginary part of the coefficient array has zero norm, False otherwise.
- Return type:
bool
- getListOfSubFields()#
Get a list of single-component fields.
- Returns:
A list containing individual scalar fields if the function space has multiple subspaces. Otherwise, returns a list containing only this field.
- Return type:
list
Notes
If the space has multiple subspaces, each subspace is extracted as an individual field.
- setListOfSubFields(listOfFields, name=[])#
Set the field coefficients from a list of single-component fields.
- Parameters:
listOfFields (
list) – A list of Field objects representing individual subspaces.name (
string, optional) – The name of the new field that contains all the given fields as subfields.
Notes
This only works if the “listOfFields” contains fields with the correct spaces, which are subspaces of the space which which this field has been initialized.
If the space has multiple subspaces, their coefficients are mapped back.
Throws an error if the input list does not match the expected size.
- describeFunctionSpace()#
Describe the structure of the function space and its subspaces.
- Returns:
A dictionary containing: - ‘type’: str - ‘scalar’, ‘vector’, or ‘mixed’ - ‘num_subspaces’: int - Number of subspaces - ‘subspaces’: list - List of subspace descriptions for mixed spaces - ‘value_size’: int - Total number of components - ‘description’: str - Human-readable description
- Return type:
dict
Examples
For a scalar field: {‘type’: ‘scalar’, ‘num_subspaces’: 0, ‘value_size’: 1, ‘description’: ‘Single scalar space’}
For a vector field: {‘type’: ‘vector’, ‘num_subspaces’: 3, ‘value_size’: 3, ‘description’: ‘Single vector space with 3 components’}
For a mixed field: {‘type’: ‘mixed’, ‘num_subspaces’: 2, ‘subspaces’: […], ‘value_size’: 4, ‘description’: ‘Mixed space with 2 subspaces’}
- getSize()#
Get the length of the coefficient array from the underlying function.
- Returns:
length of coefficient array
- Return type:
integer
- getCoefficientArray()#
Get the coefficient array of the field.
- Returns:
A complex-valued array representing the field coefficients.
- Return type:
numpy.ndarray
- getRealCoefficientArray()#
Get the real part of the coefficient array of the field.
- Returns:
A float-valued array representing the field coefficients.
- Return type:
numpy.ndarray
- getImagCoefficientArray()#
Get the imaginary part of the coefficient array of the field.
- Returns:
A float-valued array representing the field coefficients.
- Return type:
numpy.ndarray
- setCoefficientArray(array)#
Set the coefficient array of the field.
- Parameters:
array (
numpy.ndarray) – A complex-valued array of coefficients.
- setConstantValue(value)#
Set all coefficients to a constant value.
- Parameters:
value (
complexorfloat) – The constant value to assign to all coefficients.
- getPetscVector()#
Convert the field to a PETSc vector.
- Returns:
A PETSc vector created from the coefficient array.
- Return type:
petsc4py.PETSc.Vec
- conjugate()#
Compute the complex conjugate of the field.
- setBoundaryConditions(bcs)#
Apply boundary conditions to the field.
- Parameters:
bcs (
list) – A list of Dirichlet boundary conditions to be applied.
- getGradientField()#
Compute the gradient of a scalar field as a new Field.
The gradient is obtained in a vector-valued function space with dimension equal to the geometric dimension of the mesh (plus an additional spectral dimension if present). The result is computed via a tensor-based weak form and projection.
- Returns:
A new Field instance representing the gradient of the original field.
- Return type:
- calculateL2Norm()#
Compute the L2 norm of the field (and its subfields, if mixed).
For mixed or vector-valued fields, the norm is computed by summing the contributions of each scalar subfield in the mixed decomposition.
- Returns:
The L2 norm of the field.
- Return type:
float
- getVorticityField()#
Compute the vorticity field associated with a velocity field.
For a 2D velocity field, this returns a scalar vorticity field (ω_z = ∂v/∂x − ∂u/∂y). For a 3D velocity field, a vector-valued vorticity field is constructed component-wise using the curl of the velocity.
- Returns:
Vorticity field as a scalar (2D) or vector (3D) Field. Returns None if the number of components is less than 2.
- Return type:
FieldorNone
Notes
3D behaviour is marked as experimental in the implementation and may not be fully tested.
No explicit error is raised if the field is not a velocity-type vector; it is the caller’s responsibility to ensure consistency.
- exportToH5(writer, fileName=None)#
Export the field to an HDF5 file using a FELiCS writer.
- Parameters:
writer (
object) – Writer object providing anexportFieldToH5(field, fileName)method.fileName (
str, optional) – Base name for the exported dataset or file. If None, the field’snameattribute is used.
- importData(reader, importFilePath, groupName=None)#
Import data into the field using a FELiCS reader.
This is a convenience wrapper around the reader’s
importInField()method.- Parameters:
reader (
object) – Reader object providing animportInField(field, filePath, groupName)method.importFilePath (
str) – Path to the file containing the data to be imported.groupName (
strorNone, optional) – Optional group name inside the file from which to read.
- Returns:
Field– The updated field (self).listofstr– List of variable names that were not found in the file.
- evaluateUflExpression(ufl_expression, bcs=[], restartSolver=False)#
Evaluate a UFL expression and update the field accordingly.
- Parameters:
ufl_expression (
ufl.Form) – The UFL expression to evaluate.bcs (
list, optional) – A list of boundary conditions to apply.restartSolver (
bool, optional) – Whether to restart the solver instead of reusing an existing one.
- evaluateUflTensorExpression(ufl_expression, bcs=[], restartSolver=False)#
Evaluate a tensor-based UFL expression and update the field accordingly.
- Parameters:
ufl_expression (
ufl.Form) – The UFL expression to be evaluated.bcs (
list, optional) – A list of boundary conditions to be applied.restartSolver (
bool, optional) – Whether to restart the solver instead of reusing an existing one.
Notes
This method assembles and solves a tensor-based weak form.
Uses a predefined solver if available to improve performance.
The weak form includes integration over the computational domain.
- smoothUflTensorExpression(ufl_expression, smoothFactor, bcs=[], restartSolver=False)#
Smooth a tensor-based UFL expression using a diffusion-like approach.
- Parameters:
ufl_expression (
ufl.Form) – The UFL expression to be smoothed.smoothFactor (
float) – A smoothing factor controlling the influence of the gradient term.bcs (
list, optional) – A list of boundary conditions to apply.restartSolver (
bool, optional) – Whether to restart the solver instead of reusing an existing one.
Notes
This method applies a smoothing operation by adding a gradient-based regularization term to the weak form.
It is particularly useful for regularizing noisy numerical solutions.
- smooth(smoothFactor, bcs=[], restartSolver=False)#
Smooth the field using a diffusion-like approach.
This method constructs a smoothing operator that combines an identity term with a gradient-based regularization term, and then applies it to the current field by solving the corresponding weak form.
- Parameters:
smoothFactor (
float) – Smoothing factor controlling the influence of the gradient term.bcs (
list, optional) – List of boundary conditions to apply.restartSolver (
bool, optional) – If True, rebuild the smoothing solver even if one already exists.
Notes
Particularly useful for regularizing noisy numerical solutions.
The same solver is reused between calls unless
restartSolveris True.
- plot(xlim=None, ylim=None, imaginaryPart=False)#
Plotting function for debugging purposes. This function can be used, to check if a field looks as expected and rule out e.g. import problems.
Notes
This method provides a simple visualization of the field.