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

The camera with which the Astrobj is observed. More...

#include <GyotoScreen.h>

Inheritance diagram for Gyoto::Screen:
Gyoto::SmartPointee Gyoto::Object

Classes

class  Angles
 1D specifier for an arbitrary angle coordinate set. More...
 
class  Bucket
 Class containing arbitrary 2D-points. More...
 
class  Coord1dSet
 Set of 1-d coordinates: indices or angles. More...
 
class  Coord2dSet
 Class to specify a set of points on the Screen. More...
 
class  Empty
 A dummy, empty 2D set. More...
 
class  Grid
 Class containing 2D-points organized in a grid. More...
 
class  Indices
 1D specifier for an arbitrary pixel coordinate set. More...
 
class  Range
 1D coordinated specifier for a range More...
 
class  RepeatAngle
 1D specifier for an angle that is repeated. More...
 

Public Types

enum  CoordType_e { angle, pixel }
 Enum to specify whether a coordinate set (Coord1dSet or Coord2dSet) holds pixel values or angles.
 
typedef Gyoto::SmartPointer< Gyoto::SmartPointeeSubcontractor_t(Gyoto::FactoryMessenger *, std::vector< std::string > const &)
 A subcontractor builds an object upon order from the Factory. More...
 

Public Member Functions

virtual Property const * getProperties () const
 Get list of properties. More...
 
 Screen ()
 Default constructor.
 
 Screen (const Screen &)
 Copy constructor.
 
Screenclone () const
 Cloner.
 
virtual ~Screen ()
 Destructor.
 
void setProjection (const double paln, const double inclination, const double argument)
 Set inclination etc.
 
void setProjection (const double distance, const double paln, const double inclination, const double argument)
 Set distance, inclination etc.
 
void distance (double dist)
 Set distance from observer. More...
 
void dMax (double dist)
 Set ray-tracing maximum distance. More...
 
void distance (double dist, const std::string &unit)
 Set distance from observer. More...
 
void inclination (double)
 Set inclination relative to line-of-sight. More...
 
void inclination (double, const std::string &unit)
 Set inclination relative to line-of-sight. More...
 
void PALN (double)
 Set position angle of the line of nodes.
 
void PALN (double, const std::string &unit)
 Set position angle of the line of nodes.
 
void argument (double)
 Set angle beetwen line of nodes and X axis of object.
 
void argument (double, const std::string &unit)
 Set angle beetwen line of nodes and X axis of object.
 
void spectrometer (SmartPointer< Spectrometer::Generic > spectro)
 Set Screen::spectro_.
 
SmartPointer< Spectrometer::Genericspectrometer () const
 Get Screen::spectro_.
 
void freqObs (double fo)
 Set freq_obs_. More...
 
void freqObs (double fo, const std::string &unit)
 Set freq_obs_. More...
 
double freqObs () const
 Get freq_obs_.
 
double freqObs (const std::string &unit) const
 Get freq_obs_. More...
 
void setObserverPos (const double pos[4])
 Alternative way to set projection. More...
 
void observerKind (const std::string &kind)
 
std::string observerKind () const
 
void setFourVel (const double coord[4])
 Sets the observer's 4-velocity.
 
void setScreen1 (const double coord[4])
 Sets the screen vector e1.
 
void setScreen2 (const double coord[4])
 Sets the screen vector e2.
 
void setScreen3 (const double coord[4])
 Sets the screen vector e3 (normal)
 
int coordKind () const
 Get coordinate kind. More...
 
double distance () const
 Get distance from observer. More...
 
double distance (const std::string &) const
 Get distance from observer. More...
 
double dMax () const
 Get maximum ray-tracing distance. More...
 
double inclination () const
 Get inclination relative to line-of-sight. More...
 
double inclination (const std::string &) const
 Get inclination relative to line-of-sight. More...
 
double PALN () const
 Get position angle of the line of nodes.
 
double PALN (const std::string &) const
 Get position angle of the line of nodes.
 
