API reference

Contents

API reference#

Parameters#

Tools to manipulate the physical and computational parameters. It contains variables and methods designed to be a link between XCompact3d and Python applications for pre and post-processing.

class xcompact3d_toolbox.parameters.Parameters(**kwargs: Any)#

Bases: ParametersBasicParam, ParametersNumOptions, ParametersInOutParam, ParametersScalarParam, ParametersLESModel, ParametersIbmStuff, ParametersALMParam, ParametersExtras

The physical and computational parameters are built on top of traitlets. It is a framework that lets Python classes have attributes with type checking, dynamically calculated default values, and “on change” callbacks. In this way, many of the parameters are validated regarding the type, business rules, and the range of values supported by XCompact3d. There are methods to handle the parameters file (.i3d and .prm).

The parameters are arranged in different classes, but just for organization purposes, this class inherits from all of them.

In addition, there are ipywidgets for a friendly user interface, see xcompact3d_toolbox.gui.ParametersGui.

After that, it is time to read the binary arrays produced by XCompact3d and also to write a new xdmf file, so the binary fields can be opened in any external visualization tool. See more details in xcompact3d_toolbox.parameters.ParametersExtras.dataset.

Note

This is a work in progress, not all parameters are covered yet.

__init__(*, raise_warning: bool = False, **kwargs)#

Initializes the Parameters Class.

Parameters:
  • raise_warning (bool, optional) – Raise a warning instead of an error if an invalid parameter is found. By default False.

  • **kwargs – Keyword arguments for valid attributes, like nx, re and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Examples

There are a few ways to initialize the class.

First, calling it with no arguments initializes all variables with default value:

>>> prm = xcompact3d_toolbox.Parameters()

It is possible to set any value afterwards:

>>> prm.re = 1e6
>>> prm.set(
...     iibm=0,
...     p_row=4,
...     p_col=2,
... )

Second, we can specify some values, and let the missing ones be initialized with default value:

>>> prm = x3d.Parameters(
...     filename="example.i3d",
...     itype=12,
...     nx=257,
...     ny=129,
...     nz=32,
...     xlx=15.0,
...     yly=10.0,
...     zlz=3.0,
...     nclx1=2,
...     nclxn=2,
...     ncly1=1,
...     nclyn=1,
...     nclz1=0,
...     nclzn=0,
...     re=300.0,
...     init_noise=0.0125,
...     dt=0.0025,
...     ilast=45000,
...     ioutput=200,
...     iprocessing=50,
... )

And finally, it is possible to read the parameters from the disc:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="example.i3d")

