-
Notifications
You must be signed in to change notification settings - Fork 6
/
pctpairprotonsLomaLinda.cxx
259 lines (218 loc) · 8.45 KB
/
pctpairprotonsLomaLinda.cxx
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
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
#include "pctpairprotonsLomaLinda_ggo.h"
#include <iomanip>
#include <random>
#include <rtkGgoFunctions.h>
#include <itkImage.h>
#include <itkImageIterator.h>
#include <itkImageFileWriter.h>
#include <itksys/SystemTools.hxx>
// Root includes
#include <TChain.h>
#include <TROOT.h>
#include <TVector3.h>
#define MAX_RUNS 4096
struct ParticleData
{
float wepl;
itk::Vector<float,3> position0;
itk::Vector<float,3> position1;
itk::Vector<float,3> position2;
itk::Vector<float,3> position3;
};
struct ParticleDataFinal
{
float wepl;
itk::Vector<float,3> position;
itk::Vector<float,3> direction;
float time;
};
struct ParticleInfo
{
int runID;
};
bool SetTreeBranch(TChain *tree, std::string branchName, void *add, bool mandatory=true)
{
unsigned int found = 0;
tree->SetBranchStatus(branchName.c_str(), 1, &found);
if(!found)
{
if(mandatory)
{
std::cerr << "Could not load branch "
<< branchName << std::endl;
exit(EXIT_FAILURE);
}
}
else
tree->SetBranchAddress(branchName.c_str(), add);
return found;
}
void BranchParticleToPhaseSpace(struct ParticleInfo &piInput, struct ParticleData &pdInput, TChain *tree)
{
std::cout << "In BranchParticleToPhaseSpace" << std::endl;
tree->GetListOfBranches(); // force reading of chain
SetTreeBranch(tree, "projection_angle", &piInput.runID);;
SetTreeBranch(tree, "calculated_WEPL", &pdInput.wepl);
// WARNING: X and Z are purposely swap...
SetTreeBranch(tree, "v_hit0", pdInput.position0.GetDataPointer()+1);
SetTreeBranch(tree, "t_hit0", pdInput.position0.GetDataPointer());
SetTreeBranch(tree, "u_hit0", pdInput.position0.GetDataPointer()+2);
SetTreeBranch(tree, "v_hit1", pdInput.position1.GetDataPointer()+1);
SetTreeBranch(tree, "t_hit1", pdInput.position1.GetDataPointer());
SetTreeBranch(tree, "u_hit1", pdInput.position1.GetDataPointer()+2);
SetTreeBranch(tree, "v_hit2", pdInput.position2.GetDataPointer()+1);
SetTreeBranch(tree, "t_hit2", pdInput.position2.GetDataPointer());
SetTreeBranch(tree, "u_hit2", pdInput.position2.GetDataPointer()+2);
SetTreeBranch(tree, "v_hit3", pdInput.position3.GetDataPointer()+1);
SetTreeBranch(tree, "t_hit3", pdInput.position3.GetDataPointer());
SetTreeBranch(tree, "u_hit3", pdInput.position3.GetDataPointer()+2);
}
void WritePairs(const std::vector< std::pair<ParticleDataFinal,ParticleDataFinal> > &pairs, std::string fileName)
{
itk::ImageRegion<2> region;
itk::ImageRegion<2>::SizeType size;
size[0] = 5;
size[1] = pairs.size();
region.SetSize(size);
typedef itk::Vector<float,3> PixelType;
typedef itk::Image<PixelType,2> ImageType;
ImageType::Pointer img = ImageType::New();
img->SetRegions(region);
img->Allocate();
itk::ImageRegionIterator<ImageType> it(img, region);
PixelType eet;
for(size_t i=0; i<pairs.size(); i++)
{
eet[0] = pairs[i].first.wepl;
eet[1] = pairs[i].second.wepl;
eet[2] = pairs[i].second.time - pairs[i].first.time;
it.Set( pairs[i].first.position );
++it;
it.Set( pairs[i].second.position );
++it;
it.Set( pairs[i].first.direction );
++it;
it.Set( pairs[i].second.direction );
++it;
it.Set( eet );
++it;
}
// Write
typedef itk::ImageFileWriter< ImageType > WriterType;
WriterType::Pointer writer = WriterType::New();
writer->SetFileName( fileName );
writer->SetInput( img );
TRY_AND_EXIT_ON_ITK_EXCEPTION( writer->Update() );
}
int main(int argc, char * argv[])
{
GGO(pctpairprotonsLomaLinda, args_info); //RTK macro parsing options from .ggo file (rtkMacro.h)
// Create root trees
TChain *treeIn = new TChain("recoENTRY");
treeIn->AddFile(args_info.input_arg);
std::cout << "Reading in file:" << args_info.input_arg << std::endl;
if(args_info.fmpct_flag) std::cout << "In fmpct mode with ROI radius = " << args_info.roiR_arg << " mm" << std::endl;
// Branch particles
struct ParticleInfo pi;
struct ParticleData pd;
struct ParticleDataFinal pdIn;
struct ParticleDataFinal pdOut;
BranchParticleToPhaseSpace(pi, pd, treeIn);
// Init
std::vector< std::vector< std::pair<ParticleDataFinal, ParticleDataFinal> > > pairs(MAX_RUNS);
size_t nparticulesIn = treeIn->GetEntries();
std::cout << "Number of entries = " << nparticulesIn << std::endl;
size_t iIn=0;
size_t counterPairs=0;
std::cout << iIn << " particles of input phase space processed ("
<< 100*iIn/nparticulesIn << "%)" << std::flush;
// Go over root files
while(iIn<nparticulesIn)
{
if(iIn%10000==0)
std::cout << '\r'
<< iIn << " particles of input phase space processed ("
<< 100*iIn/nparticulesIn << "%)"
<< std::flush;
treeIn->GetEntry(iIn);
#if 0
std::cout << "Input Entry = " << iIn << " with data : "
<< pd.position0[0] << " " << pd.position0[1] << " " << pd.position0[2] << " "
<< pd.position1[0] << " " << pd.position1[1] << " " << pd.position1[2] << " "
<< pd.position2[0] << " " << pd.position2[1] << " " << pd.position2[2] << " "
<< pd.position3[0] << " " << pd.position3[1] << " " << pd.position3[2] << " "
<< " " << pd.wepl << std::endl;
#endif
// Part for proton track segment - corcular ROI intersection
// http://mathworld.wolfram.com/Circle-LineIntersection.html
bool interceptionFlag = true;
if(args_info.fmpct_flag){
double radius = args_info.roiR_arg; // ROI radius in mm
double dx = pd.position2[2]-pd.position1[2];
double dy = pd.position2[0]-pd.position1[0];
double dr_sq = (dx*dx) + (dy*dy);
double D = (pd.position1[2]*pd.position2[0]) - (pd.position2[2]*pd.position1[0]);
double Delta = (radius*radius*dr_sq) - (D*D);
double randomNumber = ((double) rand() / (RAND_MAX));
//std::cout << "RAND_MAX = " << RAND_MAX << std::endl;
//std::cout << "randomNumber = " << randomNumber << std::endl;
// no interception, then discard proton with randomNumber > modF_arg
if (Delta<0) {
if (randomNumber>args_info.modF_arg) interceptionFlag = false;
}
//std::cout << "Delta, randomNumber, interceptionFlag = " << Delta << " " << randomNumber << " " << interceptionFlag <<std::endl;
}
// Convert info from ParticleData to info for ParticleDataFinal
// The inner tracker planes are used
pdIn.position = pd.position1;
pdOut.position = pd.position2;
pdIn.direction = pd.position1 - pd.position0;
pdIn.direction /= pdIn.direction.GetNorm();
pdOut.direction = pd.position3 - pd.position2;
pdOut.direction /= pdOut.direction.GetNorm();
pdIn.wepl = 0.;
pdOut.wepl = pd.wepl;
pdIn.time = 0.;
pdOut.time = 0.;
#if 0
std::cout << "Output Entry = " << iIn << " with data : "
<< pdIn.position[0] << " " << pdIn.position[1] << " " << pdIn.position[2] << " "
<< pdOut.position[0] << " " << pdOut.position[1] << " " << pdOut.position[2] << " "
<< pdIn.direction[0] << " " << pdIn.direction[1] << " " << pdIn.direction[2] << " "
<< pdOut.direction[0] << " " << pdOut.direction[1] << " " << pdOut.direction[2] << " "
<< " " << pdIn.wepl << " " << pdOut.wepl << std::endl;
#endif
if( pdOut.wepl>=-20. && pdOut.wepl<=900. && interceptionFlag)
{
pairs[args_info.runID_arg].push_back( std::pair<ParticleDataFinal,ParticleDataFinal>(pdIn, pdOut) );
/* std::cout << "Test = " << iIn << " with data : "
<< pdIn.position[0] << " " << pdIn.position[1] << " " << pdIn.position[2] << " "
<< pdOut.position[0] << " " << pdOut.position[1] << " " << pdOut.position[2] << " "
<< pdIn.direction[0] << " " << pdIn.direction[1] << " " << pdIn.direction[2] << " "
<< pdOut.direction[0] << " " << pdOut.direction[1] << " " << pdOut.direction[2] << " "
<< " " << pdIn.wepl << " " << pdOut.wepl << std::endl; */
/*if (pdIn.position[0]>-45. && pdIn.position[0]<-10. && pdIn.position[1]>-20. && pdIn.position[1]<30.)*/ counterPairs++;
}
iIn++;
}
std::cout << "\r"
<< nparticulesIn << " particles of input phase space processed ("
<< 100 << "%)"
<< std::endl
<< "Writing..."
<< std::endl;
std::cout << "number of accepted pairs = " << counterPairs << std::endl;
for(unsigned int i=0; i<MAX_RUNS; i++)
{
if(pairs[i].size())
{
std::ostringstream os;
os << itksys::SystemTools::GetFilenameWithoutLastExtension(args_info.output_arg)
<< std::setw(4) << std::setfill ('0') << i
<< itksys::SystemTools::GetFilenameLastExtension(args_info.output_arg);
std::cout << "Writing into file:" << os.str() << std::endl;
WritePairs(pairs[i], os.str());
}
}
return EXIT_SUCCESS;
}