GLAMERDOC++
Gravitational Lensing Code Library
Public Member Functions | Static Public Member Functions | Friends | List of all members
PixelMap Class Reference

Image structure that can be manipulated and exported to/from fits files. More...

#include <image_processing.h>

Public Member Functions

 PixelMap (const PixelMap &pmap, double res_ratio)
 Creates a PixelMap at a different resolution. The new counts are calculated integrating over the input pixels. No interpolation or smoothing is performed. More...
 
 PixelMap (const PixelMap &other)
 
 PixelMap (PixelMap &&other)
 
 PixelMap (const PixelMap &pmap, const double *center, std::size_t Npixels)
 Creates a new PixelMap from a square region of a PixelMap. If the region exceeds the boundaries of the original map, the new map is completed with zeros. More...
 
 PixelMap (const PixelMap &pmap, long nx, long ny, std::size_t Npixels)
 Produces a square cut-out of the input PixelMap. More...
 
 PixelMap (const double *center, std::size_t Npixels, double resolution, PixelMapUnits u=PixelMapUnits::ndef)
 make square PixelMap More...
 
 PixelMap (const double *center, std::size_t Nx, std::size_t Ny, double resolution, PixelMapUnits u=PixelMapUnits::ndef)
 make rectangular PixelMap with square pixels More...
 
 PixelMap (std::string fitsfilename, double resolution=-1, PixelMapUnits u=PixelMapUnits::ndef)
 Constructs a PixelMap reading in a fits file Infos about resolution, Npixels and center are read from the header. More...
 
PixelMapoperator= (const PixelMap &other)
 
PixelMapoperator= (PixelMap &&other)
 
void ChangeUnits (PixelMapUnits u)
 
PixelMapUnits getUnits () const
 
bool valid () const
 
std::size_t size () const
 
std::size_t getNx () const
 
std::size_t getNy () const
 
double getRangeX () const
 
double getRangeY () const
 
void const getCenter (Point_2d &c) const
 
Point_2d getCenter () const
 
double getResolution () const
 
double getRA ()
 returns right accention of center
 
double getDEC ()
 returns declination of center
 
void setRAandDec (double RAin, double DECin)
 set the coordinates of center
 
void Clean ()
 
void AddImages (ImageInfo *imageinfo, int Nimages, float rescale=1.)
 Add an image to the map. More...
 
void AddImages (std::vector< ImageInfo > &imageinfo, int Nimages, float rescale=1.)
 
void AddGridBrightness (Grid &grid)
 Add an image from a the surface brightnesses of a Grid to the PixelMap.
 
void AddGridMapBrightness (const GridMap &grid)
 Add an image from a the surface brightnesses of a GridMap to the PixelMap.
 
void AddUniformImages (ImageInfo *imageinfo, int Nimages, double value)
 Add images with uniform surface brightness set by input parameter value. More...
 
PosType AddSource (Source &source)
 
PosType AddSource (Source &source, int oversample)
 Add a source to the pixel map by oversamples the source so that oversample^2 points within each pixel are averaged.
 
void AddPointSource (const Point_2d &x, double flux)
 
void copy_in (const PixelMap &pmap)
 copy a PixelMap into this one. More...
 
void paste (const PixelMap &pmap)
 Replace overlaping pixel values with those of the input map. More...
 
void paste (const PixelMap &pmap, long nx, long ny)
 paste a PixelMap on with the lower left pixel match to [nx,ny] of this
 
void duplicate (const PixelMap &pmap)
 copy a PixelMap that must be the same without creating a new one.. More...
 
template<typename T >
PosType AddSource_parallel (T &source, int oversample)
 Adds source to map. This version breaks pixels up into blocks and does them in seporate threads.
 
void AddCurve (ImageInfo *curve, double value)
 Draws a closed curve through the points in curve->imagekist. More...
 
void AddCurve (Kist< Point > *imagekist, PosType value)
 
void AddCurve (std::vector< Point_2d > &curve, double value)
 
