Gyoto
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | Friends | List of all members
Gyoto::Metric::Python Class Reference

Metric coded in Python. More...

#include <GyotoPython.h>

Inheritance diagram for Gyoto::Metric::Python:
Gyoto::Python::Object< Gyoto::Metric::Generic > Gyoto::Metric::Generic Gyoto::Python::Base Gyoto::SmartPointee Gyoto::Object Gyoto::Hook::Teller

Public Types

typedef Gyoto::SmartPointer< Gyoto::SmartPointeeSubcontractor_t(Gyoto::FactoryMessenger *, std::vector< std::string > const &)
 A subcontractor builds an object upon order from the Factory.
 

Public Member Functions

virtual Property const * getProperties () const
 Get list of properties.
 
 Python (const Python &)
 
virtual Pythonclone () const
 Virtual copy constructor.
 
void spherical (bool)
 
bool spherical () const
 
virtual std::string module () const
 Return module_.
 
virtual void module (const std::string &)
 Set module_ and import the Python module.
 
virtual std::string inlineModule () const
 Return inline_module_.
 
virtual void inlineModule (const std::string &)
 Set inline_module_ and import the Python module.
 
virtual std::string klass () const
 Retrieve class_.
 
virtual void klass (const std::string &)
 Set class_ and instantiate the Python class.
 
virtual std::vector< double > parameters () const
 Retrieve parameters_.
 
virtual void parameters (const std::vector< double > &)
 Set parameters_ and send them to pInstance_.
 
virtual void mass (double m)
 Set mass used in unitLength()
 
void gmunu (double g[4][4], const double *x) const
 
int christoffel (double dst[4][4][4], const double *x) const
 
double getRmb () const
 
double getRms () const
 
double getSpecificAngularMomentum (double rr) const
 
double getPotential (double const pos[4], double l_cst) const
 
int isStopCondition (double const coord[8]) const
 Check whether integration should stop.
 
void circularVelocity (double const pos[4], double vel[4], double dir=1.) const
 Yield circular velocity at a given position.
 
virtual void set (std::string const &key, Value val)
 
virtual void set (Property const &p, Value val)
 
virtual void set (Property const &p, Value val, std::string const &unit)
 
virtual void set (std::string const &pname, Value val, std::string const &unit)
 Set Value (expressed in unit) of a Property.
 
virtual Value get (std::string const &key) const
 
Value get (Property const &p, std::string const &unit) const
 
Value get (Property const &p) const
 
virtual Value get (std::string const &pname, std::string const &unit) const
 Get Value of a Property, converted to unit.
 
virtual int setParameter (std::string name, std::string content, std::string unit)
 
virtual void setParameter (Gyoto::Property const &p, std::string const &name, std::string const &content, std::string const &unit)
 Set parameter by Property (and name)
 
virtual void fillElement (Gyoto::FactoryMessenger *fmp) const
 
void setParameters (Gyoto::FactoryMessenger *fmp)
 
int coordKind () const
 Get coordinate kind.
 
int getRefCount ()
 
virtual void mass (const double, const std::string &unit)
 Set mass used in unitLength()
 
double mass () const
 Get mass used in unitLength()
 
double mass (const std::string &unit) const
 Get mass used in unitLength()
 
double unitLength () const
 M * G / c^2, M is in kg, unitLength in meters.
 
double unitLength (const std::string &unit) const
 unitLength expressed in specified unit
 
double deltaMin () const
 
void deltaMin (double h1)
 
double deltaMax () const
 
virtual double deltaMax (double const pos[8], double delta_max_external) const
 
void deltaMax (double h1)
 
double deltaMaxOverR () const
 Get delta_max_over_r_.
 
void deltaMaxOverR (double t)
 Set delta_max_over_r_.
 
bool keplerian () const
 Get keplerian_.
 
void keplerian (bool)
 Set keplerian_.
 
virtual void cartesianVelocity (double const coord[8], double vel[3])
 Compute xprime, yprime and zprime from 8-coordinates.
 
virtual double SysPrimeToTdot (const double coord[4], const double v[3]) const
 Compute tdot as a function of dr/dt, dtheta/dt and dphi/dt. Everything is in geometrical units.
 
virtual void zamoVelocity (double const pos[4], double vel[4]) const
 Yield ZAMO velocity at a given position.
 
virtual void nullifyCoord (double coord[8]) const
 Set tdot (coord[4]) such that coord is light-like. Everything is in geometrical units.
 
virtual void nullifyCoord (double coord[8], double &tdot2) const
 Set tdot (coord[4]) such that coord is light-like and return other possible tdot.
 
virtual void normalizeFourVel (double coord[8]) const
 Normalize fourvelvel to -1.
 
virtual void normalizeFourVel (double const pos[4], double fourvel[4]) const
 Normalize fourvelvel to -1.
 
virtual double ScalarProd (const double pos[4], const double u1[4], const double u2[4]) const
 Scalar product.
 
double norm (const double pos[4], const double u1[4]) const
 Scalar product.
 
void multiplyFourVect (double vect[4], double a) const
 multiply vector by scalar
 
void addFourVect (double u1[4], double const u2[4]) const
 add second vector to first one
 
void projectFourVect (double const pos[4], double u1[4], double const u2[4]) const
 project u1 orthogonally to u2 at pos
 
void dualOneForm (double const IN_ARRAY1_1[4], double const IN_ARRAY1_2[4], double ARGOUT_ARRAY1[4]) const
 Computes dual 1-form Compute the dual 1-form of 4-vector.
 
virtual void observerTetrad (obskind_t obskind, double const pos[4], double fourvel[4], double screen1[4], double screen2[4], double screen3[4]) const
 Computes the orthonormal local tetrad of the observer.
 
virtual void observerTetrad (double const pos[4], double fourvel[4], double screen1[4], double screen2[4], double screen3[4]) const
 Computes the orthonormal local tetrad of the observer.
 
void GramSchmidt (double const pos[4], double u0[4], double u1[4], double u2[4], double u3[4]) const
 Apply Gram-Schmidt orthonormalization to a basis.
 
virtual double gmunu (double const x[4], int mu, int nu) const
 Metric coefficients.
 
virtual void gmunu (double ARGOUT_ARRAY2[4][4], double const IN_ARRAY1[4]) const
 Metric coefficients.
 
virtual double gmunu_up (double const x[4], int mu, int nu) const
 Metric contravariant coefficients.
 
virtual void gmunu_up (double ARGOUT_ARRAY2[4][4], const double IN_ARRAY1[4]) const
 Metric contravariant coefficients.
 
virtual void jacobian (double ARGOUT_ARRAY3[4][4][4], const double IN_ARRAY1[4]) const
 Derivatives of the metric covariant coefficients.
 
virtual void gmunu_up_and_jacobian (double ARGOUT_ARRAY2[4][4], double ARGOUT_ARRAY3[4][4][4], const double IN_ARRAY1[4]) const
 gmunu_up() and jacobian() in one go
 
virtual void computeNBeta (const double coord[4], double &NN, double beta[3]) const
 Computes lapse scalar and shift vector at coord.
 
virtual double christoffel (const double coord[4], const int alpha, const int mu, const int nu) const
 Chistoffel symbol.
 
virtual int christoffel (double dst[4][4][4], const double coord[4]) const
 Chistoffel symbol.
 
