[ITK-users] [ITK] Filter working on both itk::Image and itk::VectorImage
Lowekamp, Bradley (NIH/NLM/LHC) [C]
blowekamp at mail.nih.gov
Mon Feb 27 13:46:52 EST 2017
Hello,
Sorry I didn’t see this sooner.
You want to look at the itk::Image::Rebind class. I created it for this exact purpose. [1]
You can look at itkSmoothingRecursiveGaussianImageFilter.h for an example usage.
HTH,
Brad
[1] https://itk.org/Doxygen/html/structitk_1_1Image_1_1Rebind.html
> On Feb 27, 2017, at 6:09 AM, Cyril Mory <cyril.mory at creatis.insa-lyon.fr> wrote:
>
> Hi,
>
> Still the same topic, but a different issue:
>
> I have a filter that is templated over its image input type, usually itk::Image<T, dim>. I want to modify it so that it also works when that type is itk::VectorImage<T, dim>.
>
> This filter internally defines SliceType = itk::Image<InputImageType::InternalPixelType, dim - 1>. But when processing an itk::VectorImage, I would like this SliceType to be itk::VectorImage<InputImageType::InternalPixelType, dim - 1>.
>
> I don't think there is a way in C++ to recover the itk::Image<> or the itk::VectorImage<> "un-templated", and define a type with whole new template arguments. But is there an ITK mechanism that would solve my problem ?
>
> Best regards,
>
> Cyril
>
> On 22/02/2017 17:39, Matt McCormick wrote:
>> Hi Cyril,
>>
>> Yes, this is a good approach. Thanks for sharing your investigation.
>>
>> Matt
>>
>> On Wed, Feb 22, 2017 at 5:52 AM, Cyril Mory
>> <cyril.mory at creatis.insa-lyon.fr> wrote:
>>> Answering my own question:
>>>
>>> I found a way to do what I needed using itk::NumericTraits<T> :
>>>
>>> OutputImagePixelType pix;
>>> itk::NumericTraits<OutputImagePixelType>::SetLength(pix,
>>> this->GetVectorLength());
>>> pix = itk::NumericTraits<OutputImagePixelType>::OneValue(pix) * value;
>>>
>>> It works for both cases with the same code. However, if there is a better
>>> way (cleaner, safer, more compact, ... whatever the wayin which it is
>>> better), I'd be happy to learn about it.
>>>
>>> Best regards,
>>> Cyril
>>>
>>>
>>> On 22/02/2017 10:38, Cyril Mory wrote:
>>>> Hi ITK users,
>>>>
>>>> I am writing an itk::ImageToImageFilter<T>, and I want the filter to work
>>>> with T=itk::Image as well as with T=itk::VectorImage.
>>>>
>>>> To set/get the vectorLength of an image (which is 1 for an itk::Image, and
>>>> varies for an itk::VectorImage), I have successfully used the functions
>>>> Set/GetNumberOfComponentsPerPixel from the itk::ImageBase class, which work
>>>> for both itk::Image and itk::VectorImage.
>>>>
>>>> But I cannot find a way to initialize one of the pixels to a given value,
>>>> that would work for both cases, i.e. whether that pixel is a scalar or a
>>>> variableLengthVector. For now, the best I could do is define a function
>>>> FillPixel in my filter, like this:
>>>>
>>>> OutputImagePixelType FillPixel(OutputImageInternalPixelType value) {return
>>>> value;}
>>>>
>>>> and write a template specialization when my filter is instantiated with
>>>> itk::VectorImage<float, 3>, like this:
>>>>
>>>> template <>
>>>> itk::VariableLengthVector<float>
>>>> rtk::ConstantImageSource<itk::VectorImage<float, 3> >
>>>> ::FillPixel(float value)
>>>> {
>>>> itk::VariableLengthVector<float> vect;
>>>> vect.SetSize(this->GetVectorLength());
>>>> vect.Fill(value);
>>>> return (vect);
>>>> }
>>>>
>>>> Did I miss something ? Is there a better way of doing this in ITK ?
>>>>
>>>> Looking forward to reading you,
>>>> Cyril
>>>>
>>>> _____________________________________
>>>> 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.php
>>>>
>>>> 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://public.kitware.com/mailman/listinfo/insight-users
>>>
>>> _____________________________________
>>> 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.php
>>>
>>> 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://public.kitware.com/mailman/listinfo/insight-users
>
> _____________________________________
> 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.php
>
> 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://public.kitware.com/mailman/listinfo/insight-users
> _______________________________________________
> Community mailing list
> Community at itk.org
> http://public.kitware.com/mailman/listinfo/community
More information about the Insight-users
mailing list