pyforce.tools package

Submodules

pyforce.tools.backends module

class pyforce.tools.backends.LoopProgress(msg: str, final: float = 100)[source]

Bases: object

A class to make progress bar.

Parameters:
  • msg (str) – Message to be displayed

  • final (float, optional (Default = 100)) – Maximum value for the iterations

update(step: float, percentage: bool = False)[source]

Update message to display and clears the previous one.

Parameters:
  • step (float) – Interger or float value to add at the counter.

  • percentage (boolean, optional (Default = False)) – Indicates if the bar should be displayed in %.

class pyforce.tools.backends.norms(V: dolfinx.fem.FunctionSpace, is_H1=False, metadata_degree=4)[source]

Bases: object

A class to compute norms and inner products. \(L^2\) and \(H^1\) (semi and full are implemented for both scalar and vector fields), whereas the average and the integral are available for scalar only.

Parameters:
  • V (FunctionSpace) – Functional Space onto which the Function are defined.

  • is_H1 (boolean, optional (Default = False)) – If the function belongs to \(H^1\), the forms for the inner products and norms are computed.

H1innerProd(u: dolfinx.fem.Function, v: dolfinx.fem.Function, semi=True)[source]

Computes the \(H^1\) semi or full inner product of the functions u and v over the domain

\[\langle u, v \,\rangle_{H^1} = \int_\Omega \nabla u \cdot \nabla v\,d\Omega\]
\[(u,v)_{H^1} = \int_\Omega u\cdot v \,d\Omega + \int_\Omega \nabla u\cdot \nabla v \,d\Omega\]
Parameters:
  • u (Function) – Function belonging to the same functional space V

  • v (Function) – Function belonging to the same functional space V

  • semi (boolean, optional (Default = True)) – Indicates if the semi norm must be computed.

Returns:

value\(H^1\) inner product of the functions

Return type:

float

H1norm(u: dolfinx.fem.Function, semi=True)[source]

Computes the \(H^1\) semi or full norm of the function u over the domain

\[| u |_{H^1} = \sqrt{\int_\Omega \nabla u \cdot \nabla u\,d\Omega}\]
\[\| u \|_{H^1} = \sqrt{\int_\Omega \nabla u \cdot \nabla u\,d\Omega + \int_\Omega u \cdot u\,d\Omega}\]
Parameters:
  • u (Function) – Function belonging to the same functional space V

  • semi (boolean, optional (Default = True)) – Indicates if the semi norm must be computed.

Returns:

value\(H^1\) norm of the function

Return type:

float

L2innerProd(u: dolfinx.fem.Function, v: dolfinx.fem.Function)[source]

Computes the \(L^2\) inner product of the functions u and v over the domain

\[(u,v)_{L^2}=\int_\Omega u\cdot v \,d\Omega\]
Parameters:
  • u (Function) – Function belonging to the same functional space V

  • v (Function) – Function belonging to the same functional space V

Returns:

value\(L^2\) inner product between the functions

Return type:

float

L2norm(u: dolfinx.fem.Function)[source]

Computes the \(L^2\) norm of the function u over the domain

\[\| u\|_{L^2} = \sqrt{\int_\Omega u \cdot u\,d\Omega}\]
Parameters:

u (Function) – Function belonging to the same functional space V

Returns:

value\(L^2\) norm of the function

Return type:

float

Linftynorm(u: dolfinx.fem.Function)[source]

Computes the \(L^\infty\) norm of a given function u over the domain

\[\| u \|_{L^\infty}=\max\limits_\Omega |u|\]
Parameters:

u (Function) – Function belonging to the same functional space V

Returns:

value\(L^\infty\) norm of the function

Return type:

float

average(u: dolfinx.fem.Function)[source]

Computes the integral average of a given scalar function u over the domain

\[\langle u \rangle = \frac{1}{|\Omega|}\int_\Omega u \,d\Omega\]
Parameters:

u (Function) – Function belonging to the same functional space V (it must be a scalar!)

Returns:

ave_value – Average over the domain

Return type:

float

integral(u: dolfinx.fem.Function)[source]

Computes the integral of a given scalar function u over the domain

\[\int_\Omega u \,d\Omega\]
Parameters:

u (Function) – Function belonging to the same functional space V (it must be a scalar!)

Returns:

val – Integral over the domain

