Reference on kuibit.cactus_multipoles

This module provides access to data saved by the Multipole thorn.

There are multiple classes defined in this module:

  • :py:class`~.MultipolesDir` interfaces with SimDir and organizes the data according to the variable available. This is a dictionary-like object with keys the variable names.

  • :py:class`~.MultipoleAllDets` takes all the files that correspond to a given variable and organize them according to the extraction radius.

  • :py:class`~.MultipoleOneDet` represents one single extraction radius, for a single variable. It is a dictionary-like object with keys the multipolar numbers and values the multipolar decomposition represented as represented as TimeSeries objects.

These are hierarchical classes, one containing the others, so one typically ends up with a series of brackets or dots to access the actual data. For example, if sim is a SimDir, sim.multipoles['rho_b'][100][2,2] is (2,2) decomposition of rho_b at radius 100 represented as TimeSeries.

class kuibit.cactus_multipoles.MultipoleAllDets(data, l_min=0)[source]

This class collects available surfaces with multipole data.

It is a dictionary-like object with keys the spherical surface radius, and values MultipoleOneDet object. Iteration is supported, sorted by ascending radius. You can iterate over all the radii and all the available l and m with a nested loop.

Variables:
  • radii – Available surface radii.

  • r_outer – Radius of the outermost detector.

  • l_min – l smaller than l_min are dropped.

  • outermost – Outermost detector.

  • available_lm – Available components as tuple (l,m).

  • available_l – List of available l.

  • available_m – List of available m.

Constructor.

Parameters:

data (list of tuples) – List of tuples with (multipole_l, multipole_m, extraction_radius, [timeseries]), where [timeseries] is a list of the TimeSeries associated.

copy()[source]

Return a deep copy.

Returns:

Deep copy of self.

Return type:

MultipoleAllDets

has_detector(mult_l, mult_m, dist)[source]

Check if a given multipole component extracted at a given distance is available.

Parameters:
  • mult_l (int) – Multipole component l.

  • mult_m (int) – Multipole component m.

  • dist (float) – Distance of the detector.

Returns:

If available or not.

Return type:

bool

keys()[source]

Return available extraction radii.

Returns:

Available extraction radii.

Return type:

dict_keys

class kuibit.cactus_multipoles.MultipoleOneDet(dist, data, l_min=0)[source]

This class collects multipole components of a specific variable a given spherical surface.

Multipoles are tightly connected with gravitational waves, so morally a sphere where multipoles are computed is a “detector”. Hence, the name of the class.

MultipoleOneDet is a dictionary-like object with components as the tuples (l,m) and values the corresponding multipolar decomposition as a TimeSeries object. Alternatively, this can also be called directly multipoleonedet(l,m). Iteration is supported and yields tuples (l, m, data), which can be used to loop through all the multipoles available.

The reason we allow for l_min is to remove those files that are not necessary when considering gravitational waves (l<2).

Variables:
  • dist (float) – Radius of the sphere.

  • radius – Radius of the sphere.

  • l_min – l smaller than l_min are dropped.

  • available_l – Available l values.

  • available_m – Available m values.

  • available_m – set

  • available_lm – Available (l, m) values.

  • available_lm – set of tuples

  • missing_lm – Missing (l, m) values to have all from l_min to l_max.

  • missing_lm – set of tuples

Constructor.

Parameters:
  • dist (float) – Radius of the spherical surface.

  • data (list of tuple (l, m, timeseries)) – List of tuples with the two multipolar numbers and the data as TimeSeries.

Variables:

l_min – l smaller than l_min are dropped.

copy()[source]

Return a deep copy.

Returns:

Deep copy of self.

Return type:

MultipoleOneDet

crop(init=None, end=None)[source]

Remove all the data before init and after end.

If init or end are not specified, do not crop on that side.

cropped(init=None, end=None)[source]

Return a copy where the data is cropped in the provided interval.

If init or end are not specified, do not crop on that side.

Returns:

Copy of self with data cropped.

Return type:

MultipoleOneDet

keys()[source]

Return available multipolar numbers.

Returns:

Available multipolar numbers.

Return type:

dict_keys

total_function_on_available_lm(function, *args, l_max=None, **kwargs)[source]

Evaluate function on each multipole and accumulate the result.

total_function_on_available_lm will call function with the following arguments:

function(timeseries, mult_l, mult_m, dist, *args, **kwargs)

If function does not need some paramters, it should use take the *args argument to ignore the additional paramters that are always passed (l, m, r).

Values of l larger than l_max are ignored.

This method is used to compute quantities like the total power in gravitational waves.

function can take additional paramters passed directly from total_function_on_available_lm (e.g. pcut for FFI).

Params function:

Function that has to be applied on each multipole.

Returns:

Sum of function applied to each monopole

Return type:

return type of function

class kuibit.cactus_multipoles.MultipolesDir(sd)[source]

This class provides acces to various types of multipole data in a given simulation directory.

This class is like a dictionary, you can access its values using the brackets operator, with values that are MultipoleAllDets. These contain the full multipolar description for all the available radii. Files are lazily loaded. If both HDF5 and ASCII are present, HDF5 are preferred. There’s no attempt to combine the two. Alternatively, you can access variables with get or with fields.var_name.

Constructor.

Parameters:

sd (SimDir) – Simulation directory.

get(key, default=None)[source]

Return a the multipolar data for the given variable if available, else return the default value.

Parameters:
  • key (str) – Requested variable.

  • default (any) – Returned value if key is not available.

Returns:

Collection of all the multipolar data for the given variable.

Return type:

MultipoleAllDets

keys()[source]

Return available variables with multipolar data.

Returns:

Available variables that have multipolar data.

Return type:

dict_keys