FELiCS.Equation.EquationCollection#
Module Contents#
Classes#
Constructs variational formulations for equations and manages boundary conditions. |
Attributes#
- FELiCS.Equation.EquationCollection.logger#
- class FELiCS.Equation.EquationCollection.EquationCollectionClass(param, FEMSpaces, mean, mesh)#
Constructs variational formulations for equations and manages boundary conditions.
This class generates variational formulations for various matrices, including: - A (real and imaginary parts of the linear operator) - B (real and imaginary parts of the weight/time-derivative matrix) - B_forcing (resolvent forcing norm) - B_response (resolvent response norm)
It also manages test and trial functions, boundary conditions, and mesh-based coordinate transformations. The design supports multiple physical models (momentum, mass, energy, species) and enables configuration through an external parameter object.
Initialize the EquationCollectionClass object
- Parameters:
param (
object) – Parameters for the simulation, including boundary conditions, equations, and I/O settings.FEMSpaces (
object) – Collection of Finite Element Method spaces for the mixed formulation.mean (
object) – Mean flow fields used in variational formulations.mesh (
object) – Mesh object defining the spatial domain and coordinate system.
- Variables:
x (
ufl.SpatialCoordinate) – Spatial coordinates of the computational domain.ds (
ufl.Measure) – Measure for integration over boundary subdomains.all_ds (
ufl.Measure) – Combined measure over all boundary subdomains.n_BC (
ufl.FacetNormal) – Outward-pointing unit normal on the domain boundary.n (
Tensor) – Tensor form of the boundary normal vector in the chosen coordinate system.BCs (
list) – List of Dirichlet boundary condition objects.trialFunctionsFEM (
ufl.TrialFunctions) – Mixed trial functions used in the variational problem.testFunctionsFEM (
ufl.TestFunctions) – Mixed test functions used in the variational problem.X (
listofTensor) – List of structured test functions with tensor-awareness.R (
ufl.Coefficient) – Radial coordinate value for cylindrical systems or unity for Cartesian.equationList (
list) – List of equation objects that define the weak formulations.resolventResponseIndices (
numpy.ndarray) – Indices selecting solution components used in the resolvent response norm.resolventForcingIndices (
numpy.ndarray) – Indices selecting solution components used in the resolvent forcing norm.
Notes
The class dynamically loads and constructs equation objects depending on the configuration provided in param.Case.Equations. It supports tensor-based formulations and multiple coordinate systems.
Examples
>>> eq_collection = EquationCollectionClass(param, FEMSpaces, mean, mesh) >>> A = eq_collection.getLinearOperator(mean) >>> B = eq_collection.getWeightMatrix(mean)
- getLinearOperator(meanFlow)#
Construct the linear operator matrix for the equation system.
This method assembles the linear part of the weak formulation by iterating through all equations in the equation list and adding their linear expressions.
- Parameters:
meanFlow (
object) – Mean flow properties used in constructing the linear operator- Returns:
Assembled linear operator matrix with boundary conditions applied
- Return type:
petsc4py.PETSc.Mat
- getWeightMatrix(meanFlow)#
Construct the weight matrix representing the time derivative term.
Assembles the matrix by summing contributions from all equations. Boundary conditions are included for resolvent analysis mode only.
- Parameters:
meanFlow (
object) – Mean flow properties used in constructing the weight matrix.- Returns:
Assembled weight matrix.
- Return type:
petsc4py.PETSc.Mat
Notes
For Resolvent analysis, boundary conditions are applied
For other analysis modes, no boundary conditions are applied to avoid computational issues during eigenvalue problem solving
- getFEMWeightMatrix()#
Construct the full Finite Element Method (FEM) weight matrix.
This method creates a weight matrix by performing inner products of test and trial functions across all function spaces in the mixed function space.
- Returns:
Assembled FEM weight matrix without boundary conditions
- Return type:
petsc4py.PETSc.Mat
Notes
Handles both scalar and vector function spaces
Uses conjugate of test functions for matrix construction
Includes a workaround for mesh object compatibility with different versions of DOLFINx
- getFEMDiffusionMatrix(diffusionFactor, sponge=None)#
Construct the Finite Element Method (FEM) diffusion matrix.
This method creates a diffusion matrix by computing the weak form of the diffusion operator across all function spaces in the mixed function space.
- Parameters:
diffusionFactor (
float) – Coefficient multiplying the diffusion termssponge (
floatorNone, optional) – Additional damping term to be applied to the matrix (default is None)
- Returns:
Assembled FEM diffusion matrix with boundary conditions applied
- Return type:
petsc4py.PETSc.Mat
Notes
Handles both scalar and vector function spaces
Computes second-order derivatives in x and y directions
Currently a quick implementation; a more comprehensive tensor framework is needed for future improvements
- getFullRHS(func)#
Construct the full FEM right-hand side (RHS) vector.
Forms the RHS by integrating each test function against the corresponding part of the provided function.
- Parameters:
func (
dolfinx.fem.Function) – Mixed function representing the source term or solution.- Returns:
Assembled RHS vector.
- Return type:
petsc4py.PETSc.Vec
- getNonlinearExpression(meanFlow, setBC=True)#
Construct the nonlinear vector expression for the system.
Assembles the nonlinear contributions from all equations.
- Parameters:
meanFlow (
object) – Mean flow field used in the nonlinear formulation.setBC (
bool, optional) – Whether to apply boundary conditions (default is True).
- Returns:
Assembled nonlinear vector.
- Return type:
petsc4py.PETSc.Vec
- getBilinearOperator(meanFlow)#
Construct the bilinear operator matrix.
Assembles all bilinear terms across the system equations.
- Parameters:
meanFlow (
object) – Mean flow field used in the bilinear formulation.- Returns:
Assembled bilinear operator matrix with boundary conditions applied.
- Return type:
petsc4py.PETSc.Mat
- getForcingForInputOutput(meanFlow)#
Construct the forcing vector for input-output analysis.
Creates ufl object with the linear equation system.
- Parameters:
meanFlow (
object) – Mean flow properties used to construct the forcing vector.- Returns:
Assembled complex forcing vector.
- Return type:
petsc4py.PETSc.Vec
- getResolventNorm_response(meanFlow)#
Construct the resolvent response norm matrix.
- Parameters:
meanFlow (
object) – Mean flow field. Currently not used directly here, but kept for interface consistency with other assembly routines.- Returns:
PETSc matrix representing the response norm, assembled from the UFL form stored in
self.response_vf.- Return type:
petsc4py.PETSc.Mat
- getResolventNorm_forcing(meanFlow)#
Construct the resolvent forcing norm matrix.
- Parameters:
meanFlow (
object) – Mean flow field. Currently not used directly here, but kept for interface consistency with other assembly routines.- Returns:
PETSc matrix representing the forcing norm, assembled from the UFL form stored in
self.forcing_vf.- Return type:
petsc4py.PETSc.Mat
- getResolventWeighting_FEM(meanFlow)#
Assemble the FEM weighting matrix for resolvent analysis.
The underlying UFL form for the weighting is constructed beforehand by
computeResolventFEMWeights()and stored inself.fem_weighting. This method only assembles the corresponding PETSc matrix.- Parameters:
meanFlow (
object) – Mean flow field. Currently not used directly here, but included for API compatibility with other assembly routines.- Returns:
PETSc matrix with FEM weights based on the previously defined variable projections.
- Return type:
petsc4py.PETSc.Mat
- computeResolventNorms(X, param, mean, fluc)#
Compute energy norms for resolvent forcing and response.
- Parameters:
X (
listofTensor) – Test functions with tensor awareness.param (
object) – Simulation configuration parameters.mean (
object) – Mean flow fields.fluc (
object) – Fluctuation fields object.
Notes
Only ‘TKE’ and ‘Chu’ norm types are currently implemented.
- computeResolventFEMWeights(X, param, mean, fluc)#
Compute FEM weighting matrix for resolvent input/output scaling.
- Parameters:
X (
listofTensor) – Test functions.param (
object) – Parameter configuration.mean (
object) – Mean flow fields.fluc (
object) – Fluctuation fields.
- getRestrictorMatResponse()#
Construct the response restrictor matrix for the resolvent analysis.
This matrix is diagonal and applies spatial restriction (if defined) to the response vector based on the responseDomain field in the mean flow.
- Returns:
PETSc matrix with diagonal entries corresponding to the spatial restriction mask.
- Return type:
petsc4py.PETSc.Mat
Notes
Provides a quadratic matrix, with the size of the solution space (VMixed).
Has the response restrictor values, given with the mean field, on the diagonal.
- getRestrictorMatForcing()#
Construct the forcing restrictor matrix for the resolvent analysis.
This matrix is diagonal and applies spatial restriction (if defined) to the forcing vector based on the
forcingDomainfield in the mean flow. If the forcing domain is zero everywhere, the matrix reduces to the identity.- Returns:
PETSc matrix with diagonal entries corresponding to the spatial restriction mask.
- Return type:
petsc4py.PETSc.Mat
Notes
The matrix is quadratic with the size of the solution space (
VMixed).The diagonal entries are given by the forcing restrictor values provided with the mean field (or ones everywhere if no spatial restriction is prescribed).
- getShrinkerMatResponse()#
Construct a shrinker matrix for the resolvent response vector.
- Returns:
PETSc matrix projecting full response to a reduced subspace.
- Return type:
petsc4py.PETSc.Mat
- getShrinkerMatForcing()#
Creates a (possibly rectangular) matrix which serves the purpose to “shrink” the forcing vector to the requested size (e.g. only consider the velocity components when using the forcing norm “TKE”).
- Returns:
PETSc matrix projecting full forcing vector to reduced subspace.
- Return type:
petsc4py.PETSc.Mat