go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkImageRandomSamplerSparseMask.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * http://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18#ifndef itkImageRandomSamplerSparseMask_h
19#define itkImageRandomSamplerSparseMask_h
20
22#include "itkMersenneTwisterRandomVariateGenerator.h"
23#include "itkImageFullSampler.h"
24
25namespace itk
26{
37template <class TInputImage>
38class ITK_TEMPLATE_EXPORT ImageRandomSamplerSparseMask : public ImageRandomSamplerBase<TInputImage>
39{
40public:
42
46 using Pointer = SmartPointer<Self>;
47 using ConstPointer = SmartPointer<const Self>;
48
50 itkNewMacro(Self);
51
54
56 using typename Superclass::DataObjectPointer;
59 using typename Superclass::InputImageType;
60 using typename Superclass::InputImagePointer;
64 using typename Superclass::ImageSampleType;
67 using typename Superclass::MaskType;
68
70 itkStaticConstMacro(InputImageDimension, unsigned int, Superclass::InputImageDimension);
71
73 using InputImageIndexType = typename InputImageType::IndexType;
74 using InputImagePointType = typename InputImageType::PointType;
75
77 using RandomGeneratorType = itk::Statistics::MersenneTwisterRandomVariateGenerator;
78 using RandomGeneratorPointer = typename RandomGeneratorType::Pointer;
79
80protected:
83
87 ~ImageRandomSamplerSparseMask() override = default;
88
90 void
91 PrintSelf(std::ostream & os, Indent indent) const override;
92
94 void
95 GenerateData() override;
96
97 RandomGeneratorPointer m_RandomGenerator{ RandomGeneratorType::GetInstance() };
98 InternalFullSamplerPointer m_InternalFullSampler{ InternalFullSamplerType::New() };
99
100private:
101 struct UserData
102 {
104
105 const std::vector<ImageSampleType> & m_AllValidSamples;
106 const std::vector<size_t> & m_RandomIndices;
107 std::vector<ImageSampleType> & m_Samples;
108 };
109
110 std::vector<size_t> m_RandomIndices{};
111
112 static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
113 ThreaderCallback(void * arg);
114};
115
116} // end namespace itk
117
118#ifndef ITK_MANUAL_INSTANTIATION
119# include "itkImageRandomSamplerSparseMask.hxx"
120#endif
121
122#endif // end #ifndef itkImageRandomSamplerSparseMask_h
Samples all voxels in the InputImageRegion.
SmartPointer< Self > Pointer
This class is a base class for any image sampler that randomly picks samples.
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::RegionType InputImageRegionType
typename ImageSampleContainerType::Pointer ImageSampleContainerPointer
typename InputImageType::Pointer InputImagePointer
ImageMaskSpatialObject< Self::InputImageDimension > MaskType
typename InputImageType::PixelType InputImagePixelType
Samples randomly some voxels of an image.
typename RandomGeneratorType::Pointer RandomGeneratorPointer
itk::Statistics::MersenneTwisterRandomVariateGenerator RandomGeneratorType
typename InputImageType::PointType InputImagePointType
void PrintSelf(std::ostream &os, Indent indent) const override
ITK_DISALLOW_COPY_AND_MOVE(ImageRandomSamplerSparseMask)
~ImageRandomSamplerSparseMask() override=default
typename InternalFullSamplerType::Pointer InternalFullSamplerPointer
typename InputImageType::IndexType InputImageIndexType
itkStaticConstMacro(InputImageDimension, unsigned int, Superclass::InputImageDimension)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreaderCallback(void *arg)
A class that defines an image sample, which is the coordinates of a point and its value.
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
const std::vector< ImageSampleType > & m_AllValidSamples


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo