18 #ifndef NDARRAY_UTILS_H
19 #define NDARRAY_UTILS_H
25 #include <Eigen/Dense>
28 #include <unordered_set>
48 std::unordered_set<int64_t>
getLabels(ptr<const NDArray> in);
87 bool comparable(
const NDArray* left,
const NDArray* right,
88 bool* elL = NULL,
bool* elR = NULL);
107 ptr<NDArray>
derivative(ptr<const NDArray> in);
117 ptr<NDArray>
derivative(ptr<const NDArray> in,
size_t dir);
127 ptr<NDArray>
dilate(ptr<NDArray> in,
size_t reps);
137 ptr<NDArray>
erode(ptr<NDArray> in,
size_t reps);
176 void shiftImageFFT(ptr<NDArray> inout,
size_t dim,
double dist,
209 void shearImageFFT(ptr<NDArray> inout,
size_t dim,
size_t len,
double* dist,
231 double Rx,
double Ry,
double Rz,
double stepx=1,
double stepy=1,
double stepz=1);
247 int shearTest(
double Rx,
double Ry,
double Rz,
double sx=0,
double sy=0,
double sz=0);
260 ptr<NDArray>
linearRotate(
double rx,
double ry,
double rz,
261 ptr<const NDArray> in);
275 void rotateImageKern(ptr<NDArray> inout,
double rx,
double ry,
double rz);
318 ptr<NDArray>
pseudoPolar(ptr<const NDArray> in,
size_t prdim);
330 std::vector<ptr<NDArray>>
pseudoPolar(ptr<const NDArray> in);
356 void fillCircle(ptr<NDArray> inout,
double radius,
double alpha);
393 ptr<NDArray>
concat(
const vector<ptr<NDArray>>& images,
size_t dir);
405 ptr<NDArray>
concatElevate(
const vector<ptr<NDArray>>& images);
417 ptr<NDArray>
sobelEdge(ptr<const NDArray> img);
430 ptr<NDArray>
collapseSum(ptr<const NDArray> img,
size_t dim,
bool doabs=
false);
469 ptr<NDArray>
binarize(ptr<const NDArray> in,
double t);
479 ptr<NDArray>
threshold(ptr<const NDArray> in,
double t);
515 double mse(ptr<const NDArray> a, ptr<const NDArray> b,
516 ptr<const NDArray> mask = NULL);
528 double dice(ptr<const NDArray> a, ptr<const NDArray> b, ptr<const NDArray> mask);
540 double corr(ptr<const NDArray> a, ptr<const NDArray> b,
541 ptr<const NDArray> mask = NULL);
557 double information(ptr<const NDArray> a, ptr<const NDArray> b,
559 ptr<const NDArray> mask = NULL);
596 ptr<NDArray>
varianceT(ptr<const NDArray> img);
609 ptr<NDArray>
createRandLabels(ptr<const NDArray> example,
size_t nlabels,
double sd);
622 #endif //NDARRAY_UTILS_H
void shearImageKern(ptr< NDArray > inout, size_t dim, size_t len, double *dist, double(*kern)(double, double)=npl::lanczosKern)
Performs a shear on the image where the sheared dimension (dim) will be shifted depending on the inde...
void normalizeTS(ptr< NDArray > inout)
Normalize timeseries' so that their mean is 1 and standard deviation is 0.
ptr< NDArray > standardize(ptr< const NDArray > img)
Standardizes image distribution (makes it mean = 0, variance = 1). This computation is done in place ...
ptr< NDArray > collapseSum(ptr< const NDArray > img, size_t dim, bool doabs=false)
Creates a new image with the specified dimension collapsed and the values in each output point set to...
int shearTest(double Rx, double Ry, double Rz, double sx=0, double sy=0, double sz=0)
Tests shear results. If there is a solution (error is not NAN), then these should result small errors...
ptr< NDArray > threshold(ptr< const NDArray > in, double t)
Thresholds the image, changing everything below t to 0.
ptr< NDArray > relabelConnected(ptr< NDArray > input)
Performs relabeling based on connected component using the two pass algorithm.
ptr< NDArray > medianFilter(ptr< const NDArray > img)
median filter
std::unordered_set< int64_t > getLabels(ptr< const NDArray > in)
Produces a std::unordered_set of labels within a labelmap.
int rotateImageShearFFT(ptr< MRImage > inout, double rx, double ry, double rz, double(*window)(double, double)=npl::sincWindow)
Rotates an image around the center using shear decomposition followed by FFT-based shearing...
double dice(ptr< const NDArray > a, ptr< const NDArray > b, ptr< const NDArray > mask)
Computes the dice coefficent between two labelmap images. They should be identically gridded...
ptr< NDArray > binarize(ptr< const NDArray > in, double t)
Binarize an image to 0 or 1 based on whether it is greater than or less than the given threshold (t) ...
ptr< NDArray > dilate(ptr< NDArray > in, size_t reps)
Dilate an binary array repeatedly.
ptr< NDArray > pseudoPolarZoom(ptr< const NDArray > inimg, size_t prdim)
Computes the pseudopolar-gridded fourier transform on the input image, with prdim as the pseudo-radiu...
void standardizeIP(ptr< NDArray > img)
Standardizes image distribution (makes it mean = 0, variance = 1). This computation is done in place ...
int rotateImageShearKern(ptr< MRImage > inout, double rx, double ry, double rz, double(*kern)(double, double)=npl::lanczosKern)
Rotates an image around the center using shear decomposition followed by kernel-based shearing...
double sincWindow(double x, double a)
Sinc function centered at 0, with radius a, range should be = 2a.
double corr(ptr< const NDArray > a, ptr< const NDArray > b, ptr< const NDArray > mask=NULL)
Computes the correlation between two images. They should be identically gridded.
double mse(ptr< const NDArray > a, ptr< const NDArray > b, ptr< const NDArray > mask=NULL)
Computes the mean-squared-error between two images. They should be identically gridded.
ptr< NDArray > concat(const vector< ptr< NDArray >> &images, size_t dir)
Concatinates image in the direction specified by dir. So if dir is 0, and two images, sized [32, 32, 34] and [12, 32, 34] were passed in the input vector, then the output would be [44, 32, 34].
void gaussianSmooth1D(ptr< MRImage > inout, size_t dim, double stddev)
Gaussian smooths an image in 1 direction. Questionable whether it works. Seems to shift image...
ptr< NDArray > derivative(ptr< const NDArray > in)
Computes the derivative of the image. Computes all directional derivatives of the input image and the...
ptr< NDArray > pseudoPolar(ptr< const NDArray > in, size_t prdim)
Computes the pseudopolar-gridded fourier transform on the input image, with prdim as the pseudo-radiu...
double otsuThresh(ptr< const NDArray > in)
Computes a threshold based on OTSU.
ptr< NDArray > createRandLabels(ptr< const NDArray > example, size_t nlabels, double sd)
Create a random label map with the same size/orientation as the example image, and with nlabels uniqu...
ptr< NDArray > varianceT(ptr< const NDArray > img)
Takes the variance of the higher dimensions and returns a 3D image.
double lanczosKern(double x, double a)
Derivative of lanczos kernel with respect to x.
void rotateImageKern(ptr< NDArray > inout, double rx, double ry, double rz)
Performs a rotation using 3D intperolation kernel (lanczos)
void fillCircle(ptr< NDArray > inout, double radius, double alpha)
Sets the middle of the image += radius (in index space) to 1, everything else to 0.
void fillLinear(ptr< NDArray > inout)
Fills image with the linear index at each pixel.
void shearImageFFT(ptr< NDArray > inout, size_t dim, size_t len, double *dist, double(*window)(double, double)=npl::sincWindow)
Performs a shear on the image where the sheared dimension (dim) will be shifted depending on the inde...
void shiftImageKern(ptr< NDArray > inout, size_t dd, double dist, double(*kern)(double, double)=npl::lanczosKern)
Performs unidirectional shift in the direction of +dd, of distance (in units of pixels). Uses Lanczos interpolation.
int shearDecompose(std::list< Eigen::Matrix3d > &bestshears, double Rx, double Ry, double Rz, double stepx=1, double stepy=1, double stepz=1)
Decomposes a euler angle rotation using the rotation matrix made up of R = Rx*Ry*Rz. Note that this would be multiplying the input vector by Rz then Ry, then Rx. This does not support angles > PI/4. To do that, you should first do bulk rotation using 90 degree rotations (which requires not interpolation).
double information(ptr< const NDArray > a, ptr< const NDArray > b, int bins=100, int krad=4, Metric m=METRIC_MI, ptr< const NDArray > mask=NULL)
Computes the information based difference/similarity between two images. They should be identically g...
ptr< NDArray > histEqualize(ptr< const NDArray > in)
Computes a histogram, then spaces out intensities so that each intensity has equal volume/area in the...
void fillZero(ptr< NDArray > inout)
Fills image with the linear index at each pixel.
bool comparable(const NDArray *left, const NDArray *right, bool *elL=NULL, bool *elR=NULL)
Returns whether two NDArrays have the same dimensions, and therefore can be element-by-element compar...
void fillGaussian(ptr< NDArray > inout)
Fills image with the linear index at each pixel.
ptr< NDArray > concatElevate(const vector< ptr< NDArray >> &images)
Concatinates images/arrays. 1 Extra dimension will be added, all the lower dimensions of the images m...
ptr< NDArray > sobelEdge(ptr< const NDArray > img)
Increases the number of dimensions by 1 then places the edges in each dimension at indexes matching t...
void thresholdIP(ptr< NDArray > in, double t)
Thresholds the image, changing everything below t to 0.
ptr< NDArray > linearRotate(double rx, double ry, double rz, ptr< const NDArray > in)
Performs a rotation of the image first by rotating around z, then around y, then around x...
void shiftImageFFT(ptr< NDArray > inout, size_t dim, double dist, double(*window)(double, double)=npl::sincWindow)
Performs unidirectional shift in the direction of +dd, of distance (in units of pixels), using FFT.
void binarizeIP(ptr< NDArray > in, double t)
Binarize an image to 0 or 1 based on whether it is greater than or less than the given threshold (t)...
ptr< NDArray > erode(ptr< NDArray > in, size_t reps)
Erode an binary array repeatedly.
Metric
Information-based Metric to use.