void AddCurve (std::vector< RAY > &curve, double value)
 
void drawline (double x1[], double x2[], double value, bool add)
 simple line More...
 
void DrawLine (long x0, long x1, long y0, long y1, double value, bool add)
 line by Bresenham's line algorithm
 
void DrawLineGS (long x0, long x1, long y0, long y1, double value, bool add)
 
void drawcircle (PosType r_center[], PosType radius, PosType value)
 Draws a circle. More...
 
void drawdisk (PosType r_center[], PosType radius, PosType value, int Nstrip)
 Draws a disk. More...
 
void AddGrid (const Grid &grid, double value=1.0)
 Fills in pixels where the image plane points in the grid are located with the value given. More...
 
void AddGrid (const Grid &grid, LensingVariable val)
 Fills in pixels with the selected quantity from the grid points. More...
 
void Renormalize (double factor)
 Multiplies the whole map by a scalar factor.
 
void AddValue (std::size_t i, double value)
 Adds a value to the i-th pixel.
 
void AssignValue (std::size_t i, double value)
 Assigns a value to the i-th pixel.
 
void printASCII () const
 Print an ASCII table of all the pixel values.
 
void printASCIItoFile (std::string filename) const
 Print an ASCII table of all the pixel values.
 
void printFITS (std::string filename, bool Xflip=false, bool verbose=false)
 Output the pixel map as a fits file.
 
void printFITS (std::string filename, std::vector< std::tuple< std::string, double, std::string > > &extra_header_info, bool verbose)
 
void printFITS (std::string filename, std::vector< std::string > &headercards)
 This overides all header information and relaces it with the inputs. Meant for making a modified copy. More...
 
void smooth (double sigma)
 Smoothes a map with a Gaussian kernel of width sigma (in arcseconds)
 
double getValue (std::size_t i) const
 
double & operator[] (std::size_t i)
 
double operator[] (std::size_t i) const
 
double operator() (std::size_t i) const
 
double operator() (std::size_t i, std::size_t j) const
 
double & operator() (std::size_t i, std::size_t j)
 
PixelMapoperator+= (const PixelMap &rhs)
 Add the values of another PixelMap to this one.
 
void operator+= (double f)
 
PixelMap operator+ (const PixelMap &) const
 Add two PixelMaps.
 
PixelMapoperator-= (const PixelMap &rhs)
 Subtract the values of another PixelMap from this one.
 
PixelMap operator- (const PixelMap &) const
 Subtract two PixelMaps.
 
PixelMapoperator*= (const PixelMap &rhs)
 Multiply the values of another PixelMap by this one.
 
PixelMap operator* (const PixelMap &a) const
 Multiply two PixelMaps.
 
PixelMap operator/ (const PixelMap &a) const
 Multiply two PixelMaps.
 
PixelMapoperator*= (PosType b)
 
PixelMap operator* (PosType b) const
 
std::valarray< double > & data ()
 
bool agrees (const PixelMap &other) const
 Check whether two PixelMaps agree in their physical dimensions.
 
PosType ave () const
 return average pixel value
 
PosType sum () const
 return sum of all pixel values
 
size_t size ()
 Total number of pixels.
 
double max () const
 
double min () const
 
void FindArc (PosType &radius, PosType *xc, PosType *arc_center, PosType &arclength, PosType &width, PosType threshold)
 Find arcs in image WARNING: THIS IS UNDER CONSTRUCTION!
 
long find_index (PosType const x[], long &ix, long &iy) const
 get the index for a position, returns -1 if out of map, this version returns the 2D grid coordinates More...
 
long find_index (PosType const x[]) const
 get the index for a position, returns -1 if out of map
 
long find_index (PosType x, PosType y, long &ix, long &iy) const
 get the index for a position, returns -1 if out of map, this version returns the 2D grid coordinates More...
 
long find_index (PosType x, PosType y) const
 get the index for a position, returns -1 if out of map
 
