capytaine.bem.engines module¶
Definition of the methods to build influence matrices, using possibly some sparse structures.
- class capytaine.bem.engines.BasicMatrixEngine(*, green_function=None, linear_solver='lu_decomposition')[source]¶
Bases:
MatrixEngineDefault matrix engine.
- Features:
Caching of the last computed matrices.
Supports plane symmetries and nested plane symmetries.
Linear solver can be customized. Default is lu_decomposition with caching of the LU decomposition.
- Parameters:
green_function (AbstractGreenFunction) – the low level implementation used to compute the coefficients of the matrices.
linear_solver (str or function, optional) – Setting of the numerical solver for linear problems Ax = b. It can be set with the name of a preexisting solver (available: “lu_decomposition”, “lu_decompositon_with_overwrite” and “gmres”, the former is the default choice) or by passing directly a solver function.
- build_S_matrix(mesh1, mesh2, free_surface, water_depth, wavenumber) ndarray[source]¶
Similar to
build_matrices, but returning only \(S\)
- build_fullK_matrix(mesh1, mesh2, free_surface, water_depth, wavenumber) ndarray[source]¶
Similar to
build_matrices, but returning only full \(K\) (that is the three components of the gradient, not just the normal one)
- build_matrices(mesh1, mesh2, free_surface, water_depth, wavenumber, adjoint_double_layer=True)[source]¶
Build the influence matrices between mesh1 and mesh2.
- Parameters:
mesh1 (MeshLike or list of points) – mesh of the receiving body (where the potential is measured)
mesh2 (MeshLike) – mesh of the source body (over which the source distribution is integrated)
free_surface (float) – position of the free surface (default: \(z = 0\))
water_depth (float) – position of the sea bottom (default: \(z = -\infty\))
wavenumber (float) – wavenumber (default: 1.0)
adjoint_double_layer (bool, optional) – compute double layer for direct method (F) or adjoint double layer for indirect method (T) matrices (default: True)
- Returns:
the matrices \(S\) and \(K\)
- Return type:
tuple of matrix-like (Numpy arrays or BlockCirculantMatrix)
- green_function: AbstractGreenFunction¶
- last_computed_matrices: Tuple[ndarray | BlockDiagonalMatrix | BlockCirculantMatrix, ndarray | BlockDiagonalMatrix | BlockCirculantMatrix | LUDecomposedMatrix | LUDecomposedBlockDiagonalMatrix | LUDecomposedBlockCirculantMatrix] | None¶
- linear_solver(A: ndarray | BlockDiagonalMatrix | BlockCirculantMatrix | LUDecomposedMatrix | LUDecomposedBlockDiagonalMatrix | LUDecomposedBlockCirculantMatrix, b: ndarray) ndarray[source]¶
Solve a linear system with left-hand side A and right-hand-side b
- Parameters:
A (matrix-like) – Expected to be the second output of build_matrices
b (np.ndarray) – Vector of the correct length
- Returns:
x – Vector such that A@x = b
- Return type:
np.ndarray