bolt.lib.linear_solver module

This is the module which contains the functions of the linear solver of Bolt. It performs an FFT to map the given input onto the Fourier basis, and evolves each mode of the input independantly. It is to be noted that this module can only be applied to systems with periodic boundary conditions. Additionally, this module isn’t parallelized to run across multiple devices/nodes.

class bolt.lib.linear.linear_solver.linear_solver(physical_system)[source]

Bases: object

An instance of this class’ attributes contains methods which are used in evolving the system declared under physical system linearly. The state of the system then may be determined from the attributes of the system such as the distribution function and electromagnetic fields

Methods

RK2_timestep(dt) Evolves the physical system defined using an RK2 integrator.
RK4_timestep(dt) Evolves the physical system defined using an RK4 integrator.
RK5_timestep(dt) Evolves the physical system defined using an RK5 integrator.
compute_moments(moment_name[, f, f_hat]) Used in computing the moments of the distribution function.
dump_distribution_function(file_name) This function is used to dump distribution function to a file for later usage.This dumps the complete 5D distribution function which can be used for post-processing
dump_moments(file_name) This function is used to dump variables to a file for later usage.
get_dist_func() Returns the distribution function in the same format as the nonlinear solver
load_distribution_function(file_name) This function is used to load the distribution function from the dump file that was created by dump_distribution_function.
RK2_timestep(dt)

Evolves the physical system defined using an RK2 integrator. This method is 2nd order accurate.

Parameters:
dt: double

The timestep size.

RK4_timestep(dt)

Evolves the physical system defined using an RK4 integrator. This method is 4th order accurate.

Parameters:
dt: double

The timestep size.

RK5_timestep(dt)

Evolves the physical system defined using an RK5 integrator. This method is 5th order accurate.

Parameters:
dt: double

The timestep size.

compute_moments(moment_name, f=None, f_hat=None)

Used in computing the moments of the distribution function. The moment definitions which are passed to physical system are used in computing these moment quantities.

Parameters:
moments_name : str

Pass the moment name which needs to be computed. It must be noted that this needs to be defined by the user under moments under src and passed to the physical_system object.

f/f_hat: np.ndarray

Pass this argument as well when you want to compute the moments of the input array and not the one stored by the state vector of the object.

Examples

>> solver.compute_moments(‘density’)

Will return the density of the system at its current state.

dump_distribution_function(file_name)

This function is used to dump distribution function to a file for later usage.This dumps the complete 5D distribution function which can be used for post-processing

Parameters:
file_name : The distribution_function array will be dumped to this

provided file name.

Returns:
This function returns None. However it creates a file ‘file_name.h5’,
containing the data of the distribution function

Examples

>> solver.dump_distribution_function(‘distribution_function’)

The above statement will create a HDF5 file which contains the distribution function. The data is always stored with the key ‘distribution_function’

This can later be accessed using

>> import h5py

>> h5f = h5py.File(‘distribution_function’, ‘r’)

>> f = h5f[‘distribution_function’][:]

>> h5f.close()

Alternatively, it can also be used with the load function to resume a long-running calculation.

>> solver.load_distribution_function(‘distribution_function’)

dump_moments(file_name)

This function is used to dump variables to a file for later usage.

Parameters:
file_name : str

The variables will be dumped to this provided file name.

Returns:
This function returns None. However it creates a file ‘file_name.h5’,
containing all the moments that were defined under moments_defs in
physical_system.

Examples

>> solver.dump_variables(‘boltzmann_moments_dump’)

The above set of statements will create a HDF5 file which contains the all the moments which have been defined in the physical_system object. The data is always stored with the key ‘moments’ inside the HDF5 file. Suppose ‘density’ and ‘energy’ are two these moments, and are declared the first and second in the moment_exponents object:

These variables can then be accessed from the file using

>> import h5py

>> h5f = h5py.File(‘boltzmann_moments_dump.h5’, ‘r’)

>> rho = h5f[‘moments’][:][:, :, 0]

>> energy = h5f[‘moments’][:][:, :, 1]

>> h5f.close()

However, in the case of multiple species, the following holds:

>> rho_species_1 = h5f[‘moments’][:][:, :, 0]

>> rho_species_2 = h5f[‘moments’][:][:, :, 1]

>> energy_species_1 = h5f[‘moments’][:][:, :, 2]

>> energy_species_2 = h5f[‘moments’][:][:, :, 3]

get_dist_func()[source]

Returns the distribution function in the same format as the nonlinear solver

load_distribution_function(file_name)

This function is used to load the distribution function from the dump file that was created by dump_distribution_function.

Parameters:
file_name : The distribution_function array will be loaded from this

provided file name.

Examples

>> solver.load_distribution_function(‘distribution_function’)

The above statemant will load the distribution function data stored in the file distribution_function.h5 into self.f