double argument () const
 Get angle between line of nodes and X axis of object.
 
double argument (const std::string &) const
 Get angle between line of nodes and X axis of object.
 
SmartPointer< Metric::Genericmetric () const
 Get Screen::gg_.
 
void metric (SmartPointer< Metric::Generic > gg)
 Set Screen::gg_.
 
double time () const
 Get observing date in seconds.
 
double time (const std::string &) const
 Get observing date in seconds.
 
void time (double, const std::string &)
 Set observing date in specified unit.
 
void time (double)
 Set observing date in seconds.
 
double fieldOfView () const
 Get Screen::fov_ in radians.
 
double fieldOfView (std::string const &unit) const
 Get Screen::fov_ in specified unit.
 
void fieldOfView (double)
 Set Screen::fov_ in radians.
 
void fieldOfView (double, const std::string &unit)
 Set Screen::fov_ in specified unit.
 
double azimuthalFieldOfView () const
 Get Screen::azimuthal_fov_.
 
void azimuthalFieldOfView (double ff)
 Set Screen::azimuthal_fov_.
 
void dangle1 (double)
 Set increment to first position angle.
 
void dangle1 (double, const std::string &unit)
 Set increment to first position angle in specified unit.
 
double dangle1 () const
 Get increment to first position angle.
 
double dangle1 (std::string const &unit) const
 Get increment to first position angle in specified unit.
 
void dangle2 (double)
 Set increment to second position angle.
 
void dangle2 (double, const std::string &unit)
 Set increment to second position angle in specified unit.
 
double dangle2 () const
 Get increment to second position angle.
 
double dangle2 (std::string const &unit) const
 Get increment to second position angle in specified unit.
 
void anglekind (int)
 Set Screen::anglekind_.
 
void anglekind (std::string const &)
 
std::string anglekind () const
 
size_t resolution () const
 Get Screen::npix_.
 
void resolution (size_t)
 Set Screen::npix_.
 
void mask (double const *const mm, size_t resolution=0)
 Set mask_ from array. More...
 
double const * mask () const
 Retrieve const pointer to mask_.
 
void maskFile (std::string const &fname)
 
std::string maskFile () const
 
void fitsReadMask (std::string const &fname)
 Read mask_ from FITS file.
 
void fitsWriteMask (std::string const &fname)
 Save mask_ from FITS file.
 
bool operator() (size_t, size_t)
 Whether this pixel should be ray-traced. More...
 
void getObserverPos (double dest[4]) const
 4-Position of the observer relative to the metric More...
 
void getFourVel (double dest[4]) const
 Get copy of Screen::fourvel_. More...
 
void fourVel (std::vector< double > const &)
 
std::vector< double > fourVel () const
 
void screenVector1 (std::vector< double > const &)
 
std::vector< double > screenVector1 () const
 
void screenVector2 (std::vector< double > const &)
 
std::vector< double > screenVector2 () const
 
void screenVector3 (std::vector< double > const &)
 
std::vector< double > screenVector3 () const
 
void getScreen1 (double dest[4]) const
 Get copy of Screen::screen1_. More...
 
void getScreen2 (double dest[4]) const
 Get copy of Screen::screen2_. More...
 
void getScreen3 (double dest[4]) const
 Get copy of Screen::screen3_. More...
 
void getRayCoord (double x, double y, double dest[8]) const
 Get 8-coordinate of Photon hitting screen from a given direction. More...
 
void getRayTriad (double coord[8], double Ephi[4], double Etheta[4]) const
 Get polarization triad. More...
 
void getRayCoord (const size_t i, const size_t j, double dest[8]) const
 Get 8-coordinate of Photon hitting screen pixel. More...
 
void coordToSky (const double pos[4], double dest[3], bool geometrical=false) const
 Convert metric 4-position to sky 3-position. More...
 
void skyToCoord (const double sky[3], double dest[4], bool geometrical=false) const
 Convert sky 3-position to metric 4-position. More...
 
void coordToXYZ (const double pos[4], double dest[3]) const
 Convert 4-position to 3-cartesian coordinates.
 