virtual int myrk4 (Worldline *line, state_t const &coord, double h, state_t &res) const
 RK4 integrator.
 
virtual int myrk4 (Worldline *line, const double coord[8], double h, double res[8]) const =delete
 Obsolete, update your code.
 
virtual int myrk4_adaptive (Gyoto::Worldline *line, state_t const &coord, double lastnorm, double normref, state_t &coordnew, double h0, double &h1, double deltamax=GYOTO_DEFAULT_DELTA_MAX) const
 RK4 integrator with adaptive step.
 
virtual int myrk4_adaptive (Gyoto::Worldline *line, const double coord[8], double lastnorm, double normref, double coordnew[8], double h0, double &h1, double deltamax=GYOTO_DEFAULT_DELTA_MAX) const =delete
 Obsolete, update your code.
 
virtual int diff (state_t const &x, state_t &dxdt, double mass) const
 F function such as dx/dt=F(x,cst)
 
virtual int diff (state_t const &x, state_t &dxdt) const =delete
 Obsolete, update your code.
 
virtual int diff (const double y[8], double res[8]) const =delete
 
virtual int diff31 (state_t const &x, state_t &dxdt, double mass) const
 F function such as dx/dt=F(x,cst) for 3+1 case.
 
virtual void setParticleProperties (Gyoto::Worldline *line, double const coord[8]) const
 Set Metric-specific constants of motion. Used e.g. in KerrBL.
 
void incRefCount ()
 Increment the reference counter. Warning: Don't mess with the counter.
 
int decRefCount ()
 Decrement the reference counter and return current value. Warning: Don't mess with the counter.
 
virtual bool isThreadSafe () const
 Whether this class is thread-safe.
 
Property const * property (std::string const pname) const
 Find property by name.
 
virtual void fillProperty (Gyoto::FactoryMessenger *fmp, Property const &p) const
 Output a single Property to XML.
 
std::string describeProperty (Gyoto::Property const &p) const
 Format desrciption for a property.
 
void help () const
 Print (to stdout) some help on this class.
 
virtual std::string kind () const
 Get kind_.
 
virtual void hook (Listener *listener)
 Start listening.
 
virtual void unhook (Listener *listener)
 Stop listening.
 
virtual bool hasPythonProperty (std::string const &key) const
 
virtual void setPythonProperty (std::string const &key, Value val)
 
virtual Value getPythonProperty (std::string const &key) const
 
virtual int pythonPropertyType (std::string const &key) const
 

Public Attributes

 GYOTO_OBJECT_THREAD_SAFETY
 

Static Public Attributes

static GYOTO_OBJECT Property const properties []
 

Protected Member Functions

void coordKind (int coordkind)
 Set coordkind_.
 
virtual void kind (const std::string)
 Set kind_.
 
virtual void tellListeners ()
 Call tell() on each hooked Listener.
 

Protected Attributes

double delta_min_
 Minimum integration step for the adaptive integrator.
 
double delta_max_
 Maximum integration step for the adaptive integrator.
 
double delta_max_over_r_
 Numerical tuning parameter.
 
bool keplerian_
 1 if circularVelocity should return the Newtonian Keplerian velocity, in r^-3/2
 
std::string kind_
 The "kind" that is output in the XML entity.
 
std::vector< std::string > plugins_
 The plug-ins that needs to be loaded to access this instance's class.
 
std::string module_
 Name of the Python module that holds the class.
 
std::string inline_module_
 Python source code for module that holds the class.
 
std::string class_
 Name of the Python class that we want to expose.
 
std::vector< double > parameters_
 Parameters that this class needs.
 
PyObject * pModule_
 Reference to the python module once it has been loaded.
 
PyObject * pInstance_
 Reference to the python instance once it has been instantiated.
 
PyObject * pProperties_
 Reference to the properties member.
 
PyObject * pSet_
 Reference to the (optional) Set method.
 
PyObject * pGet_
 Reference to the (optional) Get method.
 

Private Attributes

PyObject * pGmunu_
 Reference to the gmunu method.
 
PyObject * pChristoffel_
 Reference to the christoffel method.
 
PyObject * pGetRmb_
 Reference to the getRmb method.
 
PyObject * pGetRms_
 Reference to the getRms method.
 
PyObject * pGetSpecificAngularMomentum_
 Reference to the method getSpecificAngularMomentum.
 
PyObject * pGetPotential_
 Reference to the getPotential method.
 
PyObject * pIsStopCondition_
 Reference to the isStopCondition method.
 
PyObject * pCircularVelocity_
 Reference to the circularVelocity method.
 
double mass_
 Mass yielding geometrical unit (in kg).
 
int coordkind_
 Kind of coordinates (cartesian-like, spherical-like, unspecified)
 
int __defaultfeatures
 Whether some virtual methods are implemented.
 
int refCount
 Reference counter.
 
pthread_mutex_t mutex_
 A mutex.
 
ListenerItemlisteners_
 Linked list of Listener items.
 

Friends

class Gyoto::SmartPointer< Gyoto::Metric::Python >
 

Detailed Description

Metric coded in Python.

Loader for Python Metric classes. It interfaces with a Python class which must implement at least the gmunu and christoffel methods.

Other methods are optional: getRmb, getRms, getSpecificAngularMomentum, getPotential.

Use <Cartesian> or </Spherical> to select the coordinate system kind.

Sample XML file:

<?xml version="1.0" encoding="UTF-8" standalone="no"?>
<Scenery>
The goal of this example is to exhibit using a Metric and a Spectrum
coded in Python. It is functionaly equivalent to
example-fixed-star-minkowski-cartesian.xml in the main Gyoto example
directory.
An optically thin blob centered on the origin of the coordinate
system, in flat space-time. Computation uses the Cartesian
coordinates.
<Metric kind = "Python">
<Mass unit="sunmass"> 4e6 </Mass>
<Cartesian/>
<Module>gyoto_sample_metrics</Module>
<Class>Minkowski</Class>
</Metric>
<Screen>
<Distance unit="kpc"> 8 </Distance>
<Time unit="yr"> 30e3 </Time>
<FieldOfView unit="microas"> 150 </FieldOfView>
In UTF-8 locales, "microas" may be written "µas".
<Inclination unit="degree"> 90 </Inclination>
In UTF-8 locales, "degree" may be written "°".
<PALN> 0 </PALN>
<Argument> 0 </Argument>
<Resolution> 32 </Resolution>
</Screen>
<Astrobj kind = "FixedStar">
<Radius> 12 </Radius>
<Position> 0 0 0 </Position>
<UseGenericImpact/>
<Spectrum kind="Python">
<Module>gyoto_sample_spectra</Module>
<Class>PowerLaw</Class>
<Constant>0.001</Constant>
<Exponent>0.</Exponent>
</Spectrum>
<Opacity kind="Python">
<Module>gyoto_sample_spectra</Module>
<Class>PowerLaw</Class>
<Constant>0.001</Constant>
<Exponent>0.</Exponent>
</Opacity>
<OpticallyThin/>
</Astrobj>
<Delta> 1e0 </Delta>
<MinimumTime> 0. </MinimumTime>
<Quantities>
Intensity
</Quantities>
One can also specify a unit (if Gyoto was compiled with --with-udunits):
Intensity[mJy/pix²]
Intensity[mJy/µas²]
Intensity[J.m-2.s-1.sr-1.Hz-1]
Intensity[erg.cm-2.s-1.sr-1.Hz-1]
Intensity[mJy.sr-1]
Intensity[Jy.sr-1]
</Scenery>