void find_position (PosType x[], std::size_t const index) const
 get the index for a position, returns -1 if out of map
 
void find_position (PosType x[], std::size_t const ix, std::size_t const iy) const
 get the index for a position, returns -1 if out of map
 
PixelMap rotate (PosType theta, PosType scale=1)
 rotate and scale the image while keeping pixels, resoluiton More...
 
PosType linear_interpolate (PosType x[])
 interpolate to point x[] More...
 
void drawgrid (int N, PosType value)
 draw a grid on the image that divides the each demension into N cells More...
 
void drawPoints (std::vector< Point * > points, PosType size, PosType value)
 
void drawPoints (std::vector< Point > points, PosType size, PosType value)
 
void drawCurve (std::vector< Point * > points, PosType value)
 
void drawCurve (std::vector< Point > points, PosType value)
 
void drawPoints (std::vector< Point_2d > points, PosType size, PosType value)
 
void drawCurve (std::vector< Point_2d > points, PosType value)
 
void drawSquare (PosType p1[], PosType p2[], PosType value)
 Draw a rectangle. More...
 
void drawBox (PosType p1[], PosType p2[], PosType value, int Nstrip)
 Draws a box (filling the inside with horizontal lines, starting from the top)
 
void PowerSpectrum (std::vector< PosType > &power_spectrum, std::vector< PosType > &lvec, bool overwrite=true)
 Find the power spectrum of the map. More...
 
void PowerSpectrum (std::vector< PosType > &power_spectrum, const std::vector< PosType > &lbins, std::vector< PosType > &lave, bool overwrite=true)
 Find the power spectrum of the map. More...
 
void AdaptiveSmooth (PosType value)
 
void find_contour (double level, std::vector< std::vector< Point_2d > > &points, std::vector< bool > &hits_edge) const
 returns a vector of contour curves
 
void find_islands_holes (double level, std::vector< std::vector< size_t > > &points) const
 
void lens_definition (double min_sn_per_image, double pixel_threshold, int &Nimages, double &total_sig_noise_source, std::vector< size_t > &maxima_indexes, std::vector< std::vector< size_t > > &image_points, bool &lens_TF, double &level, size_t &n_pix_in_source, bool verbose=false)
 
std::vector< size_t > maxima (double minlevel) const
 find maxima that are above minlevel More...
 
int count_islands (std::vector< size_t > &pixel_index) const
 For a list of pixel indexes this will count and separated islands that are not connected. More...
 
size_t threshold (std::list< size_t > &pixel_index, PosType value)
 get a list of pixels above value
 
void flipY ()
 reflects the image about the horizontal mid-line
 
void flipX ()
 reflects the image about the vertical mid-line
 
void doubleFlip ()
 rotate the image by 180deg or equivalently reflect it through the origin
 
void recenter (PosType newcenter[2])
 recenter the map without changing anything else. More...
 
void recenter (Point_2d newcenter)
 
void convolve (PixelMap &kernel, long center_x=0, long center_y=0)
 convolve the image with a kernel. More...
 
PixelMap downsize (int n)
 Creates a PixelMap with a lower resolution. The value of the pixels are added for the new pixels. If n does not go into the orginial number of pixels evenly the right (top) redge is dropped. More...
 
void addheader (std::string label, long value, std::string comment)
 add a heaader keyword that will appear in fits output
 
void addheader (std::string label, size_t value, std::string comment)
 
void addheader (std::string label, float value, std::string comment)
 
void addheader (std::string label, double value, std::string comment)
 
void addheader (std::string label, std::string value, std::string comment)
 

Static Public Member Functions

static void swap (PixelMap &, PixelMap &)
 

Friends

void swap (PixelMap &, PixelMap &)
 

Detailed Description

Image structure that can be manipulated and exported to/from fits files.

Constructor & Destructor Documentation

◆ PixelMap() [1/6]

PixelMap::PixelMap ( const PixelMap pmap,
double  res_ratio 
)

