forked from SCIInstitute/ShapeWorks
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathImage.cpp
1001 lines (800 loc) · 31.1 KB
/
Image.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include "Image.h"
#include <itkAntiAliasBinaryImageFilter.h>
#include <itkBinaryFillholeImageFilter.h>
#include <itkBinaryThresholdImageFilter.h>
#include <itkChangeInformationImageFilter.h>
#include <itkConnectedComponentImageFilter.h>
#include <itkConstantPadImageFilter.h>
#include <itkCurvatureFlowImageFilter.h>
#include <itkDiscreteGaussianImageFilter.h>
#include <itkExtractImageFilter.h>
#include <itkGDCMImageIO.h>
#include <itkGDCMSeriesFileNames.h>
#include <itkGradientMagnitudeImageFilter.h>
#include <itkImageDuplicator.h>
#include <itkImageFileReader.h>
#include <itkImageFileWriter.h>
#include <itkImageSeriesReader.h>
#include <itkImageToVTKImageFilter.h>
#include <itkIntensityWindowingImageFilter.h>
#include <itkMultiplyImageFilter.h>
#include <itkNearestNeighborInterpolateImageFunction.h>
#include <itkOrientImageFilter.h>
#include <itkRegionOfInterestImageFilter.h>
#include <itkReinitializeLevelSetImageFilter.h>
#include <itkRelabelComponentImageFilter.h>
#include <itkResampleImageFilter.h>
#include <itkScalableAffineTransform.h>
#include <itkSigmoidImageFilter.h>
#include <itkTestingComparisonImageFilter.h>
#include <itkThresholdImageFilter.h>
#include <itkVTKImageExport.h>
#include <itkVTKImageToImageFilter.h>
#include <vtkContourFilter.h>
#include <vtkImageCast.h>
#include <vtkImageData.h>
#include <vtkImageImport.h>
#include <cmath>
#include <exception>
#include "Exception.h"
#include "MeshUtils.h"
#include "ShapeworksUtils.h"
#include "itkTPGACLevelSetImageFilter.h" // actually a shapeworks class, not itk
namespace shapeworks {
Image::Image(const Dims dims) : itk_image_(ImageType::New()) {
ImageType::RegionType region;
region.SetSize(dims);
region.SetIndex(Coord({0, 0, 0}));
itk_image_->SetRegions(region);
itk_image_->Allocate();
}
Image::Image(const vtkSmartPointer<vtkImageData> vtkImage) {
// ensure input image data is PixelType (note: it'll either be float or double)
auto cast = vtkSmartPointer<vtkImageCast>::New();
cast->SetInputData(vtkImage);
if (typeid(PixelType) == typeid(float)) {
cast->SetOutputScalarTypeToFloat();
}
cast->Update();
using FilterType = itk::VTKImageToImageFilter<Image::ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(cast->GetOutput());
filter->Update();
this->itk_image_ = cloneData(filter->GetOutput());
}
vtkSmartPointer<vtkImageData> Image::getVTKImage() const {
using connectorType = itk::ImageToVTKImageFilter<Image::ImageType>;
connectorType::Pointer connector = connectorType::New();
connector->SetInput(this->itk_image_);
connector->Update();
return connector->GetOutput();
}
Image::ImageType::Pointer Image::cloneData(const Image::ImageType::Pointer image) {
using DuplicatorType = itk::ImageDuplicator<ImageType>;
DuplicatorType::Pointer duplicator = DuplicatorType::New();
duplicator->SetInputImage(image);
duplicator->Update();
return duplicator->GetOutput();
}
Image& Image::operator=(const Image& img) {
this->itk_image_ = Image::cloneData(img.itk_image_);
return *this;
}
Image& Image::operator=(Image&& img) {
this->itk_image_ = nullptr; // make sure to free existing image by setting it to nullptr (works b/c it's a smart ptr)
this->itk_image_.Swap(img.itk_image_);
return *this;
}
Image::ImageType::Pointer Image::read(const std::string& pathname) {
if (pathname.empty()) {
throw std::invalid_argument("Empty pathname");
}
if (ShapeworksUtils::is_directory(pathname)) return readDICOMImage(pathname);
using ReaderType = itk::ImageFileReader<ImageType>;
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(pathname);
try {
reader->Update();
} catch (itk::ExceptionObject& exp) {
throw shapeworks_exception(std::string(exp.what()));
}
// reorient the image to RAI if it's not already
ImageType::Pointer img = reader->GetOutput();
if (itk::SpatialOrientationAdapter().FromDirectionCosines(img->GetDirection()) !=
itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI) {
using Orienter = itk::OrientImageFilter<ImageType, ImageType>;
Orienter::Pointer orienter = Orienter::New();
orienter->UseImageDirectionOn();
// set orientation to RAI
orienter->SetDesiredCoordinateOrientation(itk::SpatialOrientation::ITK_COORDINATE_ORIENTATION_RAI);
orienter->SetInput(img);
orienter->Update();
img = orienter->GetOutput();
}
return img;
}
Image& Image::operator-() {
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
iter.Set(-iter.Value());
++iter;
}
return *this;
}
Image Image::operator+(const Image& other) const {
Image ret(*this);
ret += other;
return ret;
}
Image& Image::operator+=(const Image& other) {
if (dims() != other.dims()) {
throw std::invalid_argument("images must have same logical dims");
}
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
itk::ImageRegionIteratorWithIndex<ImageType> otherIter(other.itk_image_,
other.itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd() && !otherIter.IsAtEnd()) {
iter.Set(iter.Value() + otherIter.Value());
++iter;
++otherIter;
}
return *this;
}
Image Image::operator-(const Image& other) const {
Image ret(*this);
ret -= other;
return ret;
}
Image& Image::operator-=(const Image& other) {
if (dims() != other.dims()) {
throw std::invalid_argument("images must have same logical dims");
}
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
itk::ImageRegionIteratorWithIndex<ImageType> otherIter(other.itk_image_,
other.itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd() && !otherIter.IsAtEnd()) {
iter.Set(iter.Value() - otherIter.Value());
++iter;
++otherIter;
}
return *this;
}
Image Image::operator+(const PixelType x) const {
Image ret(*this);
ret += x;
return ret;
}
Image& Image::operator+=(const PixelType x) {
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
iter.Set(iter.Value() + x);
++iter;
}
return *this;
}
Image Image::operator-(const PixelType x) const {
Image ret(*this);
ret -= x;
return ret;
}
Image& Image::operator-=(const PixelType x) {
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
iter.Set(iter.Value() - x);
++iter;
}
return *this;
}
Image Image::operator*(const Image& other) const {
using FilterType = itk::MultiplyImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput1(this->itk_image_);
filter->SetInput2(other.itk_image_);
filter->Update();
return Image(filter->GetOutput());
}
Image Image::operator*(const PixelType x) const {
Image ret(*this);
ret *= x;
return ret;
}
Image& Image::operator*=(const PixelType x) {
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
iter.Set(iter.Value() * x);
++iter;
}
return *this;
}
Image Image::operator/(const PixelType x) const {
Image ret(*this);
ret /= x;
return ret;
}
Image& Image::operator/=(const PixelType x) {
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
iter.Set(iter.Value() / x);
++iter;
}
return *this;
}
template <>
Image operator*(const Image& img, const double x) {
return img.operator*(x);
}
template <>
Image operator/(const Image& img, const double x) {
return img.operator/(x);
}
template <>
Image& operator*=(Image& img, const double x) {
return img.operator*=(x);
}
template <>
Image& operator/=(Image& img, const double x) {
return img.operator/=(x);
}
Image::ImageType::Pointer Image::readDICOMImage(const std::string& pathname) {
if (pathname.empty()) {
throw std::invalid_argument("Empty pathname");
}
using ReaderType = itk::ImageSeriesReader<ImageType>;
using ImageIOType = itk::GDCMImageIO;
using InputNamesGeneratorType = itk::GDCMSeriesFileNames;
ImageIOType::Pointer gdcm_io = ImageIOType::New();
InputNamesGeneratorType::Pointer input_names = InputNamesGeneratorType::New();
input_names->SetInputDirectory(pathname);
const ReaderType::FileNamesContainer& filenames = input_names->GetInputFileNames();
ReaderType::Pointer reader = ReaderType::New();
reader->SetImageIO(gdcm_io);
reader->SetFileNames(filenames);
try {
reader->Update();
} catch (itk::ExceptionObject& exp) {
throw std::invalid_argument("Failed to read DICOM from " + pathname + "(" + std::string(exp.what()) + ")");
}
return reader->GetOutput();
}
Image& Image::write(const std::string& filename, bool compressed) {
if (!this->itk_image_) {
throw std::invalid_argument("Image invalid");
}
if (filename.empty()) {
throw std::invalid_argument("Empty pathname");
}
using WriterType = itk::ImageFileWriter<ImageType>;
WriterType::Pointer writer = WriterType::New();
writer->SetInput(this->itk_image_);
writer->SetFileName(filename);
writer->SetUseCompression(compressed);
writer->Update();
return *this;
}
Image& Image::antialias(unsigned iterations, double maxRMSErr, int layers) {
if (layers < 0) {
throw std::invalid_argument("layers must be >= 0");
}
using FilterType = itk::AntiAliasBinaryImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetMaximumRMSError(maxRMSErr);
filter->SetNumberOfIterations(iterations);
if (layers) {
filter->SetNumberOfLayers(layers);
}
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::recenter() { return setOrigin(negate(size() / 2.0)); }
Image& Image::resample(const TransformPtr transform, const Point3 origin, Dims dims, const Vector3 spacing,
const ImageType::DirectionType direction, Image::InterpolationType interp) {
using FilterType = itk::ResampleImageFilter<ImageType, ImageType>;
FilterType::Pointer resampler = FilterType::New();
switch (interp) {
case Linear:
// linear interpolation is the default
break;
case NearestNeighbor: {
using InterpolatorType = itk::NearestNeighborInterpolateImageFunction<ImageType, double>;
InterpolatorType::Pointer interpolator = InterpolatorType::New();
resampler->SetInterpolator(interpolator);
break;
}
default:
throw std::invalid_argument("Unknown Image::InterpolationType");
}
resampler->SetInput(this->itk_image_);
resampler->SetTransform(transform ? transform : IdentityTransform::New());
resampler->SetOutputOrigin(origin);
resampler->SetOutputSpacing(spacing);
resampler->SetSize(dims);
resampler->SetOutputDirection(direction);
resampler->Update();
this->itk_image_ = resampler->GetOutput();
return *this;
}
Image& Image::resample(const Vector3& spacing, Image::InterpolationType interp) {
// compute logical dimensions that keep all image data for this spacing
Dims inputDims(this->dims());
Vector3 inputSpacing(this->spacing());
Dims dims({static_cast<unsigned>(std::floor(inputDims[0] * inputSpacing[0] / spacing[0])),
static_cast<unsigned>(std::floor(inputDims[1] * inputSpacing[1] / spacing[1])),
static_cast<unsigned>(std::floor(inputDims[2] * inputSpacing[2] / spacing[2]))});
Point3 new_origin = origin() + toPoint(0.5 * (spacing - inputSpacing)); // O' += 0.5 * (p' - p)
return resample(IdentityTransform::New(), new_origin, dims, spacing, coordsys(), interp);
}
Image& Image::resample(double isoSpacing, Image::InterpolationType interp) {
return resample(makeVector({isoSpacing, isoSpacing, isoSpacing}), interp);
}
Image& Image::resize(Dims dims, Image::InterpolationType interp) {
// use existing dims for any that are unspecified
Dims inputDims(this->dims());
if (dims[0] == 0) dims[0] = inputDims[0];
if (dims[1] == 0) dims[1] = inputDims[1];
if (dims[2] == 0) dims[2] = inputDims[2];
// compute new spacing so physical image size remains the same
Vector3 inputSpacing(spacing());
Vector3 spacing(makeVector({inputSpacing[0] * inputDims[0] / dims[0], inputSpacing[1] * inputDims[1] / dims[1],
inputSpacing[2] * inputDims[2] / dims[2]}));
return resample(IdentityTransform::New(), origin(), dims, spacing, coordsys(), interp);
}
bool Image::compare(const Image& other, bool verifyall, double tolerance, double precision) const {
if (tolerance > 1 || tolerance < 0) {
throw std::invalid_argument("tolerance value must be between 0 and 1 (inclusive)");
}
// we use the region of interest filter here with the full region because our
// incoming image may be the output of an ExtractImageFilter or PadImageFilter
// which modify indices and leave the origin intact. These will not compare
// properly against a saved NRRD file because the act of saving the image to
// NRRD and back in will cause the origin (and indices) to be reset.
using RegionFilterType = itk::RegionOfInterestImageFilter<ImageType, ImageType>;
RegionFilterType::Pointer region_filter = RegionFilterType::New();
region_filter->SetInput(this->itk_image_);
region_filter->SetRegionOfInterest(this->itk_image_->GetLargestPossibleRegion());
region_filter->UpdateLargestPossibleRegion();
ImageType::Pointer itk_image = region_filter->GetOutput();
// perform the same to the other image
RegionFilterType::Pointer region_filter2 = RegionFilterType::New();
region_filter2->SetInput(other.itk_image_);
region_filter2->SetRegionOfInterest(other.itk_image_->GetLargestPossibleRegion());
region_filter2->UpdateLargestPossibleRegion();
ImageType::Pointer other_itk_image = region_filter2->GetOutput();
using DiffType = itk::Testing::ComparisonImageFilter<ImageType, ImageType>;
DiffType::Pointer diff = DiffType::New();
diff->SetValidInput(other_itk_image);
diff->SetTestInput(itk_image);
diff->SetDifferenceThreshold(precision);
diff->SetToleranceRadius(0);
diff->SetVerifyInputInformation(verifyall);
try {
diff->UpdateLargestPossibleRegion();
} catch (itk::ExceptionObject& exp) {
return false;
}
auto numberOfPixelsWithDifferences = diff->GetNumberOfPixelsWithDifferences();
Dims dim = dims();
auto allowedPixelDifference = tolerance * dim[0] * dim[1] * dim[2];
if (numberOfPixelsWithDifferences > allowedPixelDifference) {
return false;
}
return true;
}
Image& Image::pad(int padding, PixelType value) { return this->pad(padding, padding, padding, value); }
Image& Image::pad(int padx, int pady, int padz, PixelType value) {
ImageType::SizeType lowerExtendRegion;
lowerExtendRegion[0] = padx;
lowerExtendRegion[1] = pady;
lowerExtendRegion[2] = padz;
ImageType::SizeType upperExtendRegion;
upperExtendRegion[0] = padx;
upperExtendRegion[1] = pady;
upperExtendRegion[2] = padz;
return this->pad(lowerExtendRegion, upperExtendRegion, value);
}
Image& Image::pad(IndexRegion& region, PixelType value) {
auto bbox = logicalBoundingBox();
// compute positive numbers to pad in each direction
ImageType::SizeType lowerExtendRegion;
lowerExtendRegion[0] = std::max(Coord::IndexValueType(0), -region.min[0]);
lowerExtendRegion[1] = std::max(Coord::IndexValueType(0), -region.min[1]);
lowerExtendRegion[2] = std::max(Coord::IndexValueType(0), -region.min[2]);
ImageType::SizeType upperExtendRegion;
upperExtendRegion[0] = std::max(Coord::IndexValueType(0), region.max[0] - bbox.max[0]);
upperExtendRegion[1] = std::max(Coord::IndexValueType(0), region.max[1] - bbox.max[1]);
upperExtendRegion[2] = std::max(Coord::IndexValueType(0), region.max[2] - bbox.max[2]);
return this->pad(lowerExtendRegion, upperExtendRegion, value);
}
Image& Image::pad(Dims lowerExtendRegion, Dims upperExtendRegion, PixelType value) {
using FilterType = itk::ConstantPadImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetPadLowerBound(lowerExtendRegion);
filter->SetPadUpperBound(upperExtendRegion);
filter->SetConstant(value);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::translate(const Vector3& v) {
AffineTransformPtr xform(AffineTransform::New());
xform->Translate(-v); // negate v because ITK applies transformations backwards.
return applyTransform(xform);
}
Image& Image::scale(const Vector3& s) {
if (s[0] == 0 || s[1] == 0 || s[2] == 0) {
throw std::invalid_argument("Invalid scale point");
}
auto origOrigin(origin()); // scale centered at origin, so temporarily set origin to be the center
recenter();
AffineTransformPtr xform(AffineTransform::New());
xform->Scale(invertValue(Vector(s))); // invert scale ratio because ITK applies transformations backwards.
applyTransform(xform);
setOrigin(origOrigin); // restore origin
return *this;
}
Image& Image::rotate(const double angle, Axis axis) {
switch (axis) {
case X:
return rotate(angle, makeVector({1.0, 0.0, 0.0}));
case Y:
return rotate(angle, makeVector({0.0, 1.0, 0.0}));
case Z:
return rotate(angle, makeVector({0.0, 0.0, 1.0}));
default:
throw std::invalid_argument("Unknown axis.");
}
}
Image& Image::rotate(const double angle, const Vector3& axis) {
if (!axis_is_valid(axis)) {
throw std::invalid_argument("Invalid axis");
}
auto origOrigin(origin()); // rotation is around origin, so temporarily set origin to be the center
recenter();
AffineTransformPtr xform(AffineTransform::New());
xform->Rotate3D(axis, -angle); // negate angle because ITK applies transformations backwards.
applyTransform(xform);
setOrigin(origOrigin); // restore origin
return *this;
}
Image& Image::applyTransform(const TransformPtr transform, Image::InterpolationType interp) {
return applyTransform(transform, origin(), dims(), spacing(), coordsys(), interp);
}
Image& Image::applyTransform(const TransformPtr transform, const Point3 origin, const Dims dims, const Vector3 spacing,
const ImageType::DirectionType coordsys, Image::InterpolationType interp) {
return resample(transform, origin, dims, spacing, coordsys, interp);
}
Image& Image::extractLabel(const PixelType label) {
binarize(label - std::numeric_limits<PixelType>::epsilon(), label + std::numeric_limits<PixelType>::epsilon());
return *this;
}
Image& Image::closeHoles(const PixelType foreground) {
using FilterType = itk::BinaryFillholeImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetForegroundValue(foreground + std::numeric_limits<PixelType>::epsilon());
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::binarize(PixelType minVal, PixelType maxVal, PixelType innerVal, PixelType outerVal) {
using FilterType = itk::BinaryThresholdImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetLowerThreshold(minVal + std::numeric_limits<PixelType>::epsilon());
filter->SetUpperThreshold(maxVal);
filter->SetInsideValue(innerVal);
filter->SetOutsideValue(outerVal);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::computeDT(PixelType isoValue) {
using FilterType = itk::ReinitializeLevelSetImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->NarrowBandingOff();
filter->SetLevelSetValue(isoValue);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::applyCurvatureFilter(unsigned iterations) {
if (iterations < 0) {
throw std::invalid_argument("iterations must be >= 0");
}
using FilterType = itk::CurvatureFlowImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetTimeStep(0.0625);
filter->SetNumberOfIterations(iterations);
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::applyGradientFilter() {
using FilterType = itk::GradientMagnitudeImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::applySigmoidFilter(double alpha, double beta) {
using FilterType = itk::SigmoidImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetAlpha(alpha);
filter->SetBeta(beta);
filter->SetOutputMinimum(0.0);
filter->SetOutputMaximum(1.0);
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::applyTPLevelSetFilter(const Image& featureImage, double scaling) {
if (!featureImage.itk_image_) {
throw std::invalid_argument("Invalid feature image");
}
using FilterType =
itk::TPGACLevelSetImageFilter<ImageType, ImageType>; // TODO: this is no longer part of ITK and should be updated
FilterType::Pointer filter = FilterType::New();
filter->SetPropagationScaling(scaling);
filter->SetCurvatureScaling(1.0);
filter->SetAdvectionScaling(1.0);
filter->SetMaximumRMSError(0.0);
filter->SetNumberOfIterations(20);
filter->SetInput(this->itk_image_);
filter->SetFeatureImage(featureImage.itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::topologyPreservingSmooth(float scaling, float sigmoidAlpha, float sigmoidBeta) {
Image featureImage(*this);
featureImage.applyGradientFilter().applySigmoidFilter(sigmoidAlpha, sigmoidBeta);
return applyTPLevelSetFilter(featureImage, scaling);
}
Image& Image::applyIntensityFilter(double minVal, double maxVal) {
using FilterType = itk::IntensityWindowingImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetWindowMinimum(minVal);
filter->SetWindowMaximum(maxVal);
filter->SetOutputMinimum(0.0);
filter->SetOutputMaximum(255.0);
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::gaussianBlur(double sigma) {
using BlurType = itk::DiscreteGaussianImageFilter<ImageType, ImageType>;
BlurType::Pointer blur = BlurType::New();
blur->SetInput(this->itk_image_);
blur->SetVariance(sigma * sigma);
blur->Update();
this->itk_image_ = blur->GetOutput();
return *this;
}
Image& Image::crop(PhysicalRegion region, const int padding) {
region.shrink(physicalBoundingBox()); // clip region to fit inside image
if (!region.valid()) {
throw std::invalid_argument("Invalid region specified (it may lie outside physical bounds of image).");
}
using FilterType = itk::RegionOfInterestImageFilter<ImageType, ImageType>;
FilterType::Pointer filter = FilterType::New();
IndexRegion indexRegion(physicalToLogical(region));
indexRegion.pad(padding);
filter->SetRegionOfInterest(ImageType::RegionType(indexRegion.min, indexRegion.size()));
filter->SetInput(this->itk_image_);
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::clip(const Plane plane, const PixelType val) {
if (!axis_is_valid(getNormal(plane))) {
throw std::invalid_argument("Invalid clipping plane (zero length normal)");
}
itk::ImageRegionIteratorWithIndex<ImageType> iter(this->itk_image_, itk_image_->GetLargestPossibleRegion());
while (!iter.IsAtEnd()) {
Vector pq(logicalToPhysical(iter.GetIndex()) - getOrigin(plane));
// if n dot pq is < 0, point q is on the back side of the plane.
if (getNormal(plane) * pq < 0.0) iter.Set(val);
++iter;
}
return *this;
}
Image& Image::reflect(const Axis& axis) {
if (!axis_is_valid(axis)) {
throw std::invalid_argument("Invalid axis");
}
Vector scale(makeVector({1, 1, 1}));
scale[axis] = -1;
using ScalableTransform = itk::ScalableAffineTransform<double, 3>;
ScalableTransform::Pointer xform(ScalableTransform::New());
xform->SetScale(scale);
Point3 currentOrigin(origin());
recenter().applyTransform(xform).setOrigin(currentOrigin);
return *this;
}
Image& Image::setOrigin(Point3 origin) {
using FilterType = itk::ChangeInformationImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetOutputOrigin(origin);
filter->ChangeOriginOn();
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::setSpacing(Vector3 spacing) {
if (spacing[0] <= 0 || spacing[1] <= 0 || spacing[2] <= 0) {
throw std::invalid_argument("Spacing cannot b <= 0");
}
using FilterType = itk::ChangeInformationImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetOutputSpacing(spacing);
filter->ChangeSpacingOn();
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::setCoordsys(ImageType::DirectionType coordsys) {
using FilterType = itk::ChangeInformationImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->SetOutputDirection(coordsys);
filter->ChangeDirectionOn();
filter->Update();
this->itk_image_ = filter->GetOutput();
return *this;
}
Image& Image::isolate() {
typedef itk::Image<unsigned char, 3> IsolateType;
typedef itk::CastImageFilter<ImageType, IsolateType> ToIntType;
ToIntType::Pointer filter = ToIntType::New();
filter->SetInput(this->itk_image_);
filter->Update();
// Find the connected components in this image.
auto cc = itk::ConnectedComponentImageFilter<IsolateType, IsolateType>::New();
cc->SetInput(filter->GetOutput());
cc->FullyConnectedOn();
cc->Update();
auto relabel = itk::RelabelComponentImageFilter<IsolateType, IsolateType>::New();
relabel->SetInput(cc->GetOutput());
relabel->SortByObjectSizeOn();
relabel->Update();
auto thresh = itk::ThresholdImageFilter<IsolateType>::New();
thresh->SetInput(relabel->GetOutput());
thresh->SetOutsideValue(0);
thresh->ThresholdBelow(0);
thresh->ThresholdAbove(1);
thresh->Update();
auto cast = itk::CastImageFilter<IsolateType, ImageType>::New();
cast->SetInput(thresh->GetOutput());
cast->Update();
this->itk_image_ = cast->GetOutput();
return *this;
}
Point3 Image::centerOfMass(PixelType minVal, PixelType maxVal) const {
itk::ImageRegionIteratorWithIndex<ImageType> imageIt(this->itk_image_, itk_image_->GetLargestPossibleRegion());
int numPixels = 0;
Point3 com({0.0, 0.0, 0.0});
while (!imageIt.IsAtEnd()) {
PixelType val = imageIt.Get();
if (val > minVal && val <= maxVal) {
numPixels++;
com += itk_image_->TransformIndexToPhysicalPoint<double>(imageIt.GetIndex());
}
++imageIt;
}
if (numPixels > 0)
com /= static_cast<double>(numPixels);
else
com = center(); // an image with no mass still has a center
return com;
}
Image::StatsPtr Image::statsFilter() {
using FilterType = itk::StatisticsImageFilter<ImageType>;
FilterType::Pointer filter = FilterType::New();
filter->SetInput(this->itk_image_);
filter->UpdateLargestPossibleRegion();
return filter;
}
Image::PixelType Image::min() {
StatsPtr filter = statsFilter();
return filter->GetMinimum();
}
Image::PixelType Image::max() {
StatsPtr filter = statsFilter();
return filter->GetMaximum();
}
Image::PixelType Image::mean() {
StatsPtr filter = statsFilter();
return filter->GetMean();
}
Image::PixelType Image::std() {
StatsPtr filter = statsFilter();
return sqrt(filter->GetVariance());
}
IndexRegion Image::logicalBoundingBox() const {
IndexRegion region(Coord({0, 0, 0}), toCoord(dims() - Dims({1, 1, 1})));
return region;
}
PhysicalRegion Image::physicalBoundingBox() const {
PhysicalRegion region(origin(), origin() + dotProduct(toVector(dims() - Dims({1, 1, 1})), spacing()));
return region;
}
PhysicalRegion Image::physicalBoundingBox(PixelType isoValue) const {
PhysicalRegion bbox;
itk::ImageRegionIteratorWithIndex<ImageType> imageIterator(itk_image_, itk_image_->GetLargestPossibleRegion());
while (!imageIterator.IsAtEnd()) {
PixelType val = imageIterator.Get();
if (val >= isoValue) bbox.expand(logicalToPhysical(imageIterator.GetIndex()));
++imageIterator;
}
return bbox;
}
PhysicalRegion Image::logicalToPhysical(const IndexRegion region) const {
return PhysicalRegion(logicalToPhysical(region.min), logicalToPhysical(region.max));
}
IndexRegion Image::physicalToLogical(const PhysicalRegion region) const {
return IndexRegion(physicalToLogical(region.min), physicalToLogical(region.max));
}
Point3 Image::logicalToPhysical(const Coord& v) const {
Point3 value;
itk_image_->TransformIndexToPhysicalPoint(v, value);
return value;
}
Coord Image::physicalToLogical(const Point3& p) const { return itk_image_->TransformPhysicalPointToIndex(p); }
Image::ImageIterator Image::iterator() {
ImageIterator iter(this->itk_image_, itk_image_->GetRequestedRegion());
return iter;
}
Mesh Image::toMesh(PixelType isoValue) const {
auto vtkImage = getVTKImage();
auto targetContour = vtkSmartPointer<vtkContourFilter>::New();
targetContour->SetInputData(vtkImage);
targetContour->SetValue(0, isoValue);
targetContour->Update();
return Mesh(targetContour->GetOutput());
}
Image::PixelType Image::evaluate(Point p) {
if (!interpolator_) {
interpolator_ = InterpolatorType::New();
interpolator_->SetInputImage(itk_image_);
}
return interpolator_->Evaluate(p);
}
TransformPtr Image::createCenterOfMassTransform() {
AffineTransformPtr xform(AffineTransform::New());
xform->Translate(-(center() - centerOfMass())); // ITK translations go in a counterintuitive direction
return xform;
}
TransformPtr Image::createRigidRegistrationTransform(const Image& target_dt, float isoValue, unsigned iterations) {
if (!target_dt.itk_image_) {
throw std::invalid_argument("Invalid target. Expected distance transform image");
}
Mesh sourceContour = toMesh(isoValue);
Mesh targetContour = target_dt.toMesh(isoValue);
try {
auto mat = MeshUtils::createICPTransform(sourceContour, targetContour, Mesh::Rigid, iterations);
return shapeworks::createTransform(ShapeworksUtils::getMatrix(mat), ShapeworksUtils::getOffset(mat));
} catch (std::invalid_argument) {
std::cerr << "failed to create ICP transform.\n";
if (sourceContour.numPoints() == 0) {
std::cerr << "\tspecified isoValue (" << isoValue << ") results in an empty mesh for source\n";
}
if (targetContour.numPoints() == 0) {
std::cerr << "\tspecified isoValue (" << isoValue << ") results in an empty mesh for target\n";
}
}
return AffineTransform::New();
}
std::ostream& operator<<(std::ostream& os, const Image& img) {
return os << "{\n\tdims: " << img.dims() << ",\n\torigin: " << img.origin() << ",\n\tsize: " << img.size()
<< ",\n\tspacing: " << img.spacing() << "\n}";
}