fluidsim.solvers.ns3d.solver

NS3D solver (fluidsim.solvers.ns3d.solver)

class fluidsim.solvers.ns3d.solver.Simul(params)[source]

Bases: SimulBasePseudoSpectral

Pseudo-spectral solver 3D incompressible Navier-Stokes equations.

Notes

This class is dedicated to solve with a pseudo-spectral method the incompressible Navier-Stokes equations (possibly with hyper-viscosity):

\[\partial_t \textbf{v} + \textbf{v} \cdot \boldsymbol{\nabla} \textbf{v} = - \boldsymbol{\nabla} p - \nu_\alpha (-\Delta)^\alpha \textbf{v},\]

where \(\textbf{v}\) is the non-divergent velocity (\(\boldsymbol{\nabla} \cdot \textbf{v} = 0\)), \(p\) is the pressure, \(\Delta\) is the 3D Laplacian operator.

In Fourier space, these equations can be written as:

\[\partial_t \hat{\textbf{v}} = \hat N(\textbf{v}) + L \hat{\textbf{v}},\]

where

\[\hat N(\textbf{v}) = -P_\perp \widehat{\textbf{v} \cdot \boldsymbol{\nabla} \textbf{v}},\]
\[L = - \nu_\alpha |\textbf{k}|^{2\alpha},\]

with \(P_\perp = (1 - \hat{\textbf{e}}_\textbf{k} \hat{\textbf{e}}_\textbf{k} \cdot)\) the operator projection on the plane perpendicular to the wave number \(\textbf{k}\). Since the flow is incompressible (\(\textbf{k} \cdot \textbf{v} = 0\)), the effect of the pressure term is taken into account with the operator \(P_\perp\).

In practice, it is more efficient to use the relation

\[\textbf{v} \cdot \boldsymbol{\nabla} \textbf{v} = \boldsymbol{\nabla} |\textbf{v}|^2/2 - \textbf{v} \times \boldsymbol{\omega},\]

with \(\boldsymbol{\omega} = \boldsymbol{\nabla} \times \textbf{v}\) the vorticity, and to compute the nonlinear term as

\[\hat N(\textbf{v}) = P_\perp \widehat{\textbf{v} \times \boldsymbol{\omega}}.\]
InfoSolver

alias of InfoSolverNS3D

static _complete_params_with_default(params)[source]

This static method is used to complete the params container.

tendencies_nonlin(state_spect=None, old=None)[source]

Compute the nonlinear tendencies.

This function has to be overridden in a child class.

Returns:
tendencies_fftfluidsim.base.setofvariables.SetOfVariables

An array containing only zeros.

Default parameters

<params NEW_DIR_RESULTS="True" ONLY_COARSE_OPER="False" f="None"
        no_vz_kz0="False" nu_2="0.0" nu_4="0.0" nu_8="0.0" nu_m4="0.0"
        projection="None" short_name_type_run="">
  <oper Lx="6.283185307179586" Ly="6.283185307179586" Lz="6.283185307179586"
        NO_SHEAR_MODES="False" coef_dealiasing="0.6666666666666666" nx="48"
        ny="48" nz="48" truncation_shape="cubic" type_fft="default"
        type_fft2d="sequential"/>  

  <time_stepping USE_CFL="True" USE_T_END="True" cfl_coef="None" deltat0="0.2"
                 deltat_max="0.2" it_end="10" max_elapsed="None" t_end="10.0"
                 type_time_scheme="RK4">
    <phaseshift_random nb_pairs="1" nb_steps_compute_new_pair="None"/>  

  </time_stepping>

  <init_fields available_types="['from_file', 'from_simul', 'in_script',
               'constant', 'noise', 'dipole']" modif_after_init="False"
               type="constant">
    <from_file path=""/>  

    <constant value="1.0"/>  

    <noise length="None" velo_max="1.0"/>  

  </init_fields>

  <forcing available_types="['in_script', 'in_script_coarse', 'pseudo_spectral',
           'proportional', 'taylor_green', 'tcrandom', 'tcrandom_anisotropic',
           'watu_coriolis']" enable="False" forcing_rate="1.0" key_forced="None"
           nkmax_forcing="5" nkmin_forcing="4" type="">
    <milestone nx_max="None">
      <objects diameter="1.0" number="2" type="cylinders"
               width_boundary_layers="0.1"/>  

      <movement type="uniform">
        <uniform speed="1.0"/>  

        <sinusoidal length="1.0" period="1.0"/>  

        <periodic_uniform length="1.0" length_acc="1.0" speed="1.0"/>  

      </movement>

    </milestone>

    <taylor_green amplitude="1.0"/>  

    <normalized constant_rate_of="None" type="2nd_degree_eq"
                which_root="minabs"/>  

    <random only_positive="False"/>  

    <tcrandom time_correlation="based_on_forcing_rate"/>  

    <tcrandom_anisotropic angle="45°" delta_angle="10°"
                          kz_negative_enable="False"/>  

    <watu_coriolis amplitude="0.05" approximate_dt="None" delta_omega_f="0.03"
                   nb_wave_makers="2" omega_f="0.3" period_forcing="1000"/>  

  </forcing>

  <output HAS_TO_SAVE="True" ONLINE_PLOT_OK="True" period_refresh_plots="1"
          sub_directory="">
    <periods_save cross_corr="0" phys_fields="0" spatial_means="0"
                  spatiotemporal_spectra="0" spect_energy_budg="0" spectra="0"
                  temporal_spectra="0"/>  

    <periods_print print_stdout="1.0"/>  

    <periods_plot phys_fields="0"/>  

    <phys_fields field_to_plot="rotz" file_with_it="False"/>  

    <spectra HAS_TO_PLOT_SAVED="False" kzkh_periodicity="0"/>  

    <spatial_means HAS_TO_PLOT_SAVED="False"/>  

    <spect_energy_budg HAS_TO_PLOT_SAVED="False"/>  

    <temporal_spectra SAVE_AS_FLOAT32="True" file_max_size="10.0"
                      probes_deltax="0.1" probes_deltay="0.1"
                      probes_deltaz="0.1" probes_region="None"/>  

    <spatiotemporal_spectra SAVE_AS_COMPLEX64="True" file_max_size="10.0"
                            probes_region="None"/>  

    <cross_corr HAS_TO_PLOT_SAVED="False" kzkh_periodicity="0"/>  

  </output>

  <preprocess enable="False" forcing_const="1.0" forcing_scale="unity"
              init_field_const="1.0" init_field_scale="unity"
              viscosity_const="1.0" viscosity_scale="enstrophy_forcing"
              viscosity_type="laplacian"/>  

