Setting up a problem¶
Defining a test matrix with xarray
¶
The shortest way to set the parameters of the potential flow problems to be solved is to build a xarray dataset.
Then the function fill_dataset
will automatically make all the simulations defined in the test matrix.
import numpy as np
import xarray as xr
import capytaine as cpt
body = ... # Set up the body and its dofs here.
test_matrix = xr.Dataset(coords={
'omega': np.linspace(0.1, 4, 40),
'wave_direction': [0, np.pi/2],
'radiating_dof': list(body.dofs),
'water_depth': [np.inf],
})
dataset = cpt.BEMSolver().fill_dataset(test_matrix, body)
The first argument of fill_dataset
is an xarray.Dataset
specifying the cases to be run. When a parameter is not specified in the
dataset, the default value is used (see next section).
The second argument is either a single FloatingBody
or a list of
FloatingBody
. In the latter case, each body in the list is studied
independently and the output dataset contains one dimension more.
fill_dataset
returns an xarray.Dataset
with the same
coordinates as its input but filled with additional output data.
See export_outputs for methods to export this dataset in various
formats.
If the coordinate theta
is added to the test matrix, the code will
compute the Kochin function for these values of \(\theta\).
The LinearPotentialFlowProblem
class¶
For a finer grain control of the problems to be solved, a list of LinearPotentialFlowProblem
should be defined.
It is defined as, e.g.:
problem = cpt.RadiationProblem(body=my_body, omega=1.0, radiating_dof="Heave")
other_problem = cpt.DiffractionProblem(body=my_body, omega=1.0, wave_direction=pi/2)
Besides the body, all the parameters are optional. The table below gives their definitions and their default values.
Parameter |
Description (unit) |
Default value |
---|---|---|
|
Position of the free surface [1] (m) |
\(0.0\) m |
|
Constant depth of water (m) |
\(\infty\) m |
|
Acceleration of gravity \(g\) (m/s²) |
\(9.81\) m/s² |
|
Water density (kg/m³) |
\(1000.0\) kg/m³ |
|
Direction of the incoming waves (rad) (for diffraction or forward speed) |
\(0.0\) rad [2] |
|
Name of radiating dof (only for radiation) |
first one found |
|
Speed of the body in the \(x\) direction (m/s) |
\(0.0\) m/s |
Warning
Unlike other software such as Nemoh, the wave direction in Capytaine is expressed in radians.
The wave height is implicitly assumed to be \(1\) m. Since all computations are linear, any wave height or motion amplitude can be retrieved by multiplying the result by the desired value.
Setting the frequency is done by passing one and only one of the following magnitude.
Parameter |
Description (unit) |
---|---|
|
Angular frequency \(\omega\) (rad/s) |
|
Frequency \(f\) (Hz) |
|
Period \(T = \frac{2\pi}{\omega}\) (s) |
|
Wavelength \(\lambda\) (m) |
|
Angular wavenumber \(k = \frac{2\pi}{\lambda}\) (rad/m) |
If no frequency is provided, a frequency omega = 1.0
rad/s is used by default.
Once the problem has been initialized, the other parameters can be retrieved as:
problem.omega
problem.freq
problem.period
problem.wavelength
problem.wavenumber
When forward speed is non zero, the encounter frequency is computed and can be retrieved as:
problem.encounter_omega
problem.encounter_freq
problem.encounter_period
problem.encounter_wavelength
problem.encounter_wavenumber
For radiation problems, setting the frequency to zero or infinity is possible. Simply pass the value 0.0 or float(‘inf’) to one of the above magnitude. The default Green function supports zero and infinite frequency in infinite depth, and infinite frequency in finite depth (but not zero frequency in finite depth).
Importing a dataset from Bemio¶
A DataFrame or a Dataset can also be created from data structures generated using the Bemio package, which reads hydrodynamic output data from NEMOH, WAMIT, and AQWA. This allows for Capytaine post-processing of hydrodynamic data generated from other BEM codes.
Bemio does not come packaged with Capytaine and needs to to be installed independently. Note that the base repository of Bemio has been archived and is only compatible with Python 2.7.x, so using a Python 3 compatible fork is recommended, available here or installed with:
pip install git+https://github.com/mancellin/bemio.git
To build the xarray dataset using Capytaine, the output files from the BEM program in
question must be read into a Bemio data_structures.ben.HydrodynamicData
class, which is
then called by assemble_dataframe or assemble_dataset. For example, to
create an xarray dataset from a WAMIT .out
file:
from bemio.io.wamit import read as read_wamit
import capytaine as cpt
bemio_data = read_wamit("myfile.out")
my_dataset = cpt.assemble_dataset(bemio_data, hydrostatics=False)
Warning
The created dataset will only contain quantities that can be directly calculated
from the values given in the original dataset. Because of this, there may be minor
differences between the variable names in an xarray dataset build with Bemio and one created
using LinearPotentialFlowResult
, even though the format will be identical. For
example, WAMIT .out
files do not contain the radii of gyration needed to calculate
the moments of inertia, so the my_dataset[‘inertia_matrix’] variable would not be included
in the above example since the rigid body mass matrix cannot be calculated.
Legacy Nemoh.cal parameters files¶
The legacy parameters files from Nemoh can be read by a dedicated function:
from capytaine.io.legacy import import_cal_file
list_of_problems = import_cal_file("path/to/Nemoh.cal")
The function returns a list of LinearPotentialFlowProblems
.
Warning
This feature is experimental.
Some of the settings in the files (such as the free surface computation or the Kochin function) are ignored for the moment.
See the example Nemoh.cal
below.
--- Environment ---
1000 ! RHO ! KG/M**3 ! Fluid specific volume
9.81 ! G ! M/S**2 ! Gravity
0 ! DEPTH ! M ! Water depth
0 0 ! Wave measurement point !! This line is ignored !!
--- Description of floating bodies ---
1 ! Number of bodies
--- Body 1 ---
Cylinder.dat ! Name of mesh file
540 300 ! Number of vertices and panels in mesh !! This line is ignored !!
6 ! Number of degrees of freedom
1 1 0 0 0 0 0 ! Surge
1 0 1 0 0 0 0 ! Sway
1 0 0 1 0 0 0 ! Heave
2 1 0 0 0 0 -7.5 ! Roll about CdG
2 0 1 0 0 0 -7.5 ! Pitch about CdG
2 0 0 1 0 0 -7.5 ! Yaw about CdG
6 ! Number of resulting generalised forces
1 1 0 0 0 0 0 ! Force in x direction !! This line is ignored !!
1 0 1 0 0 0 0 ! Force in y direction !! This line is ignored !!
1 0 0 1 0 0 0 ! Force in z direction !! This line is ignored !!
2 1 0 0 0 0 -7.5 ! Moment force in x direction about CdG !! This line is ignored !!
2 0 1 0 0 0 -7.5 ! Moment force in y direction about CdG !! This line is ignored !!
2 0 0 1 0 0 -7.5 ! Moment force in z direction about CdG !! This line is ignored !!
0 ! Number of lines of additional information !! All the additional information is ignored !!
--- Load cases to be solved ---
1 2 0.1 2.0 ! Freq type 1,2,3=[rad/s,Hz,s], Number of wave frequencies/periods, Min, and Max
1 0.0 0.0 ! Number of wave directions, Min and Max (degrees)
--- Post processing ---
!! All lines below are ignored by Capytaine !!
Note
The line setting up the frequencies is slightly different in Nemoh v2 and Nemoh v3. Both format are supported by Capytaine.
Command-line interface¶
Warning
This feature is experimental.
Capytaine comes with a command-line command capytaine
which can be used as:
$ capytaine path/to/directory/parameter.cal
The parameter file (in Nemoh.cal
format) passed as argument is read and legacy tecplot output file are written in the directory path/to/directory/results/
.
Warning
If results files already exist, they will be overwritten!
If no argument is provided to the command, the code looks for a file Nemoh.cal
in the current directory.