void computeBaseVectors ()
 Compute base vectors according to projection parameters.
 
std::ostream & print (std::ostream &) const
 Display. More...
 
std::ostream & printBaseVectors (std::ostream &) const
 Debug helper.
 
void mapPixUnit ()
 Map "pix" and "pixel" to angular pixel width in unit system. More...
 
void unmapPixUnit ()
 Unmap "pix" and "pixel" from unit system. More...
 
void fillProperty (Gyoto::FactoryMessenger *fmp, Property const &p) const
 Output a single Property to XML. More...
 
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.
 
int getRefCount ()
 Get the current number of references.
 
virtual bool isThreadSafe () const
 Whether this class is thread-safe. More...
 
void set (Property const &p, Value val)
 Set Value of a Property.
 
void set (Property const &p, Value val, std::string const &unit)
 Set Value (expressed in unit) of a Property.
 
void set (std::string const &pname, Value val)
 Set Value of a Property.
 
void set (std::string const &pname, Value val, std::string const &unit)
 Set Value (expressed in unit) of a Property.
 
Value get (Property const &p) const
 Get Value of a Property.
 
Value get (std::string const &pname) const
 Get Value of a Property.
 
Value get (Property const &p, std::string const &unit) const
 Get Value of a Property, converted to unit.
 
Value get (std::string const &pname, std::string const &unit) const
 Get Value of a Property, converted to unit.
 
Property const * property (std::string const pname) const
 Find property by name. More...
 
virtual void fillElement (Gyoto::FactoryMessenger *fmp) const
 Fill the XML element for this Object. More...
 
virtual void setParameters (Gyoto::FactoryMessenger *fmp)
 Main loop for parsing Properties from XML description. More...
 
virtual int setParameter (std::string name, std::string content, std::string unit)
 Set parameter by name. More...
 
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) More...
 
std::string describeProperty (Gyoto::Property const &p) const
 Format desrciption for a property. More...
 
void help () const
 Print (to stdout) some help on this class. More...
 

Static Public Member Functions

static SmartPointer< ScreenSubcontractor (FactoryMessenger *fmp)
 Instanciate a Screen from XML entity.
 

Public Attributes

 GYOTO_OBJECT_THREAD_SAFETY
 

Static Public Attributes

static GYOTO_OBJECT Property const properties []
 

Protected Attributes

std::string kind_
 The "kind" that is output in the XML entity. More...
 
std::vector< std::string > plugins_
 The plug-ins that needs to be loaded to access this instance's class. More...
 

Private Types

enum  anglekind_e { equatorial_angles =0, rectilinear =1, spherical_angles =2 }
 
typedef int anglekind_t
 

Private Attributes

double tobs_
 Observing date in s.
 
double fov_
 Field-of-view in rad.
 
double azimuthal_fov_
 Azimuthal field-of-view for Spherical Angles images. Maximal extent of image in the azimuthal b-angle direction.
 
size_t npix_
 Resolution in pixels.
 
double * mask_
 Mask with 0 where the ray-tracing should not be performed.
 
std::string mask_filename_
 Last read or written FITS file. More...
 
double distance_
 Distance to the observer in m.
 
double dmax_
 Maximum distance from which the photons are launched (geometrical units)
 
anglekind_t anglekind_
 Screen angles kind (0: equatorial, 1: spherical)
 
double euler_ [3]
 Euler angles. More...
 
double ex_ [3]
 Sky coordinate of base X vector.
 
double ey_ [3]
 Sky coordinate of base Y vector.
 
double ez_ [3]
 Sky coordinate of base Z vector.
 
double fourvel_ [4]
 Observer's 4-velocity.
 
double screen1_ [4]
 Screen e1 vector.
 
double screen2_ [4]
 Screen e2 vector.
 
double screen3_ [4]
 Screen e3 vector (normal)
 
double dangle1_
 Increment to first position angle of Screen; can be typically alpha if in Equatorial Angles, or a if in Spherical Angles.
 