Sample Python module:

'''Sample metrics for using with Gyoto Python plug-in
Those classes demonstrate how to use Python classes as Gyoto
Metric implementations using Gyoto's "python" plug-in. Note that
this plug-in can be renamed to whatever matches the particular
version of Python it has been built against (e.g. python3.4).
The goal is to be able to instantiate these from XML, from
Yorick... and even from Python using the gyoto extension...
Synopsis:
import gyoto.core
gyoto.core.requirePlugin("python") # or python2.7 or python3.4...
gg=gyoto.core.Metric("Python")
gg.set("Module", "gyoto_sample_metric")
gg.set("Class", "Minkowski")
Classes that aim at implementing the Gyoto::Metric::Generic
interface do so by providing the following methods:
gmunu(self, dst, pos): mandatory;
christoffel(self, dst, pos): mandatory;
getRms(self): optional
getRmb(self): optional
getSpecificAngularMomentum(self, rr): optional
getPotential(self, pos, l_cst): optional
__setattr__(self, key, value): optional, useful to react
when the C++ Metric object sets attributes:
this: if the Python extension "gyoto.core" can be imported,
it will be set to a gyoto.core.Metric instance
pointing to the C++-side instance. If the "gyoto.core"
extension cannot be loaded, this will be set to None.
spherical: when the spherical(bool t) method is called in
the C++ layer, it sets the spherical attribute in the
Python side.
mass: when the mass(double m) method is called in the C++
side, it sets the spherical attribute in the Python
side.
__setitem__(self, key, value): optional, mandatory to support the
"Parameters" Property to set arbitrary parameters for this
metric.
set(self, key, val) and get(self, key): optional, both need to be
implemented in order to support additional "Properties" with
arbitrary names. In addition, member self.properties must be set
to a dictionary listing key:datatype pairs, where both key and
datatype are strings. As of writing, only "double" is supported
as a datatype.
'''
import math
import numpy
import gyoto.core
class Minkowski:
'''Flat space metric
Implemented for both Cartesian and spherical coordinates.
'''
def __setattr__(self, key, value):
'''Set attributes.
Optional.
C++ will set several attributes. By overloading __setattr__,
one can react when that occurs, in particular to make sure `this'
knows the coordinate kind as in this example.
Attributes set by the C++ layer:
this: if the Python extension "gyoto.core" can be imported, it
will be set to a gyoto.core.Metric instance pointing to the
C++-side instance. If the "gyoto.core" extension cannot be
loaded, this will be set to None.
spherical: when the spherical(bool t) method is called in
the C++ layer, it sets the spherical attribute in the
Python side.
mass: when the mass(double m) method is called in the C++
side, it sets the spherical attribute in the Python
side.
This example initializes coordKind in the C++ side if it is
not already set, since this Minkowski class can work in
either.
'''
# First, actually store the attribute. This is what would
# happen if we did not overload __setattr__.
self.__dict__[key]=value
# Then, if key is "this", ensure this knows a valid coordKind.
if (key == "this"):
cK=value.coordKind()
if cK is gyoto.core.GYOTO_COORDKIND_UNSPECIFIED:
value.set("Spherical", False)
# We could do without this, since this will tell us later
# anyway.
else:
self.spherical = (cK is gyoto.core.GYOTO_COORDKIND_SPHERICAL)
def gmunu(self, g, x):
''' Gyoto::Metric::Generic::gmunu(double dst[4][4], const double pos[4])
Mandatory.
C++ will send two NumPy arrays.
'''
for mu in range(0, 4):
for nu in range(0, 4):
g[mu][nu]=g[nu][mu]=0
g[0][0]=-1;
if not self.spherical:
for mu in range(1, 4):
g[mu][mu]=1.
return
r=x[1]
theta=x[2]
tmp=r*math.sin(theta)
g[1][1]=1.
g[2][2]=r*r
g[3][3]=tmp*tmp
def christoffel(self, dst, x):
'''Gyoto::Metric::Generic::christoffel(double dst[4][4][4], const double pos[4])
Mandatory.
C++ will send two NumPy arrays.
'''
for alpha in range(0, 4):
for mu in range(0, 4):
for nu in range(0, 4):
dst[alpha][mu][nu]=0.
if not self.spherical:
return 0
r=x[1]
theta=x[2]
sth=math.sin(theta)
cth=math.cos(theta)
dst[1][2][2]=-r
dst[1][3][3]=-r*sth*sth
dst[2][1][2]=dst[2][2][1]= 1./r
dst[2][3][3]=-sth*cth
dst[3][1][3]=dst[3][3][1]= dst[2][1][2]
dst[3][2][3]=dst[3][3][2]= math.tan(math.pi*0.5 - x[2])
return 0
class KerrBL:
'''A Python implementation of Gyoto::Metric::KerrBL
Spin and HorizonSecurity may be set either using the Parameters
Property, as in:
<Parameters> O.5 0.01 </Parameters>
or using two distinct ad hoc "properties":
<Spin> 0.5 </Spin>
<HorizonSecurity> 0.01 </HorizonSecurity>
The various methods are essentially cut-and-paste from the C++ class.
Only gmunu and christoffel absolutely need to be implemented, but
various Astrobj require more. See below.
This is for educational and testing purposes. For all practical uses,
the C++ implementation is much faster.
'''
# This is needed to support setting Spin and HorizonSecurity by name
# in addition to the set and get methods
properties={"Spin":"double", "HorizonSecurity":"double"}
spin=0.
drhor=gyoto.core.GYOTO_KERR_HORIZON_SECURITY
rsink=2.+gyoto.core.GYOTO_KERR_HORIZON_SECURITY
def __setitem__(self, key, value):
'''Set parameters
Optional, one way to handle parameters
This is how Gyoto sends the <Parameters/> XML entity:
metric[key]=value
key=0: set spin (0.)
key=1: set drhor, thickness of sink layer around horizon (0.01)
'''
if (key==0):
self.spin = value
elif (key==1):
self.drhor = value
else:
raise IndexError
def set(self, key, value):
'''Set parameters by name
Optional, one way to handle parameters
This is how Gyoto sends custom XML entities, e.g.
<MyEntity>value</MyEntity>
metric.set("MyEntity", value)
Here:
key="Spin": set spin (0.)
key="HorizonSecurity": set drhor, thickness of sink layer
around horizon (0.01)
'''
if (key=="Spin"):
self.spin = value
elif (key=="HorizonSecurity"):
self.drhor = value
else:
raise IndexError
def get(self, key):
'''Get parameters by name
Optional, mandatory if "properties" and "set" are defined.
'''
if (key=="Spin"):
return self.spin
elif (key=="HorizonSecurity"):
return self.drhor
else:
raise IndexError
def __setattr__(self, key, value):
'''Set attributes.
Optional.
C++ will set several attributes. By overloading __setattr__,
one can react when that occurs, in particular to make sure
`this' knows the coordinate kind as in this example.
In addition, update self.rsink each time self.spin or
self.drhor change.
'''
# First, actually store the attribute. This is what would
# happen if we did not overload __setattr__.
self.__dict__[key]=value
# Then, if key is "this", ensure `this' knows a valid coordKind.
if (key == "this"):
self.this.set("Spherical", True)
elif key in ("spin", "drhor"):
if self.spin > 1:
self.rsink=drhor
else:
self.rsink=1.+math.sqrt(1.-self.spin**2)+self.drhor
def gmunu(self, g, x):
''' Gyoto::Metric::Generic::gmunu(double g[4], const double pos[4])
Mandatory.
C++ will send two NumPy arrays.
Note that the user will not be able to call this method directly but
through self.this.gmunu which has a different calling sequence:
g=self.this.gmunu(x)
'''
spin_=self.spin
a2_=spin_**2
r=x[1]
sth2=math.sin(x[2])**2
cth2=math.cos(x[2])**2
r2=r*r
sigma=r2+a2_*cth2
delta=r2-2.*r+a2_
for mu in range(0, 4):
for nu in range(0, 4):
g[mu][nu]=g[nu][mu]=0
g[0][0] = -1.+2.*r/sigma;
g[1][1] = sigma/delta;
g[2][2] = sigma;
g[3][3] = (r2+a2_+2.*r*a2_*sth2/sigma)*sth2;
g[0][3] = g[3][0] = -2*spin_*r*sth2/sigma;
def christoffel(self, dst, x):
'''Gyoto::Metric::Generic::christoffel(double dst[4][4][4], const double pos[4])
Mandatory.
C++ will send two NumPy arrays.
Like gmunu, the call will actually be:
dst=metric.gmunu(x)
where `metric' is self.this.
'''
for alpha in range(0, 4):
for mu in range(0, 4):
for nu in range(0, 4):
dst[alpha][mu][nu]=0.
spin_=self.spin
a2_=spin_**2
r=x[1]
sth=math.sin(x[2])
cth=math.cos(x[2])
sth2 = sth*sth
cth2 = cth*cth
sth4=sth2*sth2
s2th = 2.*sth*cth
c2th=cth2-sth2
s4th = 2.*s2th*c2th
s2th2= s2th*s2th
ctgth=cth/sth
r2=r*r
r4=r2*r2
r6=r4*r2;
Sigma=r2+a2_*cth2
Sigma2=Sigma*Sigma
Delta=r2-2.*r+a2_
Deltam1=1./Delta
Sigmam1=1./Sigma
Sigmam2=Sigmam1*Sigmam1
Sigmam3=Sigmam2*Sigmam1
a2cthsth=a2_*cth*sth
rSigmam1=r*Sigmam1
Deltam1Sigmam2=Deltam1*Sigmam2
r2plusa2 = r2+a2_
dst[1][1][1]=(1.-r)*Deltam1+rSigmam1;
dst[1][2][1]=dst[1][1][2]=-a2cthsth*Sigmam1;
dst[1][2][2]=-Delta*rSigmam1;
dst[1][3][3]=-Delta*sth2*(r+(a2_*(-2.*r2+Sigma)*sth2)/Sigma2)/Sigma;
dst[1][3][0]=dst[1][0][3]=spin_*Delta*(-2*r2+Sigma)*sth2*Sigmam3;
dst[1][0][0]=-Delta*(-2.*r2+Sigma)*Sigmam3;
dst[2][1][1]=a2cthsth*Deltam1*Sigmam1;
dst[2][2][1]=dst[2][1][2]=rSigmam1;
dst[2][2][2]=-a2cthsth*Sigmam1;
dst[2][3][3]=-sth*cth*Sigmam3 * (Delta*Sigma2 + 2.*r*r2plusa2*r2plusa2);
dst[2][0][3]=dst[2][3][0]=spin_*r*r2plusa2*s2th*Sigmam3;
dst[2][0][0]=-2.*a2cthsth*r*Sigmam3;
dst[3][3][1]=dst[3][1][3]=Deltam1*Sigmam2 * (r*Sigma*(Sigma-2.*r) + a2_*(Sigma-2.*r2)*sth2);
dst[3][3][2]=dst[3][2][3]=Sigmam2*ctgth * (-(Sigma+Delta)*a2_*sth2 + r2plusa2*r2plusa2);
dst[3][0][1]=dst[3][1][0]=spin_*(2.*r2-Sigma)*Deltam1Sigmam2;
dst[3][0][2]=dst[3][2][0]=-2.*spin_*r*ctgth*Sigmam2;
dst[0][3][1]=dst[0][1][3]=-spin_*sth2*Deltam1Sigmam2 * (2.*r2*r2plusa2 + Sigma*(r2-a2_));
dst[0][3][2]=dst[0][2][3]=Sigmam2*spin_*a2_*r*sth2*s2th;
dst[0][0][1]=dst[0][1][0]=(a2_+r2)*(2.*r2-Sigma)*Deltam1Sigmam2;
dst[0][0][2]=dst[0][2][0]=-a2_*r*s2th*Sigmam2;
return 0
def getRms(self):
aa=self.spin;
a2_=aa*aa
z1 = 1. + pow((1. - a2_),1./3.)*(pow((1. + aa),1./3.) + pow((1. - aa),1./3.))
z2 = pow(3.*a2_ + z1*z1,1./2.);
return (3. + z2 - pow((3. - z1)*(3. + z1 + 2.*z2),1./2.));
def getRmb(self):
spin_=self.spin
return 2.-spin_+2.*math.sqrt(1.-spin_);
def getSpecificAngularMomentum(self, rr):
aa=self.spin
sqrtr=math.sqrt(rr);
return (rr*rr-2.*aa*sqrtr+aa*aa)/(rr**1.5-2.*sqrtr+aa);
def getPotential(self, pos, l_cst):
# this is W = -ln(|u_t|) for a circular equatorial 4-velocity
# Don't try to call directly gmunu from `self',
# Call it through `this', a pointer to the C++ instance of
# Gyoto::Metric::Python
# Note that the calling sequence is not gmunu(self, gg, pos)
gg=self.this.gmunu(pos)
gtt = gg[0,0];
gtp = gg[0,3];
gpp = gg[3,3];
Omega = -(gtp + l_cst * gtt)/(gpp + l_cst * gtp) ;
W = (0.5 * math.log(abs(gtt + 2. * Omega * gtp + Omega*Omega * gpp))
- math.log(abs(gtt + Omega * gtp))) ;
return W ;
def isStopCondition(self, coord):
return coord[1] < self.rsink
def circularVelocity(self, coor, vel, d):
'''Velocity of a particle in circular motion at point
Note: this is called only if this->keplerian_ is False.
'''
sinth = math.sin(coor[2]);
coord = [coor[0], coor[1]*sinth, math.pi*0.5, coor[3]];
vel[1] = vel[2] = 0.;
#vel[3] = 1./((d*pow(coord[1], 1.5) + spin_)*sinth);
vel[3] = 1./((d*pow(coord[1], 1.5) + self.spin));
vel[0] = self.this.SysPrimeToTdot(coor, vel[1:]);
vel[3] *= vel[0];