</params>

Documentation for params

short_name_type_run: str

A short name of the simulation used to create the directory name.

NEW_DIR_RESULTS: bool

To be used only when loading a simulation. If True (default), a new directory is created to contain the results of the simulation. If False, the results of the simulation are appended in the old directory.

ONLY_COARSE_OPER: bool

To be used only when loading a simulation. If True (not default), the operator is created with a very small resolution. It is very fast but then it can not be used to process data.

nu_2: float (default = 0.)

Viscosity coefficient. Used in particular in the method fluidsim.base.solvers.pseudo_spect.SimulBasePseudoSpectral.compute_freq_diss()).

nu_8: float

Hyper-viscous coefficient of order 8. Also used in the method compute_freq_diss.

nu_4: float

Hyper-viscous coefficient of order 4.

nu_m4: float

Hypo-viscous coefficient of order -4. Hypo-viscosity affect more the large scales.

f: float (default None)

Coriolis parameter (effect of the system rotation).

no_vz_kz0: bool (default False)

If True, vz(kz=0) is 0.

projection: str (default None)

If “toroidal” or “vortical”, the solution and the equations are projected on the toroidal manifold. If “poloidal”, on the poloidal one.

Documentation for params.oper

See the documentation of fluidfft (in particular of the 3d operator class).

type_fft: str

Method for the FFT (as defined by fluidfft).

type_fft2d: str

Method for the 2d FFT.

coef_dealiasing: float

dealiasing coefficient.

nx: int

Number of points over the x-axis (last dimension in the physical space).

ny: int

Number of points over the y-axis (second dimension in the physical space).

nz: int

Number of points over the z-axis (first dimension in the physical space).

Lx, Ly and Lz: float

Length of the edges of the numerical domain.

Documentation for params.time_stepping

See fluidsim.base.time_stepping.base.

USE_T_END: bool (default True)

If True, time step until t > t_end. If False, time step until it >= it_end.

t_end: float

See documentation USE_T_END.

it_end: int

If USE_T_END is False, number of time steps.

USE_CFL: bool (default False)

If True, use an adaptive time step computed in particular with a Courant-Friedrichs-Lewy (CFL) condition.

type_time_scheme: str (default “RK4”)

Type of time scheme. Can be in (“RK2”, “RK4”).

deltat0: float (default 0.2)

If USE_CFL is False, value of the time step.

deltat_max: float (default 0.2)

Maximum value of the time step (useful when USE_CFL is True).

cfl_coef: float (default None)

If not None, clf_coef used in the CFL condition. If None, the value is choosen taking into account the time scheme.

max_elapsed: number or str (default None)

If not None, the computation stops when the elapsed time becomes larger than max_elapsed. Can be a number (in seconds) or a string (formated as “%H:%M:%S”).