double dangle2_
 Increment to second position angle of Screen; can be typically delta if in Equatorial Angles, or b if in Spherical Angles.
 
SmartPointer< Metric::Genericgg_
 The Metric in this end of the Universe.
 
SmartPointer< Spectrometer::Genericspectro_
 Gyoto::Spectrometer::Generic subclass instance used for quantities Spectrum and BinSpectrum.
 
double freq_obs_
 Frequency at which the observer observes. More...
 
obskind_t observerkind_
 What kind of observer are we considering? (At infinity, ZAMO...)
 

Friends

class Gyoto::SmartPointer< Gyoto::Screen >
 

Detailed Description

The camera with which the Astrobj is observed.

In the observer-centric point-of-view, the center of the Metric's coordinate system is positioned relatively to the observing Screen using three Euler angles and the distance (in meters). The three Euler angles are:

In addition, the Screen conveys:

The scalar FreqObs defines the observing frequency for Scenery quantity Intensity.

Likewise, a Gyoto::Spectrometer defines for which frequencies spectra are computed (when the Quantity Spectrum is requested in the Scenery).

For the sake of theoreticians, there is an alternate way of specifying the relative position of the Screen and Metric, by specifying the 4-coordinates of the Screen in the Metric's coordinate system (in that case, eerything is specified in geometrical units).

So an XML stanza for a Screen may look like that:

<Time> 1000. </Time>
<FieldOfView> 0.3141592653589793 </FieldOfView>
<Resolution> 128 </Resolution>
<Distance> 1e30 </Distance>
<PALN> 3.14159 </PALN>
<Inclination> 2.0944 </Inclination>
<Argument> -2.0944 </Argument>
<Spectrometer kind="freqlog" nsamples="10"> 17. 23. </Spectrometer>
<FreqObs> 1e20 </FreqObs>

or like that:

<Position> 1000. 1000. 0.15. 0.</Position>
<FieldOfView> 0.3141592653589793 </FieldOfView>
<Resolution> 128 </Resolution>
<Spectrometer kind="freqlog" nsamples="10"> 17. 23. </Spectrometer>
<FreqObs> 1e20 </FreqObs>

Units can be specified using the unit attribute in the XML file, for instance:

<Distance unit="kpc"> 8 </Distance>

Possible units are (with [] noting the default):

When the distance is really large and most of the ray-tracing would happen de facto in flat space, the camera is transported to a location at a reasonable distance from the metric and the images are scaled accordingly. The default value for this distance should be fine, but it can be customized using the "dmax" attribute of the "Distance" element. "dmax" is always expressed in geometrical units:

<Distance unit="kpc" dmax="1e7"> 8 </Distance>

Symptoms when dmax is too large include pixelization of the image (neighbouring photons are numerically identical) and other numerical overflows. dmax is too small when it is apparent that changing it yields projection effects. dmax must be large compared to rmax in the Astrobj and ideally, changing it by an order of magnitude should not yield significant changes in the ray-traced image.

A mask may be used to limit ray-tracing to only some portions of the field. The Scenery checks whether a mask is to be used using Screen::operator()(size_t i, size_t j). The mask can be loaded from a FITS file as a square image of doubles:

<Mask>maskfile.fits</Mask>

The mask needs to be have the same size as the Screen itself, so loading a mask also sets the resolution, and changing the resolution after setting a mask also removes the mask. The content of the Mask entity is parsed by Factory::fullPath(), so it can be an absolute path, a path relative to where the XML file is stored, or relative to the current working directory if prefixed with "`pwd`/".

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

◆ coordKind()

int Gyoto::Screen::coordKind ( ) const

Get coordinate kind.

From Screen::gg_.

◆ coordToSky()

void Gyoto::Screen::coordToSky ( const double  pos[4],
double  dest[3],
bool  geometrical = false 
) const

Convert metric 4-position to sky 3-position.

