fluidsim.solvers.ns2d.solver
NS2D solver (fluidsim.solvers.ns2d.solver
)
This module provides two classes defining the pseudo-spectral solver 2D incompressible Navier-Stokes equations (ns2d).
- class fluidsim.solvers.ns2d.solver.InfoSolverNS2D(only_root=False, **kargs)[source]
Bases:
InfoSolverPseudoSpectral
Contain the information on the solver ns2d.
- _init_root()[source]
Init. self by writing the information on the solver.
The function InfoSolverPseudoSpectral._init_root is called. We keep two classes listed by this function:
The other first-level classes for this solver are:
- class fluidsim.solvers.ns2d.solver.Simul(params)[source]
Bases:
SimulBasePseudoSpectral
Pseudo-spectral solver 2D incompressible Navier-Stokes equations.
- InfoSolver
alias of
InfoSolverNS2D
- static _complete_params_with_default(params)[source]
Complete the params container (static method).
- tendencies_nonlin(state_spect=None, old=None)[source]
Compute the nonlinear tendencies.
- Parameters:
- state_spect
fluidsim.base.setofvariables.SetOfVariables
optional
Array containing the state, i.e. the vorticity, in Fourier space. If state_spect, the variables vorticity and the velocity are computed from it, otherwise, they are taken from the global state of the simulation, self.state.
These two possibilities are used during the Runge-Kutta time-stepping.
- state_spect
- Returns:
- tendencies_fft
fluidsim.base.setofvariables.SetOfVariables
An array containing the tendencies for the vorticity.
- tendencies_fft
Notes
The 2D Navier-Stokes equation can be written as:
\[\partial_t \hat\zeta = \hat N(\zeta) - \sigma(k) \hat \zeta,\]This function compute the nonlinear term (“tendencies”) \(N(\zeta) = - \mathbf{u}\cdot \mathbf{\nabla} \zeta\).
Default parameters
<params NEW_DIR_RESULTS="True" ONLY_COARSE_OPER="False" beta="0.0" nu_2="0.0"
nu_4="0.0" nu_8="0.0" nu_m4="0.0" short_name_type_run="">
<oper Lx="8" Ly="8" NO_KY0="False" NO_SHEAR_MODES="False"
coef_dealiasing="0.6666666666666666" nx="48" ny="48"
truncation_shape="cubic" type_fft="default"/>
<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', 'jet', 'dipole']" modif_after_init="False"
type="constant">
<from_file path=""/>
<constant value="1.0"/>
<noise length="0.0" velo_max="1.0"/>
</init_fields>
<forcing available_types="['in_script', 'in_script_coarse', 'pseudo_spectral',
'proportional', 'tcrandom', 'tcrandom_anisotropic']" 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>
<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="None"
kz_negative_enable="False"/>
</forcing>
<output HAS_TO_SAVE="True" ONLINE_PLOT_OK="True" period_refresh_plots="1"
sub_directory="">
<periods_save increments="0" phys_fields="0" spatial_means="0"
spatiotemporal_spectra="0" spect_energy_budg="0" spectra="0"
spectra_multidim="0" temporal_spectra="0"/>
<periods_print print_stdout="1.0"/>
<periods_plot phys_fields="0"/>
<phys_fields field_to_plot="rot" file_with_it="False"/>
<spectra HAS_TO_PLOT_SAVED="False"/>
<spectra_multidim HAS_TO_PLOT_SAVED="False"/>
<spatial_means HAS_TO_PLOT_SAVED="False"/>
<spect_energy_budg HAS_TO_PLOT_SAVED="False"/>
<increments 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_region="None"/>
<spatiotemporal_spectra SAVE_AS_COMPLEX64="True" file_max_size="10.0"
probes_region="None"/>
</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.
Documentation for params.oper
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.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
andparams.forcing.nkmax_forcing
.
kz_negative_enable: bool
If True, modes with negative kz are also forced.
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
Documentation for params.output.spectra_multidim
Documentation for params.output.spatial_means
Documentation for params.output.spect_energy_budg
Documentation for params.output.increments
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.preprocess
Functions
|
Classes
|
Contain the information on the solver ns2d. |
|
Pseudo-spectral solver 2D incompressible Navier-Stokes equations. |