Creates a PixelMap at a different resolution. The new counts are calculated integrating over the input pixels. No interpolation or smoothing is performed.

Parameters
res_ratioresolution of map is res_ratio times the resolution of the input map

◆ PixelMap() [2/6]

PixelMap::PixelMap ( const PixelMap pmap,
const double *  center,
std::size_t  my_Npixels 
)

Creates a new PixelMap from a square region of a PixelMap. If the region exceeds the boundaries of the original map, the new map is completed with zeros.

Parameters
centerInput PixelMap (from which the stamp is taken)
my_Npixelscenter of the region to be duplicated (in rads) size of the region to be duplicated (in pixels)

◆ PixelMap() [3/6]

PixelMap::PixelMap ( const PixelMap pmap,
long  nx,
long  ny,
std::size_t  my_Npixels 
)

Produces a square cut-out of the input PixelMap.

Parameters
nxInput PixelMap (from which the stamp is taken)
nylower left pixels of in pmap
my_Npixelslower left pixels of in pmap size of the region to be duplicated (in pixels)

◆ PixelMap() [4/6]

PixelMap::PixelMap ( const double *  center,
std::size_t  Npixels,
double  resolution,
PixelMapUnits  u = PixelMapUnits::ndef 
)

make square PixelMap

Parameters
NpixelsThe location of the center of the map
resolutionNumber of pixels in one dimension of map. One dimensional range of map in whatever units the point positions are in

◆ PixelMap() [5/6]

PixelMap::PixelMap ( const double *  center,
std::size_t  Nx,
std::size_t  Ny,
double  resolution,
PixelMapUnits  u = PixelMapUnits::ndef 
)

make rectangular PixelMap with square pixels

Parameters
NxThe location of the center of the map
NyNumber of pixels in x dimension of map.
resolutionNumber of pixels in y dimension of map. One dimensional range of map in whatever units the point positions are in

◆ PixelMap() [6/6]

PixelMap::PixelMap ( std::string  fitsfilename,
double  my_res = -1,
PixelMapUnits  u = PixelMapUnits::ndef 
)

Constructs a PixelMap reading in a fits file Infos about resolution, Npixels and center are read from the header.

Parameters
fitsfilenamefile name of fits file to be read
my_resresolution (rad) of fits image if not given in fits file, use default or -1 otherwise

Member Function Documentation

◆ AddCurve()

void PixelMap::AddCurve ( ImageInfo curve,
double  value 
)

Draws a closed curve through the points in curve->imagekist.

This differs form PixelMap::AddImage() in that it draws lines between the points and takes no account of the cells that the points are in or the surface brightness. The points must be ordered already. Particularly useful for drawing the caustics that may have irregular cell sizes. The last point will be connected to the first point.

◆ AddGrid() [1/2]

void PixelMap::AddGrid ( const Grid grid,
double  value = 1.0 
)

Fills in pixels where the image plane points in the grid are located with the value given.

This is for lensing quantities and not surface brightness. If you want surface brightness use PixelMap::AddGridBrightness()

◆ AddGrid() [2/2]

void PixelMap::AddGrid ( const Grid grid,
LensingVariable  val 
)

Fills in pixels with the selected quantity from the grid points.

The grid and PixelMap do not need to be related in any way. Using this function multiple grids can be added to the same image.

This is for lensing quantities and not surface brightness. If you want surface brightness use PixelMap::AddGridBrightness()

Warning: When adding a new grid it should not overlap with any of the previously added grids.

◆ AddImages() [1/2]

void PixelMap::AddImages ( ImageInfo imageinfo,
int  Nimages,
float  rescale = 1. 
)

Add an image to the map.

If rescale==0 gives constant surface brightness, if < 0 the surface brightness is not scaled by the pixel area as for the flux (default: 1). Negative values are good for mapping some quantity independant of the pixel size

Parameters
imageinfoAn array of ImageInfo-s. There is no reason to separate images for this routine
NimagesNumber of images on input.
rescalerescales the surface brightness while leaving the image unchanged, see full notes