Member Typedef Documentation

◆ Subcontractor_t

typedef Gyoto::SmartPointer< Gyoto::SmartPointee > Gyoto::SmartPointee::Subcontractor_t(Gyoto::FactoryMessenger *, std::vector< std::string > const &)
inherited

A subcontractor builds an object upon order from the Factory.

Various classes need to provide a subcontractor to be able to instantiate themselves upon order from the Factory. A subcontractor is a function (often a static member function) which accepts a pointer to a FactoryMessenger as unique parameter, communicates with the Factory using this messenger to read an XML description of the object to build, and returns this objet. SmartPointee::Subcontractor_t* is just generic enough a typedef to cast to and from other subcontractor types: Astrobj::Subcontractor_t, Metric::Subcontractor_t, Spectrum::Subcontractor_t. A subcontractor needs to be registered using the relevant Register() function: Astrobj::Register(), Metric::Register(), Spectrum::Register().

Member Function Documentation

◆ christoffel() [1/2]

virtual double Gyoto::Metric::Generic::christoffel ( const double  coord[4],
const int  alpha,
const int  mu,
const int  nu 
) const
virtualinherited

◆ christoffel() [2/2]

virtual int Gyoto::Metric::Generic::christoffel ( double  dst[4][4][4],
const double  coord[4] 
) const
virtualinherited