It also supports the previous parameters file format (see #7):

>>> prm = xcompact3d_toolbox.Parameters(loadfile="incompact3d.prm")

Changed in version 1.2.0: The argument raise_warning changed to keyword-only.

__repr__()#

Return repr(self).

__str__()#

Representation of the parameters class, similar to the representation of the .i3d file.

from_file(filename: str | None = None, *, raise_warning: bool = False) None#

Loads the attributes from the parameters file.

It also includes support for the previous format prm (see #7).

Parameters:
  • filename (str, optional) – The filename for the parameters file. If None, it uses the filename specified in the class (default is None).

  • raise_warning (bool, optional) – Raise a warning instead of an error if an invalid parameter is found. By default False.

Raises:
  • IOError – If the file format is not i3d or prm.

  • KeyError – Exception is raised when an attributes is invalid.

Examples

>>> prm = xcompact3d_toolbox.Parameters(filename="example.i3d")
>>> prm.load()
>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.load("example.i3d")

or just:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="example.i3d")
>>> prm = xcompact3d_toolbox.Parameters(loadfile="incompact3d.prm")
from_string(string: str, *, raise_warning: bool = False) None#

Loads the attributes from a string.

Parameters:
  • filename (str, optional) – The filename for the parameters file. If None, it uses the filename specified in the class (default is None).

  • raise_warning (bool, optional) – Raise a warning instead of an error if an invalid parameter is found. By default False.

Raises:

KeyError – Exception is raised when an attributes is invalid.

get_boundary_condition(variable_name: str) dict#

This method returns the appropriate boundary parameters that are expected by the derivatives methods.

Parameters:

variable_name (str) – Variable name. The supported options are ux, uy, uz, pp and phi, otherwise the method returns the default option.

Returns:

A dict containing the boundary conditions for the variable specified.

Return type:

dict of dict

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.get_boundary_condition("ux")
{'x': {'ncl1': 1, 'ncln': 1, 'npaire': 0},
'y': {'ncl1': 1, 'ncln': 2, 'npaire': 1, 'istret': 0, 'beta': 0.75},
'z': {'ncl1': 0, 'ncln': 0, 'npaire': 1}}

It is possible to store this information as an attribute in any xarray.DataArray:

>>> DataArray.attrs["BC"] = prm.get_boundary_condition("ux")

So the correct boundary conditions will be used to compute the derivatives:

>>> DataArray.x3d.first_derivative("x")
>>> DataArray.x3d.second_derivative("x")
get_mesh(*, refined_for_ibm: bool = False) dict#

Get mesh the three-dimensional coordinate system. The coordinates are stored in a dictionary. It supports mesh refinement in y when ParametersBasicParam.istret \(\ne\) 0.

Parameters:

refined_for_ibm (bool) – If True, it returns a refined mesh as a function of ParametersIbmStuff.nraf (default is False).

Returns:

It contains the mesh point locations at three dictionary keys, for x, y and z.

Return type:

dict of numpy.ndarray

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.get_mesh()
{'x': array([0.    , 0.0625, 0.125 , 0.1875, 0.25  , 0.3125, 0.375 , 0.4375,
        0.5   , 0.5625, 0.625 , 0.6875, 0.75  , 0.8125, 0.875 , 0.9375,
        1.    ]),
 'y': array([0.    , 0.0625, 0.125 , 0.1875, 0.25  , 0.3125, 0.375 , 0.4375,
        0.5   , 0.5625, 0.625 , 0.6875, 0.75  , 0.8125, 0.875 , 0.9375,
        1.    ]),
 'z': array([0.    , 0.0625, 0.125 , 0.1875, 0.25  , 0.3125, 0.375 , 0.4375,
        0.5   , 0.5625, 0.625 , 0.6875, 0.75  , 0.8125, 0.875 , 0.9375,
        1.    ])}
load(*arg, **kwarg) None#

An alias for Parameters.from_file

set(*, raise_warning: bool = False, **kwargs) None#

Set a new value for any parameter after the initialization.

Parameters:
  • raise_warning (bool, optional) – Raise a warning instead of an error if an invalid parameter is found. By default False.

  • **kwargs – Keyword arguments for valid attributes, like nx, re and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.set(
...     iibm=0,
...     p_row=4,
...     p_col=2,
... )

Changed in version 1.2.0: The argument raise_warning changed to keyword-only.

write(filename: str | None = None) None#

Write all valid attributes to an i3d file.

An attribute is considered valid if it has a tag named group, witch assigns it to the respective namespace at the i3d file.

Parameters:

filename (str, optional) – The filename for the i3d file. If None, it uses the filename specified in the class (default is None).

Examples

>>> prm = xcompact3d_toolbox.Parameters(
...     filename="example.i3d",
...     nx=101,
...     ny=65,
...     nz=11,
...     # and so on...
... )
>>> prm.write()

or just:

>>> prm.write("example.i3d")
class xcompact3d_toolbox.parameters.ParametersALMParam(**kwargs: Any)#

Bases: HasTraits

__init__()#
class xcompact3d_toolbox.parameters.ParametersBasicParam(**kwargs: Any)#

Bases: HasTraits

__init__()#
beta#

Refinement factor in y.

Notes

Only necessary if istret \(\ne\) 0.

Type:

float

dt#

Time step \((\Delta t)\).

Type:

float

gravx#

Component of the unitary vector pointing in the gravity’s direction.

Type:

float

gravy#

Component of the unitary vector pointing in the gravity’s direction.

Type:

float

gravz#

Component of the unitary vector pointing in the gravity’s direction.

Type:

float

ifirst#

The number for the first iteration.

Type:

int

iibm#

Enables Immersed Boundary Method (IBM):

  • 0 - Off (default);

  • 1 - On with direct forcing method, i.e., it sets velocity to zero inside the solid body;

  • 2 - On with alternating forcing method, i.e, it uses Lagrangian Interpolators to define the velocity inside the body and imposes no-slip condition at the solid/fluid interface;

  • 3 - Cubic Spline Reconstruction;

Any option greater than zero activates the namespace ibmstuff, for variables like ParametersIbmStuff.nobjmax and ParametersIbmStuff.nraf.

Type:

int

iin#

Defines perturbation at the initial condition:

  • 0 - No random noise (default);

  • 1 - Random noise with amplitude of init_noise;

  • 2 - Random noise with fixed seed (important for reproducibility, development and debugging) and amplitude of init_noise;

  • 3 - Read inflow planes.

Notes

The exactly behavior may be different according to each flow configuration.

Type:

int

ilast#

The number for the last iteration.

Type:

int

ilesmod#

Enables Large-Eddy methodologies:

Type:

int

inflow_noise#

Random number amplitude at inflow boundary (where \(x=0\)).

Notes

Only necessary if nclx1 is equal to 2.

Type:

float

init_noise#

Random number amplitude at initial condition.

Notes

Only necessary if iin \(\ne\) 0. The exactly behavior may be different according to each flow configuration.

Type:

float

ipost#

Enables online postprocessing at a frequency ParametersInOutParam.iprocessing:

  • 0 - No;

  • 1 - Yes (default).

Note

The computation for each case is specified at the BC.<flow-configuration>.f90 file.

Type:

int

istret#

Controls mesh refinement in y:

  • 0 - No refinement (default);

  • 1 - Refinement at the center;

  • 2 - Both sides;

  • 3 - Just near the bottom.

Notes

See beta.

Type:

int

itype#

Sets the flow configuration, each one is specified in a different BC.<flow-configuration>.f90 file (see Xcompact3d/src), they are:

  • 0 - User configuration;

  • 1 - Turbidity Current in Lock-Release;

  • 2 - Taylor-Green Vortex;

  • 3 - Periodic Turbulent Channel;

  • 4 - Periodic Hill

  • 5 - Flow around a Cylinder;

  • 6 - Debug Schemes (for developers);

  • 7 - Mixing Layer;

  • 8 - Jet;

  • 9 - Turbulent Boundary Layer;

  • 10 - ABL;

  • 11 - Uniform;

  • 12 - Sandbox.

Type:

int

ivisu#

Enables store snapshots at a frequency ParametersInOutParam.ioutput:

  • 0 - No;

  • 1 - Yes (default).

Type:

int

nclx1#

Boundary condition for velocity field where \(x=0\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - Inflow.

Type:

int

nclxn#

Boundary condition for velocity field where \(x=L_x\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - Convective outflow.

Type:

int

ncly1#

Boundary condition for velocity field where \(y=0\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - No-slip.

Type:

int

nclyn#

Boundary condition for velocity field where \(y=L_y\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - No-slip.

Type:

int

nclz1#

Boundary condition for velocity field where \(z=0\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - No-slip.

Type:

int

nclzn#

Boundary condition for velocity field where \(z=L_z\), the options are:

  • 0 - Periodic;

  • 1 - Free-slip;

  • 2 - No-slip.

Type:

int

numscalar#

Number of scalar fraction, which can have different properties. Any option greater than zero activates the namespace ParametersScalarParam.

Type:

int

nx#

Number of mesh points.

Notes

See xcompact3d_toolbox.mesh.Coordinate.possible_grid_size for recommended grid sizes.

Type:

int

ny#

Number of mesh points.

Notes

See xcompact3d_toolbox.mesh.Coordinate.possible_grid_size for recommended grid sizes.

Type:

int

nz#

Number of mesh points.

Notes

See xcompact3d_toolbox.mesh.Coordinate.possible_grid_size for recommended grid sizes.

Type:

int

p_col#

Defines the domain decomposition for (large-scale) parallel computation.

Notes

The product p_row * p_col must be equal to the number of computational cores where XCompact3d will run. More information can be found at 2DECOMP&FFT.

p_row = p_col = 0 activates auto-tunning.

Type:

int

p_row#

Defines the domain decomposition for (large-scale) parallel computation.

Notes

The product p_row * p_col must be equal to the number of computational cores where XCompact3d will run. More information can be found at 2DECOMP&FFT.

p_row = p_col = 0 activates auto-tunning.

Type:

int

re#

Reynolds number \((Re)\).

Type:

float

xlx#

Domain size.

Type:

float

yly#

Domain size.

Type:

float

zlz#

Domain size.

Type:

float

class xcompact3d_toolbox.parameters.ParametersExtras(**kwargs: Any)#

Bases: HasTraits

Extra utilities that are not present at the parameters file, but are useful for Python applications.

__init__()#
dataset#

An object that reads and writes the raw binary files from XCompact3d on-demand.

Notes

All arrays are wrapped into Xarray objects (xarray.DataArray or xarray.Dataset), take a look at xarray’s documentation, specially, see Why xarray? Xarray has many useful methods for indexing, comparisons, reshaping and reorganizing, computations and plotting.

Consider using hvPlot to explore your data interactively, see how to plot Gridded Data.

Examples

The first step is specify the filename properties. If the simulated fields are named like ux-000.bin, they are in the default configuration, there is no need to specify filename properties. But just in case, it would be like:

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator = "-",
...     file_extension = ".bin",
...     number_of_digits = 3
... )

If the simulated fields are named like ux0000, the parameters are:

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator = "",
...     file_extension = "",
...     number_of_digits = 4
... )

Data type is defined by xcompact3d_toolbox.param:

>>> import numpy
>>> xcompact3d_toolbox.param["mytype] = numpy.float64 # if double precision
>>> xcompact3d_toolbox.param["mytype] = numpy.float32 # if single precision

Now it is possible to customize the way the dataset will be handled:

>>> prm.dataset.set(
...     data_path = "./data/",
...     drop_coords = "",
...     set_of_variables = {"ux", "uy", "uz"},
...     snapshot_step = "ioutput",
...     snapshot_counting = "ilast",
...     stack_scalar = True,
...     stack_velocity = False,
... )

Note

For convenience, data_path is set as "./data/" relative to the filename of the parameters file when creating a new instance of Parameters (i.g., if filename = "./example/input.i3d" then data_path = "./example/data/").

There are many ways to load the arrays produced by your numerical simulation, so you can choose what best suits your post-processing application. See the examples:

  • Load one array from the disc:

    >>> ux = prm.dataset.load_array("ux-0000.bin")
    
  • Load the entire time series for a given variable:

    >>> ux = prm.dataset.load_time_series("ux")
    >>> uy = prm.dataset.load_time_series("uy")
    >>> uz = prm.dataset.load_time_series("uz")
    

    or just:

    >>> ux = prm.dataset["ux"]
    >>> uy = prm.dataset["uy"]
    >>> uz = prm.dataset["uz"]
    

    You can organize them using a dataset:

    >>> dataset = xarray.Dataset()
    >>> for var in "ux uy uz".split():
    ...     dataset[var] = prm.dataset[var]
    
  • Load all variables from a given snapshot:

    >>> snapshot = prm.dataset.load_snapshot(10)
    

    or just:

    >>> snapshot = prm.dataset[10]
    
  • Loop through all snapshots, loading them one by one:

    >>> for ds in prm.dataset:
    ...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
    ...     prm.dataset.write(data = vort, file_prefix = "w3")
    
  • Loop through some snapshots, loading them one by one, with the same arguments of a classic Python range, for instance, from 0 to 100 with a step of 5:

    >>> for ds in prm.dataset(0, 101, 5):
    ...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
    ...     prm.dataset.write(data = vort, file_prefix = "w3")
    
  • Or simply load all snapshots at once (if you have enough memory):

    >>> ds = prm.dataset[:]
    

And finally, it is possible to produce a new xdmf file, so all data can be visualized on any external tool:

>>> prm.dataset.write_xdmf()
Type:

xcompact3d_toolbox.io.Dataset

dx#

Mesh resolution.

Type:

float

dy#

Mesh resolution.

Type:

float

dz#

Mesh resolution.

Type:

float

filename#

Filename for the parameters file.

Type:

str

mesh#

Mesh object.

Type:

xcompact3d_toolbox.mesh.Mesh3D

ncores#

Number of computational cores where XCompact3d will run.

Type:

int

size#

Auxiliary variable indicating the demand for storage.

Type:

str

class xcompact3d_toolbox.parameters.ParametersIbmStuff(**kwargs: Any)#

Bases: HasTraits

__init__()#
iforces#

Enables forces calculation.

Type:

int

izap#

How many points to skip for reconstruction ranging from 0 to 3, the recommended is 1.

Type:

int

nobjmax#

Maximum number of objects in any direction. It is defined automatically at gene_epsi_3d.

Type:

int

npif#

Number of Points for the reconstruction ranging from 1 to 3, the recommended is 2.

Type:

int

nraf#

Level of refinement to find the surface of the immersed object, when ParametersBasicParam.iibm is equal to 2.

Type:

int

class xcompact3d_toolbox.parameters.ParametersInOutParam(**kwargs: Any)#

Bases: HasTraits

__init__()#
icheckpoint#

Frequency for writing restart file.

Type:

int

ioutput#

Frequency for visualization (3D snapshots).

Type:

int

iprocessing#

Frequency for online postprocessing.

Type:

int

irestart#

Reads initial flow field if equals to 1.

Type:

int

nvisu#

Size for visual collection.

Type:

int

class xcompact3d_toolbox.parameters.ParametersLESModel(**kwargs: Any)#

Bases: HasTraits

__init__()#
iwall#

int:

jles#

Chooses LES model, they are:

  • 0 - No model (DNS);

  • 1 - Phys Smag;

  • 2 - Phys WALE;

  • 3 - Phys dyn. Smag;

  • 4 - iSVV.

Type:

int

maxdsmagcst#

float:

nSmag#

float:

smagcst#

float:

smagwalldamp#

int:

walecst#

float:

class xcompact3d_toolbox.parameters.ParametersNumOptions(**kwargs: Any)#

Bases: HasTraits

__init__()#
cnu#

Ratio between hyperviscosity at \(k_m=2/3\pi\) and \(k_c= \pi\).

Type:

float

ifirstder#

Scheme for first order derivative:

  • 1 - 2nd central;

  • 2 - 4th central;

  • 3 - 4th compact;

  • 4 - 6th compact (default).

Type:

int

iimplicit#

Time integration scheme:

  • 0 - Off (default);

  • 1 - With backward Euler for Y diffusion;

  • 2 - With CN for Y diffusion.

Type:

int

isecondder#

Scheme for second order derivative:

  • 1 - 2nd central;

  • 2 - 4th central;

  • 3 - 4th compact;

  • 4 - 6th compact (default);

  • 5 - Hyperviscous 6th.

Type:

int

itimescheme#

Time integration scheme:

  • 1 - Forwards Euler;

  • 2 - Adams-bashforth 2;

  • 3 - Adams-bashforth 3 (default);

  • 4 - Adams-bashforth 4 (Not Implemented);

  • 5 - Runge-kutta 3;

  • 6 - Runge-kutta 4 (Not Implemented).

Type:

int

nu0nu#

Ratio between hyperviscosity/viscosity at nu.

Type:

float

class xcompact3d_toolbox.parameters.ParametersScalarParam(**kwargs: Any)#

Bases: HasTraits

__init__()#
cp#

Initial concentration(s).

Type:

list of float

nclxS1#

Boundary condition for scalar field(s) where \(x=0\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Inflow.

Type:

int

nclxSn#

Boundary condition for scalar field(s) where \(x=L_x\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Convective outflow.

Type:

int

nclyS1#

Boundary condition for scalar field(s) where \(y=0\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Dirichlet.

Type:

int

nclySn#

Boundary condition for scalar field(s) where \(y=L_y\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Dirichlet.

Type:

int

nclzS1#

Boundary condition for scalar field(s) where \(z=0\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Dirichlet.

Type:

int

nclzSn#

Boundary condition for scalar field(s) where \(z=L_z\), the options are:

  • 0 - Periodic;

  • 1 - No-flux;

  • 2 - Dirichlet.

Type:

int

ri#

Richardson number(s).

Type:

list of float

sc#

Schmidt number(s).

Type:

list of float

scalar_lbound#

Lower scalar bound(s), for clipping methodology.

Type:

list of float

scalar_ubound#

Upper scalar bound(s), for clipping methodology.

Type:

list of float

uset#

Settling velocity(s).

Type:

list of float

Parameters - GUI#

Manipulate the physical and computational parameters, just like xcompact3d_toolbox.parameters.Parameters, but with ipywidgets.

class xcompact3d_toolbox.gui.ParametersGui(**kwargs: Any)#

Bases: Parameters

This class is derivated from xcompact3d_toolbox.parameters.Parameters, including all its features. In addition, there is a two way link between the parameters and their widgets. Control them with code and/or with the graphical user interface.

__call__(*args: str) VBox#

Returns widgets on demand.

Parameters:

*args (str) – Name(s) for the desired widget(s).

Returns:

Widgets for an user friendly interface.

Return type:

ipywidgets.VBox

Examples

>>> prm = xcompact3d_toolbox.ParametersGui()
>>> prm("nx", "xlx", "dx", "nclx1", "nclxn")
__init__(**kwargs)#

Initializes the Parameters Class.

Parameters:

**kwargs – Keyword arguments for xcompact3d_toolbox.parameters.Parameters.

Creates a two-way link between the value of an attribute and its widget. This method is called at initialization, but provides an easy way to link any new variable.

Examples

>>> prm = xcompact3d_toolbox.ParametersGui(loadfile="example.i3d")
>>> prm.link_widgets()

Mesh#

Objects to handle the coordinates and coordinate system. Note they are an attribute at xcompact3d_toolbox.parameters.ParametersExtras, so they work together with all the other parameters. They are presented here for reference.

class xcompact3d_toolbox.mesh.Coordinate(**kwargs: Any)#

Bases: HasTraits

A coordinate.

Thanks to traitlets, the attributes can be type checked, validated and also trigger “on change” callbacks. It means that:

All these functionalities aim to make a user-friendly interface, where the consistency between different coordinate parameters is ensured even when they change at runtime.

Parameters:
  • length (float) – Length of the coordinate (default is 1.0).

  • grid_size (int) – Number of mesh points (default is 17).

  • delta (float) – Mesh resolution (default is 0.0625).

  • is_periodic (bool) – Specifies if the boundary condition is periodic (True) or not (False) (default is False).

Notes

There is no need to specify both length and delta, because they are a function of each other, the missing value is automatically computed from the other.

Returns:

Coordinate

Return type:

xcompact3d_toolbox.mesh.Coordinate

__array__() ndarray#

This method makes the coordinate automatically work as a numpy like array in any function from numpy.

Returns:

A numpy array.

Return type:

numpy.ndarray

Examples

>>> from xcompact3d_toolbox.mesh import Coordinate
>>> import numpy
>>> coord = Coordinate(length=1.0, grid_size=9)
>>> numpy.sin(coord)
array([0.        , 0.12467473, 0.24740396, 0.36627253, 0.47942554,
       0.58509727, 0.68163876, 0.7675435 , 0.84147098])
>>> numpy.cos(coord)
array([1.        , 0.99219767, 0.96891242, 0.93050762, 0.87758256,
        0.81096312, 0.73168887, 0.64099686, 0.54030231])
__init__(**kwargs)#

Initializes the Coordinate class.

Parameters:

**kwargs – Keyword arguments for attributes, like grid_size, length and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Examples

>>> from xcompact3d_toolbox.mesh import Coordinate
>>> coord = Coordinate(length=1.0, grid_size=9, is_periodic=False)
__len__()#

Make the coordinate work with the Python function len.

Returns:

Coordinate size (grid_size)

Return type:

int

Examples

>>> from xcompact3d_toolbox.mesh import Coordinate
>>> coord = Coordinate(grid_size=9)
>>> len(coord)
9
__repr__()#

Return repr(self).

property possible_grid_size: list#

Possible values for grid size.

Due to restrictions at the FFT library, they must be equal to:

\[n = 2^{1+a} \times 3^b \times 5^c,\]

if the coordinate is periodic, and:

\[n = 2^{1+a} \times 3^b \times 5^c + 1,\]

otherwise, where \(a\), \(b\) and \(c\) are non negative integers.

Additionally, the derivative’s stencil imposes that \(n \ge 8\) if periodic and \(n \ge 9\) otherwise.

Returns:

Possible values for grid size

Return type:

list

Notes

There is no upper limit, as long as the restrictions are satisfied.

Examples

>>> from xcompact3d_toolbox.mesh import Coordinate
>>> coordinate(is_periodic=True).possible_grid_size
[8, 10, 12, 16, 18, 20, 24, 30, 32, 36, 40, 48, 50, 54, 60, 64, 72, 80,
 90, 96, 100, 108, 120, 128, 144, 150, 160, 162, 180, 192, 200, 216, 240,
 250, 256, 270, 288, 300, 320, 324, 360, 384, 400, 432, 450, 480, 486,
 500, 512, 540, 576, 600, 640, 648, 720, 750, 768, 800, 810, 864, 900,
 960, 972, 1000, 1024, 1080, 1152, 1200, 1250, 1280, 1296, 1350, 1440,
 1458, 1500, 1536, 1600, 1620, 1728, 1800, 1920, 1944, 2000, 2048, 2160,
 2250, 2304, 2400, 2430, 2500, 2560, 2592, 2700, 2880, 2916, 3000, 3072,
 3200, 3240, 3456, 3600, 3750, 3840, 3888, 4000, 4050, 4096, 4320, 4374,
 4500, 4608, 4800, 4860, 5000, 5120, 5184, 5400, 5760, 5832, 6000, 6144,
 6250, 6400, 6480, 6750, 6912, 7200, 7290, 7500, 7680, 7776, 8000, 8100,
 8192, 8640, 8748, 9000]
>>> coordinate(is_periodic=False).possible_grid_size
[9, 11, 13, 17, 19, 21, 25, 31, 33, 37, 41, 49, 51, 55, 61, 65, 73, 81,
 91, 97, 101, 109, 121, 129, 145, 151, 161, 163, 181, 193, 201, 217, 241,
 251, 257, 271, 289, 301, 321, 325, 361, 385, 401, 433, 451, 481, 487,
 501, 513, 541, 577, 601, 641, 649, 721, 751, 769, 801, 811, 865, 901,
 961, 973, 1001, 1025, 1081, 1153, 1201, 1251, 1281, 1297, 1351, 1441,
 1459, 1501, 1537, 1601, 1621, 1729, 1801, 1921, 1945, 2001, 2049, 2161,
 2251, 2305, 2401, 2431, 2501, 2561, 2593, 2701, 2881, 2917, 3001, 3073,
 3201, 3241, 3457, 3601, 3751, 3841, 3889, 4001, 4051, 4097, 4321, 4375,
 4501, 4609, 4801, 4861, 5001, 5121, 5185, 5401, 5761, 5833, 6001, 6145,
 6251, 6401, 6481, 6751, 6913, 7201, 7291, 7501, 7681, 7777, 8001, 8101,
 8193, 8641, 8749, 9001]
set(**kwargs) None#

Set a new value for any parameter after the initialization.

Parameters:

**kwargs – Keyword arguments for attributes, like grid_size, length and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Examples

>>> from xcompact3d_toolbox.mesh import Coordinate
>>> coord = Coordinate()
>>> coord.set(length=1.0, grid_size=9, is_periodic=False)
property size: int#

An alias for grid_size.

Returns:

Grid size

Return type:

int

property vector: ndarray#

Construct a vector with numpy.linspace and return it.

Returns:

Numpy array

Return type:

numpy.ndarray

class xcompact3d_toolbox.mesh.Istret(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)#

Bases: IntEnum

Mesh refinement type.

Added in version 1.2.0.

__format__(format_spec, /)#

Convert to a string according to format_spec.

__new__(value)#
class xcompact3d_toolbox.mesh.Mesh3D(**kwargs: Any)#

Bases: HasTraits

A three-dimensional coordinate system

Parameters:

Notes

mesh is in fact an attribute of xcompact3d_toolbox.parameters.ParametersExtras, so there is no need to initialize it manually for most of the common use cases. The features of each coordinate are copled by a two-way link with their corresponding values at the Parameters class. For instance, the length of each of them is copled to xlx, yly and zlz, grid size to nx, ny and nz and so on.

Returns:

Coordinate system

Return type:

xcompact3d_toolbox.mesh.Mesh3D

__init__(**kwargs)#

Initializes the 3DMesh class.

Parameters:

**kwargs – Keyword arguments for each coordinate (x, y and z), containing a dict with the parameters for them, like grid_size, length and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid coordinate.

Examples

>>> from xcompact3d_toolbox.mesh import Mesh3D
>>> mesh = Mesh3D(
...     x=dict(length=4.0, grid_size=65, is_periodic=False),
...     y=dict(length=1.0, grid_size=17, is_periodic=False, istret=0),
...     z=dict(length=1.0, grid_size=16, is_periodic=True),
... )
__len__()#

Make the coordinate work with the Python function len.

Returns:

Mesh size is calculated by multiplying the size of the three coordinates

Return type:

int

__repr__()#

Return repr(self).

copy()#

Return a copy of the Mesh3D object.

drop(*args) dict#

Get the coordinates in a dictionary, where the keys are their names and the values are their vectors. It is possible to drop any of the coordinates in case they are needed to process planes. For instance:

  • Drop x if working with yz planes;

  • Drop y if working with xz planes;

  • Drop z if working with xy planes.

Parameters:

*args (str or list of str) – Name of the coordinate(s) to be dropped

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Returns:

A dict containing the desired coordinates

Return type:

dict of numpy.ndarray

get() dict#

Get the three coordinates in a dictionary, where the keys are their names (x, y and z) and the values are their vectors.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Returns:

A dict containing the coordinates

Return type:

dict of numpy.ndarray

Notes

It is an alias for Mesh3D.drop(None).

set(**kwargs) None#

Set new values for any of the coordinates after the initialization.

Parameters:

**kwargs – Keyword arguments for each coordinate (x, y and z), containing a dict with the parameters for them, like grid_size, length and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid attribute.

Examples

>>> from xcompact3d_toolbox.mesh import Mesh3D
>>> mesh = Mesh3D()
>>> mesh.set(
...     x=dict(length=4.0, grid_size=65, is_periodic=False),
...     y=dict(length=1.0, grid_size=17, is_periodic=False, istret=0),
...     z=dict(length=1.0, grid_size=16, is_periodic=True),
... )
property size#

Mesh size

Returns:

Mesh size is calculated by multiplying the size of the three coordinates

Return type:

int

class xcompact3d_toolbox.mesh.StretchedCoordinate(**kwargs: Any)#

Bases: Coordinate

Another coordinate, as a subclass of Coordinate. It includes parameters and methods to handle stretched coordinates, which is employed by XCompact3d at the vertical dimension y.

Parameters:
  • length (float) – Length of the coordinate (default is 1.0).

  • grid_size (int) – Number of mesh points (default is 17).

  • delta (float) – Mesh resolution (default is 0.0625).

  • is_periodic (bool) – Specifies if the boundary condition is periodic (True) or not (False) (default is False).

  • istret (int) –

    Type of mesh refinement:

    • 0 - No refinement (default);

    • 1 - Refinement at the center;

    • 2 - Both sides;

    • 3 - Just near the bottom.

  • beta (float) – Refinement parameter.

Notes

There is no need to specify both length and delta, because they are a function of each other, the missing value is automatically computed from the other.

Returns:

Stretched coordinate

Return type:

xcompact3d_toolbox.mesh.StretchedCoordinate

__array__()#

This method makes the coordinate automatically work as a numpy like array in any function from numpy.

Returns:

A numpy array.

Return type:

numpy.ndarray

Examples

>>> from xcompact3d_toolbox.mesh import StretchedCoordinate
>>> import numpy
>>> coord = StretchedCoordinate(length=1.0, grid_size=9)
>>> numpy.sin(coord)
array([0.        , 0.12467473, 0.24740396, 0.36627253, 0.47942554,
       0.58509727, 0.68163876, 0.7675435 , 0.84147098])
>>> numpy.cos(coord)
array([1.        , 0.99219767, 0.96891242, 0.93050762, 0.87758256,
        0.81096312, 0.73168887, 0.64099686, 0.54030231])
__repr__()#

Return repr(self).

Reading and writing files#

Useful objects to read and write the binary fields produced by XCompact3d.

class xcompact3d_toolbox.io.Dataset(**kwargs: Any)#

Bases: HasTraits

An object that reads and writes the raw binary files from XCompact3d on-demand.

Parameters:
  • data_path (str) – The path to the folder where the binary fields are located (default is "./data/").

  • note:: (..) – The default "./data/" is relative to the path to the parameters file when initialized from xcompact3d_toolbox.parameters.ParametersExtras.

  • drop_coords (str) – If working with two-dimensional planes, specify which of the coordinates should be dropped, i.e., "x", "y" or "z", or leave it empty for 3D fields (default is "").

  • filename_properties (FilenameProperties) – Specifies filename properties for the binary files, like the separator, file extension and number of digits.

  • set_of_variables (set) – The methods in this class will try to find all variables per snapshot, use this parameter to work with just a few specified variables if you need to speedup your application (default is an empty set).

  • snapshot_counting (str) – The parameter that controls the number of timesteps used to produce the datasets (default is "ilast").

  • snapshot_step (str) – The parameter that controls the number of timesteps between each snapshot, it is often "ioutput" or "iprocessing" (default is "ioutput").

  • stack_scalar (bool) – When True, the scalar fields will be stacked in a new coordinate n, otherwise returns one array per scalar fraction (default is False).

  • stack_velocity (bool) – When True, the velocity will be stacked in a new coordinate i , otherwise returns one array per velocity component (default is False).

Notes

__call__(*args) Iterator[xr.Dataset]#

Yields selected snapshots, so the application can iterate over them, loading one by one, with the same arguments of a classic Python range.

Parameters:

*args (int) – Same arguments used for a range.

Yields:

xarray.Dataset – Dataset containing the arrays loaded from the disc with the appropriate dimensions, coordinates and attributes.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

Iterate over some snapshots, loading them one by one, with the same arguments of a classic Python range, for instance, from 0 to 100 with a step of 5:

>>> for ds in prm.dataset(0, 101, 5):
...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
...     prm.dataset.write(data=vort, file_prefix="w3")
__getitem__(arg: int | slice | str) DataArray | Dataset#

Get specified items from the disc.

Note

Make sure to have enough memory to load many files at the same time.

Parameters:

arg (int or slice or str) –

Specifies the items to load from the disc, depending on the type of the argument:

Returns:

Xarray objects containing values loaded from the disc with the appropriate dimensions, coordinates and attributes.

Return type:

xarray.Dataset or xarray.DataArray

Raises:

TypeError – Raises type error if arg is not an integer, string or slice

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     drop_coords="z",
...     stack_scalar=True,
...     stack_velocity=True,
... )
  • Load the entire time series for a given variable:

    >>> ux = prm.dataset["ux"]
    >>> uy = prm.dataset["uy"]
    >>> uz = prm.dataset["uz"]
    

    or organize them using a dataset:

    >>> dataset = xarray.Dataset()
    >>> for var in "ux uy uz".split():
    ...     dataset[var] = prm.dataset[var]
    
  • Load all variables from a given snapshot:

    >>> snapshot = prm.dataset[10]
    
  • Load many snapshots at once with a slice, for instance, from 0 to 100 with a step of 10:

    >>> snapshots = prm.dataset[0:101:10]
    
  • Or simply load all snapshots at once (if you have enough memory):

    >>> snapshots = prm.dataset[:]
    
__init__(**kwargs)#

Initializes the Dataset class.

Parameters:
  • filename_properties (dict, optional) – Keyword arguments for FilenamePropertie, like separator, file_extension and so on.

  • **kwargs – Keyword arguments for the parameters, like data_path, drop_coords and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid parameter.

Returns:

An object to read and write the raw binary files from XCompact3d on-demand.

Return type:

Dataset

__iter__()#

Yields all the snapshots, so the application can iterate over them.

Yields:

xarray.Dataset – Dataset containing the arrays loaded from the disc with the appropriate dimensions, coordinates and attributes.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

Iterate over all snapshots, loading them one by one:

>>> for ds in prm.dataset:
...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
...     prm.dataset.write(data=vort, file_prefix="w3")
__len__() int#

Make the dataset work with the Python function len.

Returns:

Total of snapshots as a function of snapshot_counting and snapshot_step.

Return type:

int

__repr__()#

Return repr(self).

load_array(filename: str, attrs: dict | None = None, *, add_time: bool = True) type[DataArray]#

This method reads a binary field from XCompact3d with numpy.fromfile and wraps it into a xarray.DataArray with the appropriate dimensions, coordinates and attributes.

Parameters:
  • filename (str) – Name of the file.

  • add_time (bool, optional) – Add time as a coordinate (default is True).

  • attrs (dict_like, optional) – Attributes to assign to the new instance xarray.DataArray.

Returns:

Data array containing values loaded from the disc.

Return type:

xarray.DataArray

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

Load one array from the disc:

>>> ux = prm.dataset.load_array("ux-000.bin")

Changed in version 1.2.0: The argument add_time changed to keyword-only.

load_snapshot(numerical_identifier: int, list_of_variables: list | None = None, stack_scalar: bool | None = None, stack_velocity: bool | None = None, *, add_time: bool = True) type[Dataset]#

Load the variables for a given snapshot.

Parameters:
  • numerical_identifier (int) – The number of the snapshot.

  • list_of_variables (list, optional) – List of variables to be loaded, if None, it uses Dataset.set_of_variables, if Dataset.set_of_variables is empty, it automatically loads all arrays from this snapshot, by default None.

  • add_time (bool, optional) – Add time as a coordinate, by default True.

  • stack_scalar (bool, optional) – When true, the scalar fields will be stacked in a new coordinate n, otherwise returns one array per scalar fraction. If none, it uses Dataset.stack_scalar, by default None.

  • stack_velocity (bool, optional) – When true, the velocity will be stacked in a new coordinate i, otherwise returns one array per velocity component. If none, it uses Dataset.stack_velocity, by default None.

Returns:

Dataset containing the arrays loaded from the disc with the appropriate dimensions, coordinates and attributes.

Return type:

xarray.Dataset

Raises:

IOError – Raises IO error if it does not find any variable for this snapshot.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

Load all variables from a given snapshot:

>>> snapshot = prm.dataset.load_snapshot(10)

or just

>>> snapshot = prm.dataset[10]
load_time_series(array_prefix: str) type[DataArray]#

Load the entire time series for a given variable.

Note

Make sure to have enough memory to load all files at the same time.

Parameters:

array_prefix (str) – Name of the variable, for instance ux, uy, uz, pp, phi1.

Returns:

DataArray containing the time series loaded from the disc, with the appropriate dimensions, coordinates and attributes.

Return type:

xarray.DataArray

Raises:

IOError – Raises IO error if it does not find any snapshot for this variable.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

Load the entire time series for a given variable:

>>> ux = prm.dataset.load_time_series("ux")
>>> uy = prm.dataset.load_time_series("uy")
>>> uz = prm.dataset.load_time_series("uz")

or just:

>>> ux = prm.dataset["ux"]
>>> uy = prm.dataset["uy"]
>>> uz = prm.dataset["uz"]

You can organize them using a dataset:

>>> dataset = xarray.Dataset()
>>> for var in "ux uy uz".split():
...     dataset[var] = prm.dataset[var]
load_wind_turbine_data(file_pattern: str | None = None) type[Dataset]#

Load the data produced by wind turbine simulations.

Note

This feature is experimental

Parameters:

file_pattern (str, optional) – A filename pattern used to locate the files with glob.iglob. If None, it is obtained from datapath, i.g., if datapath = ./examples/Wind-Turbine/data then file_pattern = ./examples/Wind-Turbine/data/../*.perf. By default None.

Returns:

A dataset with all variables as a function of the time

Return type:

xarray.Dataset

Examples

>>> prm = xcompact3d_toolbox.Parameters(loadfile="NREL-5MW.i3d")
>>> ds = prm.dataset.load_wind_turbine_data()
>>> ds
<xarray.Dataset>
Dimensions:          (t: 21)
Coordinates:
* t                (t) float64 0.0 2.0 4.0 6.0 8.0 ... 34.0 36.0 38.0 40.0
Data variables: (12/14)
    Number of Revs   (t) float64 0.0 0.4635 0.9281 1.347 ... 7.374 7.778 8.181
    GeneratorSpeed   (t) float64 0.0 149.3 133.6 123.2 ... 122.9 122.9 123.0
    GeneratorTorque  (t) float64 0.0 2.972e+04 3.83e+04 ... 4.31e+04 4.309e+04
    BladePitch1      (t) float64 0.0 12.0 14.21 13.21 ... 11.44 11.44 11.44
    BladePitch2      (t) float64 0.0 12.0 14.21 13.21 ... 11.44 11.44 11.44
    BladePitch3      (t) float64 0.0 12.0 14.21 13.21 ... 11.44 11.44 11.44
    ...               ...
    Ux               (t) float64 0.0 15.0 15.0 15.0 15.0 ... 15.0 15.0 15.0 15.0
    Uy               (t) float64 0.0 -1.562e-05 2.541e-05 ... 7.28e-07 7.683e-07
    Uz               (t) float64 0.0 1.55e-06 ... -1.828e-06 2.721e-06
    Thrust           (t) float64 9.39e+05 1.826e+05 ... 4.084e+05 4.066e+05
    Torque           (t) float64 8.78e+06 1.268e+06 ... 4.231e+06 4.203e+06
    Power            (t) float64 1.112e+07 1.952e+06 ... 5.362e+06 5.328e+06
set(**kwargs)#

Set new values for any of the properties after the initialization.

Parameters:
  • filename_properties (dict, optional) – Keyword arguments for FilenameProperties, like separator, file_extension and so on.

  • **kwargs – Keyword arguments for the parameters, like data_path, drop_coords and so on.

Raises:

KeyError – Exception is raised when an Keyword arguments is not a valid parameter.

Examples

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )
write(data: DataArray | Dataset, file_prefix: str | None = None)#

Write an array or dataset to raw binary files on the disc, in the same order that Xcompact3d would do, so they can be easily read with 2DECOMP.

In order to avoid overwriting any relevant field, only the variables in a dataset with an attribute called file_name will be written.

Coordinates are properly aligned before writing.

If n is a valid coordinate (for scalar fractions) in the array, one numerated binary file will be written for each scalar field.

If i is a valid coordinate (for velocity components) in the array, one binary file will be written for each of them (x, y or z).

If t is a valid coordinate (for time) in the array, one numerated binary file will be written for each available time.

Parameters:
Raises:

IOError – Raises IO error data is not an xarray.DataArray or xarray.Dataset.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )
  • From a dataset, write only the variables with the attribute file_name, notice that ux and uy will not be overwritten because them do not have the attribute file_name:

    >>> for ds in prm.dataset:
    ...     ds["vort"] = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative(
    ...         "y"
    ...     )
    ...     ds["vort"].attrs["file_name"] = "vorticity"
    ...     prm.dataset.write(ds)
    
  • Write an array:

    >>> for ds in prm.dataset:
    ...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
    ...     vort.attrs["file_name"] = "vorticity"
    ...     prm.dataset.write(vort)
    

    or

    >>> for ds in prm.dataset:
    ...     vort = ds.uy.x3d.first_derivative("x") - ds.ux.x3d.first_derivative("y")
    ...     prm.dataset.write(data=vort, file_prefix="vorticity")
    

Note

It is not recommended to load the arrays with add_time = False when planning to write the results in a time series (e.g., vort-000.bin, vort-001.bin, vort-002.bin, …)

write_xdmf(xdmf_name: str = 'snapshots.xdmf', *, float_precision: int | None = None) None#

Write the xdmf file, so the results from the simulation and its postprocessing can be opened in an external visualization tool, like Paraview.

Make sure to set all the parameters in this object properly.

If set_of_objects is empty, the files are obtained automatically with glob.glob.

Parameters:
  • xdmf_name (str, optional) – Filename for the xdmf file, by default “snapshots.xdmf”

  • float_precision (int, optional) – Number of digits for the float precision on the output file, by default None. If None, it uses 8 for single precision and 16 for double precision arrays.

Raises:

IOError – Raises IO error if it does not find any file for this simulation.

Examples

Initial setup:

>>> prm = xcompact3d_toolbox.Parameters(loadfile="input.i3d")
>>> prm.dataset.set(
...     filename_properties=dict(
...         separator="-", file_extension=".bin", number_of_digits=3
...     ),
...     stack_scalar=True,
...     stack_velocity=True,
... )

It is possible to produce a new xdmf file, so all data can be visualized on any external tool:

>>> prm.dataset.write_xdmf()

Changed in version 1.2.0: Added argument float_precision.

class xcompact3d_toolbox.io.FilenameProperties(**kwargs: Any)#

Bases: HasTraits

Filename properties are important to guarantee consistency for input/output operations. This class makes xcompact3d-toolbox work with different types of file names for the binary fields produced from the numerical simulations and their pre/postprocessing.

Parameters:
  • separator (str) – The string used as separator between the name of the variable and its numeration, it can be an empty string (default is "-").

  • file_extension (str) – The file extension that identify the raw binary files from XCompact3d, it can be an empty string (default is ".bin").

  • number_of_digits (int) – The number of numerical digits used to identify the time series (default is 3).

  • scalar_num_of_digits (int) – The number of numerical digits used to identify each scalar field (default is 1).

Notes

FilenameProperties is in fact an attribute of xcompact3d_toolbox.io.Dataset, so there is no need to initialize it manually for most of the common use cases.

__init__(**kwargs)#

Initializes the object.

Parameters:

**kwargs – Keyword arguments for the parameters, like separator, file_extension and so on.

Raises:

KeyError – Raises an error if the user tries to set an invalid parameter.

Returns:

Filename properties

Return type:

xcompact3d_toolbox.io.FilenameProperties

__repr__()#

Return repr(self).

get_filename_for_binary(prefix: str, counter: int, data_path='') str#

Get the filename for an array.

Parameters:
  • prefix (str) – Name of the array.

  • counter (int or str) – The number that identifies this array, it can be the string "*" if the filename is going to be used with glob.glob.

  • data_path (str) – Path to the folder where the data is stored.

Returns:

The filename

Return type:

str

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator="-", file_extension=".bin", number_of_digits=3
... )
>>> prm.dataset.filename_properties.get_filename_for_binary("ux", 10)
'ux-010.bin'
>>> prm.dataset.filename_properties.get_filename_for_binary("ux", "*")
'ux-???.bin'
>>> prm.dataset.filename_properties.set(
...     separator="", file_extension="", number_of_digits=4
... )
>>> prm.dataset.filename_properties.get_filename_for_binary("ux", 10)
'ux0010'
>>> prm.dataset.filename_properties.get_filename_for_binary("ux", "*")
'ux????'
get_info_from_filename(filename: str) tuple[int, str]#

Get information from the name of a binary file.

Parameters:

filename (str) – The name of the array.

Returns:

A tuple with the name of the array and the number that identifies it.

Return type:

tuple[int, str]

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator="-", file_extension=".bin", number_of_digits=3
... )
>>> prm.dataset.filename_properties.get_info_from_filename("ux-010.bin")
(10, 'ux')
>>> prm.dataset.filename_properties.set(
...     separator="", file_extension="", number_of_digits=4
... )
>>> prm.dataset.filename_properties.get_info_from_filename("ux0010")
(10, 'ux')
get_name_from_filename(filename: str) str#

Same as get_info_from_filename, but just returns the name.

get_num_from_filename(filename: str) int#

Same as get_info_from_filename, but just returns the number.

set(**kwargs) None#

Set new values for any of the properties after the initialization.

Parameters:

**kwargs – Keyword arguments for parameters, like separator, file_extension and so on.

Raises:

KeyError – Raises an error if the user tries to set an invalid parameter.

Examples

If the simulated fields are named like ux-000.bin, they are in the default configuration, there is no need to specify filename properties. But just in case, it would be like:

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator="-", file_extension=".bin", number_of_digits=3
... )

If the simulated fields are named like ux0000, the parameters are:

>>> prm = xcompact3d_toolbox.Parameters()
>>> prm.dataset.filename_properties.set(
...     separator="", file_extension="", number_of_digits=4
... )

Computation and Plotting#

The data structure is provided by xarray, that introduces labels in the form of dimensions, coordinates and attributes on top of raw NumPy-like arrays, which allows for a more intuitive, more concise, and less error-prone developer experience.

See xarray’s User Guide for a complete overview about its data structures and built-in functions for indexing, selecting, computing, plotting and much more. It integrates tightly with dask for parallel computing.

Consider using hvPlot to explore your data interactively, see how to plot Gridded Data.

Xcompact3d-toolbox adds extra functions on top of xarray.DataArray and xarray.Dataset, all the details are described below.

class xcompact3d_toolbox.array.X3dDataArray(data_array)#

An accessor with extra utilities for xarray.DataArray.

cumtrapz(dim)#

See cumulative_trapezoid.

Deprecated since version 1.2.0: SciPy 1.12.0 deprecated scipy.integrate.cumtrapz in favor of scipy.integrate.cumulative_trapezoid. Since it is running under the hood, it makes sense to rename the method on xcompact3d-toolbox as well.

cumulative_trapezoid(dim)#

Cumulatively integrate xarray.DataArray in direction dim using the composite trapezoidal rule. It is a wrapper for scipy.integrate.cumulative_trapezoid. Initial value is defined to zero.

Parameters:

dim (str) – Coordinate used for the integration.

Returns:

Integrated

Return type:

xarray.DataArray

Examples

>>> da.x3d.cumulative_trapezoid("t")

Added in version 1.2.0.

first_derivative(dim)#

Compute first derivative with the 4th order accurate centered scheme.

It is fully functional with all boundary conditions available on XCompact3d and stretched mesh in the vertical direction (y). The attribute BC is used to store Boundary Condition information in a dictionary (see examples), default is ncl1 = ncln = 2 and npaire = 1.

Parameters:

dim (str) – Coordinate used for the derivative.

Returns:

differentiated

Return type:

xarray.DataArray

Examples

>>> da.attrs['BC'] = {
...     'x': {
...         'ncl1': 1,
...         'ncln': 1,
...         'npaire': 0
...     },
...     'y': {
...         'ncl1': 2,
...         'ncln': 1,
...         'npaire': 1,
...         'istret': 0,
...         'beta': 1.0
...     },
...     'z': {
...         'ncl1': 0,
...         'ncln': 0,
...         'npaire': 1
... }
>>> da.x3d.first_derivative("x")

or just:

>>> prm = xcompact3d_toolbox.Parameters()
>>> da.attrs["BC"] = prm.get_boundary_condition("ux")
>>> da.x3d.first_derivative("x")
pencil_decomp(*args)#

Coerce the data array into dask array.

It applies chunk=-1 for all coordinates listed in args, which means no decomposition, and 'auto' to the others, resulting in a pencil decomposition for parallel evaluation.

For customized chunks adjust, see xarray.DataArray.chunk.

Parameters:

arg (str or sequence of str) – Dimension(s) to apply no decomposition.

Returns:

chunked

Return type:

xarray.DataArray

Raises:

ValueError – args must be valid dimensions in the data array

Examples

>>> da.x3d.pencil_decomp("x")  # Pencil decomposition
>>> da.x3d.pencil_decomp("t")
>>> da.x3d.pencil_decomp("y", "z")  # Slab decomposition
second_derivative(dim)#

Compute second derivative with the 4th order accurate centered scheme.

It is fully functional with all boundary conditions available on Xcompact3d and stretched mesh in y direction. The attribute BC is used to store Boundary Condition information in a dictionary (see examples), default is ncl1 = ncln = 2 and npaire = 1.

Parameters:

dim (str) – Coordinate used for the derivative.

Returns:

differentiated

Return type:

xarray.DataArray

Examples

>>> da.attrs['BC'] = {
...     'x': {
...         'ncl1': 1,
...         'ncln': 1,
...         'npaire': 0
...     },
...     'y': {
...         'ncl1': 2,
...         'ncln': 1,
...         'npaire': 1
...         'istret': 0,
...         'beta': 1.0
...     },
...     'z': {
...         'ncl1': 0,
...         'ncln': 0,
...         'npaire': 1
... }
>>> da.x3d.second_derivative("x")

or just:

>>> prm = xcompact3d_toolbox.Parameters()
>>> da.attrs["BC"] = prm.get_boundary_condition("ux")
>>> da.x3d.second_derivative("x")
simps(*args)#

See simpson.

Deprecated since version 1.2.0: SciPy 1.12.0 deprecated scipy.integrate.simps in favor of scipy.integrate.simpson. Since it is running under the hood, it makes sense to rename the method on xcompact3d-toolbox as well.

simpson(*args)#

Integrate xarray.DataArray in direction(s) args using the composite Simpson’s rule. It is a wrapper for scipy.integrate.simpson.

Parameters:

arg (str or sequence of str) – Dimension(s) to compute integration.

Returns:

Integrated

Return type:

xarray.DataArray

Raises:

ValueError – args must be valid dimensions in the data array

Examples

>>> da.x3d.simpson("x")
>>> da.x3d.simpson("t")
>>> da.x3d.simpson("x", "y", "z")

Added in version 1.2.0.

class xcompact3d_toolbox.array.X3dDataset(data_set)#

An accessor with extra utilities for xarray.Dataset.

cumtrapz(dim)#

See cumulative_trapezoid.

Deprecated since version 1.2.0: SciPy 1.12.0 deprecated scipy.integrate.cumtrapz in favor of scipy.integrate.cumulative_trapezoid. Since it is running under the hood, it makes sense to rename the method on xcompact3d-toolbox as well.

cumulative_trapezoid(dim)#

Cumulatively integrate all arrays in this dataset in direction dim using the composite trapezoidal rule. It is a wrapper for scipy.integrate.cumulative_trapezoid. Initial value is defined to zero.

Parameters:

dim (str) – Coordinate used for the integration.

Returns:

Integrated

Return type:

xarray.Dataset

Examples

>>> ds.x3d.cumulative_trapezoid("t")

Added in version 1.2.0.

pencil_decomp(*args)#

Coerce all arrays in this dataset into dask arrays.

It applies chunk=-1 for all coordinates listed in args, which means no decomposition, and "auto" to the others, resulting in a pencil decomposition for parallel evaluation.

For customized chunks adjust, see xarray.Dataset.chunk.

Parameters:

arg (str or sequence of str) – Dimension(s) to apply no decomposition.

Returns:

chunked

Return type:

xarray.Dataset

Raises:

ValueError – args must be valid dimensions in the dataset

Examples

>>> ds.x3d.pencil_decomp("x")  # Pencil decomposition
>>> ds.x3d.pencil_decomp("t")
>>> ds.x3d.pencil_decomp("y", "z")  # Slab decomposition
simps(*args)#

See simpson.

Deprecated since version 1.2.0: SciPy 1.12.0 deprecated scipy.integrate.simps in favor of scipy.integrate.simpson. Since it is running under the hood, it makes sense to rename the method on xcompact3d-toolbox as well.

simpson(*args)#

Integrate all arrays in this dataset in direction(s) args using the composite Simpson’s rule. It is a wrapper for scipy.integrate.simpson.

Parameters:

arg (str or sequence of str) – Dimension(s) to compute integration.

Returns:

Integrated

Return type:

xarray.Dataset

Raises:

ValueError – args must be valid dimensions in the dataset

Examples

>>> ds.x3d.simpson("x")
>>> ds.x3d.simpson("t")
>>> ds.x3d.simpson("x", "y", "z")

Added in version 1.2.0.

Sandbox#

The new Sandbox Flow Configuration (itype = 12) aims to break many of the barriers to entry in a Navier-Stokes solver. The idea is to easily provide everything that XCompact3d needs from a Python Jupyter Notebook, like initial conditions, solid geometry, boundary conditions, and the parameters. For students in computational fluid dynamics, it provides a direct hands-on experience and a safe place for practicing and learning, while for advanced users and code developers, it works as a rapid prototyping tool. For more details, see:

exception xcompact3d_toolbox.sandbox.DimensionNotFoundError(dim)#

Raised when a dimension is not found in the DataArray.

class xcompact3d_toolbox.sandbox.Geometry(data_array: DataArray)#

An accessor with some standard geometries for xarray.DataArray. Use them in combination with the arrays initialized at xcompact3d_toolbox.sandbox.init_epsi and the new xcompact3d_toolbox.genepsi.gene_epsi_3d.

ahmed_body(*, scale: float = 1.0, angle: float = 45.0, wheels: bool = False, remp: bool = True, **kwargs) DataArray#

Draw an Ahmed body.

Parameters:
  • scale (float) – Ahmed body’s scale (the default is 1).

  • angle (float) – Ahmed body’s angle at the back, in degrees (the default is 45).

  • wheel (bool) – Draw “wheels” if True (the default is False).

  • remp (bool) – Adds the geometry to the xarray.DataArray if True and removes it if False (the default is True).

  • **kwargs (float) – Ahmed body’s center.

Returns:

Array with(out) the Ahmed body

Return type:

xarray.DataArray

Raises:

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.ahmed_body(x=2)

Changed in version 1.2.0: All arguments changed to keyword-only.

box(*, remp: bool = True, **kwargs) DataArray#

Draw a box.

Parameters:
  • remp (bool) – Adds the geometry to the xarray.DataArray if True and removes it if False (the default is True).

  • **kwargs (tuple of float) – Box’s boundaries.

Returns:

Array with(out) the box

Return type:

xarray.DataArray

Raises:

KeyError – Boundaries coordinates must be valid dimensions

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.box(x=(2,5), y=(0,1))

Changed in version 1.2.0: All arguments changed to keyword-only.

cylinder(*, radius: float = 0.5, axis: str = 'z', height: float | None = None, remp: bool = True, **kwargs) DataArray#

Draw a cylinder.

Parameters:
  • radius (float) – Cylinder’s radius (the default is 0.5).

  • axis (str) – Cylinder’s axis (the default is "z").

  • height (float or None) – Cylinder’s height (the default is None), if None, it will take the entire axis, otherwise \(\pm h/2\) is considered from the center.

  • remp (bool) – Adds the geometry to the xarray.DataArray if True and removes it if False (the default is True).

  • **kwargs (float) – Cylinder’s center point.

Returns:

Array with(out) the cylinder

Return type:

xarray.DataArray

Raises:

KeyError – Center coordinates must be valid dimensions

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.cylinder(x=4.0, y=5.0)

Changed in version 1.2.0: All arguments changed to keyword-only.

from_stl(*, filename: str | None = None, stl_mesh: Mesh | None = None, origin: dict | None = None, rotate: dict | None = None, scale: float | None = None, user_tol: float = 6.283185307179586, remp: bool = True) DataArray#

Load a STL file and compute if the nodes of the computational mesh are inside or outside the object. In this way, the customized geometry can be used at the flow solver.

The methodology is based on the work of:

  • Jacobson, A., Kavan, L., & Sorkine-Hornung, O. (2013). Robust inside-outside segmentation using generalized winding numbers. ACM Transactions on Graphics (TOG), 32(4), 1-12.

The Python implementation is an adaptation from inside-3d-mesh (licensed under the MIT License), by @marmakoide.

To maximize the performance here at the toolbox, from_stl is powered by Numba, that translates Python functions to optimized machine code at runtime. This method is compatible with Dask for parallel computation. In addition, just the subdomain near the object is tested, to save computational time.

Note

The precision of the method is influenced by the complexity of the STL mesh, there is no guarantee it will work for all geometries. This feature is experimental, its interface may change in future releases.

Parameters:
  • filename (str, optional) – Filename of the STL file to be loaded and included in the cartesian domain, by default None

  • scale (float, optional) – This parameters can be used to scale up the object when greater than one and scale it down when smaller than one, by default None

  • rotate (dict, optional) – Rotate the object, including keyword arguments that are expected by stl.mesh.Mesh.rotate, like axis, theta and point. For more details, see numpy-stl’s documentation. By default None

  • origin (dict, optional) – Specify the location of the origin point for the geometry. It is considered as the minimum value computed from all points in the object for each coordinate, after scaling and rotating them. The keys of the dictionary are the coordinate names (x, y and z) and the values are the origin on that coordinate. For missing keys, the value is assumed as zero. By default None

  • stl_mesh (stl.mesh.Mesh, optional) – For very customizable control over the 3D object, you can provide it directly. Note that none of the arguments above are applied in this case. For more details about how to create and modify the geometry, see numpy-stl’s documentation. By default None

  • user_tol (float, optional) – Control the tolerance used to compute if a mesh node is inside or outside the object. Values smaller than the default may reduce the number of false negatives. By default \(2\pi\)

  • remp (bool, optional) – Add the geometry to the xarray.DataArray if True and removes it if False, by default True

Returns:

Array with(out) the customized geometry

Return type:

xarray.DataArray

Raises:

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm, dask=True)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.from_stl(
...         filename="My_file.stl",
...         scale=1.0,
...         rotate=dict(axis=[0, 0.5, 0], theta=math.radians(90)),
...         origin=dict(x=2, y=1, z=0),
...     )

Changed in version 1.2.0: All arguments changed to keyword-only.

mirror(dim: str = 'x') DataArray#

Mirror the \(\epsilon\) array with respect to the central plane in the direction dim.

Parameters:

dim (str) – Reference for the mirror (the default is x).

Returns:

Mirrored array

Return type:

xarray.DataArray

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.cylinder(x=4, y=5).geo.mirror("x")
sphere(*, radius: float = 0.5, remp: bool = True, **kwargs) DataArray#

Draw a sphere.

Parameters:
  • radius (float) – Sphere’s radius (the default is 0.5).

  • remp (bool) – Adds the geometry to the xarray.DataArray if True and removes it if False (the default is True).

  • **kwargs (float) – Sphere’s center.

Returns:

Array with(out) the sphere

Return type:

xarray.DataArray

Raises:

KeyError – Center coordinates must be valid dimensions

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.sphere(x=1, y=1, z=1)

Changed in version 1.2.0: All arguments changed to keyword-only.

square(*, length: float = 1.0, thickness: float = 0.1, remp: bool = True, **kwargs) DataArray#

Draw a squared frame.

Parameters:
  • length (float) – Frame’s external length (the default is 1).

  • thickness (float) – Frames’s thickness (the default is 0.1).

  • remp (bool) – Adds the geometry to the xarray.DataArray if True and removes it if False (the default is True).

  • **kwargs (float) – Frames’s center.

Returns:

Array with(out) the squared frame

Return type:

xarray.DataArray

Raises:

KeyError – Center coordinates must be valid dimensions

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)
>>> for key in epsi.keys():
>>>     epsi[key] = epsi[key].geo.square(x=5, y=2, z=1)

Changed in version 1.2.0: All arguments changed to keyword-only.

xcompact3d_toolbox.sandbox.init_dataset(prm: Parameters) xr.Dataset#

This function initializes a xarray.Dataset including all variables that should be provided to XCompact3d and the sandbox flow configuration, according to the computational and physical parameters.

Parameters:

prm (xcompact3d_toolbox.parameters.Parameters) – Contains the computational and physical parameters.

Returns:

Each variable is initialized with np.zeros(dtype=xcompact3d_toolbox.param["mytype"]) and wrapped into a xarray.Dataset with the proper size, dimensions, coordinates and attributes, check them for more details. The variables are:

  • bxx1, bxy1, bxz1 - Inflow boundary condition for ux, uy and uz, respectively (if nclx1 = 2);

  • noise_mod_x1 - for random noise modulation at inflow boundary condition (if nclx1 = 2);

  • bxphi1 - Inflow boundary condition for scalar field(s) (if nclx1 = 2 and numscalar > 0);

  • byphi1 - Bottom boundary condition for scalar field(s) (if nclyS1 = 2, numscalar > 0 and uset = 0);

  • byphin - Top boundary condition for scalar field(s) (if nclySn = 2, numscalar > 0 and uset = 0);

  • ux, uy, uz - Initial condition for velocity field;

  • phi - Initial condition for scalar field(s) (if numscalar > 0);

  • vol_frc - Integral operator employed for flow rate control in case of periodicity in x direction (nclx1 = 0 and nclxn = 0). Xcompact3d will compute the volumetric integration as I = sum(vol_frc * ux) and them will correct streamwise velocity as ux = ux / I, so, set vol_frc properly.

Return type:

xarray.Dataset

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> dataset = xcompact3d_toolbox.init_dataset(prm)
>>> #
>>> # Code here your customized flow configuration
>>> #
>>> prm.dataset.write(dataset)  # write the files to the disc
xcompact3d_toolbox.sandbox.init_epsi(prm: Parameters, *, dask: bool = False) dict[str, xr.DataArray]#

Initializes the \(\epsilon\) arrays that define the solid geometry for the Immersed Boundary Method.

Parameters:
  • prm (xcompact3d_toolbox.parameters.Parameters) – Contains the computational and physical parameters.

  • dask (bool) – Defines the lazy parallel execution with dask arrays. See xcompact3d_toolbox.array.x3d.pencil_decomp().

Returns:

A dictionary containing the epsi(s) array(s):

  • epsi (nx, ny, nz) if iibm != 0;

  • xepsi (nxraf, ny, nz) if iibm = 2;

  • yepsi (nx, nyraf, nz) if iibm = 2;

  • zepsi (nx, ny, nzraf) if iibm = 2.

Each one initialized with np.zeros(dtype=np.bool) and wrapped into a xarray.DataArray with the proper size, dimensions and coordinates. They are used to define the object(s) that is(are) going to be inserted into the cartesian domain by the Immersed Boundary Method (IBM). They should be set to one (True) at the solid points and stay zero (False) at the fluid points, some standard geometries are provided by the accessor xcompact3d_toolbox.sandbox.Geometry.

Return type:

dict of xarray.DataArray

Examples

>>> prm = xcompact3d_toolbox.Parameters()
>>> epsi = xcompact3d_toolbox.init_epsi(prm)

Changed in version 1.2.0: The argument dask changed to keyword-only.

Genepsi#

This module generates all the files necessary for our customized Immersed Boundary Method, based on Lagrange reconstructions. It is an adaptation to Python from the original Fortran code and methods from:

  • Gautier R., Laizet S. & Lamballais E., 2014, A DNS study of jet control with microjets using an alternating direction forcing strategy, Int. J. of Computational Fluid Dynamics, 28, 393–410.

gene_epsi_3d is powered by Numba, it translates Python functions to optimized machine code at runtime. Numba-compiled numerical algorithms in Python can approach the speeds of C or FORTRAN.

xcompact3d_toolbox.genepsi.gene_epsi_3d(epsi_in_dict, prm)#

This function generates all the Auxiliary files necessary for our customize IBM, based on Lagrange reconstructions. The arrays can be initialized with xcompact3d_toolbox.sandbox.init_epsi(), then, some standard geometries are provided by the accessor xcompact3d_toolbox.sandbox.Geometry. Notice that you can apply our own routines for your own objects. The main outputs of the function are written to disc at the files epsilon.bin, nobjx.dat, nobjy.dat, nobjz.dat, nxifpif.dat, nyifpif.dat, nzifpif.dat, xixf.dat, yiyf.dat and zizf.dat. They will be used by Xcompact3d and the sandbox flow configuration.

Parameters:
Returns:

All computed variables are returned in a Dataset if prm.iibm >= 2, but just for reference, since all the relevant values are written to the disc.

Return type:

xarray.Dataset or None

Examples

>>> prm = x3d.Parameters()
>>> epsi = x3d.sandbox.init_epsi(prm)
>>> for key in epsi.keys():
...     epsi[key] = epsi[key].geo.cylinder(x=4, y=5)
>>> dataset = x3d.gene_epsi_3d(epsi, prm)

Remember to set the number of objects after that if prm.iibm >= 2:

>>> if prm.iibm >= 2:
...     prm.nobjmax = dataset.obj.size

Sample Data#

xcompact3d_toolbox.tutorial.open_dataset(name: str, **kws) tuple[Dataset, Parameters]#

Open a dataset from the online repository (requires internet).

If a local copy is found then always use that to avoid network traffic.

Available datasets: * "cylinder": Flow around a cylinder

Parameters: