bolt.lib.nonlinear_solver module

This is the module where the main solver object for the nonlinear solver of bolt is defined. This solver object stores the details of the system defined under physical_system, and is evolved using the methods of this module.

The solver has the option of using 2 different methods:

  • A semi-lagrangian scheme based on Cheng-Knorr(1978) which uses advective interpolation.(non-conservative)

    • The interpolation schemes available are linear and cubic spline.
  • Finite volume scheme(conservative):

    • Riemann solvers available are the local Lax-Friedrichs and 1st order upwind scheme.
    • The reconstruction schemes available are minmod, PPM, and WENO5
class bolt.lib.nonlinear.nonlinear_solver.nonlinear_solver(physical_system, performance_test_flag=False)[source]

Bases: object

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

Relevant physical information is obtained by coarse graining this system by taking moments of the distribution function. This is performed by the compute_moments() method.

Methods

compute_moments(moment_name[, f]) Used in computing the moments of the distribution function.
dump_EM_fields(file_name) This function is used to EM fields to a file for later usage.
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.
jia_timestep(dt) Advances the system using the Jia split scheme.
lie_timestep(dt) Advances the system using a lie-split scheme.
load_EM_fields(file_name) This function is used to load the EM fields from the dump file that was created by dump_EM_fields.
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.
print_performance_timings(N_iters) This function is used to check the timings of each of the functions which are used during the process of a single-timestep.
strang_timestep(dt) Advances the system using a strang-split scheme.
swss_timestep(dt) Advances the system using a SWSS-split scheme.
compute_moments(moment_name, f=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: af.Array

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’)

The above line will lookup the definition for ‘density’ and calculate the same accordingly

dump_EM_fields(file_name)

This function is used to EM fields to a file for later usage. This dumps all the EM fields quantities E1, E2, E3, B1, B2, B3 which can then be used later for post-processing

Parameters:
file_name : The EM_fields 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 EM fields.

Examples

>> solver.dump_EM_fields(‘data_EM_fields’)

The above statement will create a HDF5 file which contains the EM fields data. The data is always stored with the key ‘EM_fields’

This can later be accessed using

>> import h5py

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

>> EM_fields = h5f[‘EM_fields’][:]

>> E1 = EM_fields[:, :, 0]

>> E2 = EM_fields[:, :, 1]

>> E3 = EM_fields[:, :, 2]

>> B1 = EM_fields[:, :, 3]

>> B2 = EM_fields[:, :, 4]

>> B3 = EM_fields[:, :, 5]

>> h5f.close()

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

>> solver.load_EM_fields(‘data_EM_fields’)

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.h5’, ‘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’, ‘mom_v1_bulk’ and ‘energy’ are the 3 functions defined. Then the moments get stored in alphabetical order, ie. ‘density’, ‘energy’…:

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]

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

>> h5f.close()

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

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

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

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

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

>> mom_p1_species_1 = h5f[‘moments’][:][:, :, 4]

>> mom_p1_species_2 = h5f[‘moments’][:][:, :, 5]

jia_timestep(dt)

Advances the system using the Jia split scheme. reference:<https://www.sciencedirect.com/science/article/pii/S089571771000436X>

NOTE: This scheme is computationally expensive, and should only be used for testing/debugging

Parameters:
dt : double

Time-step size to evolve the system

lie_timestep(dt)

Advances the system using a lie-split scheme. This scheme is 1st order accurate in time.

Parameters:
dt : double

Time-step size to evolve the system

load_EM_fields(file_name)

This function is used to load the EM fields from the dump file that was created by dump_EM_fields.

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

provided file name.

Examples

>> solver.load_EM_fields(‘data_EM_fields’)

The above statemant will load the EM fields data stored in the file data_EM_fields.h5 into self.cell_centered_EM_fields

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

print_performance_timings(N_iters)

This function is used to check the timings of each of the functions which are used during the process of a single-timestep.

strang_timestep(dt)

Advances the system using a strang-split scheme. This scheme is 2nd order accurate in time.

Parameters:
dt : double

Time-step size to evolve the system

swss_timestep(dt)

Advances the system using a SWSS-split scheme. This scheme is 2nd order accurate in time.

Parameters:
dt : double

Time-step size to evolve the system