◆ circularVelocity()

void Gyoto::Metric::Python::circularVelocity ( double const  pos[4],
double  vel[4],
double  dir = 1. 
) const
virtual

Yield circular velocity at a given position.

Give the velocity of a massive particle in circular orbit at the given position projected onto the equatorial plane. Such a velocity may not exist everywhere (or anywhere) for a given metric. This method is intended to be used by Astrobj classes such as Torus or ThinDisk.

If keplerian_ is set to true, this method should return the Keplerian velcity instead (derived classes should ensure this, see KerrBL::circularVelocity() for instance).

The default implementation throws an error if keplerian_ is set to false.

Parameters
posinput: position,
veloutput: velocity,
dir1 for corotating, -1 for counterrotating.

Reimplemented from Gyoto::Metric::Generic.

◆ clone()

virtual Python * Gyoto::Metric::Python::clone ( ) const
virtual

Virtual copy constructor.

Reimplemented from Gyoto::Metric::Generic.

◆ computeNBeta()

virtual void Gyoto::Metric::Generic::computeNBeta ( const double  coord[4],
double &  NN,
double  beta[3] 
) const
virtualinherited

Computes lapse scalar and shift vector at coord.

Reimplemented in Gyoto::Metric::KerrBL, and Gyoto::Metric::NumericalMetricLorene.

◆ coordKind()

void Gyoto::Metric::Generic::coordKind ( int  coordkind)
protectedinherited

Set coordkind_.

coordkind(int coordkind) is protected because, for most Metrics, it should not be changed in runtime. Set coordinate kind

◆ deltaMax() [1/3]

double Gyoto::Metric::Generic::deltaMax ( ) const
inherited

Get delta_max_

◆ deltaMax() [2/3]

virtual double Gyoto::Metric::Generic::deltaMax ( double const  pos[8],
double  delta_max_external 
) const
virtualinherited

Get delta max at a given position

Parameters
pos4-position
[optional]delta_max_external external constraint on delta_max
Returns
the smallest value between delta_max_, delta_max_external, and R*delta_max_over_r_ where R is pos[1] in spherical coordinates and max(x1, x2, x3) in Cartesian coordinates.

◆ deltaMax() [3/3]

void Gyoto::Metric::Generic::deltaMax ( double  h1)
inherited

Set delta_max_

◆ deltaMin() [1/2]

double Gyoto::Metric::Generic::deltaMin ( ) const
inherited

Get delta_min_

◆ deltaMin() [2/2]

void Gyoto::Metric::Generic::deltaMin ( double  h1)
inherited

Set delta_min_

◆ describeProperty()

std::string Gyoto::Object::describeProperty ( Gyoto::Property const &  p) const
inherited

Format desrciption for a property.

Returns a string containing the name(s) and type of the property, as well as whether it supports unit.

◆ diff()

virtual int Gyoto::Metric::Generic::diff ( state_t const &  x,
state_t &  dxdt,
double  mass 
) const
virtualinherited

F function such as dx/dt=F(x,cst)

Reimplemented in Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::RotStar3_1, and Gyoto::Metric::Minkowski.

◆ diff31()

virtual int Gyoto::Metric::Generic::diff31 ( state_t const &  x,
state_t &  dxdt,
double  mass 
) const
virtualinherited

F function such as dx/dt=F(x,cst) for 3+1 case.

Reimplemented in Gyoto::Metric::KerrBL, and Gyoto::Metric::NumericalMetricLorene.

◆ dualOneForm()

void Gyoto::Metric::Generic::dualOneForm ( double const  IN_ARRAY1_1[4],
double const  IN_ARRAY1_2[4],
double  ARGOUT_ARRAY1[4] 
) const
inherited

Computes dual 1-form Compute the dual 1-form of 4-vector.

Parameters
IN_ARRAY1_14-position;
IN_ARRAY1_2quadrivector;
[out]ARGOUT_ARRAY1output 1-form

◆ fillElement()

virtual void Gyoto::Python::Object< Gyoto::Metric::Generic >::fillElement ( Gyoto::FactoryMessenger fmp) const
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ fillProperty()

virtual void Gyoto::Object::fillProperty ( Gyoto::FactoryMessenger fmp,
Property const &  p 
) const
virtualinherited

Output a single Property to XML.

The base implementation decides what to do based on the p.type. The format matches how setParameters() an setParameter() would interpret the XML descition.

Overriding this method should be avoided, but makes sense in some cases (for instance Screen::fillProperty() selects a different unit for Distance based on its magnitude, so that stellar sizes are expressed in solar radii while smaller sizes can be expressed in meters and larger sizes in parsecs).

Overriding implementation should fall-back on calling the implementation in the direct parent class:

class A: public Object {};
class B: public A {
using B::setParameter;
virtual void fillProperty(Gyoto::FactoryMessenger *fmp,
Property const &p) const ;
};
void B::fillProperty(Gyoto::FactoryMessenger *fmp,
Property const &p) const {
if (name=="Duff") fmp->doSomething();
else A::fillProperty(fmp, p);
}
Factory / SmartPointee::Subcontractor_t interface.
Definition GyotoFactoryMessenger.h:92
Object with properties.
Definition GyotoObject.h:152
Property that can be set and got using standard methods.
Definition GyotoProperty.h:608

Reimplemented in Gyoto::Scenery, Gyoto::Astrobj::DirectionalDisk, Gyoto::Astrobj::Disk3D, Gyoto::Astrobj::DynamicalDisk, Gyoto::Astrobj::EquatorialHotSpot, Gyoto::Astrobj::NeutronStarModelAtmosphere, Gyoto::Astrobj::PatternDisk, Gyoto::Astrobj::PolishDoughnut, Gyoto::Screen, Gyoto::Metric::Shift, Gyoto::Astrobj::Star, Gyoto::Spectrometer::Uniform, and Gyoto::Astrobj::XillverReflection.

◆ get() [1/3]

Value Gyoto::Python::Object< Gyoto::Metric::Generic >::get ( Property const &  p) const
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ get() [2/3]

Value Gyoto::Python::Object< Gyoto::Metric::Generic >::get ( Property const &  p,
std::string const &  unit 
) const
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ get() [3/3]