◆ AddImages() [2/2]

void PixelMap::AddImages ( std::vector< ImageInfo > &  imageinfo,
int  Nimages,
float  rescale = 1. 
)
Parameters
imageinfoAn array of ImageInfo-s. There is no reason to separate images for this routine
NimagesNumber of images on input.
rescalerescales the surface brightness while leaving the image unchanged, see full notes

◆ AddUniformImages()

void PixelMap::AddUniformImages ( ImageInfo imageinfo,
int  Nimages,
double  value 
)

Add images with uniform surface brightness set by input parameter value.

This does not use the surface brightnesses stored in the image points.

Parameters
imageinfoAn array of ImageInfo-s. There is no reason to separate images for this routine

◆ convolve()

void PixelMap::convolve ( PixelMap kernel,
long  center_x = 0,
long  center_y = 0 
)

convolve the image with a kernel.

It is assumed that the size of the kernel is much smaller than the image and that the kernal has the same pixel size as the image.

◆ copy_in()

void PixelMap::copy_in ( const PixelMap pmap)

copy a PixelMap into this one.

The size, resolution and center of the pixel maps are not changed and do not need to match. The input pixel map is added while conserving the area integral of the map within the area of overlaping pixels.

◆ count_islands()

int PixelMap::count_islands ( std::vector< size_t > &  pixel_index) const

For a list of pixel indexes this will count and separated islands that are not connected.

On return, 'pixel_index' is ordered into groups and the 'heads' list points to the first elemant in each group plus the end of the list so that heads[i] to heads[i+1] is a group for 0 <= i <= ngroups. The number of groups is returned which is also heads.size() - 1

◆ downsize()

PixelMap PixelMap::downsize ( int  n)

Creates a PixelMap with a lower resolution. The value of the pixels are added for the new pixels. If n does not go into the orginial number of pixels evenly the right (top) redge is dropped.

Parameters
nnumber of pixels each direction added into each new pixel

◆ drawcircle()

void PixelMap::drawcircle ( PosType  r_center[],
PosType  radius,
PosType  value 
)

Draws a circle.

Parameters
r_centercenter of circle
radiusradius of circle
valuevalue that it is set to on the map

◆ drawdisk()

void PixelMap::drawdisk ( PosType  r_center[],
PosType  radius,
PosType  value,
int  Nstrip 
)

Draws a disk.

Parameters
r_centercenter of disk
radiusradius of disk
valuevalue that it is set to on the map
Nstripnumber of lines we want

◆ drawgrid()

void PixelMap::drawgrid ( int  N,
PosType  value 
)

draw a grid on the image that divides the each demension into N cells

Draws a grid.

◆ drawline()

void PixelMap::drawline ( double  x1[],
double  x2[],
double  value,
bool  add 
)

simple line

Draws a line between two points on the image by setting the pixels equal to value.

TODO: Could be improved by detecting if the line passes through the map at all before starting. Could also be improved by making the line fatter by including neighbor points.

Parameters
x1one end point of line
x2other end point of line
valuevalue that it is set to on the map
addtrue : add value, false replace with value

◆ drawSquare()

void PixelMap::drawSquare ( PosType  p1[],
PosType  p2[],
PosType  value 
)

Draw a rectangle.

Draws a square.

◆ duplicate()

void PixelMap::duplicate ( const PixelMap pmap)

copy a PixelMap that must be the same without creating a new one..

This avoids calling a any constructor or destructor.

◆ find_index() [1/2]

long PixelMap::find_index ( PosType const  x[],
long &  ix,
long &  iy 
) const

get the index for a position, returns -1 if out of map, this version returns the 2D grid coordinates

get the index for a position, returns -1 if out of map

◆ find_index() [2/2]

long PixelMap::find_index ( PosType  x,
PosType  y,
long &  ix,
long &  iy 
) const

get the index for a position, returns -1 if out of map, this version returns the 2D grid coordinates

get the index for a position, returns -1 if out of map

