19 #ifndef rtkJosephForwardAttenuatedProjectionImageFilter_h 20 #define rtkJosephForwardAttenuatedProjectionImageFilter_h 22 #include "rtkConfiguration.h" 42 template <
class TInput,
class TCoordRepType,
class TOutput = TInput>
50 m_AttenuationRay[i] = 0;
51 m_AttenuationPixel[i] = 0;
66 return !(*
this != other);
71 const double stepLengthInVoxel,
72 const TCoordRepType weight,
76 const double w = weight * stepLengthInVoxel;
78 m_AttenuationRay[threadId] += w * (p + m_AttenuationMinusEmissionMapsPtrDiff)[i];
79 m_AttenuationPixel[threadId] += w * (p + m_AttenuationMinusEmissionMapsPtrDiff)[i];
86 m_AttenuationMinusEmissionMapsPtrDiff = pd;
91 return m_AttenuationRay;
96 return m_AttenuationPixel;
118 template <
class TInput,
class TOutput>
135 return !(*
this != other);
141 TInput ex2 = exp(-m_AttenuationRay[threadId] * stepInMM.GetNorm());
144 if (m_AttenuationPixel[threadId] > 0)
146 wf = (m_Ex1[threadId] - ex2) / m_AttenuationPixel[threadId];
150 wf = m_Ex1[threadId] * stepInMM.GetNorm();
153 m_Ex1[threadId] = ex2;
154 m_AttenuationPixel[threadId] = 0;
155 sumValue += wf * volumeValue;
161 m_AttenuationRay = attenuationRayVector;
166 m_AttenuationPixel = attenuationPixelVector;
187 template <
class TInput,
class TOutput>
204 return !(*
this != other);
209 const TInput & input,
211 const TOutput & rayCastValue,
218 output = input + rayCastValue;
219 m_Attenuation[threadId] = 0;
226 m_Attenuation = attenuationVector;
259 typename TInputImage::PixelType,
261 class TProjectedValueAccumulation =
268 TInterpolationWeightMultiplication,
269 TProjectedValueAccumulation,
279 TInterpolationWeightMultiplication,
280 TProjectedValueAccumulation,
291 static constexpr
unsigned int InputImageDimension = TInputImage::ImageDimension;
297 #ifdef itkOverrideGetNameOfClassMacro 310 GenerateInputRequestedRegion()
override;
313 BeforeThreadedGenerateData()
override;
318 VerifyInputInformation()
const override;
322 #ifndef ITK_MANUAL_INSTANTIATION 323 # include "rtkJosephForwardAttenuatedProjectionImageFilter.hxx"
typename TOutputImage::PixelType OutputPixelType
typename TInputImage::PixelType InputPixelType
void SetAttenuationVector(TInput *attenuationVector)
InterpolationWeightMultiplicationAttenuated()
constexpr vcl_size_t ITK_MAX_THREADS
Joseph forward projection.
Function to accumulate the ray casting on the projection.
TOutput * GetAttenuationPixel()
bool operator!=(const ProjectedValueAccumulationAttenuated &) const
typename TPixelType::ValueType ValueType
TInput * m_AttenuationRay
TOutput * GetAttenuationRay()
typename OutputImageType::RegionType OutputImageRegionType
bool operator==(const ProjectedValueAccumulationAttenuated &other) const
bool operator!=(const InterpolationWeightMultiplicationAttenuated &) const
void operator()(const ThreadIdType threadId, TOutput &sumValue, const TInput volumeValue, const VectorType &stepInMM)
void SetAttenuationPixelVector(TInput *attenuationPixelVector)
std::ptrdiff_t m_AttenuationMinusEmissionMapsPtrDiff
bool operator==(const ComputeAttenuationCorrection &other) const
TInput * m_AttenuationPixel
void SetAttenuationMinusEmissionMapsPtrDiff(std::ptrdiff_t pd)
Joseph forward projection.
unsigned int ThreadIdType
TOutput operator()(const ThreadIdType threadId, const double stepLengthInVoxel, const TCoordRepType weight, const TInput *p, const int i)
void SetAttenuationRayVector(TInput *attenuationRayVector)
Function to multiply the interpolation weights with the projected volume values and attenuation map...
void operator()(const ThreadIdType threadId, const TInput &input, TOutput &output, const TOutput &rayCastValue, const VectorType &, const VectorType &, const VectorType &, const VectorType &, const VectorType &)
bool operator==(const InterpolationWeightMultiplicationAttenuated &other) const
Function to compute the attenuation correction on the projection.
bool operator!=(const ComputeAttenuationCorrection &) const