virtual Value Gyoto::Python::Object< Gyoto::Metric::Generic >::get ( std::string const &  pname) const
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ getPotential()

double Gyoto::Metric::Python::getPotential ( double const  pos[4],
double  l_cst 
) const
virtual

Returns potential W=-ln(|u_t|) for a cst specific angular momentum l_cst Should be implemented in derived classes if useful If called on the base class, returns an error

Reimplemented from Gyoto::Metric::Generic.

◆ getProperties()

virtual Property const * Gyoto::Metric::Python::getProperties ( ) const
virtual

Get list of properties.

This method is declared automatically by the GYOTO_OBJECT macro and defined automatically by the GYOTO_PROPERTY_END macro.

Reimplemented from Gyoto::Metric::Generic.

◆ getRmb()

double Gyoto::Metric::Python::getRmb ( ) const
virtual

Returns the marginally bound radius Should be implemented in derived classes if useful If called on the base class, returns an error

Reimplemented from Gyoto::Metric::Generic.

◆ getRms()

double Gyoto::Metric::Python::getRms ( ) const
virtual

Returns the marginally stable (ISCO) radius Should be implemented in derived classes if useful If called on the base class, returns an error

Reimplemented from Gyoto::Metric::Generic.

◆ getSpecificAngularMomentum()

double Gyoto::Metric::Python::getSpecificAngularMomentum ( double  rr) const
virtual

Returns the specific angular momentum l=-u_phi/u_t Should be implemented in derived classes if useful If called on the base class, returns an error

Reimplemented from Gyoto::Metric::Generic.

◆ gmunu() [1/2]

virtual void Gyoto::Metric::Generic::gmunu ( double  ARGOUT_ARRAY2[4][4],
double const  IN_ARRAY1[4] 
) const
virtualinherited

Metric coefficients.

The default implementation calls double gmunu(const double * x, int mu, int nu) const.

Parameters
[out]ARGOUT_ARRAY2(g) 4x4 array to store the coeefficients
[in]IN_ARRAY1(x) 4-position at which to compute the coefficients;
Returns
Metric coefficient gμ,ν at point x

Reimplemented in Gyoto::Metric::Complex, Gyoto::Metric::Hayward, Gyoto::Metric::KerrBL, Gyoto::Metric::KerrKS, Gyoto::Metric::Minkowski, Gyoto::Metric::Shift, Gyoto::Metric::ChernSimons, Gyoto::Metric::Complex, Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::RezzollaZhidenko, Gyoto::Metric::RotStar3_1, Gyoto::Metric::SchwarzschildHarmonic, and Gyoto::Metric::Shift.

◆ gmunu() [2/2]

virtual double Gyoto::Metric::Generic::gmunu ( double const  x[4],
int  mu,
int  nu 
) const
virtualinherited

Metric coefficients.

The default implementation calls Metric:: gmunu(double g[4][4], const double * pos) const

Parameters
x4-position at which to compute the coefficient;
mu1st index of coefficient, 0≤μ≤3;
nu2nd index of coefficient, 0≤ν≤3;
Returns
Metric coefficient gμ,ν at point x

Reimplemented in Gyoto::Metric::Minkowski, Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::KerrKS, Gyoto::Metric::ChernSimons, Gyoto::Metric::Complex, Gyoto::Metric::Hayward, Gyoto::Metric::KerrBL, Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::RezzollaZhidenko, Gyoto::Metric::RezzollaZhidenko, Gyoto::Metric::RotStar3_1, Gyoto::Metric::RotStar3_1, Gyoto::Metric::SchwarzschildHarmonic, Gyoto::Metric::SchwarzschildHarmonic, and Gyoto::Metric::Shift.

◆ gmunu_up() [1/2]

virtual void Gyoto::Metric::Generic::gmunu_up ( double  ARGOUT_ARRAY2[4][4],
const double  IN_ARRAY1[4] 
) const
virtualinherited

Metric contravariant coefficients.

The default implementation inverts the covariant coefficients matrix.

Reimplemented in Gyoto::Metric::Hayward, Gyoto::Metric::KerrBL, Gyoto::Metric::KerrKS, Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::Shift, Gyoto::Metric::Shift, and Gyoto::Metric::ChernSimons.

◆ gmunu_up() [2/2]

virtual double Gyoto::Metric::Generic::gmunu_up ( double const  x[4],
int  mu,
int  nu 
) const
virtualinherited

Metric contravariant coefficients.

The default implementation calls Metric:: gmunu_up(double g[4][4], const double * pos) const

Parameters
x4-position at which to compute the coefficient;
mu1st index of coefficient, 0≤μ≤3;
nu2nd index of coefficient, 0≤ν≤3;
Returns
Metric coefficient gμ,ν at point x

Reimplemented in Gyoto::Metric::ChernSimons, Gyoto::Metric::Hayward, Gyoto::Metric::KerrBL, and Gyoto::Metric::Shift.

◆ gmunu_up_and_jacobian()

virtual void Gyoto::Metric::Generic::gmunu_up_and_jacobian ( double  ARGOUT_ARRAY2[4][4],
double  ARGOUT_ARRAY3[4][4][4],
const double  IN_ARRAY1[4] 
) const
virtualinherited

gmunu_up() and jacobian() in one go

Reimplemented in Gyoto::Metric::KerrKS.

◆ GramSchmidt()

void Gyoto::Metric::Generic::GramSchmidt ( double const  pos[4],
double  u0[4],
double  u1[4],
double  u2[4],
double  u3[4] 
) const
inherited

Apply Gram-Schmidt orthonormalization to a basis.

On input, u0 to u3 must be four non-zero norm, independent 4-vectors. On output, they will form an orthonormal basis.

Parameters
[in]posposition,
[in,out]u0basis vector
[in,out]u1basis vector
[in,out]u2basis vector
[in,out]u3basis vector

◆ help()

void Gyoto::Object::help ( ) const
inherited

Print (to stdout) some help on this class.

Describe all properties that this instance supports.

◆ hook()

virtual void Gyoto::Hook::Teller::hook ( Listener listener)
virtualinherited

Start listening.

Use from a Hook::Listener object method:

teller->hook(this)

where "this" is a Listener and "teller" is a Teller.

Use unhook() later to stop listening to a given Teller.

Parameters
listenerpointer to the new listener

◆ inlineModule() [1/2]

virtual std::string Gyoto::Metric::Python::inlineModule ( ) const
virtual

Return inline_module_.

Reimplemented from Gyoto::Python::Base.

◆ inlineModule() [2/2]

virtual void Gyoto::Metric::Python::inlineModule ( const std::string &  )
virtual

Set inline_module_ and import the Python module.

Side effects:

Reimplemented from Gyoto::Python::Base.

◆ isStopCondition()

int Gyoto::Metric::Python::isStopCondition ( double const  coord[8]) const
virtual

Check whether integration should stop.

The integrating loop will ask this the Metric through this method whether or not it is happy to continue the integration. Typically, the Metric should answer 0 when everything is fine, 1 when too close to the event horizon, inside the BH...

Parameters
coord8-coordinate vector to check.

Reimplemented from Gyoto::Metric::Generic.

◆ isThreadSafe()

virtual bool Gyoto::Object::isThreadSafe ( ) const
virtualinherited

Whether this class is thread-safe.