Parameters
[in]pos4-position in metric coordinates.
[in]dest3-position in plane of the sky: Cartesian East, North, front.
[in]geometricalif true, #dest will be in geometrical units instead of meters.

◆ 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.

◆ distance() [1/4]

void Gyoto::Screen::distance ( double  dist)

Set distance from observer.

Parameters
distDistance in meters.

◆ distance() [2/4]

void Gyoto::Screen::distance ( double  dist,
const std::string &  unit 
)

Set distance from observer.

Parameters
distthe distance expressed in the specified unit;
unitconvertible to meters

◆ distance() [3/4]

double Gyoto::Screen::distance ( ) const

Get distance from observer.

In meters.

◆ distance() [4/4]

double Gyoto::Screen::distance ( const std::string &  ) const

Get distance from observer.

In specified unit.Get distance from observer

◆ dMax() [1/2]

void Gyoto::Screen::dMax ( double  dist)

Set ray-tracing maximum distance.

Parameters
distDistance in geometrical units.

◆ dMax() [2/2]

double Gyoto::Screen::dMax ( ) const

Get maximum ray-tracing distance.

In geometrical units.

◆ fillElement()

virtual void Gyoto::Object::fillElement ( Gyoto::FactoryMessenger fmp) const
virtualinherited

Fill the XML element for this Object.

The base implementation simply calls fillProperty() for each Property defined for the Object.

Derived classes should avoid overriding fillElement(). It may make sense occasionally, e.g. to make sure that the metric is output first.

To customize how a given Property is rendered, it is better to override fillProperty().

If this method is overridden, the implementation should in general call fillElement() on the direct base.

Reimplemented in Gyoto::Scenery, Gyoto::Spectrometer::Complex, Gyoto::Astrobj::Complex, and Gyoto::Metric::Complex.

◆ fillProperty()

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

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;
Property const &p) const ;
};
void B::fillProperty(Gyoto::FactoryMessenger *fmp,
Property const &p) const {
if (name=="Duff") fmp->doSomething();
else A::fillProperty(fmp, p);
}

Reimplemented from Gyoto::Object.

◆ freqObs() [1/3]

void Gyoto::Screen::freqObs ( double  fo)

Set freq_obs_.

Parameters
fodouble: observing frequency in Hz

◆ freqObs() [2/3]

void Gyoto::Screen::freqObs ( double  fo,
const std::string &  unit 
)

Set freq_obs_.

Parameters
fodouble: observing frequency (or wavelength) in "unit"
unitstring: unit in which fo is expressed, convertible to Herz or meters or energy.

◆ freqObs() [3/3]

double Gyoto::Screen::freqObs ( const std::string &  unit) const

Get freq_obs_.

Parameters
unitstring: unit in which freq_obs_ should be returned is expressed, convertible to Herz or meters or energy.

◆ getFourVel()

void Gyoto::Screen::getFourVel ( double  dest[4]) const

Get copy of Screen::fourvel_.

Parameters
[out]fourvelpreallocated 4-element array

◆ getObserverPos()

void Gyoto::Screen::getObserverPos ( double  dest[4]) const

4-Position of the observer relative to the metric

A Screen is positioned relative to the observer with four elements: Screen::distance, Screen::inclination, Screen::paln and Screen::argument.

This function returns the position of the observer relative to the metric system in Screen::gg_, using these parameters. The output parameter is coord.

Parameters
[out]coordposition of the observer. Must be preallocated.

◆ getProperties()

virtual Property const* Gyoto::Screen::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::Object.

◆ getRayCoord() [1/2]

void Gyoto::Screen::getRayCoord ( double  x,
double  y,
double  dest[8] 
) const

Get 8-coordinate of Photon hitting screen from a given direction.

Similar to Screen::getObserverPos() but will return in addition the 4-velocity of a photon corresponding to the sky direction given by x and y.

Parameters
[in]xRA (d_alpha*cos(delta)) offset in radians;
[in]yDec offset (d_delta) in radians;
[out]destposition-velocity of the observer Photon. Preallocated.

◆ getRayCoord() [2/2]