◆ find_islands_holes()

void PixelMap::find_islands_holes ( double  level,
std::vector< std::vector< size_t > > &  points 
) const

find all the points above level divided into seprated groups

Groups with points are connected regions above level. Groups without points are regions surrounded by regions above level, i.e. holes.

◆ lens_definition()

void PixelMap::lens_definition ( double  min_sn_per_image,
double  pixel_threshold,
int &  Nimages,
double &  total_sig_noise_source,
std::vector< size_t > &  maxima_indexes,
std::vector< std::vector< size_t > > &  image_points,
bool &  lens_TF,
double &  level,
size_t &  n_pix_in_source,
bool  verbose = false 
)

This applies a definition of lesning for a resolved source based on that of Sonnenfeld et al. 2023

This is meant to be used on a signal-to-noise map of the lensed source only.

The definition does the following tests: 1) define a footprint at surface brightness level pixel_threshold 2) romove disconected region that have signal-to-noise below min_sn_per_image 3) if there are multiple images or a hole in the image lens_TF =true 4) if there is one image with no hole increase the threshold and apply 3) until the it is classified as a lens or it reaches the maximum surface brightness level

Parameters
min_sn_per_imagesignal-to-noise required for a seporate image (ex 10)
pixel_thresholdsignal-to-noise threshold that defines the footprint (ex. 2)
Nimagesthe number of images
total_sig_noise_sourcegives the total signal-to-noise of all images
maxima_indexesindex of maxima
lens_TFwhether it passes lens diffintion
levellevels on which the multiple images are defined
n_pix_in_sourcenumber of pixels in footprint

◆ linear_interpolate()

PosType PixelMap::linear_interpolate ( PosType  x[])

interpolate to point x[]

bilinear interpolation

◆ maxima()

std::vector< size_t > PixelMap::maxima ( double  minlevel) const

find maxima that are above minlevel

find the index of the pixels that are larger than all its neighbors

◆ paste()

void PixelMap::paste ( const PixelMap pmap)

Replace overlaping pixel values with those of the input map.

No attempt is made to interpolate or average the pixels of pmap. No integration is done. If the resolution of input map is higher than the reslution of the map then the pixels values will just be that of the last pixel visited while going through them.

This can be used to sew tiles together into a larger map.

◆ PowerSpectrum() [1/2]

void PixelMap::PowerSpectrum ( std::vector< PosType > &  power_spectrum,
const std::vector< PosType > &  lbins,
std::vector< PosType > &  lave,
bool  overwrite = true 
)
inline

Find the power spectrum of the map.

Parameters
power_spectrumoutput power spectrum
lbinsoutput l values of bands
laveoutput l values of bands
overwriteif false add power to existing power_spectrum (used for averaging over many fields

◆ PowerSpectrum() [2/2]

void PixelMap::PowerSpectrum ( std::vector< PosType > &  power_spectrum,
std::vector< PosType > &  lvec,
bool  overwrite = true 
)
inline

Find the power spectrum of the map.

Parameters
power_spectrumoutput power spectrum
lvecoutput l values of bands
overwriteif false add power to existing power_spectrum (used for averaging over many fields

◆ printFITS()

void PixelMap::printFITS ( std::string  filename,
std::vector< std::string > &  headercards 
)

This overides all header information and relaces it with the inputs. Meant for making a modified copy.

Parameters
filenamefile to create
headercardsheader information in cfitsio "card" format

◆ recenter() [1/2]

void PixelMap::recenter ( Point_2d  newcenter)
Parameters
newcenterin radians

◆ recenter() [2/2]

void PixelMap::recenter ( PosType  newcenter[2])

recenter the map without changing anything else.

Parameters
newcenterin radians

◆ rotate()

PixelMap PixelMap::rotate ( PosType  theta,
PosType  scale = 1 
)

rotate and scale the image while keeping pixels, resoluiton

Parameters
thetacounter-clockwise rotation (radians)
scalescale <1 shrinks it

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