ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkVideoStream.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
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 * https://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 itkVideoStream_h
19#define itkVideoStream_h
20
22#include "itkImage.h"
23
24namespace itk
25{
26
41template <typename TFrameType>
42class ITK_TEMPLATE_EXPORT VideoStream : public TemporalDataObject
43{
44public:
45 ITK_DISALLOW_COPY_AND_MOVE(VideoStream);
46
53
54 using FrameType = TFrameType;
55 using FramePointer = typename FrameType::Pointer;
56 using FrameConstPointer = typename FrameType::ConstPointer;
57 using typename Superclass::BufferType;
58
59 using SpatialRegionType = typename FrameType::RegionType;
60 using IndexType = typename FrameType::IndexType;
61 using PixelType = typename FrameType::PixelType;
62 using PointType = typename FrameType::PointType;
63 using SpacingType = typename FrameType::SpacingType;
64 using SizeType = typename FrameType::SizeType;
65 using DirectionType = typename FrameType::DirectionType;
66 using NumberOfComponentsPerPixelType = unsigned int;
67
69 using SpatialRegionMapType = typename std::map<SizeValueType, SpatialRegionType>;
70 using PointMapType = typename std::map<SizeValueType, PointType>;
71 using DirectionMapType = typename std::map<SizeValueType, DirectionType>;
72 using SpacingMapType = typename std::map<SizeValueType, SpacingType>;
73 using NumberOfComponentsPerPixelMapType = typename std::map<SizeValueType, NumberOfComponentsPerPixelType>;
74
76 static constexpr unsigned int FrameDimension = FrameType::ImageDimension;
77 static unsigned int
79 {
80 return FrameType::ImageDimension;
81 }
82
83 itkNewMacro(Self);
84
86 itkOverrideGetNameOfClassMacro(VideoStream);
87
89 void
90 SetMinimumBufferSize(SizeValueType minimumNumberOfFrames);
91
98 void
100
103 BufferType *
105 {
106 return reinterpret_cast<BufferType *>(m_DataObjectBuffer.GetPointer());
107 }
108 const BufferType *
110 {
111 return reinterpret_cast<BufferType *>(m_DataObjectBuffer.GetPointer());
112 }
113
115 void
117
125 void
130
131 const SpatialRegionMapType &
136 void
141
142 const SpatialRegionMapType &
147 void
152
153 const SpacingMapType &
155 {
156 return m_SpacingCache;
157 }
158 void
160 {
161 m_SpacingCache = map;
162 }
163
164 const PointMapType &
166 {
167 return m_OriginCache;
168 }
169 void
171 {
172 m_OriginCache = map;
173 }
174
175 const DirectionMapType &
177 {
178 return m_DirectionCache;
179 }
180 void
185
186 const NumberOfComponentsPerPixelMapType &
191 void
196
198 void
200
206 FrameType *
208 const FrameType *
209 GetFrame(SizeValueType frameNumber) const;
212 void
214
215 const SpatialRegionType &
217
219 void
221
222 const SpatialRegionType &
224
226 void
228
229 const SpatialRegionType &
231
233 void
235
236 const SpacingType &
237 GetFrameSpacing(SizeValueType frameNumber) const;
238
240 void
242
243 const PointType &
244 GetFrameOrigin(SizeValueType frameNumber) const;
245
247 void
249
250 const DirectionType &
252
254 void
256
259
263 void
265
269 void
271
275 void
277
281 void
283
287 void
289
293 void
295
299 void
301
334 void
336
341 void
342 Graft(const DataObject * data) override;
343
344protected:
345 VideoStream() = default;
346 ~VideoStream() override = default;
347
348 void
349 PrintSelf(std::ostream & os, Indent indent) const override
350 {
351 Superclass::Print(os, indent);
352 }
353
360
367
368}; // end class VideoStream
369
370} // end namespace itk
371
372#ifndef ITK_MANUAL_INSTANTIATION
373# include "itkVideoStream.hxx"
374#endif
375
376#endif
Base class for all data objects in ITK.
Control indentation during Print() invocation.
Definition itkIndent.h:50
void Print(std::ostream &os, Indent indent=0) const
Implements transparent reference counting.
RingBuffer< DataObject > BufferType
BufferType::Pointer m_DataObjectBuffer
TFrameType FrameType
const SpacingMapType & GetSpacingCache() const
WeakPointer< const Self > ConstWeakPointer
void GetNumberOfComponentsPerPixelCache(NumberOfComponentsPerPixelMapType map)
RingBuffer< DataObject > BufferType
BufferType * GetFrameBuffer()
void PrintSelf(std::ostream &os, Indent indent) const override
void SetFrameLargestPossibleSpatialRegion(SizeValueType frameNumber, SpatialRegionType region)
void SetAllBufferedSpatialRegions(SpatialRegionType region)
const DirectionType & GetFrameDirection(SizeValueType frameNumber) const
VideoStream()=default
void SetAllRequestedSpatialRegions(SpatialRegionType region)
void SetFrameRequestedSpatialRegion(SizeValueType frameNumber, SpatialRegionType region)
typename FrameType::RegionType SpatialRegionType
typename std::map< SizeValueType, NumberOfComponentsPerPixelType > NumberOfComponentsPerPixelMapType
void SetSpacingCache(SpacingMapType map)
~VideoStream() override=default
typename FrameType::Pointer FramePointer
const BufferType * GetFrameBuffer() const
void SetAllFramesDirection(DirectionType direction)
void SetAllLargestPossibleSpatialRegions(SpatialRegionType region)
const SpatialRegionType & GetFrameRequestedSpatialRegion(SizeValueType frameNumber) const
typename std::map< SizeValueType, DirectionType > DirectionMapType
typename std::map< SizeValueType, SpacingType > SpacingMapType
typename FrameType::PixelType PixelType
void InitializeEmptyFrames()
const SpatialRegionType & GetFrameBufferedSpatialRegion(SizeValueType frameNumber) const
const DirectionMapType & GetDirectionCache() const
void SetMinimumBufferSize(SizeValueType minimumNumberOfFrames)
const NumberOfComponentsPerPixelMapType & GetNumberOfComponentsPerPixelCache() const
SmartPointer< const Self > ConstPointer
void SetBufferedSpatialRegionCache(SpatialRegionMapType map)
const SpatialRegionType & GetFrameLargestPossibleSpatialRegion(SizeValueType frameNumber) const
const PointType & GetFrameOrigin(SizeValueType frameNumber) const
void SetDirectionCache(DirectionMapType map)
FrameType * GetFrame(SizeValueType frameNumber)
void SetFrameBufferedSpatialRegion(SizeValueType frameNumber, SpatialRegionType region)
typename FrameType::IndexType IndexType
typename FrameType::DirectionType DirectionType
SpatialRegionMapType m_RequestedSpatialRegionCache
void SetFrameBuffer(BufferType *buffer)
const FrameType * GetFrame(SizeValueType frameNumber) const
typename FrameType::SpacingType SpacingType
typename FrameType::PointType PointType
TemporalDataObject Superclass
void SetOriginCache(PointMapType map)
void Graft(const DataObject *data) override
const PointMapType & GetOriginCache() const
typename std::map< SizeValueType, PointType > PointMapType
const NumberOfComponentsPerPixelType & GetFrameNumberOfComponentsPerPixel(SizeValueType frameNumber) const
NumberOfComponentsPerPixelMapType m_NumberOfComponentsPerPixelCache
unsigned int NumberOfComponentsPerPixelType
void SetAllFramesSpacing(SpacingType spacing)
SpatialRegionMapType m_BufferedSpatialRegionCache
static constexpr unsigned int FrameDimension
typename std::map< SizeValueType, SpatialRegionType > SpatialRegionMapType
void SetLargestPossibleSpatialRegionCache(SpatialRegionMapType map)
void SetFrameOrigin(SizeValueType frameNumber, PointType origin)
const SpatialRegionMapType & GetRequestedSpatialRegionCache() const
void SetAllFramesNumberOfComponentsPerPixel(NumberOfComponentsPerPixelType n)
typename FrameType::ConstPointer FrameConstPointer
void SetAllFramesOrigin(PointType origin)
SmartPointer< Self > Pointer
void SetRequestedSpatialRegionCache(SpatialRegionMapType map)
typename FrameType::SizeType SizeType
static unsigned int GetFrameDimension()
void SetFrame(SizeValueType frameNumber, FramePointer frame)
const SpatialRegionMapType & GetBufferedSpatialRegionCache() const
void SetFrameSpacing(SizeValueType frameNumber, SpacingType spacing)
const SpatialRegionMapType & GetLargestPossibleSpatialRegionCache() const
const SpacingType & GetFrameSpacing(SizeValueType frameNumber) const
void SetFrameDirection(SizeValueType frameNumber, DirectionType direction)
SpatialRegionMapType m_LargestPossibleSpatialRegionCache
void SetFrameNumberOfComponentsPerPixel(SizeValueType frameNumber, unsigned int n)
Implements a weak reference to an object.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition itkIntTypes.h:86