base-metric module

This module contains the class, that defines a general spacetime. All other metric classes derive from it. It has methods, that can be used as utility functions. Users are recommended to inherit this class to create user-defined metric classes.

Docstring for base_metric.py module

This module defines the BaseMetric class, which is the base class for all metrics in EinsteinPy. This class contains several utilities, that are used in einsteinpy.metric to define classes for vacuum solutions. Users are advised to inherit this class, while defining their own metric classes. Two parameters to note are briefly described below:

  • metric_cov: User-supplied function, defining the covariant form of the metric tensor. Users need to supply just this to completely determine the metric tensor, as the contravariant form is calculated and accessed through a predefined method, metric_contravariant().

  • perturbation: User-supplied function, defining a perturbation to the metric. Currently, no checks are performed to ascertain the physicality of the resulting perturbed metric. Read the documentation on metric_covariant() below, to learn more.

Also, note that, users should call metric_covariant to access the perturbed, covariant form of the metric. For unperturbed underlying metric, users should call metric_cov, which returns the metric, that they had supplied.

class einsteinpy.metric.base_metric.BaseMetric(coords, M, a=<Quantity 0.>, Q=<Quantity 0. C>, name='Base Metric', metric_cov=None, christoffels=None, f_vec=None, perturbation=None)[source]

Bases: object

For defining a general Metric

Constructor

Parameters
  • coords (*) – Coordinate system, in which Metric is to be represented

  • M (Quantity) – Mass of gravitating body, e.g. Black Hole

  • a (Quantity) – Dimensionless Spin Parameter Defaults to 0

  • Q (Quantity) – Charge on gravitating body, e.g. Black Hole Defaults to 0 C

  • name (str, optional) – Name of the Metric Tensor. Defaults to "Base Metric"

  • metric_cov (callable, optional) – Function, defining Covariant Metric Tensor It should return a real-valued tensor (2D Array), at supplied coordinates Defaults to None Consult pre-defined classes for function definition

  • christoffels (callable, optional) – Function, defining Christoffel Symbols It should return a real-valued (4,4,4) array, at supplied coordinates Defaults to None Consult pre-defined classes for function definition

  • f_vec (callable, optional) – Function, defining RHS of Geodesic Equation It should return a real-valued (8) vector, at supplied coordinates Defaults to None Consult pre-defined classes for function definition

  • perturbation (callable, optional) – Function, defining Covariant Perturbation Tensor It should return a real-valued tensor (2D Array), at supplied coordinates Defaults to None Function definition similar to metric_cov

static sigma(r, theta, M, a)[source]

Returns the value of r**2 + alpha**2 * cos(theta)**2 Specific to Boyer-Lindquist coordinates Applies to Kerr Geometry

Parameters
Returns

The value of r**2 + alpha**2 * cos(theta)**2

Return type

float

static delta(r, M, a, Q=0)[source]

Returns the value of r**2 - r_s * r + alpha**2 + r_Q**2 Specific to Boyer-Lindquist coordinates Applies to Kerr & Kerr-Newman Geometries

Parameters
Returns

The value of r**2 - r_s * r + alpha**2 + r_Q**2

Return type

float

static rho(r, theta, M, a)[source]

Returns the value of sqrt(r**2 + alpha**2 * cos(theta)**2) == sqrt(sigma) Specific to Boyer-Lindquist coordinates Applies to Kerr-Newman Geometry

Parameters
Returns

The value of sqrt(r**2 + alpha**2 * cos(theta)**2) == sqrt(sigma)

Return type

float

static schwarzschild_radius(M)[source]

Returns Schwarzschild Radius

Parameters

M (float or Quantity) – Mass of gravitating body

Returns

Schwarzschild Radius for a given mass

Return type

float

static alpha(M, a)[source]

Returns Rotational Length Parameter (alpha) that is used in the Metric. Following equations are relevant: alpha = J / Mc a = Jc / GM**2 alpha = GMa / c**2 where, ‘a’ is the dimensionless Spin Parameter (0 <= a <= 1)

Parameters
Returns

Rotational Length Parameter

Return type

float

Raises

ValueError – If a is not between 0 and 1

singularities()[source]

Returns the Singularities of the Metric Depends on the choice of Coordinate Systems Applies to Kerr and Kerr-Newman Geometries

Returns

Dictionary of singularities in the geometry { "inner_ergosphere": function(theta), "inner_horizon": float, "outer_horizon": float, "outer_ergosphere": function(theta) }

Return type

dict

Raises

CoordinateError – If einsteinpy.metric.* does not have the metric in the coordinate system, the metric object has been instantiated with

metric_covariant(x_vec)[source]

Returns Covariant Metric Tensor Adds Kerr-Schild (Linear) Perturbation to metric, if perturbation is defined in Metric object Currently, this does not consider Gauge Fixing or any physical checks for the perturbation matrix. Please exercise caution while using perturbation.

Parameters

x_vec (array_like) – Position 4-Vector

Returns

Covariant Metric Tensor Numpy array of shape (4,4)

Return type

ndarray

metric_contravariant(x_vec)[source]

Returns Contravariant Metric Tensor Adds Kerr-Schild (Linear) Perturbation to metric, if perturbation is not None in Metric object Currently, this does not consider Gauge Fixing or any physical checks for the perturbation matrix. Please exercise caution while using perturbation.

Parameters

x_vec (array_like) – Position 4-Vector

Returns

Contravariant Metric Tensor

Return type

ndarray

calculate_trajectory(start_lambda=0.0, end_lambda=10.0, stop_on_singularity=True, OdeMethodKwargs={'stepsize': 0.001}, return_cartesian=False)[source]

Deprecated in 0.4.0. Please use einsteinpy.Geodesic.

Calculate trajectory in manifold according to geodesic equation