void Gyoto::Screen::getRayCoord ( const size_t  i,
const size_t  j,
double  dest[8] 
) const

Get 8-coordinate of Photon hitting screen pixel.

Similar to Screen::getObserverPos() but will return in addition the 4-velocity of a photon corresponding to the sky direction given by x and y.

Parameters
[in]i,jpixel coordinates
[out]destposition-velocity of the Photon. Preallocated.

◆ getRayTriad()

void Gyoto::Screen::getRayTriad ( double  coord[8],
double  Ephi[4],
double  Etheta[4] 
) const

Get polarization triad.

Computes the polarization triad (Ephi,Etheta,k) from the tangent to photon geodesic k.

Parameters
[in]coordPhoton coordinate + tangent vector (double[8])
[out]Ephifirst polarisation direction. Preallocated. Default: NULL.
[out]Ethetasecond polarisation direction. Preallocated. Default: NULL.

◆ getScreen1()

void Gyoto::Screen::getScreen1 ( double  dest[4]) const

Get copy of Screen::screen1_.

Parameters
[out]destpreallocated 4-element array

◆ getScreen2()

void Gyoto::Screen::getScreen2 ( double  dest[4]) const

Get copy of Screen::screen2_.

Parameters
[out]destpreallocated 4-element array

◆ getScreen3()

void Gyoto::Screen::getScreen3 ( double  dest[4]) const

Get copy of Screen::screen3_.

Parameters
[out]destpreallocated 4-element array

◆ help()

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

Print (to stdout) some help on this class.

Describe all properties that this instance supports.

◆ inclination() [1/4]

void Gyoto::Screen::inclination ( double  )

Set inclination relative to line-of-sight.

Inclination of z-axis relative to line-of-sight, or inclination of equatorial plane relative to plane of the sky, in radians

◆ inclination() [2/4]

void Gyoto::Screen::inclination ( double  ,
const std::string &  unit 
)

Set inclination relative to line-of-sight.

Inclination of z-axis relative to line-of-sight, or inclination of equatorial plane relative to plane of the sky, in specified unit.

◆ inclination() [3/4]

double Gyoto::Screen::inclination ( ) const

Get inclination relative to line-of-sight.

Inclination of z-axis relative to line-of-sight, or inclination of equatorial plane relative to plane of the sky, in radians.

◆ inclination() [4/4]

double Gyoto::Screen::inclination ( const std::string &  ) const

Get inclination relative to line-of-sight.

Inclination of z-axis relative to line-of-sight, or inclination of equatorial plane relative to plane of the sky, in specified unit.

◆ 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.

◆ mapPixUnit()

void Gyoto::Screen::mapPixUnit ( )

Map "pix" and "pixel" to angular pixel width in unit system.

"pix" or "pixel" can then be used in units.

There is only one unit system in Gyoto: "pix" can therefore be registered only for one Screen at a time. See Gyoto::Units.

The unit must later be unmapped with unmapPixUnit().

◆ mask()

void Gyoto::Screen::mask ( double const *const  mm,
size_t  resolution = 0 
)

Set mask_ from array.

mm will be copied. mm must be a square resolution x resolution array. If mm==NULL, just deallocate mask_.

◆ operator()()

bool Gyoto::Screen::operator() ( size_t  ,
size_t   
)

Whether this pixel should be ray-traced.

If mask_ is not set, always true. Else, true for non-zero cells in mask_.

◆ print()

std::ostream& Gyoto::Screen::print ( std::ostream &  ) const

Display.

Debug helper

◆ 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.

◆ setObserverPos()

void Gyoto::Screen::setObserverPos ( const double  pos[4])

Alternative way to set projection.

Beware : paln can not be set this way, setting later other parameters change the observer's coordinates. For observationnal ray-tracing purposes, prefer setProjection().

Parameters
[in]posposition of observer in Screen's coordinate system. Content is copied.

◆ setParameter() [1/2]

virtual int Gyoto::Object::setParameter ( std::string  name,
std::string  content,
std::string  unit 
)
virtualinherited

