forked from SimonRit/PCT
-
Notifications
You must be signed in to change notification settings - Fork 0
/
pctProtonPairsToBackProjection.h
148 lines (111 loc) · 4.88 KB
/
pctProtonPairsToBackProjection.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
#ifndef __pctProtonPairsToBackProjection_h
#define __pctProtonPairsToBackProjection_h
#include "rtkConfiguration.h"
#include "pctBetheBlochFunctor.h"
#include <rtkQuadricShape.h>
#include <rtkThreeDCircularProjectionGeometry.h>
#include <itkInPlaceImageFilter.h>
#include <mutex>
namespace pct
{
template <class TInputImage, class TOutputImage>
class ITK_EXPORT ProtonPairsToBackProjection :
public itk::InPlaceImageFilter<TInputImage,TOutputImage>
{
public:
/** Standard class typedefs. */
typedef ProtonPairsToBackProjection Self;
typedef itk::InPlaceImageFilter<TInputImage,TOutputImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
typedef std::vector<std::string> FileNamesContainer;
typedef itk::Vector<float, 3> ProtonPairsPixelType;
typedef itk::Image<ProtonPairsPixelType,2> ProtonPairsImageType;
typedef ProtonPairsImageType::Pointer ProtonPairsImagePointer;
typedef typename itk::Image< unsigned int,
TInputImage::ImageDimension> CountImageType;
typedef typename CountImageType::Pointer CountImagePointer;
typedef TOutputImage OutputImageType;
typedef typename OutputImageType::Pointer OutputImagePointer;
typedef typename OutputImageType::RegionType OutputImageRegionType;
typedef rtk::QuadricShape RQIType;
typedef rtk::ThreeDCircularProjectionGeometry GeometryType;
typedef typename GeometryType::Pointer GeometryPointer;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Run-time type information (and related methods). */
itkTypeMacro(ProtonPairsToBackProjection, itk::InPlaceImageFilter);
/** Set the vector of strings that contains the file names. Files
* are processed in sequential order. */
void SetProtonPairsFileNames (const FileNamesContainer &name)
{
if ( m_ProtonPairsFileNames != name)
{
m_ProtonPairsFileNames = name;
this->Modified();
}
}
const FileNamesContainer & GetProtonPairsFileNames() const
{
return m_ProtonPairsFileNames;
}
/** Get/Set the most likely path type. Can be "schulte" or "polynomial" */
itkGetMacro(MostLikelyPathType, std::string);
itkSetMacro(MostLikelyPathType, std::string);
itkGetMacro(MostLikelyPathPolynomialDegree, int);
itkSetMacro(MostLikelyPathPolynomialDegree, int);
/** Get/Set the boundaries of the object. */
itkGetMacro(QuadricIn, RQIType::Pointer);
itkSetMacro(QuadricIn, RQIType::Pointer);
itkGetMacro(QuadricOut, RQIType::Pointer);
itkSetMacro(QuadricOut, RQIType::Pointer);
/** Get/Set the count of proton pairs per pixel. */
itkGetMacro(Counts, CountImagePointer);
itkSetMacro(Counts, CountImagePointer);
/** Get/Set the ionization potential used in the Bethe-Bloch equation. */
itkGetMacro(IonizationPotential, double);
itkSetMacro(IonizationPotential, double);
/** Get / Set the object pointer to projection geometry */
itkGetMacro(Geometry, GeometryPointer);
itkSetMacro(Geometry, GeometryPointer);
/** Get / Set the boolean desactivating rotation to bin in coordinate orientation. Default is off. */
itkGetMacro(DisableRotation, bool);
itkSetMacro(DisableRotation, bool);
itkBooleanMacro(DisableRotation);
protected:
ProtonPairsToBackProjection();
virtual ~ProtonPairsToBackProjection() {}
virtual void BeforeThreadedGenerateData() ITK_OVERRIDE;
virtual void GenerateData() ITK_OVERRIDE;
virtual void AfterThreadedGenerateData() ITK_OVERRIDE;
/** The two inputs should not be in the same space so there is nothing
* to verify. */
virtual void VerifyInputInformation() const ITK_OVERRIDE {}
private:
ProtonPairsToBackProjection(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
/** A list of filenames to be processed. */
FileNamesContainer m_ProtonPairsFileNames;
std::string m_MostLikelyPathType;
int m_MostLikelyPathPolynomialDegree;
/** Count event in each thread */
CountImagePointer m_Counts;
/** The two quadric functions defining the object support. */
RQIType::Pointer m_QuadricIn;
RQIType::Pointer m_QuadricOut;
/** Ionization potential used in the Bethe Bloch equation */
double m_IonizationPotential;
/** The functor to convert energy loss to attenuation */
Functor::IntegratedBetheBlochProtonStoppingPowerInverse<float, double> *m_ConvFunc;
ProtonPairsImageType::Pointer m_ProtonPairs;
/** RTK geometry object */
GeometryPointer m_Geometry;
/** Disable rotation to bin in coordinate orientation. Default is off. */
bool m_DisableRotation = false;
std::mutex m_Mutex;
};
} // end namespace pct
#ifndef ITK_MANUAL_INSTANTIATION
#include "pctProtonPairsToBackProjection.txx"
#endif
#endif