mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
imageUtils.hpp
Go to the documentation of this file.
1/** \file imageUtils.hpp
2 * \author Jared R. Males
3 * \brief Header for the image processing utilities
4 * \ingroup image_processing_files
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2020 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#ifndef improc_imageUtils_hpp
28#define improc_imageUtils_hpp
29
30#include <cstdint>
31#include <cmath>
32
33#include "imageTransforms.hpp"
34
35namespace mx
36{
37namespace improc
38{
39
40/// Get the mxlib standard center coordinate of an image
41/**
42 * \returns 0.5*(rows_cols-1)
43 */
44template <typename realT = double>
45realT imCen( int rows_cols /**< [in] the number of rows or columns in the image */ )
46{
47 return 0.5 * ( 1.0 * rows_cols - 1 );
48}
49
50/// Get the mxlib standard center x-coordinate of an image
51/**
52 * \returns 0.5*(rows-1)
53 */
54template <typename realT = double, typename imT>
55realT imCenX( const imT &im /**< [in] the image to find the center of */ )
56{
57 return imCen<realT>( im.rows() );
58}
59
60/// Get the mxlib standard center y-coordinate of an image
61/**
62 * \returns 0.5*(cols-1)
63 */
64template <typename realT = double, typename imT>
65realT imCenY( const imT &im /**< [in] the image to find the center of */ )
66{
67 return imCen<realT>( im.cols() );
68}
69
70/** \ingroup image_utils
71 *@{
72 */
73
74template <typename T>
75constexpr T invalidNumber()
76{
77 return -3e38;
78}
79
80/// Check if the number is nan, using several different methods
81inline bool IsNan( float value )
82{
83 return ( ( ( ( *(uint32_t *)&value ) & 0x7fffffff ) > 0x7f800000 ) || ( value == invalidNumber<float>() ) ||
84 !std::isfinite( value ) );
85}
86
87/// Reflect pixel coordinates across the given center pixel.
88/**
89 */
90template <typename realT>
91int reflectImageCoords( int &x1, ///< [out] the reflected x coordinate
92 int &y1, ///< [out] the reflected y coordinate
93 int x0, ///< [in] the input x coordinate
94 int y0, ///< [in] the input y coordinate
95 realT xc, ///< [in] the center pixel x coordinate
96 realT yc ///< [in] the center pixel y coordinate
97)
98{
99 x1 = xc - ( x0 - xc );
100 y1 = yc - ( y0 - yc );
101
102 return 0;
103}
104
105/// Zero any NaNs in an image
106/**
107 * \overload
108 */
109template <class imageT, typename valueT>
110void zeroNaNs( imageT &im, ///< [in.out] image which will have any NaN pixels set to zero
111 valueT val ///< [in] [optional] The value to set NaN pixels too. Default is 0.
112)
113{
114 for( int c = 0; c < im.cols(); ++c )
115 {
116 for( int r = 0; r < im.rows(); ++r )
117 {
118 if( IsNan( im( r, c ) ) )
119 {
120 im( r, c ) = val;
121 }
122 }
123 }
124}
125
126/// Zero any NaNs in an image
127template <class imageT>
128void zeroNaNs( imageT &im /**< [in.out] image which will have any NaN pixels set to zero */ )
129{
130 typename imageT::Scalar zero = 0;
131 zeroNaNs<imageT, typename imageT::Scalar>( im, zero );
132}
133
134/// Zero any NaNs in an image cube
135/** This version fills in a mask with 1s where there were nans, 0s elsewhere.
136 *
137 */
138template <class cubeT, class maskCubeT>
139void zeroNaNCube( cubeT &imc, /**< [in.out] cube which will have any NaN pixels set to zero */
140 maskCubeT *mask /**< [out] a 1/0 mask with 1 indicating which pixels where nan */
141)
142{
143 if( mask )
144 {
145 mask->resize( imc.rows(), imc.cols(), imc.planes() );
146 mask->setZero();
147 }
148
149 for( int p = 0; p < imc.planes(); ++p )
150 {
151 for( int c = 0; c < imc.cols(); ++c )
152 {
153 for( int r = 0; r < imc.rows(); ++r )
154 {
155 if( IsNan( imc.image( p )( r, c ) ) )
156 {
157 imc.image( p )( r, c ) = 0;
158 if( mask )
159 {
160 ( *mask ).image( p )( r, c ) = 1;
161 }
162 }
163 }
164 }
165 }
166}
167
168/// Zero any NaNs in an image cube
169template <class cubeT>
170void zeroNaNCube( cubeT &imc /**< [in.out] cube which will have any NaN pixels set to zero */ )
171{
172 return zeroNaNCube<cubeT, cubeT>( imc, nullptr );
173}
174
175/// Calculate the mean value of an image
176/**
177 *
178 * \returns the mean of the input image
179 */
180template <class calcT, class imageT>
181calcT imageMean( imageT &im /**< [in] the image of which to calculate the mean*/ )
182{
183 return static_cast<calcT>( im.sum() ) / ( im.rows() * im.cols() );
184}
185
186/// Calculate the mean value of an image over a mask
187/**
188 *
189 * \returns the mean of the input image in the masked region
190 */
191template <class calcT, class imageT, class maskT>
192calcT imageMean( imageT &im, /**< [in] the image of which to calculate the mean*/
193 const maskT &mask /**< [in] a 1/0 mask where 1 defines the good pixels */ )
194{
195 return static_cast<calcT>( ( im * mask ).sum() ) / ( mask.sum() );
196}
197
198/// Calculate the variance of an image given its mean
199/**
200 *
201 * \returns the variance of the input image
202 */
203template <typename calcT, class imageT>
204calcT imageVariance( imageT &im /**< [in] the image of which to calculate the variance*/,
205 calcT mn /**< [in] the mean value of the image w.r.t. which to calcualate the variance */
206)
207{
208 return ( im.template cast<calcT>() - mn ).square().sum() / ( im.rows() * im.cols() );
209}
210
211/// Calculate the variance of an image given its mean
212/**
213 *
214 * \returns the variance of the input image
215 */
216template <typename calcT, class imageT, class maskT>
217calcT imageVariance( imageT &im /**< [in] the image of which to calculate the variance*/,
218 calcT mn, /**< [in] the mean value of the image w.r.t. which to calcualate the variance */
219 const maskT &mask /**< [in] a 1/0 mask where 1 defines the good pixels */
220)
221{
222 return ( im.template cast<calcT>() * mask - mn ).square().sum() / ( mask.sum() );
223}
224
225/// Calculate the center of light of an image
226/** Note that the sum of the image should be > 0.
227 *
228 */
229template <typename imageT>
230int imageCenterOfLight( typename imageT::Scalar &x, ///< [out] the x coordinate of the center of light [pixels]
231 typename imageT::Scalar &y, ///< [out] the y coordinate of hte center of light [pixels]
232 const imageT &im ///< [in] the image to centroid
233)
234{
235 x = 0;
236 y = 0;
237
238 typename imageT::Scalar sum = im.sum();
239
240 if( sum == 0 )
241 {
242 x = 0.5 * ( im.rows() - 1.0 );
243 y = 0.5 * ( im.cols() - 1.0 );
244 return 0;
245 }
246
247 for( int j = 0; j < im.cols(); ++j )
248 {
249 for( int i = 0; i < im.rows(); ++i )
250 {
251 x += ( i + 1 ) * im( i, j );
252 y += ( j + 1 ) * im( i, j );
253 }
254 }
255
256 x = x / sum - 1;
257 y = y / sum - 1;
258
259 return 0;
260}
261
262/// Find the maximum in an image at sub-pixel resolution by interpolation
263/** Uses imageMagnify() to expand the image to the desired scale. Because of the
264 * scaling used in imageMagnify, the desired scale may not be exact. As a result
265 * the actual scale is returned in scale_x and scale_y.
266 *
267 */
268template <typename floatT, typename imageT, typename magImageT, typename transformT>
269int imageMaxInterp( floatT &x, ///< [out] the x-position of the maximum, in pixels of the input image
270 floatT &y, ///< [out] the y-position of the maximum, in pixels of the input image
271 floatT &scale_x, ///< [in.out] the desired scale or resolution, in pixels < 1, in the x direction.
272 ///< On output contains the actual scale calculated.
273 floatT &scale_y, ///< [in.out] the desired scale or resolution, in pixels < 1, in the y direction.
274 ///< On output contains the actual scale calculated.
275 magImageT &magIm, ///< [in] the magnified image. This is used as working memory, will be resized.
276 const imageT &im, ///< [in] the image to find the maximum of
277 transformT trans ///< [in] the transform to use for interpolation
278)
279{
280 floatT magSize_x = ceil( ( im.rows() - 1.0 ) / scale_x ) + 1;
281 floatT magSize_y = ceil( ( im.cols() - 1.0 ) / scale_y ) + 1;
282
283 floatT mag_x = ( (floatT)magSize_x - 1.0 ) / ( (floatT)im.rows() - 1.0 );
284 floatT mag_y = ( (floatT)magSize_y - 1.0 ) / ( (floatT)im.cols() - 1.0 );
285
286 scale_x = 1.0 / mag_x;
287 scale_y = 1.0 / mag_y;
288
289 magIm.resize( magSize_x, magSize_y );
290
291 imageMagnify( magIm, im, trans );
292
293 int ix, iy;
294 magIm.maxCoeff( &ix, &iy );
295 x = ix * scale_x;
296 y = iy * scale_y;
297
298 return 0;
299}
300
301/// Find the maximum in an image at sub-pixel resolution by cubic convolution interpolation
302/** Uses imageMagnify() to expand the image to the desired scale. Because of the
303 * scaling used in imageMagnify, the desired scale may not be exact. As a result
304 * the actual scale is returned in scale_x and scale_y.
305 *
306 */
307template <typename floatT, typename imageT, typename magImageT>
308int imageMaxInterp( floatT &x, ///< [out] the x-position of the maximum, in pixels of the input image
309 floatT &y, ///< [out] the y-position of the maximum, in pixels of the input image
310 floatT &scale_x, ///< [in.out] the desired scale or resolution, in pixels < 1, in the x direction.
311 ///< On output contains the actual scale calculated.
312 floatT &scale_y, ///< [in.out] the desired scale or resolution, in pixels < 1, in the y direction.
313 ///< On output contains the actual scale calculated.
314 magImageT &magIm, ///< [in] the magnified image. This is used as working memory, will be resized.
315 const imageT &im ///< [in] the image to find the maximum of
316)
317{
318 return imageMaxInterp( x, y, scale_x, scale_y, magIm, im, cubicConvolTransform<typename imageT::Scalar>() );
319}
320
321/// Combine two images, each with their own mask defining good pixels.
322/** The combined image is made up of the pixels in im1 where mask1 is 1, and the pixels of im2 where mask2 is 1.
323 * If a pixel in both mask1 and mask2 has a value of 1, that pixel in the combo is the average of im1 and im2.
324 * All other pixels are set to 0 in the combined image.
325 *
326 * Separate template types are used for each argument to allow references, etc.
327 *
328 * \tparam imageT the eigen-like array type of the combined image
329 * \tparam imageT1 the eigen-like array type of image 1
330 * \tparam imageT2 the eigen-like array type of mask 1
331 * \tparam imageT3 the eigen-like array type of image 2
332 * \tparam imageT4 the eigen-like array type of mask 2
333 */
334template <typename imageT, typename imageT1, typename imageT2, typename imageT3, typename imageT4>
335void combine2ImagesMasked( imageT &combo, ///< [out] the combined image. will be resized.
336 const imageT1 &im1, ///< [in] the first image
337 const imageT2 &mask1, ///< [in] the mask for the first image
338 const imageT3 &im2, ///< [in] the second image
339 const imageT4 &mask2 ///< [in] the mask for the second image
340)
341{
342 combo.resize( im1.rows(), im2.cols() );
343
344 for( int c = 0; c < combo.cols(); ++c )
345 {
346 for( int r = 0; r < combo.rows(); ++r )
347 {
348 if( mask1( r, c ) == 1 && mask2( r, c ) == 0 )
349 combo( r, c ) = im1( r, c );
350 else if( mask2( r, c ) == 1 & mask1( r, c ) == 0 )
351 combo( r, c ) = im2( r, c );
352 else if( mask1( r, c ) == 1 && mask2( r, c ) == 1 )
353 combo( r, c ) = 0.5 * ( im1( r, c ) + im2( r, c ) );
354 else
355 combo( r, c ) = 0;
356 }
357 }
358}
359
360template <typename eigenT, typename eigenTin>
361void removeRowsAndCols( eigenT &out, const eigenTin &in, int st, int w )
362{
363
364 out.resize( in.rows() - w, in.cols() - w );
365
366 out.topLeftCorner( st, st ) = in.topLeftCorner( st, st );
367
368 out.bottomLeftCorner( in.rows() - ( st + w ), st ) = in.bottomLeftCorner( in.rows() - ( st + w ), st );
369
370 out.topRightCorner( st, in.cols() - ( st + w ) ) = in.topRightCorner( st, in.cols() - ( st + w ) );
371
372 out.bottomRightCorner( in.rows() - ( st + w ), in.cols() - ( st + w ) ) =
373 in.bottomRightCorner( in.rows() - ( st + w ), in.cols() - ( st + w ) );
374}
375
376template <typename eigenT, typename eigenTin>
377void removeRows( eigenT &out, const eigenTin &in, int st, int w )
378{
379
380 out.resize( in.rows() - w, in.cols() );
381
382 out.topLeftCorner( st, in.cols() ) = in.topLeftCorner( st, in.cols() );
383
384 out.bottomLeftCorner( in.rows() - ( st + w ), in.cols() ) =
385 in.bottomLeftCorner( in.rows() - ( st + w ), in.cols() );
386}
387
388template <typename eigenT, typename eigenTin>
389void removeCols( eigenT &out, const eigenTin &in, int st, int w )
390{
391
392 out.resize( in.rows(), in.cols() - w );
393
394 out.topLeftCorner( in.rows(), st ) = in.topLeftCorner( in.rows(), st );
395
396 out.topRightCorner( in.rows(), in.cols() - ( st + w ) ) = in.topRightCorner( in.rows(), in.cols() - ( st + w ) );
397}
398
399/// Copy one image to another, with no transformation
400/** This is merely memcpy
401 *
402 * \returns dest
403 */
404void *imcpy( void *dest, ///< [out] the address of the first pixel in the destination image
405 void *src, ///< [in] the address of the first pixel in the source image
406 size_t width, ///< [in] the width in pixels of size szof
407 size_t height, ///< [in] the height in pixels of size szof
408 size_t szof ///< [in] the size in bytes of a one pixel
409);
410
411/// Copy one image to another, flipping up-down
412/** This is a reversed row-by-row memcpy
413 *
414 * \returns dest
415 */
416void *imcpy_flipUD( void *dest, ///< [out] the address of the first pixel in the destination image
417 void *src, ///< [in] the address of the first pixel in the source image
418 size_t width, ///< [in] the width in pixels of size szof
419 size_t height, ///< [in] the height in pixels of size szof
420 size_t szof ///< [in] the size in bytes of a one pixel
421);
422
423/// Copy one image to another, flipping left-right
424/**
425 *
426 * \returns dest
427 */
428void *imcpy_flipLR( void *dest, ///< [out] the address of the first pixel in the destination image
429 void *src, ///< [in] the address of the first pixel in the source image
430 size_t width, ///< [in] the width in pixels of size szof
431 size_t height, ///< [in] the height in pixels of size szof
432 size_t szof ///< [in] the size in bytes of a one pixel
433);
434
435/// Copy one image to another, flipping up-down and left-right
436/**
437 *
438 * \returns dest
439 */
440void *imcpy_flipUDLR( void *dest, ///< [out] the address of the first pixel in the destination image
441 void *src, ///< [in] the address of the first pixel in the source image
442 size_t width, ///< [in] the width in pixels of size szof
443 size_t height, ///< [in] the height in pixels of size szof
444 size_t szof ///< [in] the size in bytes of a one pixel
445);
446
447///@}
448
449} // namespace improc
450} // namespace mx
451
452#endif // improc_imageUtils_hpp
void imageMagnify(arrOutT &transim, const arrInT &im, transformT trans)
Magnify an image.
void * imcpy(void *dest, void *src, size_t width, size_t height, size_t szof)
Copy one image to another, with no transformation.
int imageMaxInterp(floatT &x, floatT &y, floatT &scale_x, floatT &scale_y, magImageT &magIm, const imageT &im, transformT trans)
Find the maximum in an image at sub-pixel resolution by interpolation.
calcT imageVariance(imageT &im, calcT mn)
Calculate the variance of an image given its mean.
int imageCenterOfLight(typename imageT::Scalar &x, typename imageT::Scalar &y, const imageT &im)
Calculate the center of light of an image.
void combine2ImagesMasked(imageT &combo, const imageT1 &im1, const imageT2 &mask1, const imageT3 &im2, const imageT4 &mask2)
Combine two images, each with their own mask defining good pixels.
bool IsNan(float value)
Check if the number is nan, using several different methods.
calcT imageMean(imageT &im)
Calculate the mean value of an image.
void zeroNaNs(imageT &im, valueT val)
Zero any NaNs in an image.
void * imcpy_flipLR(void *dest, void *src, size_t width, size_t height, size_t szof)
Copy one image to another, flipping left-right.
void zeroNaNCube(cubeT &imc, maskCubeT *mask)
Zero any NaNs in an image cube.
void * imcpy_flipUDLR(void *dest, void *src, size_t width, size_t height, size_t szof)
Copy one image to another, flipping up-down and left-right.
void * imcpy_flipUD(void *dest, void *src, size_t width, size_t height, size_t szof)
Copy one image to another, flipping up-down.
int reflectImageCoords(int &x1, int &y1, int x0, int y0, realT xc, realT yc)
Reflect pixel coordinates across the given center pixel.
Image interpolation and transformation.
realT imCen(int rows_cols)
Get the mxlib standard center coordinate of an image.
realT imCenX(const imT &im)
Get the mxlib standard center x-coordinate of an image.
realT imCenY(const imT &im)
Get the mxlib standard center y-coordinate of an image.
The mxlib c++ namespace.
Definition mxError.hpp:106
Transformation by cubic convolution interpolation.