Return type:

float

pyforce.tools.functions_list module

class pyforce.tools.functions_list.FunctionsList(function_space: dolfinx.fem.FunctionSpace | None = None, dofs: int | None = None)[source]

Bases: object

A class wrapping a list of functions. They are stored as a list of np.ndarray.

Parameters:
  • function_space (FunctionSpace, optional) – Functional Space onto which the Function are defined. If not provided, dofs must be specified.

  • dofs (int, optional) – Degrees of freedom of the functions \(\mathcal{N}_h\). Required if function_space is not provided.

append(dofs_fun: ndarray) None[source]

Extend the current list by adding a new function. The dolfinx.fem.Function element is stored in a list as a numpy array, to avoid problems when the number of elements becomes large. The input can be either a np.ndarray or a dolfinx.fem.Function, in the latter case it is mapped to np.ndarray.

Parameters:

dofs_fun (np.ndarray) – Functions to be appended.

clear() None[source]

Clear the storage.

copy()[source]

Defining the copy of the _list of elements

delete(idx: int) None[source]

Delete a single element in position idx.

Parameters:

idx (int) – Integers indicating the position inside the _list to delete.

lin_combine(vec: ndarray, use_numpy=True) ndarray[source]

Linearly combine functions (a_i = vec[i]) in the list \(\{\phi_i\}_{i=0}^N\):

\[\sum_{i=0}^N a_i\cdot \phi_i\]

given N = len(vec) and \(\mathbf{a}\in\mathbb{R}^N\).

Parameters:
  • vec (np.ndarray) – Iterable containing the coefficients of the linear combination.

  • use_numpy (boolean, optional (Default=True)) – If True the functions are treated as np.ndarray, otherwise the formulation in dolfinx is used.

Returns:

combination – Function object storing the result of the linear combination

Return type:

np.ndarray or Function

map(idx: int) dolfinx.fem.Function[source]

Mapping the element in position idx into a dolfinx.fem.Function.

Parameters:

idx (int) – Integers indicating the position inside the _list.

Returns:

eval_fun – Evaluated dofs into a Function

Return type:

Function

return_matrix() ndarray[source]

Returns the list of arrays as a matrix \(S\in\mathbb{R}^{\mathcal{N}_h\times N_s}\).

sort(order: list) None[source]

Sorting the list according to order - iterable of indices.

Parameters:

order (list) – List of indices for the sorting.

pyforce.tools.functions_list.train_test_split(params: list, fun_list: FunctionsList, test_size: float = 0.33, random_state: int = 42)[source]

This function can be used to create a train and test using sklearn capabilities.

Parameters:
  • params (list) – List of parameters to be split.

  • fun_list (FunctionsList) – Object containing the functions as a list of arrays to convert.

  • test_size (float) – DimensionFunctional space of the functions (the dofs should be compliant).

  • random_state (int, optional (Default = 42)) – Random seed for the splitting algorithm.

Returns:

  • train_params (list) – List of the train parameters.

  • test_params (list) – List of the test parameters.

  • train_fun (list) – List of the train functions.

  • test_fun (list) – List of the test functions.

pyforce.tools.plotting module

pyforce.tools.plotting.grids(fun: dolfinx.fem.Function, varname='u', log_plot: bool = False, mag_plot: bool = True, **kwargs)[source]

Creates a PyVista grid from a dolfinx.fem.Function and returns the warped or glyph representation.

Parameters:
  • fun (dolfinx.fem.Function) – The function representing the field to be visualized.

  • varname (str) – The name to assign to the data in the PyVista grid (default is ‘u’).

  • log_plot (bool) – If True, apply a logarithmic plot to scalar data (default is False).

  • mag_plot (bool) – If True, creates a vector warp of the grid. If False, uses glyphs (default is True).

  • **kwargs – Additional keyword arguments passed to the function.

Returns:

grid – The resulting PyVista grid, which can be visualized using PyVista plotting functions.

Return type:

pyvista.UnstructuredGrid

pyforce.tools.plotting.plot_function(fun: dolfinx.fem.Function, filename: str | None = None, format: str = 'png', varname: str = 'u', clim=None, colormap=<matplotlib.colors.LinearSegmentedColormap object>, resolution=[1080, 720], zoom=1.0, title=None, **kwargs)[source]

Python function to plot a scalar field.