Return True if this object is thread-safe, i.e. if an instance and its clone can be used in parallel threads (in the context of Scenery::raytrace()). Known objects which are not thread-safe include Lorene metrics and everything from the Python plug-in.

The default implementation considers that the class itself is thread safe and recurses into the declared properties to check whether they are safe too. Classes that abide to the Object/Property paradigm and are themselves thread-safe have nothing special to do.

Objects that clone children in their copy constructor that are not declared as properties must take these children into account.

Classes that are never thread-safe must declare it. It acn be easily done using GYOTO_OBJECT_THREAD_SAFETY in the class declaration and GYOTO_PROPERTY_THREAD_UNSAFE in the class definition.

◆ jacobian()

virtual void Gyoto::Metric::Generic::jacobian ( double  ARGOUT_ARRAY3[4][4][4],
const double  IN_ARRAY1[4] 
) const
virtualinherited

Derivatives of the metric covariant coefficients.

The default implementation evaluates them numerically. The gmunu matrix is assumed to be symmetrical but no other assumptions are made at the moment.

Reimplemented in Gyoto::Metric::Complex, Gyoto::Metric::NumericalMetricLorene, Gyoto::Metric::Shift, and Gyoto::Metric::KerrKS.

◆ kind() [1/2]

virtual std::string Gyoto::Object::kind ( ) const
virtualinherited

Get kind_.

Reimplemented in Gyoto::Spectrometer::Uniform.

◆ kind() [2/2]

virtual void Gyoto::Object::kind ( const std::string  )
protectedvirtualinherited

Set kind_.

kind(const std::string) is protected because, for most Objects, it should not be changed in runtime. Set kind_

◆ klass() [1/2]

virtual std::string Gyoto::Metric::Python::klass ( ) const
virtual

Retrieve class_.

Reimplemented from Gyoto::Python::Base.

◆ klass() [2/2]

virtual void Gyoto::Metric::Python::klass ( const std::string &  c)
virtual

Set class_ and instantiate the Python class.

Sets pInstance_.

This generic implementation takes care of the common ground, but does not set 'this' or call parameters(parameters_). Therefore, all the derived classes should reimplement this method and at least call Python::Base::klass(c) and parameters(parameters_). Between the two is the right moment to check that the Python class implements the required API and to cache PyObject* pointers to class methods.

Reimplemented from Gyoto::Python::Base.

◆ mass()

virtual void Gyoto::Metric::Python::mass ( double  )
virtual

Set mass used in unitLength()

Reimplemented from Gyoto::Metric::Generic.

◆ module() [1/2]

virtual std::string Gyoto::Metric::Python::module ( ) const
virtual

Return module_.

Reimplemented from Gyoto::Python::Base.

◆ module() [2/2]

virtual void Gyoto::Metric::Python::module ( const std::string &  )
virtual

Set module_ and import the Python module.

Side effects:

Reimplemented from Gyoto::Python::Base.

◆ myrk4()

virtual int Gyoto::Metric::Generic::myrk4 ( Worldline line,
state_t const &  coord,
double  h,
state_t &  res 
) const
virtualinherited

RK4 integrator.

Reimplemented in Gyoto::Metric::KerrBL.

◆ myrk4_adaptive()

virtual int Gyoto::Metric::Generic::myrk4_adaptive ( Gyoto::Worldline line,
state_t const &  coord,
double  lastnorm,
double  normref,
state_t &  coordnew,
double  h0,
double &  h1,
double  deltamax = GYOTO_DEFAULT_DELTA_MAX 
) const
virtualinherited

RK4 integrator with adaptive step.

Reimplemented in Gyoto::Metric::KerrBL, and Gyoto::Metric::RotStar3_1.

◆ norm()

double Gyoto::Metric::Generic::norm ( const double  pos[4],
const double  u1[4] 
) const
inherited

Scalar product.

Compute the norm of the quadrivector u1 in this Metric, at point pos expressed in coordinate system sys.

Parameters
[in]pos4-position;
[in]u1quadrivector;
Returns
||u1||

◆ normalizeFourVel() [1/2]

virtual void Gyoto::Metric::Generic::normalizeFourVel ( double const  pos[4],
double  fourvel[4] 
) const
virtualinherited

Normalize fourvelvel to -1.

First computes threevel as xiprime=xidot/x0dot for i in {1, 2, 3}, then computes x0dot using SyPrimeToTdot, then computes again xidot as xidot=xiprime*x0dot.

Parameters
[in]pos4-position;
[in,out]fourvel4-velocity, will be renormalized.

◆ normalizeFourVel() [2/2]

virtual void Gyoto::Metric::Generic::normalizeFourVel ( double  coord[8]) const
virtualinherited

Normalize fourvelvel to -1.

First computes threevel as xiprime=xidot/x0dot for i in {1, 2, 3}, then computes x0dot using SyPrimeToTdot, then computes again xidot as xidot=xiprime*x0dot.

Parameters
[in,out]coord8-position, coord[4-7] will be set according to the other elements;

◆ nullifyCoord() [1/2]

virtual void Gyoto::Metric::Generic::nullifyCoord ( double  coord[8]) const
virtualinherited

Set tdot (coord[4]) such that coord is light-like. Everything is in geometrical units.

Set coord[4] so that the 4-velocity coord[4:7] is lightlike, i.e. of norm 0. There may be up to two solutions. coord[4] is set to the hightest. The lowest can be retrieved using nullifyCoord(double coord[8], double& tdot2) const. Everything is expressed in geometrical units.

Parameters
[in,out]coord8-position, coord[4] will be set according to the other elements;

Reimplemented in Gyoto::Metric::KerrBL.

◆ nullifyCoord() [2/2]

virtual void Gyoto::Metric::Generic::nullifyCoord ( double  coord[8],
double &  tdot2 
) const
virtualinherited

Set tdot (coord[4]) such that coord is light-like and return other possible tdot.

Set coord[4] so that the 4-velocity coord[4:7] is lightlike, i.e. of norm 0. There may be up to two solutions. coord[4] is set to the hightest. The lowest can be retrieved in tdot2. Everything is expressed in geometrical units.

Parameters
[in,out]coord8-position, coord[4] will be set according to the other elements;
[out]tdot2will be set to the smallest solution

Reimplemented in Gyoto::Metric::KerrBL.

◆ observerTetrad() [1/2]

virtual void Gyoto::Metric::Generic::observerTetrad ( double const  pos[4],
double  fourvel[4],
double  screen1[4],
double  screen2[4],
double  screen3[4] 
) const
virtualinherited

Computes the orthonormal local tetrad of the observer.

Parameters
[in]posposition,
[in]fourvelobserver 4-velocity (norm -1)
[out]screen1first vector in the screen plane
[out]screen2second vector in the screen plane
[out]screen3vector normal to the screen

Reimplemented in Gyoto::Metric::KerrBL.

◆ observerTetrad() [2/2]

virtual void Gyoto::Metric::Generic::observerTetrad ( obskind_t  obskind,
double const  pos[4],
double  fourvel[4],
double  screen1[4],
double  screen2[4],
double  screen3[4] 
) const
virtualinherited

Computes the orthonormal local tetrad of the observer.

