EinsteinPy is an open-source pure Python package dedicated to the study of problems arising
in General Relativity and gravitational physics. Using EinsteinPy, it is possible to approach
problems symbolically as well as numerically.
On the symbolic side, EinsteinPy provides a robust API, which allows users to access several
predefined metrics or to define custom metrics and perform symbolic calculations on them.
Computation of quantities, such as terms of metric tensors, Christoffel symbols, Riemann
Curvature tensor, Ricci tensor, stress-energy tensor and more are all supported and extensible,
since the symbolic modules are built on top of SymPy.
On the numerical side, EinsteinPy provides tools to calculate and visualize geodesics, metric
singularities and hypersurface embeddings in certain spacetimes. We hope to extend the package to
include more features in the future.
EinsteinPy works on recent versions of Python and is released under the MIT license,
hence allowing commercial use of the library.
Most of the general relativity packages currently available in the
gravitational physics research community demand knowledge of
paid-for software / computer algebra systems or programming experience
in languages like C/C++. Many of the people working in the field of
gravitational physics come from a theoretical background and have
little programming experience. As such, they find using these libraries
difficult. EinsteinPy is motivated by these problems and aims to provide
a high level of abstraction, that enables anyone to do general relativity
calculations without being bogged down by the implementation details or
obscure errors. Since the library is written in Python, it is easy to
use and understand, which also makes it a perfect fit for teaching purposes.
Moreover, EinsteinPy being FOSS
also helps interested users understand the code and modify it to suit
their needs, which in turn grows the project and positively adds to the
research community.
EinsteinPy is developed by an open community. Hence, all contributions are more than
welcome! Visit our GitHub to view EinsteinPy’s source code and to report errors or suggest
improvements. For information on contributing code to the project, head to the developer guide.
Please join our Element channel or Gitter chat room for general discussions
and further queries. Release announcements take place on our mailing list.
Please use the links below or on the navigation bar to explore some examples or the API
documentation. If you want to view the documentation for a specific version, please use the
pop-up menu on the bottom-right corner of the page.
This page provides a high-level overview of some of the concepts in
General Relativity. For more detailed treatments, see the references
at the end of this page.
Einstein’s Field Equations (EFE) relate local spacetime curvature
with local energy and momentum. In short, these equations determine the metric tensor
of a spacetime, given the arrangement of stress-energy. The EFE are given as follows:
Here, \(R_{\mu\nu}\) is the Ricci tensor, \(R\) is the
scalar curvature (trace of Ricci tensor), \(g_{\mu\nu}\)
is the metric tensor, \(\Lambda\) is the cosmological constant and
lastly, \(T_{\mu\nu}\) denotes the stress-energy tensor.
All other variables hold their usual meaning. If we introduce the
Einstein tensor \(G_{\mu\nu} = R_{\mu\nu} - \frac{1}{2}Rg_{\mu\nu}\),
then the EFE can be written as:
These equations form a ten-component tensor equation, which collectively
denotes a system of coupled non-linear partial differential equations. These are
usually intractable to approach analytically. However, under certain conditions,
the EFE can be simplified and solved to yield exact solutions.
For example, if the Einstein tensor is assumed to vanish, which implies that the
stress-energy tensor also vanishes, then the EFE reduces to the vacuum Einstein
equations, whose solutions are called vacuum solutions. The Minkowski,
Schwarzschild and Kerr spacetimes are some examples of vacuum solutions.
Similarly, if an electromagnetic field is assumed to be present and also
the only source of non-gravitational energy, then the EFE reduces to the
source-free Maxwell-Einstein equations, whose exact solutions are termed
electrovacuum solutions. The Kerr-Newman and Reissner-Nordström spacetimes are
examples of electrovacuum solutions.
The metric tensor denotes a solution to the EFE. As such, it is a fundamental
entity in general relativity, that captures the geometry of spacetime. It is
a symmetric, indefinite rank-2 tensor, which can be represented by a matrix.
It can be thought of as characterizing the differential line element for
a given geometry:
To get the matrix representation, the coefficients, \(g_{\mu\nu}\), in the above equation
are substituted into the corresponding places designated by the indices,
\(\mu\) and \(\nu\), in a matrix. For example, the metric tensor in the spherical-polar
coordinate system can be written as:
Note that the off-diagonal components are 0, since we are using an
orthogonal basis. However, it is not always the case and in general,
\(g_{\mu\nu} \ne 0\), for \(\mu \ne \nu\).
The metric tensor is also used to define the operations of lowering
and raising indices. For this, we also need the inverse metric tensor, \(g^{\mu\nu}\)
(or contravariant metric tensor to be precise), which is defined using the following identity:
Imagine a bug traveling on a sheet of paper folded into a cone. The
bug can’t see up and down. So, it lives in a 2D world. But it still
experiences the curvature of the surface (space), as, after a long journey,
it would return to the starting position.
Mathematically, the curvature of space is characterized by the rank-4
Riemann Curvature tensor, whose contraction gives the rank-2 Ricci
tensor. Taking the trace of the Ricci tensor yields the rank-0 Ricci
Scalar or Scalar Curvature.
Straight lines are used to describe the shortest path between two points in Euclidean
space. Geodesics extend this notion of “shortest path” to curved spaces. In GR, geodesics
represent the shortest path between two points in a possibly curved spacetime. They are
completely characterized by the metric tensor and its derivatives, resulting in a set
of coupled ordinary differential equations:
Here, \(x\) denotes the geodesic and the derivatives are taken with respect to \(s\), an
affine parameter that uniquely parameterizes \(x\). \(\Gamma^{\mu}_{\alpha\beta}\) denotes
Christoffel symbols of the second kind, which is essentially an array of partial derivatives of the
metric tensor, computed with respect to the coordinate basis. It can be thought of as a “connection”
between the derivates of nearby points along the geodesic and is given as:
The Riemann Curvature tensor encapsulates the idea of curvature in GR. It can be written in
a condensed notation using the Christoffel symbols and their derivaties:
A space with zero curvature implies that the Riemann tensor is zero and vice-versa. Since,
we are dealing with tensors, i.e. objects that operate independently of basis choice, this statement
holds for any coordinate system, i.e. \(R = 0\) in one coordinate system implies \(R = 0\)
and by extension, zero curvature in all coordinate systems.
The Ricci tensor is another geometrical object in GR that is related to curvature. It can be obtained
by contracting the first and third indices of the Riemann tensor:
\[R_{\mu\nu} = R^{\rho}_{\mu\rho\nu}\]
\(R_{\mu\nu}\) can be thought of as quantifying the deformation of a shape as it is
translated along a given geodesic. The trace of the Ricci tensor gives the Scalar Curvature, \(R\):
\[R = g^{\mu\nu}R_{\mu\nu}\]
\(R\) relates the volume of infinitesimal geodesic balls in curved space to that in Euclidean space.
With this, we end our short and superficial look into some of the basic quantities that are used to
characterize the structure of spacetime in General Relativity. Readers, who are interested in gaining a deeper
understanding, are strongly recommended to peruse the resources listed below.
You can install the most recent stable version of EinsteinPy in the
following ways. At the moment, using pip is recommended, since
the package on conda-forge is currently a version behind PyPI.
We are working on updating it.
If you want to contribute to EinsteinPy, please see the Developer Guide
for instructions on how to set up your development environment and start
contributing.
Now that you have EinsteinPy installed, you can try out some of the
examples listed in the Examples section. Below, we demonstrate a simple
example of plotting a precessing timelike geodesic in Schwarzschild spacetime.
This page gives a brief overview of some of the features of EinsteinPy. For more complete hands-on tutorials, please refer to the Jupyter notebooks in Examples.
EinsteinPy provides a way to define the background geometry for relativistic dynamics. This geometry has a central operating quantity known as the Metric tensor, that encapsulates all the geometrical and topological information about the 4D spacetime.
EinsteinPy provides a BaseMetric class, that has various utility functions and a proper template, that can be used to define custom Metric classes. All predefined classes in metric derive from this class.
The central quantities required to simulate the trajectory of a particle in a gravitational field are the metric derivatives, that can be succinctly written using Christoffel Symbols. EinsteinPy provides an easy-to-use interface to calculate these symbols for the predefined metrics. To perform calculations for general metrics, see the symbolic module.
BaseMetric also provides support for f_vec and perturbation, where f_vec corresponds to the RHS of the geodesic equation and perturbation is a linear Kerr-Schild Perturbation that can be defined on the underlying metric. Note that EinsteinPy does not perform physical checks on perturbation currently. Users should exercise caution while using it.
EinsteinPy provides an intuitive interface for calculating timelike and nulllike geodesics for vacuum solutions. Below, we calculate the orbit of a precessing particle in Schwarzschild spacetime.
First, we import all the relevant modules and classes and define the initial conditions for our test particle, such as its initial 3-position in spherical coordinates and corresponding 3-momentum (or 3-velocity, given that we are working in geometric units with \(M = 1\)).
importnumpyasnpfromeinsteinpy.geodesicimportTimelikefromeinsteinpy.plotting.geodesicimportStaticGeodesicPlotterposition=[40.,np.pi/2,0.]momentum=[0.,0.,3.83405]a=0.# Spin = 0 for a Schwarzschild black holesteps=5500# Number of steps to be taken by the solverdelta=1# Step size
After defining our initial conditions, we can use Timelike to create a Geodesic object, that automatically calculates the trajectory.
EinsteinPy currently supports 3 coordinate systems, namely, Cartesian, Spherical and Boyer-Lindquist. The coordinates module provides a way to convert between these coordinate systems. Below, we show how to convert 4-positions and velocities (defined alongside positions) between Cartesian and Boyer-Lindquist coordinates.
importnumpyasnpfromastropyimportunitsasufromeinsteinpy.coordinatesimportBoyerLindquistDifferential,CartesianDifferential,Cartesian,BoyerLindquistM=1e20*u.kga=0.5*u.one# 4-Positiont=1.*u.sx=0.2*u.kmy=0.15*u.kmz=0.*u.km_4pos_cart=Cartesian(t,x,y,z)# The keyword arguments, M & a are required for conversion to & from Boyer-Lindquist coordinates_4pos_bl=_4pos_cart.to_bl(M=M,a=a)print(_4pos_bl)cartsn_pos=_4pos_bl.to_cartesian(M=M,a=a)print(cartsn_pos)# With position & velocityv_x=150*u.km/u.sv_y=250*u.km/u.sv_z=0.*u.km/u.spos_vel_cart=CartesianDifferential(t,x,y,z,v_x,v_y,v_z)# Converting to Boyer-Lindquist coordinatespos_vel_bl=pos_vel_cart.bl_differential(M=M,a=a)print(pos_vel_bl)# Converting back to Cartesian coordinatespos_vel_cart=pos_vel_bl.cartesian_differential(M=M,a=a)print(pos_vel_cart)
You can also pass a einsteinpy.metric.* object to the differential object to set v_t. For usage without units, see the functions in einsteinpy.coordinates.util.
EinsteinPy also supports a robust symbolic module that allows users to access several predefined metrics, or to define custom metrics and perform symbolic calculations. A short example with a predefined metric is shown below.
fromsympyimportsimplifyfromeinsteinpy.symbolicimportSchwarzschild,ChristoffelSymbols,EinsteinTensorm=Schwarzschild()ch=ChristoffelSymbols.from_metric(m)G1=EinsteinTensor.from_metric(m)print(f"ch(1, 2, k) = {simplify(ch[1,2,:])}")# k is an index in [0, 1, 2, 3]print(G1.arr)
Here, we have listed some Jupyter notebooks that exemplify the usage of various features offered by EinsteinPy.
We have broadly divided them into two categories: Symbolic and Numerical computations. To run these notebooks
on the web, please use the binder link provided below.
syms=sympy.symbols("t chi theta phi")t,ch,th,ph=symsm=sympy.diag(-1,cos(t)**2,cos(t)**2*sinh(ch)**2,cos(t)**2*sinh(ch)**2*sin(th)**2).tolist()metric=MetricTensor(m,syms)
Calculating the Einstein Tensor (with both indices covariant)¶
Importing some of the predefined tensors. All the metrics are comprehensively listed in EinsteinPy documentation.¶
[1]:
importsympyfromsympyimportsimplifyfromeinsteinpy.symbolicimportRicciScalarfromeinsteinpy.symbolic.predefinedimportSchwarzschild,DeSitter,AntiDeSitter,Minkowski,findsympy.init_printing()# for pretty printing
The curavture is -12 which is in-line with the theoretical results
Symbolically Understanding Christoffel Symbol and Riemann Curvature Tensor using EinsteinPy¶
[1]:
importsympyfromeinsteinpy.symbolicimportMetricTensor,ChristoffelSymbols,RiemannCurvatureTensorsympy.init_printing()# enables the best printing available in an environment
Defining the metric tensor for 3d spherical coordinates¶
[2]:
syms=sympy.symbols('r theta phi')# define the metric for 3d spherical coordinatesmetric=[[0foriinrange(3)]foriinrange(3)]metric[0][0]=1metric[1][1]=syms[0]**2metric[2][2]=(syms[0]**2)*(sympy.sin(syms[1])**2)# creating metric objectm_obj=MetricTensor(metric,syms)m_obj.tensor()
Calculating the christoffel symbols for Schwarzschild Spacetime Metric¶
The expressions are unsimplified
[7]:
syms=sympy.symbols("t r theta phi")G,M,c,a=sympy.symbols("G M c a")# using metric values of schwarschild space-time# a is schwarzschild radiuslist2d=[[0foriinrange(4)]foriinrange(4)]list2d[0][0]=1-(a/syms[1])list2d[1][1]=-1/((1-(a/syms[1]))*(c**2))list2d[2][2]=-1*(syms[1]**2)/(c**2)list2d[3][3]=-1*(syms[1]**2)*(sympy.sin(syms[2])**2)/(c**2)sch=MetricTensor(list2d,syms)sch.tensor()
syms=sympy.symbols("t chi theta phi")t,ch,th,ph=symsm=sympy.diag(-1,cos(t)**2,cos(t)**2*sinh(ch)**2,cos(t)**2*sinh(ch)**2*sin(th)**2).tolist()metric=MetricTensor(m,syms)
Calculating the Weyl Tensor (with all indices covariant)¶
These notebooks demonstrate the features of the symbolic module in EinsteinPy.
You can use this module to access several pre-defined metrics or define custom ones,
and calculate quantities such as the Christoffel symbols, Ricci tensor, Ricci scalar,
Riemann curvature tensor, Einstein tensor and Weyl tensor. Since this module is based on
SymPy, you can also use the lambdify() function to convert EinsteinPy’s symbolic
expressions into NumPy-compatible functions. Indicial operations are also possible using this
module, as shown in the notebook on contravariant and covariant indices.
%matplotlib widget
# Use %matplotlib ipympl instead if you are working in Visual Studio Code.# You may need to `pip install ipympl` first.sgpl=StaticGeodesicPlotter()sgpl.animate(geod,interval=1)sgpl.show()
The plotted embedding has initial Schwarzschild radial coordinate to be greater than schwarzschild radius but the embedding can be defined for coordinates greater than 9m/4. The Schwarzschild spacetime is a static spacetime and thus the embeddings can be obtained by considering fermat’s surfaces of stationary time coordinate and thus this surface also represent the spacial geometry on which light rays would trace their paths along geodesics of this surface (spacially)!
Shadow cast by a thin emission disk around a Schwarzschild black hole¶
geod=Timelike(metric="Schwarzschild",metric_params=(a,),position=position,momentum=momentum,steps=15543,# As close as we can get before the integration becomes highly unstabledelta=0.0005,return_cartesian=True,omega=0.01,# Small omega values lead to more stable integrationsuppress_warnings=True,# Uncomment to view the tolerance warning)geod
Note that GeodesicPlotter automatically switches between “Static” and “Interactive” plots. Since we are in a Jupyter Notebook or Interactive Environment, it uses the “Interactive” backend.
Data type cannot be displayed: application/vnd.plotly.v1+json
[17]:
gpl.clear()# In Interactive mode, `clear()` must be called before drawing another plot, to avoid overlapgpl.plot2D(geod,coordinates=(1,2),color="green",title="Top/Bottom View")# "top" / "bottom" viewgpl.show()
Data type cannot be displayed: application/vnd.plotly.v1+json
Data type cannot be displayed: application/vnd.plotly.v1+json
Each interactive plot instance (gpl) has a fig attribute that is a plotly.graph_objs._figure.Figure object. So, one can easily make modifications to the plots as they would to any plotly figure. For instance, we can update the title to clarify that the plot above is showing a face-on view.
# Metric or Black Hole parameters - Mass, M and Spin Parameter, aM=4e30*u.kga1=0.4*u.onea2=0.9*u.one# Extremal Kerr Black Hole# Coordinate object to initialize metric with# Note that, for this example# the coordinate values below are irrelevantbl=BoyerLindquistDifferential(t=0.*u.s,r=1e3*u.m,theta=np.pi/2*u.rad,phi=np.pi*u.rad,v_r=0.*u.m/u.s,v_th=0.*u.rad/u.s,v_p=0.*u.rad/u.s,)# Defining two Kerr Black Holes, one with a higher spin parameterkerr1=Kerr(coords=bl,M=M,a=a1)kerr2=Kerr(coords=bl,M=M,a=a2)# Getting the list of singularitiessing_dict1=kerr1.singularities()sing_dict2=kerr2.singularities()# Let's check the contents of the dicts# 'ergosphere' entries should be functionsprint(sing_dict1,sing_dict2,sep="\n\n")
{'inner_ergosphere': <function BaseMetric.singularities.<locals>._in_ergo at 0x7f733c8b0f70>, 'inner_horizon': 247.98878315867296, 'outer_horizon': 5692.939432136259, 'outer_ergosphere': <function BaseMetric.singularities.<locals>._out_ergo at 0x7f730f324af0>}
{'inner_ergosphere': <function BaseMetric.singularities.<locals>._in_ergo at 0x7f733c062310>, 'inner_horizon': 1675.668821582463, 'outer_horizon': 4265.259393712469, 'outer_ergosphere': <function BaseMetric.singularities.<locals>._out_ergo at 0x7f733c0620d0>}
# Sampling Polar Angle for plotting in Polar Coordinatestheta=np.linspace(0,2*np.pi,100)# Ergospheres# These are functionsEi1,Eo1=sing_dict1["inner_ergosphere"],sing_dict1["outer_ergosphere"]Ei2,Eo2=sing_dict2["inner_ergosphere"],sing_dict2["outer_ergosphere"]# Creating lists of points on Ergospheres for different polar angles, for both black holesEi1_list,Eo1_list=Ei1(theta),Eo1(theta)Ei2_list,Eo2_list=Ei2(theta),Eo2(theta)# For Black Hole 1 (a = 0.4)Xei1=Ei1_list*np.sin(theta)Yei1=Ei1_list*np.cos(theta)Xeo1=Eo1_list*np.sin(theta)Yeo1=Eo1_list*np.cos(theta)# For Black Hole 2 (a = 0.9)Xei2=Ei2_list*np.sin(theta)Yei2=Ei2_list*np.cos(theta)Xeo2=Eo2_list*np.sin(theta)Yeo2=Eo2_list*np.cos(theta)# Event HorizonsHi1,Ho1=sing_dict1["inner_horizon"],sing_dict1["outer_horizon"]Hi2,Ho2=sing_dict2["inner_horizon"],sing_dict2["outer_horizon"]# For Black Hole 1 (a = 0.4)Xhi1=Hi1*np.sin(theta)Yhi1=Hi1*np.cos(theta)Xho1=Ho1*np.sin(theta)Yho1=Ho1*np.cos(theta)# For Black Hole 2 (a = 0.9)Xhi2=Hi2*np.sin(theta)Yhi2=Hi2*np.cos(theta)Xho2=Ho2*np.sin(theta)Yho2=Ho2*np.cos(theta)
fig,(ax1,ax2)=plt.subplots(1,2,figsize=(15,7.5))ax1.fill(Xei1,Yei1,'b',Xeo1,Yeo1,'r',Xhi1,Yhi1,'b',Xho1,Yho1,'r',alpha=0.3)ax1.set_title(f"$M = {M}, a = {a1}$",fontsize=18)ax1.set_xlabel("X",fontsize=18)ax1.set_ylabel("Y",fontsize=18)ax1.set_xlim([-6100,6100])ax1.set_ylim([-6100,6100])ax2.fill(Xei2,Yei2,'b',Xeo2,Yeo2,'r',Xhi2,Yhi2,'b',Xho2,Yho2,'r',alpha=0.3)ax2.set_title(f"$M = {M}, a = {a2}$",fontsize=18)ax2.set_xlabel("X",fontsize=18)ax2.set_ylabel("Y",fontsize=18)ax2.set_xlim([-6100,6100])ax2.set_ylim([-6100,6100])
[4]:
(-6100.0, 6100.0)
The surfaces are clearly visible in the plots. Going radially inward, we have Outer Ergosphere, Outer Event Horizon, Inner Event Horizon and Inner Ergosphere. We can also observe the following:
As \(a \to 1\) (its maximum attainable value), the individual singularities become prominent.
As \(a \to 0\), some singularities appear to fade away, leaving us with a single surface, that is the Event Horizon of a Schwarzschild black hole.
Note that, we are working in M-Units (\(G = c = M = 1\)). Also, setting momentum’s \(\phi\)-component to negative, implies an initial retrograde trajectory.
[7]:
position=[2.5,np.pi/2,0.]momentum=[0.,0.,-2.]a=0.99steps=7440# As close as we can get before the integration becomes highly unstabledelta=0.0005omega=0.01suppress_warnings=True
Here, omega, the coupling between the hamiltonian flows, needs to be decreased in order to decrease numerical errors and increase integration stability. Reference: https://arxiv.org/abs/2010.02237.
Also, suppress_warnings has been set to True, as the error would grow exponentially, very close to the black hole.
sgpl=StaticGeodesicPlotter(bh_colors=("red","blue"))sgpl.plot2D(geod,coordinates=(1,2),figsize=(10,10),color="indigo")# Plot X vs Ysgpl.show()
As can be seen in the plot above, the photon’s trajectory is reversed, due to frame-dragging effects, so that, it moves in the direction of the black hole’s spin, before eventually falling into the black hole.
Also, the last few steps seem to have a larger delta, but that is simply because of huge numerical errors, as the particle has crossed the Event Horizon.
Visualizing Precession in Schwarzschild Spacetime¶
Note that, we are working in M-Units (\(G = c = M = 1\)). Also, since the Schwarzschild spacetime has spherical symmetry, the values of the angular components do not affect the end result (We can always rotate our coordinate system to bring the geodesic in the equatorial plane). Hence, we set \(\theta = \pi / 2\) (equatorial plane), with initial \(p_\theta = 0\), which implies, that the geodesic should stay in the equatorial plane.
Apsidal Precession is easily observed through the plot, above, and as expected, the geodesic is confined to the equatorial plane. We can visualize this better, with a few 2D plots.
[5]:
sgpl.clear()sgpl.plot2D(geod,coordinates=(1,2))# Plot X & Ysgpl.show()
[16]:
sgpl.clear()sgpl.plot2D(geod,coordinates=(1,3))# Plot X & Zsgpl.show()
[17]:
sgpl.clear()sgpl.plot2D(geod,coordinates=(2,3))# Plot Y & Zsgpl.show()
Let’s see, how the coordinates change, with \(\lambda\) (Affine Parameter).
[18]:
sgpl.clear()sgpl.parametric_plot(geod,colors=("cyan","magenta","lime"))# Plot X, Y, Z vs Lambda (Affine Parameter)sgpl.show()
The lag between \(X1\) (\(x\)-component) and \(X2\) (\(y\)-component), in the parametric plot, reaffirms the results from above. Also, \(X3 \approx 0\) (\(z\)-component) in this plot, which is expected, as the geodesic never leaves the equatorial plane.
These notebooks show how to use the various numerical methods offered by EinsteinPy.
Using these, you can compute and plot geodesics for manifolds described by (electro-)vacuum
solutions to Einstein’s field equations. You can also calculate and view manifold singularities
in these cases or calculate black hole shadows and hypersurface embeddings for certain spacetimes.
EinsteinPy’s plotting module also enables making animations like these:
Welcome to the API documentation of EinsteinPy. Please navigate through the given modules to get to
know the API of the classes and methods. If you find anything missing, please open an issue in the repo .
This module uses Forward Mode Automatic Differentiation
to calculate metric derivatives. Currently, integrators of
orders 2, 4, 6 and 8 have been implemented.
Geodesic Integrator, based on [1].
This module uses Forward Mode Automatic Differentiation
to calculate metric derivatives to machine precision
leading to stable simulations.
References
Constructor
Parameters:
metric (callable) – Metric Function. Currently, these metrics are supported:
1. Schwarzschild
2. Kerr
3. KerrNewman
metric_params (array_like) – Tuple of parameters to pass to the metric
E.g., (a,) for Kerr
q0 (array_like) – Initial 4-Position
p0 (array_like) – Initial 4-Momentum
time_like (bool, optional) – Determines type of Geodesic
True for Time-like geodesics
False for Null-like geodesics
Defaults to True
steps (int) – Number of integration steps
Defaults to 50
delta (float) – Initial integration step-size
Defaults to 0.5
rtol (float) – Relative Tolerance
Defaults to 1e-2
atol (float) – Absolute Tolerance
Defaults to 1e-2
omega (float) – Coupling between Hamiltonian Flows
Smaller values imply smaller integration error, but too
small values can make the equation of motion non-integrable.
For non-capture trajectories, omega=1.0 is recommended.
For trajectories, that either lead to a capture or a grazing
geodesic, a decreased value of 0.01 or less is recommended.
Defaults to 1.0
suppress_warnings (bool) – Whether to suppress warnings during simulation
Warnings are shown for every step, where numerical errors
exceed specified tolerance (controlled by rtol and atol)
Defaults to False
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.
classeinsteinpy.metric.base_metric.BaseMetric(coords, M, a=<Quantity0.>, Q=<Quantity0.C>, name='BaseMetric', metric_cov=None, christoffels=None, f_vec=None, perturbation=None)[source]¶
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 0C
name (str, optional) – Name of the Metric Tensor. Defaults to "BaseMetric"
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
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)
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)}
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)
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.
The Reissner-Nordström metric in spherical coordinates
A static solution to the Einstein-Maxwell field equations,
which corresponds to the gravitational field of a charged,
non-rotating, spherically symmetric body of mass M.
Parameters:
c (Basic or int or float) – Any value to assign to speed of light. Defaults to c.
G (Basic or int or float) – Any value to assign to the Newton’s (or gravitational) constant. Defaults to G.
eps_0 (Basic or int or float) – Any value to assign to the electric constant or permittivity of free space. Defaults to eps_0.
sch (Basic or int or float) – Any value to assign to Schwarzschild Radius of the central object.
Defaults to r_s.
Q (Basic or int or float) – Any value to assign to eletric charge of the central object.
Defaults to Q.
Exact gravitational wave solution without diffraction.
Class. Quantum Grav., 16:L75-78, 1999.
D. Kramer.
An exact solution describing an axisymmetric gravitational wave propagating in the
z-direction in closed form. This solution to Einstein’s vacuum field equations has the
remarkable property that the curvature invariants decrease monotonically with increasing radial
distance from the axis and vanish at infinity. The solution is regular at the symmetry axis.
Parameters:
C (Basic or int or float) – Constant for Bessel metric, the choice of the constant is not really relavent for details see the paper. Defaults to ‘C’.
syms (tuple or list) – List of crucial symbols dentoting time-axis and/or spacial axis.
For example, in case of 4D space-time, the arrangement would look like [t, x1, x2, x3].
config (str) – Configuration of contravariant and covariant indices in tensor.
‘u’ for upper and ‘l’ for lower indices. Defaults to ‘ll’.
parent_metric (MetricTensor or None) – Metric Tensor for some particular space-time which is associated with this Tensor.
variables (tuple or list or set) – List of symbols used in expressing the tensor,
other than symbols associated with denoting the space-time axis.
Calculates in real-time if left blank.
functions (tuple or list or set) – List of symbolic functions used in epressing the tensor.
Calculates in real-time if left blank.
name (str or None) – Name of the Tensor. Defaults to “GenericTensor”.
Raises:
TypeError – Raised when arr is not a list, sympy array or numpy array.
TypeError – Raised when config is not of type str or contains characters other than ‘l’ or ‘u’
TypeError – Raised when arguments syms, variables, functions have data type other than list, tuple or set.
TypeError – Raised when argument parent_metric does not belong to MetricTensor class and isn’t None.
ValueError – Raised when argument syms does not agree with shape of argument arr
Returns lambdified function of symbolic tensors.
This means that the returned functions can accept numerical values and return numerical quantities.
Parameters:
*args – The variable number of arguments accept sympy symbols.
The returned function accepts arguments in same order as initially defined in *args.
Uses sympy symbols from class attributes syms and variables (in the same order) if no *args is passed
Leaving *args empty is recommended.
Returns:
tuple – arguments to be passed in the returned function in exact order.
function – Lambdified function which accepts and returns numerical quantities.
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘u’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Defaults to None.
This module contains the class for defining a Metric of an arbitrary spacetime, symbolically. Note that, only the coordinate symbols are to
be supplied to the syms parameter, while arr takes the metric (as a SymPy array), which may contain several constants. The symbols in
syms define the basis to perform several operations on the metric, such as symbolic differentiation. Symbols, for the constants in the
metric, should be defined independently and used directly in the specification of arr. Please check the metric definitions in
einsteinpy.symbolic.predefined for examples of doing this.
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘lll’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if not initialized with ‘from_metric’. Defaults to None.
Returns:
New tensor with new configuration. Defaults to ‘lll’
Get Riemann Tensor calculated from Christoffel Symbols.
Reimann Tensor is given as:
\[R^{t}{}_{s r n}=\Gamma^{t}{}_{s n, r} - \Gamma^{t }{}_{s r, n } +
\Gamma^{p}{}_{s n}\Gamma^{t}{}_{p r} - \Gamma^{p}{}_{s r}\Gamma^{t}{}_{p n}\]
Parameters:
chris (ChristoffelSymbols) – Christoffel Symbols from which Riemann Curvature Tensor to be calculated
parent_metric (MetricTensor or None) – Corresponding Metric for the Riemann Tensor.
None if it should inherit the Parent Metric of Christoffel Symbols.
Defaults to None.
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘llll’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if not initialized with ‘from_metric’. Defaults to None.
Returns:
New tensor with new configuration. Configuration defaults to ‘llll’
This module contains the basic classes for obtaining Ricci Tensor and Ricci Scalar
related to a Metric belonging to any arbitrary space-time symbolically:
parent_metric (MetricTensor or None) – Corresponding Metric for the Ricci Tensor.
None if it should inherit the Parent Metric of Riemann Tensor.
Defaults to None.
parent_metric (MetricTensor or None) – Corresponding Metric for the Ricci Tensor.
None if it should inherit the Parent Metric of Christoffel Symbols.
Defaults to None.
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘ul’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if somehow does not have a parent metric. Defaults to None.
Returns:
New tensor with new configuration. Defaults to ‘ul’
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘ul’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if somehow does not have a parent metric. Defaults to None.
Returns:
New tensor with new configuration. Defaults to ‘ul’
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘ul’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if somehow does not have a parent metric. Defaults to None.
Returns:
New tensor with new configuration. Defaults to ‘ul’
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘llll’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if not initialized with ‘from_metric’. Defaults to None.
Returns:
New tensor with new configuration. Configuration defaults to ‘llll’
Changes the index configuration(contravariant/covariant)
Parameters:
newconfig (str) – Specify the new configuration. Defaults to ‘ul’
metric (MetricTensor or None) – Parent metric tensor for changing indices.
Already assumes the value of the metric tensor from which it was initialized if passed with None.
Compulsory if not initialized with ‘from_metric’. Defaults to None.
Returns:
New tensor with new configuration. Configuration defaults to ‘ul’
Constructor
Initialize mass and embedding initial radial coordinate in appropiate units
in order to render the plots of the surface in finite grid. The initial r
is taken to be just greater than schwarzschild radius but it is important
to note that the embedding breaks at r < 9m/4.
Obtain the same values as of the get_values function but reshapes them to obtain
values for all points on the solid of revolution about Z axis (as the
embedding is symmetric in angular coordinates).
Parameters:
alpha (float) – scaling factor to obtain the step size for incrementing r
Returns:
(~numpy.array of X, ~numpy.array of Y, ~numpy.array of Z) values in cartesian coordinates
obtained after applying solid of revolution
This module stores the common utilities of the project, such as a DualNumber
class, that defines dual numbers, used for Forward Mode Auto Differentiation.
Numbers of the form, \(a + b\epsilon\), where
\(\epsilon^2 = 0\) and \(\epsilon \ne 0\).
Their addition and multiplication properties make them
suitable for Automatic Differentiation (AD).
EinsteinPy uses AD for solving Geodesics in arbitrary spacetimes.
This module is based on [1].
This module defines the BaseError class which is the base class for all custom Errors in EinsteinPy, and the CoordinateError class, which is a child class used for raising exceptions when the geometry does not support the supplied coordinate system.
Class for automatically switching between Matplotlib and Plotly depending on platform used.
Constructor
Parameters:
ax (Axes) – Matplotlib Axes object
To be deprecated in Version 0.5.0
Since Version 0.4.0, StaticGeodesicPlotter
automatically creates a new Axes Object.
Defaults to None
bh_colors (tuple, optional) – 2-Tuple, containing hexcodes (Strings) for the colors,
used for the Black Hole Event Horizon (Outer) and Ergosphere (Outer)
Defaults to ("#000","#FFC")
draw_ergosphere (bool, optional) – Whether to draw the ergosphere
Defaults to True
Important Bodies.
Contains some predefined bodies of the Solar System:
* Sun (☉)
* Earth (♁)
* Moon (☾)
* Mercury (☿)
* Venus (♀)
* Mars (♂)
* Jupiter (♃)
* Saturn (♄)
* Uranus (⛢)
* Neptune (♆)
* Pluto (♇)
and a way to define new bodies (Body class).
Data references can be found in constant
Base Class for defining Geodesics
Working in Geometrized Units (M-Units),
with \(c = G = M = k_e = 1\)
Constructor
Parameters:
metric (str) – Name of the metric. Currently, these metrics are supported:
1. Schwarzschild
2. Kerr
3. KerrNewman
metric_params (array_like) – Tuple of parameters to pass to the metric
E.g., (a,) for Kerr
position (array_like) – 3-Position
4-Position is initialized by taking t=0.0
momentum (array_like) – 3-Momentum
4-Momentum is calculated automatically,
considering the value of time_like
time_like (bool, optional) – Determines type of Geodesic
True for Time-like geodesics
False for Null-like geodesics
Defaults to True
return_cartesian (bool, optional) – Whether to return calculated positions in Cartesian Coordinates
This only affects the coordinates. Momenta are dimensionless
quantities, and are returned in Spherical Polar Coordinates.
Defaults to True
kwargs (dict) – Keyword parameters for the Geodesic Integrator
See ‘Other Parameters’ below.
steps (int) – Number of integration steps
Defaults to 50
delta (float) – Initial integration step-size
Defaults to 0.5
rtol (float) – Relative Tolerance
Defaults to 1e-2
atol (float) – Absolute Tolerance
Defaults to 1e-2
omega (float) – Coupling between Hamiltonian Flows
Smaller values imply smaller integration error, but too
small values can make the equation of motion non-integrable.
For non-capture trajectories, omega=1.0 is recommended.
For trajectories, that either lead to a capture or a grazing
geodesic, a decreased value of 0.01 or less is recommended.
Defaults to 1.0
suppress_warnings (bool) – Whether to suppress warnings during simulation
Warnings are shown for every step, where numerical errors
exceed specified tolerance (controlled by rtol and atol)
Defaults to False
omega (float) – Coupling between Hamiltonian Flows
Smaller values imply smaller integration error, but too
small values can make the equation of motion non-integrable.
For non-capture trajectories, omega=1.0 is recommended.
For trajectories, that either lead to a capture or a grazing
geodesic, a decreased value of 0.01 or less is recommended.
Defaults to 1.0
suppress_warnings (bool) – Whether to suppress warnings during simulation
Warnings are shown for every step, where numerical errors
exceed specified tolerance (controlled by rtol and atol)
Defaults to False
metric (str) – Name of the metric. Currently, these metrics are supported:
1. Schwarzschild
2. Kerr
3. KerrNewman
metric_params (array_like) – Tuple of parameters to pass to the metric
E.g., (a,) for Kerr
position (array_like) – 3-Position
4-Position is initialized by taking t=0.0
momentum (array_like) – 3-Momentum
4-Momentum is calculated automatically,
considering the value of time_like
return_cartesian (bool, optional) – Whether to return calculated positions in Cartesian Coordinates
This only affects the coordinates. The momenta dimensionless
quantities, and are returned in Spherical Polar Coordinates.
Defaults to True
kwargs (dict) – Keyword parameters for the Geodesic Integrator
See ‘Other Parameters’ below.
steps (int) – Number of integration steps
Defaults to 50
delta (float) – Initial integration step-size
Defaults to 0.5
rtol (float) – Relative Tolerance
Defaults to 1e-2
atol (float) – Absolute Tolerance
Defaults to 1e-2
omega (float) – Coupling between Hamiltonian Flows
Smaller values imply smaller integration error, but too
small values can make the equation of motion non-integrable.
For non-capture trajectories, omega=1.0 is recommended.
For trajectories, that either lead to a capture or a grazing
geodesic, a decreased value of 0.01 or less is recommended.
Defaults to 1.0
suppress_warnings (bool) – Whether to suppress warnings during simulation
Warnings are shown for every step, where numerical errors
exceed specified tolerance (controlled by rtol and atol)
Defaults to False
metric (str) – Name of the metric. Currently, these metrics are supported:
1. Schwarzschild
2. Kerr
3. KerrNewman
metric_params (array_like) – Tuple of parameters to pass to the metric
E.g., (a,) for Kerr
position (array_like) – 3-Position
4-Position is initialized by taking t=0.0
momentum (array_like) – 3-Momentum
4-Momentum is calculated automatically,
considering the value of time_like
return_cartesian (bool, optional) – Whether to return calculated positions in Cartesian Coordinates
This only affects the coordinates. The momenta dimensionless
quantities, and are returned in Spherical Polar Coordinates.
Defaults to True
kwargs (dict) – Keyword parameters for the Geodesic Integrator
See ‘Other Parameters’ below.
steps (int) – Number of integration steps
Defaults to 50
delta (float) – Initial integration step-size
Defaults to 0.5
rtol (float) – Relative Tolerance
Defaults to 1e-2
atol (float) – Absolute Tolerance
Defaults to 1e-2
omega (float) – Coupling between Hamiltonian Flows
Smaller values imply smaller integration error, but too
small values can make the equation of motion non-integrable.
For non-capture trajectories, omega=1.0 is recommended.
For trajectories, that either lead to a capture or a grazing
geodesic, a decreased value of 0.01 or less is recommended.
Defaults to 1.0
suppress_warnings (bool) – Whether to suppress warnings during simulation
Warnings are shown for every step, where numerical errors
exceed specified tolerance (controlled by rtol and atol)
Defaults to False
An example to showcase the usage of the various modules in einsteinpy.
Here, we assume a Schwarzschild spacetime and obtain a test particle orbit, that
shows apsidal precession.
Returns:
geod – Timelike Geodesic, defining test particle trajectory
As with any open-source project, we rely on the community for implementing
new ideas and fixing old bugs. You can contribute to EinsteinPy in many ways,
for example, by reporting errors in code or documentation, writing new code,
improving documentation, or even by writing tutorials or blog posts. If you
are interested in contributing to EinsteinPy, please go through the following
sections.
For starters, we recommend checking out the “good first issue” tag on
the issue tracker. Those issues should be relatively easy to fix
and would require minimal knowledge about the library. However, you
will need to possess some knowledge about how git works. git
is a decentralized version control system that preserves codebase history,
keeps track of changes and allows for collaboration. If you are new to
git and version control, try following the GitHub Quickstart tutorial.
If you already know all this and would like to contribute, then
that’s awesome! But before coding away, please make sure your
intended addition or change is in line with the project’s scope and goals.
You can do this by chatting with us on Element or Gitter, or mailing us.
All code changes & additions should be accompanied by altered or new
unit-tests, so that the code coverage increases or stays the same. Automated
services will ensure that your code works across the various platforms that
EinsteinPy supports.
Note that EinsteinPy is a Python 3-only library. So, you need to have Python 3
installed in your system and you need
to be familiar with it. Then, you can follow these steps to set up a
platform-agnostic development environment and start contributing to EinsteinPy:
Create a virtual environment and activate it using the following commands:
$python-mvenv<your-venv-name>
$source<your-venv-name>/bin/activate# Linux/macOS
$./<your-venv-name>/Scripts/activate# Windows
Learn more about Python virtual environments here.
After activating the environment, install EinsteinPy in editable or development mode like so:
$pipinstall-e./einsteinpy/[dev]
[dev] ensures that all the dependencies required for development and
testing are installed, while the -e or --editable flag ensures that any changes you make to the code
take effect immediately (after you save the changes).
Create and switch to a new branch like so:
$gitcheckout-b<your-branch-name>
Usually, the branch name should be kept similar to the issue/topic you are working on.
Make changes to the code and add unit-tests (if applicable). Ensure that the code is formatted properly. See the Code Linting & Testing section below.
After you are done, commit your changes and push them to your forked repository.
Finally, open a Pull Request (PR) to the main branch of the einsteinpy/einsteinpy repository. You can also request a review from specific maintainers of the project. The maintainers will review your PR and might ask you to make some changes. If there are any changes required, you can make them in the same branch and push them to your forked repository. The PR will automatically update with the new changes.
When your PR is approved and ready to be merged, we will ask you to update the project CHANGELOG and if needed, we might also ask you to squash your commits.
If you are facing issues with any of the above steps, please feel free to ask for help on Element or Gitter.
After you push your code and open a PR to the main branch, your code will be
tested automatically by the continuous integration (CI) tools. If they encounter any errors,
you will be notified by CI with the details about the error, which
you can use to make necessary changes to fix the error. These errors usually come
in the form of code quality errors, failed tests and decreased coverage.
Code quality is an important aspect of any software project. It ensures that
the code is readable, maintainable, and bug-free. The quality of a piece of code
depends on its formatting, style and complexity. To maintain consistent code
formatting and style, we use black, isort and mypy. For convenience,
we have set up tox, which runs all these in a single, short command:
$cd./einsteinpy/
$tox-ereformat
If your PR is failing quality checks, executing the above commands should fix it.
Also, to ensure that the code remains understandable and maintainable over time,
we enforce Cyclomatic Complexity (CC) checks on the codebase using CodeFactor/CodeClimate.
CC is a measure of the complexity of a program. The lower the CC, the easier it is to
understand and maintain the code. If your PR is running into CodeFactor or CodeClimate
errors, you will have to refactor your code, so that its CC is below a certain threshold.
To check the CC of your code locally, you can use radon:
Since you have made changes to the codebase, it is likely that some of the
unit-tests that were previously passing will now fail. CI will alert you to these
failures on your PR with the error details. You can use these details to
fix the errors, commit & push the changes and make the tests pass. In case you want to
debug the test errors locally, you can do so by running the following command:
Note that this can take a while. If the error is isolated to only one file or a few files,
you can choose to test only that or those file(s). For example, to execute all the
einsteinpy.metric-related tests, you can run:
This command also reports the overall code coverage after all the tests are run. Code coverage
is a measure of how much of the codebase is covered by the tests. It is a good
practice to have a high code coverage, so that the tests can catch any bugs that
might be introduced in the codebase. We use codecov to track the code coverage
of the project. You can see the current code coverage of the project
here. If your PR is failing coverage
checks, you will have to add more tests to increase the coverage.
After you have implemented your bugfix or your shiny new feature, you should also
add some documentation, so that users, maintainers and future contributors can
understand how to use or make changes to your code. All of EinsteinPy’s non-API
documentation is stored in text files under docs/source.
If you think anything can be improved there, please edit the files and open a PR.
The API docs comprising the docstrings of the Python code follow
numpydoc guidelines.
If you come across any inconsistency or opportunity for improvement, feel free to
edit the docstrings and submit a PR.
We use Sphinx to generate the overall API + non-API documentation, which is then
hosted here, courtesy of Read the Docs. To
build the docs locally, you can run the following commands:
This should open the built documentation website in a web browser at
http://127.0.0.1:8000. If you want hot-reloading, you can use
sphinx-autobuild instead of sphinx-build, after installing it using
pipinstallsphinx-autobuild. This will automatically rebuild the docs and
refresh the browser tab whenever you make changes to the source files.
In addition to the usual documentation, the GitHub Wiki for EinsteinPy is open to everybody. Please feel free to add
new content there.
Great job 🎊🎊. Your PR just got approved & merged. Now how do you ensure your local
main branch is up-to-date with upstream/main? Like so:
$gitremoteaddupstreamhttps://github.com/einsteinpy/einsteinpy.git# Set up upstream
$gitcheckoutmain
$gitfetchupstream
$gitmergeupstream/main
$gitbranch-d<your-branch-name># Delete your branch
$gitpushoriginmain
And now, you are all set to create another branch and contribute further to EinsteinPy. Have fun coding 😀!
Not only can things break at any given time, but different people also have different
use cases for the project. If you find anything that doesn’t work as expected
or have suggestions, please refer to the issue tracker on GitHub.
This major release brings a lot of improvements to the numerical module of the project. COVID19 is proving
to be very difficult for India and we are trying to cope up. Please forgive us for any issues you faced with
the v0.3.1 and the documentation.
This major release would bring some very important improvements. This release fixes a very crucial
bug with sympy. Fixes coordinate conversions so they don’t fail on edge cases anymore.
EinsteinPy now uses GitHub Actions for macOS builds. Big changes to the plotting module.
The release comes for the paper of EinsteinPy. The release marks the beginning of Google Summer of Code 2020.
The release also brings a new rays module, which will form the base for null geodesics in future releases.
This minor release would bring improvements and new feature additions to the already existing symbolic calculations module along
with performance boosts of order of 15x.
This release brings a lot of new features for the EinsteinPy Users.
A better API, intuitive structure and easy coordinates handling! This major release
comes before Python in Astronomy 2019 workshop and brings a lots of cool stuff.
Part of this release is sponsored by ESA/ESTEC - Adv. Concepts & Studies Office
(European Space Agency), through Summer of Code in Space (SOCIS) 2019 program.
This is a short-term supported version and will be supported only until December 2019.
For any feature request, write a mail to developers@einsteinpy.org describing what you need.
The old StaticGeodesicPlotter has been renamed to
einsteinpy.plotting.senile.StaticGeodesicPlotter, please adjust
your imports accordingly
The old ScatterGeodesicPlotter has been renamed to
einsteinpy.plotting.senile.ScatterGeodesicPlotter, please adjust
your imports accordingly.
einsteinpy.metric.Schwarzschild,
einsteinpy.metric.Kerr, and
einsteinpy.metric.KerrNewman now have different signatures for
class methods, and they now explicitly support einsteinpy.coordinates
coordinate objects. Check out the notebooks and their respective documentation.
The old coordinates conversion in einsteinpy.utils has been deprecated.
The old symbolic module in einsteinpy.utils has been moved to
einsteinpy.symbolic.
This is a major first release for world’s first actively maintained python library
for General Relativity and Numerical methods. This major release just comes before
the Annual AstroMeet of IIT Mandi, AstraX. This will be a short term support version
and will be supported only until late 2019.
The community of participants in EinsteinPy is made up of members from around
the globe with a diverse set of skills, personalities and experiences. It is
through these differences that our community experiences success and continued
growth. We expect everyone in our community to follow these guidelines when
interacting with others both inside and outside of our community. Our goal is
to keep ours a positive, inclusive, successful and growing community.
As members of the community,
we pledge to treat all people with respect and provide a harassment- and
bullying-free environment, regardless of sex, sexual orientation and/or
gender identity, disability, physical appearance, body size, race,
nationality, ethnicity, and religion. In particular, sexual language and
imagery, sexist, racist, or otherwise exclusionary jokes are not appropriate.
we pledge to respect the work of others by recognizing
acknowledgment/citation requests of original authors. As authors, we pledge
to be explicit about how we want our work to be cited or acknowledged.
we pledge to welcome those interested in joining the community, and realize
that including people with a variety of opinions and backgrounds will only
serve to enrich our community. In particular, discussions relating to
pros/cons of various technologies, programming languages, and so on are
welcome, but these should be done with respect, taking proactive measures to
ensure that all participants are heard and feel confident that they can
freely express their opinions.
we pledge to welcome questions and answer them respectfully, paying
particular attention to those new to the community. We pledge to provide
respectful criticisms and feedback in forums, especially in discussion
threads resulting from code contributions.
we pledge to be conscientious of the perceptions of the wider community and
to respond to criticism respectfully. We will strive to model behaviors that
encourage productive debate and disagreement, both within our community and
where we are criticized. We will treat those outside our community with the
same respect as people within our community.
we pledge to help the entire community follow the code of conduct, and to
not remain silent when we see violations of the code of conduct. We will take
action when members of our community violate this code such as contacting
shreyas@einsteinpy.org (all emails sent to this address will be treated with
the strictest confidence) or talking privately with the person.
This code of conduct applies to all community situations online and offline,
including mailing lists, forums, social media, conferences, meetings,
associated social events, and one-to-one interactions.
Parts of this code of conduct have been adapted from the Astropy code of
conduct.
The Code of Conduct for the EinsteinPy community is licensed under a Creative
Commons Attribution 4.0 International License. We encourage other communities
related to ours to use or adapt this code as they see fit.