ITK
4.2.0
Insight Segmentation and Registration Toolkit
Main Page
Related Pages
Modules
Namespaces
Classes
Files
Examples
File List
File Members
ITK
Modules
Core
Transform
include
itkMatrixOffsetTransformBase.h
Go to the documentation of this file.
1
/*=========================================================================
2
*
3
* Copyright Insight Software Consortium
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 __itkMatrixOffsetTransformBase_h
19
#define __itkMatrixOffsetTransformBase_h
20
21
#include <iostream>
22
23
#include "
itkMatrix.h
"
24
#include "
itkTransform.h
"
25
#include "
itkMacro.h
"
26
27
namespace
itk
28
{
70
template
<
71
class
TScalarType = double,
// Data type for scalars
72
unsigned
int
NInputDimensions = 3,
// Number of dimensions in the input space
73
unsigned
int
NOutputDimensions = 3>
74
// Number of dimensions in the output space
75
class
MatrixOffsetTransformBase
:
76
public
Transform
<TScalarType, NInputDimensions, NOutputDimensions>
77
{
78
public
:
80
typedef
MatrixOffsetTransformBase
Self
;
81
typedef
Transform
<TScalarType,
82
NInputDimensions,
83
NOutputDimensions>
Superclass
;
84
85
typedef
SmartPointer<Self>
Pointer
;
86
typedef
SmartPointer<const Self>
ConstPointer
;
87
89
itkTypeMacro(
MatrixOffsetTransformBase
,
Transform
);
90
92
itkNewMacro(
Self
);
93
95
itkStaticConstMacro(
InputSpaceDimension
,
unsigned
int
, NInputDimensions);
96
itkStaticConstMacro(
OutputSpaceDimension
,
unsigned
int
, NOutputDimensions);
97
itkStaticConstMacro(
ParametersDimension
,
unsigned
int
,
98
NOutputDimensions * ( NInputDimensions + 1 ) );
100
102
typedef
typename
Superclass::ParametersType
ParametersType
;
103
typedef
typename
Superclass::ParametersValueType
ParametersValueType
;
104
106
typedef
typename
Superclass::JacobianType
JacobianType
;
107
109
typedef
typename
Superclass::ScalarType
ScalarType
;
110
112
typedef
Vector
<TScalarType,
113
itkGetStaticConstMacro(
InputSpaceDimension
)>
InputVectorType
;
114
typedef
Vector
<TScalarType,
115
itkGetStaticConstMacro(
OutputSpaceDimension
)>
OutputVectorType
;
116
typedef
typename
OutputVectorType::ValueType
OutputVectorValueType
;
118
120
typedef
CovariantVector
<TScalarType,
121
itkGetStaticConstMacro(
InputSpaceDimension
)>
122
InputCovariantVectorType
;
123
typedef
CovariantVector
<TScalarType,
124
itkGetStaticConstMacro(
OutputSpaceDimension
)>
125
OutputCovariantVectorType
;
127
128
typedef
typename
Superclass::InputVectorPixelType
InputVectorPixelType
;
129
typedef
typename
Superclass::OutputVectorPixelType
OutputVectorPixelType
;
130
132
typedef
typename
Superclass::InputDiffusionTensor3DType
133
InputDiffusionTensor3DType
;
134
typedef
typename
Superclass::OutputDiffusionTensor3DType
135
OutputDiffusionTensor3DType
;
136
138
typedef
typename
Superclass::InputSymmetricSecondRankTensorType
139
InputSymmetricSecondRankTensorType
;
140
typedef
typename
Superclass::OutputSymmetricSecondRankTensorType
141
OutputSymmetricSecondRankTensorType
;
142
143
typedef
CovariantVector<TScalarType, InputDiffusionTensor3DType::Dimension>
144
InputTensorEigenVectorType
;
145
147
typedef
vnl_vector_fixed<TScalarType,
148
itkGetStaticConstMacro(
InputSpaceDimension
)>
149
InputVnlVectorType
;
150
typedef
vnl_vector_fixed<TScalarType,
151
itkGetStaticConstMacro(
OutputSpaceDimension
)>
152
OutputVnlVectorType
;
154
156
typedef
Point
<TScalarType,
157
itkGetStaticConstMacro(
InputSpaceDimension
)>
158
InputPointType
;
159
typedef
typename
InputPointType::ValueType
InputPointValueType
;
160
typedef
Point
<TScalarType,
161
itkGetStaticConstMacro(
OutputSpaceDimension
)>
162
OutputPointType
;
163
typedef
typename
OutputPointType::ValueType
OutputPointValueType
;
165
167
typedef
Matrix
<TScalarType, itkGetStaticConstMacro(
OutputSpaceDimension
),
168
itkGetStaticConstMacro(
InputSpaceDimension
)>
169
MatrixType
;
170
typedef
typename
MatrixType::ValueType
MatrixValueType
;
172
174
typedef
Matrix
<TScalarType, itkGetStaticConstMacro(
InputSpaceDimension
),
175
itkGetStaticConstMacro(
OutputSpaceDimension
)>
176
InverseMatrixType
;
177
178
typedef
InputPointType
CenterType
;
179
180
typedef
OutputVectorType
OffsetType
;
181
typedef
typename
OffsetType::ValueType
OffsetValueType
;
182
183
typedef
OutputVectorType
TranslationType
;
184
185
typedef
typename
TranslationType::ValueType
TranslationValueType
;
186
189
typedef
typename
Superclass::InverseTransformBaseType
InverseTransformBaseType
;
190
typedef
typename
InverseTransformBaseType::Pointer
InverseTransformBasePointer
;
191
195
virtual
void
SetIdentity
(
void
);
196
208
virtual
void
SetMatrix
(
const
MatrixType
& matrix)
209
{
210
m_Matrix
= matrix; this->
ComputeOffset
();
211
this->
ComputeMatrixParameters
();
212
m_MatrixMTime
.
Modified
(); this->
Modified
();
return
;
213
}
215
223
virtual
const
MatrixType
&
GetMatrix
()
const
224
{
225
return
m_Matrix
;
226
}
227
236
void
SetOffset
(
const
OutputVectorType
& offset)
237
{
238
m_Offset
= offset; this->
ComputeTranslation
();
239
this->
Modified
();
return
;
240
}
242
248
const
OutputVectorType
&
GetOffset
(
void
)
const
249
{
250
return
m_Offset
;
251
}
252
275
void
SetCenter
(
const
InputPointType
& center)
276
{
277
m_Center
= center; this->
ComputeOffset
();
278
this->
Modified
();
return
;
279
}
281
288
const
InputPointType
&
GetCenter
()
const
289
{
290
return
m_Center
;
291
}
292
299
void
SetTranslation
(
const
OutputVectorType
& translation)
300
{
301
m_Translation
= translation; this->
ComputeOffset
();
302
this->
Modified
();
return
;
303
}
305
312
const
OutputVectorType
&
GetTranslation
(
void
)
const
313
{
314
return
m_Translation
;
315
}
316
321
void
SetParameters
(
const
ParametersType
& parameters);
322
324
const
ParametersType
&
GetParameters
(
void
)
const
;
325
327
virtual
void
SetFixedParameters
(
const
ParametersType
&);
328
330
virtual
const
ParametersType
&
GetFixedParameters
(
void
)
const
;
331
343
void
Compose
(
const
Self
*other,
bool
pre = 0);
344
353
OutputPointType
TransformPoint
(
const
InputPointType
& point)
const
;
354
355
using
Superclass::TransformVector
;
356
357
OutputVectorType
TransformVector
(
const
InputVectorType
& vector)
const
;
358
359
OutputVnlVectorType
TransformVector
(
const
InputVnlVectorType
& vector)
const
;
360
361
OutputVectorPixelType
TransformVector
(
const
InputVectorPixelType
& vector)
const
;
362
363
using
Superclass::TransformCovariantVector
;
364
365
OutputCovariantVectorType
TransformCovariantVector
(
const
InputCovariantVectorType
& vector)
const
;
366
367
OutputVectorPixelType
TransformCovariantVector
(
const
InputVectorPixelType
& vector)
const
;
368
369
using
Superclass::TransformDiffusionTensor3D
;
370
371
OutputDiffusionTensor3DType
TransformDiffusionTensor3D
(
const
InputDiffusionTensor3DType
& tensor)
const
;
372
373
OutputVectorPixelType
TransformDiffusionTensor3D
(
const
InputVectorPixelType
& tensor )
const
;
374
375
using
Superclass::TransformSymmetricSecondRankTensor
;
376
OutputSymmetricSecondRankTensorType
TransformSymmetricSecondRankTensor
(
const
InputSymmetricSecondRankTensorType
& tensor )
const
;
377
378
OutputVectorPixelType
TransformSymmetricSecondRankTensor
(
const
InputVectorPixelType
& tensor )
const
;
379
389
virtual
void
ComputeJacobianWithRespectToParameters
(
const
InputPointType
& x,
JacobianType
& j)
const
;
390
394
virtual
void
ComputeJacobianWithRespectToPosition
(
const
InputPointType
& x,
JacobianType
& jac)
const
;
395
399
virtual
void
ComputeInverseJacobianWithRespectToPosition
(
const
InputPointType
& x,
JacobianType
& jac)
const
;
400
419
bool
GetInverse
(
Self
*inverse)
const
;
421
423
virtual
InverseTransformBasePointer
GetInverseTransform
()
const
;
424
430
virtual
bool
IsLinear
()
const
431
{
432
return
true
;
433
}
434
435
#if !defined(ITK_LEGACY_REMOVE)
436
public
:
437
#else
438
protected
:
439
#endif
440
442
const
InverseMatrixType
&
GetInverseMatrix
(
void
)
const
;
443
protected
:
444
452
MatrixOffsetTransformBase
(
const
MatrixType
& matrix,
const
OutputVectorType
& offset);
453
MatrixOffsetTransformBase
(
unsigned
int
paramDims);
454
MatrixOffsetTransformBase
();
456
458
virtual
~MatrixOffsetTransformBase
();
459
461
void
PrintSelf
(std::ostream & s,
Indent
indent)
const
;
462
463
464
const
InverseMatrixType
&
GetVarInverseMatrix
(
void
)
const
465
{
466
return
m_InverseMatrix
;
467
}
468
void
SetVarInverseMatrix
(
const
InverseMatrixType
& matrix)
const
469
{
470
m_InverseMatrix
= matrix;
m_InverseMatrixMTime
.
Modified
();
471
}
472
bool
InverseMatrixIsOld
(
void
)
const
473
{
474
if
(
m_MatrixMTime
!=
m_InverseMatrixMTime
)
475
{
476
return
true
;
477
}
478
else
479
{
480
return
false
;
481
}
482
}
483
484
virtual
void
ComputeMatrixParameters
(
void
);
485
486
virtual
void
ComputeMatrix
(
void
);
487
488
void
SetVarMatrix
(
const
MatrixType
& matrix)
489
{
490
m_Matrix
= matrix;
m_MatrixMTime
.
Modified
();
491
}
492
493
virtual
void
ComputeTranslation
(
void
);
494
495
void
SetVarTranslation
(
const
OutputVectorType
& translation)
496
{
497
m_Translation
= translation;
498
}
499
500
virtual
void
ComputeOffset
(
void
);
501
502
void
SetVarOffset
(
const
OutputVectorType
& offset)
503
{
504
m_Offset
= offset;
505
}
506
507
void
SetVarCenter
(
const
InputPointType
& center)
508
{
509
m_Center
= center;
510
}
511
private
:
512
513
MatrixOffsetTransformBase
(
const
Self
& other);
514
const
Self
&
operator=
(
const
Self
&);
515
516
MatrixType
m_Matrix
;
// Matrix of the transformation
517
OutputVectorType
m_Offset
;
// Offset of the transformation
518
mutable
InverseMatrixType
m_InverseMatrix
;
// Inverse of the matrix
519
mutable
bool
m_Singular
;
// Is m_Inverse singular?
520
521
InputPointType
m_Center
;
522
OutputVectorType
m_Translation
;
523
525
TimeStamp
m_MatrixMTime
;
526
mutable
TimeStamp
m_InverseMatrixMTime
;
527
};
// class MatrixOffsetTransformBase
528
}
// namespace itk
529
530
// Define instantiation macro for this template.
531
#define ITK_TEMPLATE_MatrixOffsetTransformBase(_, EXPORT, TypeX, TypeY) \
532
namespace itk \
533
{ \
534
_( 3 ( class EXPORT MatrixOffsetTransformBase<ITK_TEMPLATE_3 TypeX> ) ) \
535
namespace Templates \
536
{ \
537
typedef MatrixOffsetTransformBase<ITK_TEMPLATE_3 TypeX> MatrixOffsetTransformBase##TypeY; \
538
} \
539
}
540
541
#if ITK_TEMPLATE_EXPLICIT
542
// template < class TScalarType, unsigned int NInputDimensions, unsigned int
543
// NOutputDimensions>
544
// const unsigned int itk::MatrixOffsetTransformBase< TScalarType,
545
// NInputDimensions, NOutputDimensions >::InputSpaceDimension;
546
// template < class TScalarType, unsigned int NInputDimensions, unsigned int
547
// NOutputDimensions>
548
// const unsigned int itk::MatrixOffsetTransformBase< TScalarType,
549
// NInputDimensions, NOutputDimensions >::OutputSpaceDimension;
550
// template < class TScalarType, unsigned int NInputDimensions, unsigned int
551
// NOutputDimensions>
552
// const unsigned int itk::MatrixOffsetTransformBase< TScalarType,
553
// NInputDimensions, NOutputDimensions >::ParametersDimension;
554
#include "Templates/itkMatrixOffsetTransformBase+-.h"
555
#endif
556
557
#if ITK_TEMPLATE_TXX
558
#include "itkMatrixOffsetTransformBase.hxx"
559
#endif
560
561
#endif
/* __itkMatrixOffsetTransformBase_h */
562
Generated on Tue Jul 10 2012 23:36:08 for ITK by
1.8.1