[Insight-users] exclude pixel values (background) from metrics

Luis Ibanez luis.ibanez at kitware.com
Wed Mar 17 18:44:36 EDT 2010


Hi Darren,

Excellent,

In that case you should be able to eliminate the background
by using:

intensityThreshold = 1;
metric->SetFixedImageSamplesIntensityThreshold( intensityThreshold );

If you want to be 100% sure, it may be worth running a
BinaryThresholdImageFilter first and saving the output
for further visualization. That will provide a good sanity
check.


     Regards,


            Luis


------------------------------------------------------------------------------------------
On Wed, Mar 17, 2010 at 5:22 PM, Darren Weber
<darren.weber.lists at gmail.com> wrote:
>
> Hi Luis et al.,
>
> I've created the input images by exporting from Photoshop with a defined
> background value (RGB: 0 0 0).  It's these values that I want to exclude
> from the metric in the image registration method.
>
> The min color is confirmed by the identify program of ImageMagick:
>
>   Colorspace: RGB
>   Depth: 16-bit
>   Channel depth:
>     gray: 16-bit
>   Channel statistics:
>     Gray:
>       min: 0 (0)
>       max: 46335 (0.707027)
>       mean: 25664.6 (0.391617)
>       standard deviation: 15257.7 (0.232818)
>       kurtosis: -0.857968
>       skewness: -0.771191
>
>
> TIA,
> Darren
>
>
>
> On Tue, Mar 16, 2010 at 10:25 AM, Luis Ibanez <luis.ibanez at kitware.com>
> wrote:
>>
>> Hi Darren,
>>
>> 1) As described in the ITK Software Guide
>>
>>         http://www.itk.org/ItkSoftwareGuide.pdf
>>
>>   in the Chapter 7 "Reading and Writing Images",
>>   in section 7.1, pdf-page: 296:
>>
>>  ITK takes the pixel data from the input file and casts it to
>>  the PixelType that you used for instantiating your image
>>  at compilation time. The casting is performed using
>>  C-Language rules, and not truncation, overflow or underflow
>>  detection, nor scaling is applied.
>>
>>  So, if your TIFF file contains 16bit integers values per pixel
>>  and you instantiate the image using a "float" pixel type, then
>>  the values of the integers will be casted (without any rescaling)
>>  in to floats.
>>
>>
>> 2)  If you have a use for the original "16bits" data, then it may
>>     be worth reading the file in an image of the exact same type
>>     and then using a CastImageFilter to convert it to float.
>>
>>     On the other hand, if you only need the float image for
>>     performing the registration, then reading directly into a
>>     float image is actually saving you some memory.
>>
>> 3)   "Black" pixels is quite an ambiguous measure, since,
>>     depending on the image viewer that you are using, it may
>>     well be that what you "see" as black, is actually a level
>>     of intensity with values around 20,000 (for example).
>>
>>     What you really need is a sort of histogram that will allow
>>     you to identify the background of the image,
>>
>>     or
>>
>>     A visualization method that allows you to identify the
>>     intensity level of the background.
>>
>>     You may want to consider:
>>
>>       Insight/Code/Algorithms/
>>                        itkOtsuThresholdImageCalculator.h
>>
>>     and the Insight Journal papers:
>>
>>     "Kappa Sigma Clipping"
>>     http://www.insight-journal.org/browse/publication/132
>>     http://hdl.handle.net/1926/367
>>
>>     "Robust Automatic Threshold Selection"
>>     http://www.insight-journal.org/browse/publication/134
>>     http://hdl.handle.net/1926/370
>>
>>     as options for estimating a suitable threshold.
>>
>>
>>
>>
>>   Regards,
>>
>>
>>         Luis
>>
>>
>>
>> ----------------------------------------------------------------
>> On Mon, Mar 15, 2010 at 3:09 PM, Darren Weber
>> <darren.weber.lists at gmail.com> wrote:
>> >
>> > Hi Luis et al.,
>> >
>> > Now the ITK library is compiled and installed with the additional
>> > features
>> > available.  Now this code will compile:
>> >
>> > typedef itk::MeanSquaresImageToImageMetric<iBWImgType,iBWImgType>
>> > MetricType;
>> > MetricType::Pointer metric = MetricType::New();
>> > metric->SetUseFixedImageSamplesIntensityThreshold( true );
>> > metric->SetFixedImageSamplesIntensityThreshold( intensityThreshold );
>> > //metric->SetUseSequentialSampling( true );
>> > //metric->SetUseAllPixels( true );
>> >
>> >
>> > Now a question arises about the range of values for the
>> > "intensityThreshold".  The input pixel and image typedefs are:
>> >
>> > const unsigned short dimImg = 2;
>> > typedef float iPixType;
>> > typedef itk::Image< iPixType, dimImg > iBWImgType;
>> > typedef itk::ImageFileReader< iBWImgType > bwImageReaderType;
>> >
>> >
>> > The image reader is pulling in .tif images, which ImageMagick identifies
>> > as
>> > (see verbose output in attachment):
>> >
>> > $  identify testdata/bw/section0003_w1.tif
>> > testdata/bw/section0003_w1.tif TIFF 549x539 549x539+0+0 16-bit Grayscale
>> > DirectClass 615KB 0.010u 0:00.030
>> >
>> > When this image is read into the iPixType, I assume the range of these
>> > "unsigned short" values is mapped or cast into a range of "float"
>> > values.
>> > Is the new range of values from 0.0 to 1.0, is that correct for the
>> > "float"
>> > pixel type in ITK?  How does the ITK image reader convert or cast the
>> > pixel
>> > values to float?  (The reason the input images are specified as "float"
>> > rather than "unsigned short" is to have them input to the registration
>> > method as floats - should this be done with an explicit cast filter?)
>> >
>> > The input images have a black background (RGB value: 0 0 0).  It would
>> > be
>> > great to threshold the MeanSquaresImageToImageMetric to exclude all the
>> > black pixels.  In the range of values for the input pixel data that are
>> > an
>> > unsigned short (0 to 65535), any pixel above 0 would be
>> >
>> > unsigned short intensityThreshold = 1;
>> >
>> > However, the intensityThreshold is declared and initialized as:
>> >
>> > iBWImgType::PixelType intensityThreshold = 0.0;
>> >
>> > Is this a reasonable value for this float intensityThreshold:
>> >
>> > iBWImgType::PixelType intensityThreshold = 1.0 / USHRT_MAX;
>> >
>> >
>> > TIA,
>> > Darren
>> >
>> >
>> >
>> >
>> > On Sat, Mar 6, 2010 at 12:15 PM, Luis Ibanez <luis.ibanez at kitware.com>
>> > wrote:
>> >>
>> >> Hi Darren,
>> >>
>> >> The feature that you are referring to, is only available
>> >> in the new version of the metrics that are in the
>> >> Code/Review directory.
>> >>
>> >> In order to use this classes you must reconfigure your
>> >> binary build of ITK, by rerunning CMake and turning
>> >> ON the CMake variables:
>> >>
>> >>   *  ITK_USE_REVIEW
>> >>   *  OPTIMIZED_REGISTRATION_METHODS
>> >>
>> >> Then you can rebuild your application.
>> >>
>> >>
>> >>     Regards,
>> >>
>> >>
>> >>            Luis
>> >>
>> >>
>> >>
>> >>
>> >> -----------------------------------------------------------------------------------
>> >> On Thu, Mar 4, 2010 at 5:56 PM, Darren Weber
>> >> <darren.weber.lists at gmail.com> wrote:
>> >> >
>> >> > In an image registration framework, the metric determines image
>> >> > "matches".
>> >> > Let's suppose we have large images that contain a majority of
>> >> > background
>> >> > color (or transparency).  In a metric like the
>> >> > itk::MeanSquaresImageToImageMetric, the majority of the
>> >> > sqaured-diff-values
>> >> > would be zero and the remainder of the "relevant" pixel differences
>> >> > may
>> >> > be
>> >> > "swamped" or "lost" in the division by a large N.
>> >> > Also, it appears the default metric sampling is a random-sample of
>> >> > the
>> >> > image
>> >> > pixels (this might be modified with methods like
>> >> > "SetUseAllPixels(bool)").
>> >> >  In this case, it is possible, perhaps likely, that few of the
>> >> > "relevant"
>> >> > pixels will be evaluated in the metric.
>> >> > See:
>> >> >
>> >> >
>> >> > http://www.itk.org/Doxygen316/html/classitk_1_1ImageToImageMetric.html#42b876134388099afbf34b14faf83cdb
>> >> > This documentation suggests there are also options to define:
>> >> > IntensityThreshold, Masks, and SequentialSampling
>> >> > In my reading of the doxy page (see next link), there is a method
>> >> > called
>> >> > "SetFixedImageSamplesIntensityThreshold" that might be useful to
>> >> > define
>> >> > a
>> >> > background value (like zero) to be excluded from the metric.  Is that
>> >> > right?
>> >> >
>> >> >
>> >> > http://www.itk.org/Doxygen316/html/classitk_1_1MeanSquaresImageToImageMetric-members.html
>> >> > What is the logic of this "threshold" method?  From the
>> >> > documentation,
>> >> > it
>> >> > appears to be a minima value.  So if the background intensity is
>> >> > zero,
>> >> > the
>> >> > threshold value of 1 would include all pixels > 0 (note the threshold
>> >> > value
>> >> > must be the same type as the fixed image PixelType).
>> >> > There is no equivalent method for the moving image,
>> >> > like "SetMovingImageSamplesIntensityThreshold".  Is this unnecessary
>> >> > because
>> >> > the metric only scans the pixels of the moving image that fall within
>> >> > the
>> >> > "domain" of fixed image after the transform & interpolation?
>> >> > I could not locate this method in the class hierarchy - is it
>> >> > available
>> >> > (inherited) in all metrics?
>> >> > I'm confused because these methods do not appear to be available
>> >> > to itk::MeanSquaresImageToImageMetric
>> >> > e.g.:  CODE:
>> >> >     typedef itk::MeanSquaresImageToImageMetric< iBWImgType,
>> >> > iBWImgType >
>> >> > MetricType;
>> >> >     MetricType::Pointer metric = MetricType::New();
>> >> >     metric->SetFixedImageSamplesIntensityThreshold( 0.0 );
>> >> >     metric->SetUseAllPixels( true );
>> >> > e.g.:  COMPILATION:
>> >> > itkImageRigid2DCoregistration.cxx: In function ‘int main(int,
>> >> > char**)’:
>> >> > itkImageRigid2DCoregistration.cxx:442: error: ‘class
>> >> > itk::MeanSquaresImageToImageMetric<iBWImgType, iBWImgType>’ has no
>> >> > member
>> >> > named ‘SetFixedImageSamplesIntensityThreshold’
>> >> > itkImageRigid2DCoregistration.cxx:443: error: ‘class
>> >> > itk::MeanSquaresImageToImageMetric<iBWImgType, iBWImgType>’ has no
>> >> > member
>> >> > named ‘SetUseAllPixels’
>> >> > gmake[2]: ***
>> >> >
>> >> >
>> >> > [CMakeFiles/itkImageRigid2DCoregistration.dir/itkImageRigid2DCoregistration.cxx.o]
>> >> > Error 1
>> >> >
>> >> >
>> >> >
>> >> > TIA,
>> >> > Darren
>> >> >
>> >> > _____________________________________
>> >> > Powered by www.kitware.com
>> >> >
>> >> > Visit other Kitware open-source projects at
>> >> > http://www.kitware.com/opensource/opensource.html
>> >> >
>> >> > Kitware offers ITK Training Courses, for more information visit:
>> >> > http://www.kitware.com/products/protraining.html
>> >> >
>> >> > Please keep messages on-topic and check the ITK FAQ at:
>> >> > http://www.itk.org/Wiki/ITK_FAQ
>> >> >
>> >> > Follow this link to subscribe/unsubscribe:
>> >> > http://www.itk.org/mailman/listinfo/insight-users
>> >> >
>> >> >
>> >
>> >
>
>


More information about the Insight-users mailing list