Parameters:
  • fun (Function) – Field to plot.

  • varname (str, optional (Default = 'u')) – Name of the variable.

  • filename (str (Default = None)) – If None, the plot is shown; otherwise this is the name of the file to save.

  • clim (optional (Default = None)) – Colorbar limit, if None the mininum and maximum of fun are computed.

  • colormap (optional (Default = jet)) – Colormap for the plot

  • resolution (list, optional (Default = [1080, 720])) – Resolution of the image

  • zoom (float (Default = 1)) – Zoom of the plot.

  • title (str (Default = None)) – If None no title is displayed.

  • **kwargs – Additional keyword arguments passed to the function.

pyforce.tools.timer module

class pyforce.tools.timer.Timer[source]

Bases: object

start()[source]

Start a new timer

stop()[source]

Stop the timer, and report the elapsed time

exception pyforce.tools.timer.TimerError[source]

Bases: Exception

A custom exception used to report errors in use of Timer class

pyforce.tools.write_read module

pyforce.tools.write_read.ImportH5(V: dolfinx.fem.FunctionSpace, filename: str, var_name: str, verbose=False)[source]

This function can be used to load from xdmf/h5 files scalar or vector list of functions.

Parameters:
  • V (FunctionSpace) – Functional space in which the functions live.

  • varname (str) – Name of the variable.

  • filename (str) – Name of the file to read as xdmf/h5 file.

  • verbose (boolean, (Default = False)) – If True, printing is enabled

Returns:

  • snap (FunctionsList) – Imported list of functions

  • params_np (list) – Imported list of parameters

class pyforce.tools.write_read.ImportOF(path: str, extract_dofs=False)[source]

Bases: object

A class used to import data from OpenFOAM and interpolate them onto dolfinx mesh.

The fluidfoam library (see https://github.com/fluiddyn/fluidfoam) is exploited for the import process
  • Supported OpenFoam Versions : 2.4.0, 4.1 to 9 (org), v1712plus to v2212plus (com)

Parameters:
  • path (str) – Path of the OpenFOAM case directory.

  • extract_dofs (boolean, optional (Default=False)) – If true, the dofs \(x,\,y,\,z\) of the OpenFOAM mesh are imported.

foam_to_dolfinx(V: dolfinx.fem.FunctionSpace, snaps: list, variables: list, cut_value=0.0, verbose=None)[source]

Converting the OpenFOAM data imported into FunctionsList based on dolfinx structures. The snapshots are projected onto a suitable functional space.

Parameters:
  • V (FunctionSpace) – Functional Space onto which data are projected.

  • snaps (list) – List of numpy.ndarray of size (ndofs, ndim) with the OpenFOAM snapshots.

  • variables (list) – Lisr of spatial variables, either 2D (e.g., [‘x’, ‘y’], [‘x’, ‘z’], [‘y’, ‘z’]) or 3D (e.g., [‘x’, ‘y’, ‘z’]).

  • cut_value (float, (Default = 0.)) – For 2D meshes, one dimension is cut and this is the associated value.

  • verbose (str, (Default = None)) – If not None, printing is enabled with input message.

Returns:

snap_dolfinx – Project list of functions in dolfinx.

Return type:

FunctionsList

import_field(var_name: str, vector=False, verbose=True)[source]

Importing all time instances (skipping zero folder) from OpenFOAM directory.

Parameters:
  • var_name (str) – Name of the field to import.

  • vector (boolean, (Default: False)) – Labelling if the field is a vector.

  • verbose (boolean, (Default = True)) – If True, printing is enabled

Returns:

  • field (list) – Imported list of functions (each element is a numpy.ndarray), sorted in time.

  • time_instants (list) – Sorted list of time.

pyforce.tools.write_read.StoreFunctionsList(domain, snap: FunctionsList, var_name: str, filename: str, order=None)[source]

This function can be used to save in xdmf/h5 files scalar or vector list of functions.

Parameters:
  • domain (dolfinx.mesh) – Mesh containing the geometry.

  • snap (FunctionsList) – List of functions to save.

  • varname (str) – Name of the variable.

  • filename (str) – Name of the file to save as xdmf/h5 file.

  • order (optional (Default = None)) – It must be passed as a list of integers containing the ordered indeces.

Module contents

pyforce/tools.

Basic tools for the pyforce library.