Set parameter by name.

This function is used when parsing an XML description, if no Property of this name is found. Overriding implementation should fall-back on calling the direct's parent implementation:

class A: public Object {};
class B: public A {
using B::setParameter;
virtual int setParameter(std::string name,
std::string content,
std::string unit);
};
int B::setParameter(std::string name,
std::string content,
std::string unit) {
if (name=="Duff") doSomething(content, unit);
else return A::setParameter(name, content, unit);
return 0; // name was known
}
Parameters
nameXML name of the parameter (XML entity). This may have a path component, e.g. "Astrobj::Radius", in which case a property named "Astrobj" will be sought in the current object, and setParameter will be called recusrively on this Astrobj with Radius as name.
contentstring representation of the value
unitstring representation of the unit
Returns
0 if this parameter is known, 1 if it is not.

Reimplemented in Gyoto::Astrobj::Star, Gyoto::Metric::RotStar3_1, Gyoto::Metric::KerrKS, and Gyoto::Astrobj::EquatorialHotSpot.

◆ setParameter() [2/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.

◆ setParameters()

virtual void Gyoto::Object::setParameters ( Gyoto::FactoryMessenger fmp)
virtualinherited

Main loop for parsing Properties from XML description.

This function queries the FactoryMessenger for elements to parse, and tries to matche each element to a Property to set it accordingly.

Any class that tries to be buildable from XML must supply a subcontractor (for base classes such as Metric, Astrobj, Spectrum and Spectrometer, it is done as a template that must be specialized for each class).

This subcontractor typically looks somewhat like this:

SmartPointer<Metric::Generic>
Gyoto::Metric::MyKind::Subcontractor(FactoryMessenger* fmp) {
SmartPointer<MyKind> gg = new MyKind();
gg -> setParameters(fmp);
return gg;
}
 Although this is discouraged, it is possible to override the
 following functions to customize how XML entities are parsed:
   - setParameters() if low-level access to the
     FactoryMessenger is required;
   - setParameter(std::string name,
                  std::string content,
                  std::string unit)
     to interpret an entity that does not match a Property
     (e.g. alternative name);
   - setParameter(Gyoto::Property const &p,
                  std::string const &name,
                  std::string const &content,
                  std::string const &unit)
     to change how a Property is interpreted.

Reimplemented in Gyoto::Astrobj::Generic, Gyoto::Photon, Gyoto::Astrobj::Star, Gyoto::Spectrometer::Complex, Gyoto::Spectrometer::Uniform, Gyoto::Astrobj::Complex, Gyoto::Astrobj::OscilTorus, Gyoto::Astrobj::EquatorialHotSpot, Gyoto::Metric::Complex, and Gyoto::Metric::Shift.

◆ skyToCoord()

void Gyoto::Screen::skyToCoord ( const double  sky[3],
double  dest[4],
bool  geometrical = false 
) const

Convert sky 3-position to metric 4-position.

Parameters
[in]sky3-position in plane of the sky.
[in]dest4-position in metric coordinates (dest[0] is not modified).
[in]geometricalset to true if #sky is in geometrical units instead of meters.

◆ unmapPixUnit()

void Gyoto::Screen::unmapPixUnit ( )

Unmap "pix" and "pixel" from unit system.

See also mapPixUnit().

Member Data Documentation

◆ euler_

double Gyoto::Screen::euler_[3]
private

Euler angles.

The angles are position angle of the line of nodes (North of East), inclination (0 = face-on), argument of X axis. We use the z-x-z convention. See http://en.wikipedia.org/wiki/Euler_angles

◆ freq_obs_

double Gyoto::Screen::freq_obs_
private

Frequency at which the observer observes.

For the quantity Intensity

◆ 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>

is written.

◆ mask_filename_

std::string Gyoto::Screen::mask_filename_
private

Last read or written FITS file.

Used when saving to XML: if the mask was saved or loaded from FITS file, output this file name in the XML.

◆ 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: