#include <vtkImageT2StarMap.h>
Inheritance diagram for vtkImageT2StarMap:


This filter must have 3 or greater inputs and produces an output with 6 components. The inputs must all have the same extents, they must have one scalar component, and the same scalar type. Each input is an echo from a multi-echo T2*-weighted MRI sequence which will be used to calculate the T2* time in seconds. The output will have the same extents as the inputs, but will be type double. The first component of the output will be the T2* time in seconds, the second component will be the constant A as shown in the function below, and the second to last component will be the r-squared value from the curve fitting. The third and fourth components will vary depending on the type of T2* fitting selected. For VTK_TTWOSTAR_LINEAR and VTK_TTWOSTAR_NON_LINEAR the third and fourth components will be zero. For VTK_TTWOSTAR_NON_LINEAR_WITH_CONSTANT the third component will be the value C as shown below and the fourth component will be zero. For VTK_TTWOSTAR_NON_LINEAR_SINC_MOD_THREE_PARAM the third component will be zero and the fourth component will be the fitted gradient map (Gb). For VTK_TTWOSTAR_NON_LINEAR_SINC_MOD_THREE_PARAM_WITH_CONSTANT the third component will be c as shown below and the fourth component will be the fitted gradient map. The last component will be the error improvement from the VTK_TTWOSTAR_NON_LINEAR fittings or zero for the VTK_TTWOSTAR_LINEAR fitting.
NOTE: The nonlinear fitting functions are only enabled if the levmar library is used. If the levmar library is used, then this file must be released under the GNU General Public License.
VTK_TTWOSTAR_LINEAR (Linear least squares): S = Ae^-(TE/T2*)
VTK_TTWOSTAR_NON_LINEAR (Non-linear least squares using Levenberg-Marquardt): S = Ae^-(TE/T2)
VTK_TTWOSTAR_NON_LINEAR_WITH_CONSTANT (Non-linear least squares using Levenberg-Marquardt): S = Ae^-(TE/T2) + C
VTK_TTWOSTAR_NON_LINEAR_SINC_MOD_THREE_PARAM (Non-linear least squares using Levenberg-Marquardt): S = A*|sinc(gamma/2*deltaB*TE)|*e^-(TE/T2)
VTK_TTWOSTAR_NON_LINEAR_SINC_MOD_THREE_PARAM_WITH_CONSTANT (Non-linear least squares using Levenberg-Marquardt): S = A*|sinc(gamma/2*deltaB*TE)|*e^-(TE/T2) + C
The echo times are set via the SetEchoTimes function. If any of the VTK_TTWOSTAR_NON_LINEAR_WITH_CONSTANT are selected then no single echo time may appear more than once in the echo array. If SmartT2StarFiltering is set then the filter will set T2* times that are greater than or equal to the value specified by MaxT2StarTime to MaxT2StarTime, which is set to 10 seconds by default. If PerformR2StarMapping is set to On then the T2* output will be the inverse of T2* or R2* in Hz (1/s). Gamma (the gyromagnetic ratio) is set using
IMPORTANT NOTE: When any of the sinc modulated functions are set the last input after the echo images must be the calculated gradient (Gb) image. This class cannot perform this type of fitting without it and will exit with an error if it has not been added as an input.
Definition at line 118 of file vtkImageT2StarMap.h.
| vtkImageT2StarMap::vtkImageT2StarMap | ( | ) | [protected] |
| vtkImageT2StarMap::~vtkImageT2StarMap | ( | ) | [protected] |
| static vtkImageT2StarMap* vtkImageT2StarMap::New | ( | ) | [static] |
Reimplemented from vtkAlgorithm.
| vtkImageT2StarMap::vtkTypeRevisionMacro | ( | vtkImageT2StarMap | , | |
| vtkImageMultipleInputFilter | ||||
| ) |
| void vtkImageT2StarMap::PrintSelf | ( | ostream & | os, | |
| vtkIndent | indent | |||
| ) | [virtual] |
Reimplemented from vtkImageMultipleInputFilter.
| void vtkImageT2StarMap::SetEchoTimes | ( | vtkDoubleArray * | _arg | ) |
Set the echo times for the T2* map calculation.
| vtkDoubleArray* vtkImageT2StarMap::GetEchoTimes | ( | ) |
Get the echo times for the T2* map calculation.
| virtual void vtkImageT2StarMap::SetT2StarFitType | ( | int | ) | [virtual] |
Set/Get the T2* fitting type.
| virtual int vtkImageT2StarMap::GetT2StarFitType | ( | ) | [virtual] |
Set/Get the T2* fitting type.
| void vtkImageT2StarMap::SetT2StarFitTypeToLinear | ( | ) | [inline] |
Set/Get the T2* fitting type.
Definition at line 135 of file vtkImageT2StarMap.h.
References VTK_TTWOSTAR_LINEAR.
| void vtkImageT2StarMap::SetT2StarFitTypeToNonLinear | ( | ) |
Set/Get the T2* fitting type.
| void vtkImageT2StarMap::SetT2StarFitTypeToNonLinearWithConstant | ( | ) |
Set/Get the T2* fitting type.
| void vtkImageT2StarMap::SetT2StarFitTypeToNonLinearSincModThreeParam | ( | ) |
Set/Get the T2* fitting type.
| void vtkImageT2StarMap::SetT2StarFitTypeToNonLinearSincModThreeParamWithConstant | ( | ) |
Set/Get the T2* fitting type.
| virtual void vtkImageT2StarMap::SetSmartT2StarFiltering | ( | int | ) | [virtual] |
Set/Get whether to use smart T2 filtering.
| virtual int vtkImageT2StarMap::GetSmartT2StarFiltering | ( | ) | [virtual] |
Set/Get whether to use smart T2 filtering.
| virtual void vtkImageT2StarMap::SmartT2StarFilteringOn | ( | ) | [virtual] |
Set/Get whether to use smart T2 filtering.
| virtual void vtkImageT2StarMap::SmartT2StarFilteringOff | ( | ) | [virtual] |
Set/Get whether to use smart T2 filtering.
| void vtkImageT2StarMap::SetMaxT2StarTime | ( | double | time | ) |
Set smart T2 filtering max T2 time.
| double vtkImageT2StarMap::GetMaxT2StarTime | ( | ) |
Set smart T2 filtering max T2 time.
| virtual void vtkImageT2StarMap::SetPerformR2StarMapping | ( | int | ) | [virtual] |
Set/Get whether to perform R2 mapping.
| virtual int vtkImageT2StarMap::GetPerformR2StarMapping | ( | ) | [virtual] |
Set/Get whether to perform R2 mapping.
| virtual void vtkImageT2StarMap::PerformR2StarMappingOn | ( | ) | [virtual] |
Set/Get whether to perform R2 mapping.
| virtual void vtkImageT2StarMap::PerformR2StarMappingOff | ( | ) | [virtual] |
Set/Get whether to perform R2 mapping.
| virtual void vtkImageT2StarMap::SetGamma | ( | double | ) | [virtual] |
Set/Get the gyromagnetic ratio (Gamma) used in the sinc modulated fitting.
| virtual double vtkImageT2StarMap::GetGamma | ( | ) | [virtual] |
Set/Get the gyromagnetic ratio (Gamma) used in the sinc modulated fitting.
| void vtkImageT2StarMap::ExecuteInformation | ( | vtkImageData ** | inputs, | |
| vtkImageData * | output | |||
| ) | [protected, virtual] |
Reimplemented from vtkImageMultipleInputFilter.
| void vtkImageT2StarMap::ExecuteInformation | ( | ) | [inline, protected, virtual] |
Reimplemented from vtkImageMultipleInputFilter.
Definition at line 176 of file vtkImageT2StarMap.h.
References vtkImageMultipleInputFilter::ExecuteInformation().
| void vtkImageT2StarMap::ThreadedExecute | ( | vtkImageData ** | inDatas, | |
| vtkImageData * | outData, | |||
| int | extent[6], | |||
| int | id | |||
| ) | [protected, virtual] |
Reimplemented from vtkImageMultipleInputFilter.
| void vtkImageT2StarMap::InitOutput | ( | int | outExt[6], | |
| vtkImageData * | outData | |||
| ) | [protected] |
vtkDoubleArray* vtkImageT2StarMap::EchoTimes [protected] |
Definition at line 183 of file vtkImageT2StarMap.h.
int vtkImageT2StarMap::T2StarFitType [protected] |
Definition at line 184 of file vtkImageT2StarMap.h.
int vtkImageT2StarMap::SmartT2StarFiltering [protected] |
Definition at line 185 of file vtkImageT2StarMap.h.
double vtkImageT2StarMap::MaxT2StarTime [protected] |
Definition at line 186 of file vtkImageT2StarMap.h.
int vtkImageT2StarMap::PerformR2StarMapping [protected] |
Definition at line 187 of file vtkImageT2StarMap.h.
double vtkImageT2StarMap::Gamma [protected] |
Definition at line 188 of file vtkImageT2StarMap.h.
1.4.7