Documentation for params.time_stepping.phaseshift_random

Documentation for params.init_fields

See fluidsim.base.init_fields.

type: str (default constant)

Name of the initialization method.

available_types: list

Actually not a parameter; just a hint to set type.

modif_after_init: bool (default False)

Used internally when reloading some simulations.

Documentation for params.init_fields.from_file

path: str

Documentation for params.init_fields.constant

value: float (default 1.)

Documentation for params.init_fields.noise

velo_max: float (default 1.)

Maximum velocity.

length: float (default 0.)

The smallest (cutoff) scale in the noise.

Documentation for params.forcing

How the forcing is normalized

Documentation for params.forcing.milestone

Documentation for params.forcing.milestone.objects
Documentation for params.forcing.milestone.movement
Documentation for params.forcing.milestone.movement.uniform
Documentation for params.forcing.milestone.movement.sinusoidal
Documentation for params.forcing.milestone.movement.periodic_uniform

Documentation for params.forcing.taylor_green

Documentation for params.forcing.normalized

Documentation for params.forcing.random

Documentation for params.forcing.tcrandom

Documentation for params.forcing.tcrandom_anisotropic

See fluidsim.base.forcing.anisotropic.

angle: float or str

Angle between the wavevector and the horizontal characterising the forcing region.

delta_angle: float or None

Control the shape of the forcing region in k-space. If None, (khmin, khmax, kvmin, kvmax) are computed from the angle, params.forcing.nkmin_forcing and params.forcing.nkmax_forcing.

kz_negative_enable: bool

If True, modes with negative kz are also forced.

Documentation for params.forcing.watu_coriolis

Documentation for params.output

See fluidsim.output.base

ONLINE_PLOT_OK: bool (default: True)

If True, the online plots are enabled.

period_refresh_plots: float (default: 1)

Period of refreshment of the online plots.

HAS_TO_SAVE: bool (default: True)

If False, nothing new is saved in the directory of the simulation.

sub_directory: str (default: “”)

A name of a subdirectory where the directory of the simulation is saved.

Documentation for params.output.periods_save

Periods (float, in equation time) to set when the specific outputs are saved.

Documentation for params.output.periods_print

Periods (float, in equation time) to set when the printing specific outputs are called.

Documentation for params.output.periods_plot

Periods (float, in equation time) to set when the plots of the specific outputs are called.

Documentation for params.output.phys_fields

Documentation for params.output.spectra

HAS_TO_PLOT_SAVED : bool (False)

If True, some curves can be plotted during the run.

kzkh_periodicity : int (0)

Periodicity of saving of (kz, kh) spectra (compared to standard spectra).

Documentation for params.output.spatial_means

Documentation for params.output.spect_energy_budg

HAS_TO_PLOT_SAVED : bool (False)

If True, some curves can be plotted during the run.

Documentation for params.output.temporal_spectra

probes_deltax, probes_deltay and probes_deltaz: float (default: 0.1)

Probes spacing in the x, y and z directions (in params.oper.Li unit).

probes_region: tuple (default:None)

Boundaries of the region in the simulation domain were probes are set.

probes_region = (xmin, xmax, ymin, ymax, zmin, zmax), in params.oper.Lx unit.

If None, set to the whole simulation domain.

file_max_size: float (default: 10.0)

Maximum size of one time series file, in megabytes.

SAVE_AS_FLOAT32: bool (default: True)

If set to False, probes data is saved as float64.

Warning : saving as float64 reduces digital noise at high frequencies, but double the size of the output!

Documentation for params.output.spatiotemporal_spectra

probes_region: int tuple (default:None)

Boundaries of the region to record in the spectral domain.

probes_region = (ikxmax, ikymax, ikzmax), in nondimensional units (mode indices).

The resulting ranges for each wavevectors are :

ikx in [0 , ikxmax]

iky in [-ikymax+1 , ikymax]

ikz in [-ikzmax+1 , ikzmax]

If None, set all ikmax = 4.

file_max_size: float (default: 10.0)

Maximum size of one time series file, in megabytes.

SAVE_AS_COMPLEX64: bool (default: True)

If set to False, probes data is saved as complex128.

Warning : saving as complex128 reduces digital noise at high frequency, but doubles the size of the output!

Documentation for params.output.cross_corr

HAS_TO_PLOT_SAVED : bool (False)

If True, some curves can be plotted during the run.

kzkh_periodicity : int (0)

Periodicity of saving of (kz, kh) spectra (compared to standard spectra).

Documentation for params.preprocess

Classes

InfoSolverNS3D([only_root])

Simul(params)

Pseudo-spectral solver 3D incompressible Navier-Stokes equations.