Parameters
obskindinput: kind of observer (eg: "ZAMO","KeplerianObserver"...)
posinput: position,
fourveloutput: observer 4-velocity (norm -1)
screen1output: first vector in the screen plane
screen2output: second vector in the screen plane
screen3output: vector normal to the screen

Reimplemented in Gyoto::Metric::Minkowski.

◆ parameters() [1/2]

virtual std::vector< double > Gyoto::Metric::Python::parameters ( ) const
virtual

Retrieve parameters_.

Reimplemented from Gyoto::Python::Base.

◆ parameters() [2/2]

virtual void Gyoto::Metric::Python::parameters ( const std::vector< double > &  )
virtual

Set parameters_ and send them to pInstance_.

The parameters are sent to the class instance using the setitem method with numerical keys.

Reimplemented from Gyoto::Python::Base.

◆ property()

Property const * Gyoto::Object::property ( std::string const  pname) const
inherited

Find property by name.

Look into the Property list for a Property whose name (or name_false, for a boolean Property) is pname. Return a const pointer to the first such property found, or NULL if none is found.

◆ ScalarProd()

virtual double Gyoto::Metric::Generic::ScalarProd ( const double  pos[4],
const double  u1[4],
const double  u2[4] 
) const
virtualinherited

Scalar product.

Compute the scalarproduct of the two quadrivectors u1 and u2 in this Metric, at point pos expressed in coordinate system sys.

Parameters
pos4-position;
u11st quadrivector;
u22nd quadrivector;
Returns
u1*u2

Reimplemented in Gyoto::Metric::Hayward, Gyoto::Metric::KerrBL, and Gyoto::Metric::RotStar3_1.

◆ set() [1/3]

virtual void Gyoto::Python::Object< Gyoto::Metric::Generic >::set ( Property const &  p,
Value  val 
)
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ set() [2/3]

virtual void Gyoto::Python::Object< Gyoto::Metric::Generic >::set ( Property const &  p,
Value  val,
std::string const &  unit 
)
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ set() [3/3]

virtual void Gyoto::Python::Object< Gyoto::Metric::Generic >::set ( std::string const &  pname,
Value  val 
)
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ setParameter() [1/2]

virtual void Gyoto::Object::setParameter ( Gyoto::Property const &  p,
std::string const &  name,
std::string const &  content,
std::string const &  unit 
)
virtualinherited

Set parameter by Property (and name)

This function is used when parsing an XML description, if Property (p) of this name is found (i.e. either p.name or p.name_false is equal to name). Implementation should fall-back on calling the direct's parent implementation:

class A: public Object {};
class B: public A {
using B::setParameter;
virtual void setParameter(Gyoto::Property const &p,
std::string name,
std::string content,
std::string unit);
};
void B::setParameter(Gyoto::Property const &p,
std::string name,
std::string content,
std::string unit) {
if (name=="Duff") doSomething(content, unit);
else A::setParameter(p, name, content, unit);
}
Parameters
pProperty that matches name (p.name == name or p.name_false == name)
nameXML name of the parameter (XML entity)
contentstring representation of the value
unitstring representation of the unit

Reimplemented in Gyoto::Astrobj::PolishDoughnut.

◆ setParameter() [2/2]

virtual int Gyoto::Python::Object< Gyoto::Metric::Generic >::setParameter ( std::string  name,
std::string  content,
std::string  unit 
)
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ setParameters()

void Gyoto::Python::Object< Gyoto::Metric::Generic >::setParameters ( Gyoto::FactoryMessenger fmp)
inlinevirtualinherited

Reimplemented from Gyoto::Object.

◆ SysPrimeToTdot()

virtual double Gyoto::Metric::Generic::SysPrimeToTdot ( const double  coord[4],
const double  v[3] 
) const
virtualinherited

Compute tdot as a function of dr/dt, dtheta/dt and dphi/dt. Everything is in geometrical units.

Parameters
coord4-position (geometrical units);
v3-velocity dx1/dx0, dx2/dx0, dx3/dx0;
Returns
tdot = dx0/dtau.

◆ tellListeners()

virtual void Gyoto::Hook::Teller::tellListeners ( )
protectedvirtualinherited

Call tell() on each hooked Listener.

Whenever a Teller mutates, it should warn any Listener hooked to it using tellListeners().

◆ unhook()

virtual void Gyoto::Hook::Teller::unhook ( Listener listener)
virtualinherited

Stop listening.

Use from a Hook::Listener object method:

teller->unhook(this)

where "this" is a Listener, "teller" is a Teller, and "this" has called teller->hook(this) previously.

Parameters
listenerpointer to the listener

◆ unitLength()

double Gyoto::Metric::Generic::unitLength ( ) const
inherited

M * G / c^2, M is in kg, unitLength in meters.

Metrics implementations are free to express lengths and distances in whatever unit they see fit (presumably most often geometrical units). This function returns this unit in SI (meters).

◆ zamoVelocity()

virtual void Gyoto::Metric::Generic::zamoVelocity ( double const  pos[4],
double  vel[4] 
) const
virtualinherited

Yield ZAMO velocity at a given position.

Give the velocity of a zero angular momentul observer (whatever is closest to "at rest"). The default implementation simply projects (1, 0, 0, 0) othogonally along ephi and normalizes it, thus ensuring that vel is orthogonal to ephi.

Parameters
posinput: position,
veloutput: velocity,

Reimplemented in Gyoto::Metric::KerrBL.

Member Data Documentation

◆ __defaultfeatures

int Gyoto::Metric::Generic::__defaultfeatures
privateinherited

Whether some virtual methods are implemented.

The default implementations of some methods call one-another. This member is used internally to avoid infinite recursion.

◆ class_

std::string Gyoto::Python::Base::class_
protectedinherited

Name of the Python class that we want to expose.

Property name: Class.

◆ delta_max_over_r_

double Gyoto::Metric::Generic::delta_max_over_r_
protectedinherited

Numerical tuning parameter.

Ensure that delta (the numerical integration step) is never larger than a fraction of the distance between the current location and the center of the coordinate system.

For investigations close to the event horizon, 0.5 is usually fine. If high accuracy is needed long after deflection (weak lensing), then this must be smaller. A good test is to look at a MinDistance map for a FixedStar: it must be smooth.

◆ kind_

std::string Gyoto::Object::kind_
protectedinherited

The "kind" that is output in the XML entity.

E.g. for an Astrobj, fillElement() will ensure

<Astrobj kind="kind_" ...>...</Astrobj>
virtual std::string kind() const
Get kind_.

is written.

◆ module_

std::string Gyoto::Python::Base::module_
protectedinherited

Name of the Python module that holds the class.

For instance, if the class is implemented in toto.py, the module name is "toto". Property name: Module.

◆ mutex_

pthread_mutex_t Gyoto::SmartPointee::mutex_
privateinherited

A mutex.

When compiled with libpthread

◆ parameters_

std::vector<double> Gyoto::Python::Base::parameters_
protectedinherited

Parameters that this class needs.

A list of parameters (doubles) can be passed in the Property Parameters. They will be sent to the Python instance using setitem.

◆ plugins_

std::vector<std::string> Gyoto::Object::plugins_
protectedinherited

The plug-ins that needs to be loaded to access this instance's class.

E.g. for an Astrobj, fillElement() will ensure

<Astrobj ... plugin="plugins_">...</Astrobj>

is written.


The documentation for this class was generated from the following file: