CASToR  1.1
Tomographic Reconstruction (PET/SPECT)
 All Classes Files Functions Variables Typedefs Macros Groups Pages
oImageSpace.cc
Go to the documentation of this file.
1 
2 /*
3  Implementation of class oImageSpace
4 
5  - separators: X
6  - doxygen: X
7  - default initialization: X
8  - CASTOR_DEBUG: X
9  - CASTOR_VERBOSE: none
10 */
11 
12 
19 #include "oImageSpace.hh"
20 #include "vOptimizer.hh"
21 
22 
23 // =====================================================================
24 // ---------------------------------------------------------------------
25 // ---------------------------------------------------------------------
26 // =====================================================================
27 
29 {
30  // Initialize the member variables to their default values
31  m_loadedSensitivity = false;
32  m_loadedAnatomical = false;
33  m_loadedMask = false;
34  mp_ID = NULL;
35  m_verbose = -1;
36  m_nbBackwardImages = 1;
37  // Set all images pointers to NULL
38  m4p_image = NULL;
39  m4p_forwardImage = NULL;
40  m6p_backwardImage = NULL;
41  m5p_sensitivity = NULL;
42  m4p_attenuation = NULL;
45  m4p_outputImage = NULL;
46  mp_visitedVoxelsImage = NULL;
47  m2p_projectionImage = NULL;
48  m4p_anatomicalImage = NULL;
49  mp_maskImage = NULL;
50 }
51 
52 // =====================================================================
53 // ---------------------------------------------------------------------
54 // ---------------------------------------------------------------------
55 // =====================================================================
56 
58 
59 // =====================================================================
60 // ---------------------------------------------------------------------
61 // ---------------------------------------------------------------------
62 // =====================================================================
63 
65 {
66  // Verbose
67  if (m_verbose>=3) Cout("oImageSpace::InstantiateImage() -> Initialize to 0."<< endl);
68 
69  // For each temporal dimensions, the number of related images is recovered from the number of
70  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
71  // are initialized (standard reconstruction), this is equivalent to the standard calls to
72  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
73 
74  // First pointer is for the number of time basis functions
75  m4p_image = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
76  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
77  {
78  // Second pointer is for the number of respiratory basis functions
79  m4p_image[tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
80  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
81  {
82  // Third pointer is for the number of cardiac basis functions
83  m4p_image[tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
84  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
85  {
86  // Fourth pointer is for the 3D space
87  m4p_image[tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
88  // Initialize to 0.
89  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++) m4p_image[tbf][rbf][cbf][v] = 0.;
90  }
91  }
92  }
93 }
94 
95 // =====================================================================
96 // ---------------------------------------------------------------------
97 // ---------------------------------------------------------------------
98 // =====================================================================
99 
101 {
102  // Check the first pointer
103  if (m4p_image)
104  {
105  // Verbose
106  if (m_verbose>=3) Cout("oImageSpace::DeallocateImage() -> Free memory"<< endl);
107  // Loop on time basis functions
108  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++) if (m4p_image[tbf])
109  {
110  // Loop on respiratory basis functions
111  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++) if (m4p_image[tbf][rbf])
112  {
113  // Loop on cardiac basis functions
114  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++) if (m4p_image[tbf][rbf][cbf])
115  {
116  free(m4p_image[tbf][rbf][cbf]);
117  }
118  free(m4p_image[tbf][rbf]);
119  }
120  free(m4p_image[tbf]);
121  }
122  free(m4p_image);
123  }
124  // Reset the pointer to NULL
125  m4p_image = NULL;
126 }
127 
128 // =====================================================================
129 // ---------------------------------------------------------------------
130 // ---------------------------------------------------------------------
131 // =====================================================================
132 
134 {
135  // Verbose
136  if (m_verbose>=3) Cout("oImageSpace::InstantiateForwardImage() -> Initialize to 0."<< endl);
137 
138  // For each temporal dimensions, the number of related images is recovered from the number of
139  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
140  // are initialized (standard reconstruction), this is equivalent to the standard calls to
141  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
142 
143  // First pointer is for the number of time basis functions
144  m4p_forwardImage = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
145  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
146  {
147  // Second pointer is for the number of respiratory basis functions
148  m4p_forwardImage[tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
149  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
150  {
151  // Third pointer is for the number of cardiac basis functions
152  m4p_forwardImage[tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
153  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
154  {
155  // Fourth pointer is for the 3D space
156  m4p_forwardImage[tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
157  // Initialize to 0.
158  for (INTNB v=0; v<mp_ID->GetNbVoxXYZ(); v++) m4p_forwardImage[tbf][rbf][cbf][v] = 0.;
159  }
160  }
161  }
162 }
163 
164 // =====================================================================
165 // ---------------------------------------------------------------------
166 // ---------------------------------------------------------------------
167 // =====================================================================
168 
170 {
171  // Check the first pointer
172  if (m4p_forwardImage)
173  {
174  // Verbose
175  if (m_verbose>=3) Cout("oImageSpace::DeallocateForwardImage() -> Free memory"<< endl);
176  // Loop on time basis functions
177  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++) if (m4p_forwardImage[tbf])
178  {
179  // Loop on respiratory basis functions
180  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++) if (m4p_forwardImage[tbf][rbf])
181  {
182  // Loop on cardiac basis functions
183  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++) if (m4p_forwardImage[tbf][rbf][cbf])
184  {
185  free(m4p_forwardImage[tbf][rbf][cbf]);
186  }
187  free(m4p_forwardImage[tbf][rbf]);
188  }
189  free(m4p_forwardImage[tbf]);
190  }
191  free(m4p_forwardImage);
192  }
193  // Reset the pointer to NULL
194  m4p_forwardImage = NULL;
195 }
196 
197 // =====================================================================
198 // ---------------------------------------------------------------------
199 // ---------------------------------------------------------------------
200 // =====================================================================
201 
203 {
204  // Verbose
205  if (m_verbose>=3) Cout("oImageSpace::InstantiateBackwardImageFromDynamicBasis() -> Initialize to 0."<< endl);
206 
207  // For each temporal dimensions, the number of related images is recovered from the number of
208  // intrinsic time basis functions stored in mp_ID. If no intrinsic temporal basis functions
209  // are initialized (standard reconstruction), this is equivalent to the standard calls to
210  // GetNbTimeFrames() / GetNbRespGates() / GetNbCardGates()
211 
212  // First pointer is for the number of backward images (in case the optimizer does need more than 1)
213  m_nbBackwardImages = a_nbBackwardImages;
214  m6p_backwardImage = (FLTNB******)malloc(m_nbBackwardImages*sizeof(FLTNB*****));
215  for (int img=0; img<m_nbBackwardImages; img++)
216  {
217  // Second pointer is for the number of threads, in order to have thread-safe backward projections
218  m6p_backwardImage[img] = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
219  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
220  {
221  // The third pointer is for the number of time basis functions
222  m6p_backwardImage[img][th] = (FLTNB****)malloc(mp_ID->GetNbTimeBasisFunctions()*sizeof(FLTNB***));
223  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
224  {
225  // The fourth pointer is for the number of respiratory basis functions
226  m6p_backwardImage[img][th][tbf] = (FLTNB***)malloc(mp_ID->GetNbRespBasisFunctions()*sizeof(FLTNB**));
227  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
228  {
229  // The fifth pointer is for the number of cardiac basis functions
230  m6p_backwardImage[img][th][tbf][rbf] = (FLTNB**)malloc(mp_ID->GetNbCardBasisFunctions()*sizeof(FLTNB*));
231  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
232  {
233  // The sixth pointer is for the 3D space
234  m6p_backwardImage[img][th][tbf][rbf][cbf] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
235  }
236  }
237  }
238  }
239  }
240 }
241 
242 // =====================================================================
243 // ---------------------------------------------------------------------
244 // ---------------------------------------------------------------------
245 // =====================================================================
246 
248 {
249  // Check the first pointer
250  if (m6p_backwardImage)
251  {
252  // Verbose
253  if (m_verbose>=3) Cout("oImageSpace::DeallocateBackwardImageFromDynamicBasis() -> Free memory" << endl);
254  // Loop on backward images (in case the optimizer does need more than 1)
255  for (int img=0; img<m_nbBackwardImages; img++) if (m6p_backwardImage[img]!=NULL)
256  {
257  // Loop on threads
258  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++) if (m6p_backwardImage[img][th])
259  {
260  // Loop on time basis functions
261  for (int tbf=0 ; tbf<mp_ID->GetNbTimeBasisFunctions() ; tbf++) if (m6p_backwardImage[img][th][tbf])
262  {
263  // Loop on respiratory basis functions
264  for (int rbf=0 ; rbf<mp_ID->GetNbRespBasisFunctions() ; rbf++) if (m6p_backwardImage[img][th][tbf][rbf])
265  {
266  // Loop on cardiac basis functions
267  for (int cbf=0 ; cbf<mp_ID->GetNbCardBasisFunctions() ; cbf++) if (m6p_backwardImage[img][th][tbf][rbf][cbf])
268  {
269  free(m6p_backwardImage[img][th][tbf][rbf][cbf]);
270  }
271  free(m6p_backwardImage[img][th][tbf][rbf]);
272  }
273  free(m6p_backwardImage[img][th][tbf]);
274  }
275  free(m6p_backwardImage[img][th]);
276  }
277  free(m6p_backwardImage[img]);
278  }
279  free(m6p_backwardImage);
280  }
281  // Reset the pointer to NULL
282  m6p_backwardImage = NULL;
283 }
284 
285 // =====================================================================
286 // ---------------------------------------------------------------------
287 // ---------------------------------------------------------------------
288 // =====================================================================
289 
291 {
292  // Verbose
293  if (m_verbose>=3) Cout("oImageSpace::InstantiateBackwardImageFromDynamicBins() -> Initialize to 0."<< endl);
294  // Not related to an optimizer, so we consider only one image
295  m6p_backwardImage = (FLTNB******)malloc(1*sizeof(FLTNB*****));
296  // Projection threads
297  m6p_backwardImage[0] = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
298  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
299  {
300  // Time frames
301  m6p_backwardImage[0][th] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
302  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
303  {
304  // Respiratory gates
305  m6p_backwardImage[0][th][fr] = (FLTNB***)malloc(mp_ID->GetSensNbRespGates()*sizeof(FLTNB**));
306  for (int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
307  {
308  // Cardiac gates
309  m6p_backwardImage[0][th][fr][rg] = (FLTNB**)malloc(mp_ID->GetSensNbCardGates()*sizeof(FLTNB*));
310  for (int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
311  {
312  // Number of voxels
313  m6p_backwardImage[0][th][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
314  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++) m6p_backwardImage[0][th][fr][rg][cg][v] = 0.;
315  }
316  }
317  }
318  }
319 }
320 
321 // =====================================================================
322 // ---------------------------------------------------------------------
323 // ---------------------------------------------------------------------
324 // =====================================================================
325 
327 {
328  // Check the first pointer
329  if (m6p_backwardImage)
330  {
331  // Verbose
332  if (m_verbose>=3) Cout("oImageSpace::DeallocateBackwardImageFromDynamicBins() -> Free memory"<< endl);
333  // We assumed only one image when considering dynamic bins
334  if (m6p_backwardImage[0]!=NULL)
335  {
336  // Loop on threads
337  for (int th=0 ; th<mp_ID->GetNbThreadsForProjection() ; th++) if (m6p_backwardImage[0][th])
338  {
339  // Loop on time frames
340  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++) if (m6p_backwardImage[0][th][fr])
341  {
342  // Loop on respiratory gates
343  for (int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++) if (m6p_backwardImage[0][th][fr][rg])
344  {
345  // Loop on cardiac gates
346  for (int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++) if (m6p_backwardImage[0][th][fr][rg][cg])
347  {
348  free(m6p_backwardImage[0][th][fr][rg][cg]);
349  }
350  free(m6p_backwardImage[0][th][fr][rg]);
351  }
352  free(m6p_backwardImage[0][th][fr]);
353  }
354  free(m6p_backwardImage[0][th]);
355  }
356  free(m6p_backwardImage[0]);
357  }
358  free(m6p_backwardImage);
359  }
360  // Reset the pointer to NULL
361  m6p_backwardImage = NULL;
362 }
363 
364 // =====================================================================
365 // ---------------------------------------------------------------------
366 // ---------------------------------------------------------------------
367 // =====================================================================
368 
369 void oImageSpace::InstantiateSensitivityImage(const string& a_pathToSensitivityImage)
370 {
371  // -----------------------------------------------------------------------------------------------------
372  // Case 1: a path to a sensitivity image is provided, meaning that we will not compute the sensitivity
373  // on-the-fly. In that case, we do not take the multi-threading into account. Also note that in
374  // histogram mode, the same sensitivity image will be used in optimization algorithms for all
375  // subsets.
376  // -----------------------------------------------------------------------------------------------------
377 
378  if (!a_pathToSensitivityImage.empty())
379  {
380  // Verbose
381  if (m_verbose>=3) Cout("oImageSpace::InstantiateSensitivityImage() -> Will be loaded from '" << a_pathToSensitivityImage << "'" << endl);
382  // Allocate for all threads first
383  m5p_sensitivity = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
384  // Allocate the rest only for the first thread
385  int thread_0 = 0;
386  m5p_sensitivity[thread_0] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
387  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
388  {
389  m5p_sensitivity[thread_0][fr] = (FLTNB***)malloc(mp_ID->GetNbRespGates()*sizeof(FLTNB**));
390  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
391  {
392  m5p_sensitivity[thread_0][fr][rg] = (FLTNB**)malloc(mp_ID->GetNbCardGates()*sizeof(FLTNB*));
393  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
394  m5p_sensitivity[thread_0][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
395  }
396  }
397  // Make all thread pointers pointing to the first
398  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
399  m5p_sensitivity[th] = m5p_sensitivity[thread_0];
400  // Set the flag for loaded sensitivity to true
401  m_loadedSensitivity = true;
402  }
403 
404  // -----------------------------------------------------------------------------------------------------
405  // Case 2: no sensitivity provided, meaning we will commpute it on-the-fly. We thus need it to be
406  // thread-safe, so we allocate for all threads.
407  // -----------------------------------------------------------------------------------------------------
408 
409  else
410  {
411  // Verbose
412  if (m_verbose>=3) Cout("oImageSpace::InstantiateSensitivityImage() -> For all threads"<< endl);
413  // Allocate for all threads first
414  m5p_sensitivity = (FLTNB*****)malloc(mp_ID->GetNbThreadsForProjection()*sizeof(FLTNB****));
415  // Allocate the rest only for all threads
416  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
417  {
418  m5p_sensitivity[th] = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
419  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
420  {
421  m5p_sensitivity[th][fr] = (FLTNB***)malloc(mp_ID->GetNbRespGates()*sizeof(FLTNB**));
422  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
423  {
424  m5p_sensitivity[th][fr][rg] = (FLTNB**)malloc(mp_ID->GetNbCardGates()*sizeof(FLTNB*));
425  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
426  {
427  m5p_sensitivity[th][fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
428  // Initialize to 0.
429  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++) m5p_sensitivity[th][fr][rg][cg][v] = 0.;
430  }
431  }
432  }
433  }
434  // Set the flag for loaded sensitivity to false
435  m_loadedSensitivity = false;
436  }
437 }
438 
439 // =====================================================================
440 // ---------------------------------------------------------------------
441 // ---------------------------------------------------------------------
442 // =====================================================================
443 
445 {
446  // Check the first pointer
447  if (m5p_sensitivity)
448  {
449  // Verbose
450  if (m_verbose>=3) Cout("oImageSpace::DeallocateSensitivityImage() -> Free memory"<< endl);
451  // If loaded sensitivity then fully deallocate only the first thread pointer
453  {
454  int thread_0 = 0;
455  if (m5p_sensitivity[thread_0])
456  {
457  // Loop on time frames
458  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++) if (m5p_sensitivity[thread_0][fr])
459  {
460  // Loop on respiratory gates
461  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++) if (m5p_sensitivity[thread_0][fr][rg])
462  {
463  // Loop on cardiac gates
464  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++) if (m5p_sensitivity[thread_0][fr][rg][cg])
465  {
466  free(m5p_sensitivity[thread_0][fr][rg][cg]);
467  }
468  free(m5p_sensitivity[thread_0][fr][rg]);
469  }
470  free(m5p_sensitivity[thread_0][fr]);
471  }
472  free(m5p_sensitivity[thread_0]);
473  }
474  }
475  // Otherwise, loop on all threads
476  else
477  {
478  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++) if (m5p_sensitivity[th])
479  {
480  // Loop on time frames
481  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++) if (m5p_sensitivity[th][fr])
482  {
483  // Loop on respiratory gates
484  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++) if (m5p_sensitivity[th][fr][rg])
485  {
486  // Loop on cardiac gates
487  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++) if (m5p_sensitivity[th][fr][rg][cg])
488  {
489  free(m5p_sensitivity[th][fr][rg][cg]);
490  }
491  free(m5p_sensitivity[th][fr][rg]);
492  }
493  free(m5p_sensitivity[th][fr]);
494  }
495  free(m5p_sensitivity[th]);
496  }
497  }
498  // Free the first pointer
499  free(m5p_sensitivity);
500  }
501  // And set the pointer to NULL
502  m5p_sensitivity = NULL;
503 }
504 
505 // =====================================================================
506 // ---------------------------------------------------------------------
507 // ---------------------------------------------------------------------
508 // =====================================================================
509 
510 int oImageSpace::InitAnatomicalImage(const string& a_pathToImage)
511 {
512  // Allocate and read only if a file is provided
513  if (!a_pathToImage.empty())
514  {
515  // Verbose
516  if (m_verbose>=3) Cout("oImageSpace::InitAnatomicalImage() -> From file '" << a_pathToImage << "'"<< endl);
517  // Open file
518  ifstream image_file(a_pathToImage.c_str(), ios::in|ios::binary);
519  if (!image_file.is_open())
520  {
521  Cerr("***** oImageSpace::InitAnatomicalImage() -> Input file '" << a_pathToImage << "' is missing or corrupted !" << endl);
522  return 1;
523  }
524  // Allocate memory (if not already done)
525  if (!m_loadedAnatomical)
526  {
527  m4p_anatomicalImage = (FLTNB****)malloc(mp_ID->GetNbTimeFrames()*sizeof(FLTNB***));
528  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
529  {
530  m4p_anatomicalImage[fr] = (FLTNB***)malloc(mp_ID->GetNbRespGates()*sizeof(FLTNB**));
531  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
532  {
533  m4p_anatomicalImage[fr][rg] = (FLTNB**)malloc(mp_ID->GetNbCardGates()*sizeof(FLTNB*));
534  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
535  {
536  m4p_anatomicalImage[fr][rg][cg] = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
537  }
538  }
539  }
540  }
541  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
543  {
544  Cerr("***** oImageSpace::InitAnatomicalImage() -> Error reading interfile image '" << a_pathToImage << "' !" << endl);
545  // Read failure implies that the anatomical image will never be used, so free all the allocated memory
547  return 1;
548  }
549  // Flag as loaded
550  m_loadedAnatomical = true;
551  // Close file
552  image_file.close();
553  }
554  // End
555  return 0;
556 }
557 
558 // =====================================================================
559 // ---------------------------------------------------------------------
560 // ---------------------------------------------------------------------
561 // =====================================================================
562 
564 {
565  // Check the first pointer
567  {
568  // Verbose
569  if (m_verbose>=3) Cout("oImageSpace::DeallocateAnatomicalImage() -> Free memory" << endl);
570  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++) if (m4p_anatomicalImage[fr])
571  {
572  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++) if (m4p_anatomicalImage[fr][rg])
573  {
574  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++) if (m4p_anatomicalImage[fr][rg][cg])
575  {
576  free(m4p_anatomicalImage[fr][rg][cg]);
577  }
578  free(m4p_anatomicalImage[fr][rg]);
579  }
580  free(m4p_anatomicalImage[fr]);
581  }
582  free(m4p_anatomicalImage);
583  }
584  // Reset the pointer to NULL
585  m4p_anatomicalImage = NULL;
586 }
587 
588 // =====================================================================
589 // ---------------------------------------------------------------------
590 // ---------------------------------------------------------------------
591 // =====================================================================
592 
593 int oImageSpace::InitMaskImage(const string& a_pathToImage)
594 {
595  // Allocate and read only if a file is provided
596  if (!a_pathToImage.empty())
597  {
598  // Verbose
599  if (m_verbose>=3) Cout("oImageSpace::InitMaskImage() -> From file '" << a_pathToImage << "'"<< endl);
600  // Open file
601  ifstream image_file(a_pathToImage.c_str(), ios::in|ios::binary);
602  if (!image_file.is_open())
603  {
604  Cerr("***** oImageSpace::InitMaskImage() -> Input file '" << a_pathToImage << "' is missing or corrupted !" << endl);
605  return 1;
606  }
607  // Allocate memory (if not already done)
608  if (!m_loadedMask)
609  {
610  mp_maskImage = (FLTNB*)malloc(mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
611  }
612  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
614  {
615  Cerr("***** oImageSpace::InitMaskImage() -> Error reading interfile image '" << a_pathToImage << "' !" << endl);
616  // Read failure implies that the mask image will never be used, so free all the allocated memory
618  return 1;
619  }
620  // Flag as loaded
621  m_loadedMask = true;
622  // Close file
623  image_file.close();
624  }
625  // End
626  return 0;
627 }
628 
629 // =====================================================================
630 // ---------------------------------------------------------------------
631 // ---------------------------------------------------------------------
632 // =====================================================================
633 
635 {
636  // Check pointer
637  if (mp_maskImage)
638  {
639  // Verbose
640  if (m_verbose>=3) Cout("oImageSpace::DeallocateMaskImage() -> Free memory" << endl);
641  free(mp_maskImage);
642  }
643  // Set the pointer to NULL
644  mp_maskImage = NULL;
645 }
646 
647 
648 
649 
650 // =====================================================================
651 // ---------------------------------------------------------------------
652 // ---------------------------------------------------------------------
653 // =====================================================================
654 /*
655  \fn InstantiateOutputImage()
656  \brief Instanciate Image matrices dedicated to output writing on disk
657  \details Additionnal output image matrix is needed in the case the reconstruction uses intrinsic temporal basis functions
658  In this case, the image matrices are defined in the temporal image basis functions space, therefore requiring an additional step to recover the images in the regular image-space.
659  If no intrinsic temporal basis functions are used, m4p_outputImage just refers to the address of the image matrix containing the regular image.
660 */
662 {
663  if(m_verbose>=3) Cout("oImageSpace::InstantiateOutputImage ..."<< endl);
664 
665  // Here, we will allocate the first 3 pointers to the dynamic dimensions
667  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
668  {
669  m4p_outputImage[fr] = new FLTNB**[mp_ID->GetNbRespGates()];
670  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
671  {
672  m4p_outputImage[fr][rg] = new FLTNB*[mp_ID->GetNbCardGates()];
673  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
674  {
675  // For the last pointer (i.e. the voxels), we allocate them only if we strictly need them, which
676  // means in the case where one of the dynamic dimensions makes internal use of basis functions,
677  // otherwise, we will just make this output image pointing to the forward image which is useless
678  // at the moment of computing the output image.
679  if ( mp_ID->GetTimeStaticFlag() &&
682  {
683  // In this case, the time frames and basis functions represent the same thing, same for respiratory/cardiac
684  // gates and their equivalent basis functions
685  m4p_outputImage[fr][rg][cg] = m4p_forwardImage[fr][rg][cg];
686  }
687  else
688  {
689  // We allocate the output image
690  m4p_outputImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
691  }
692  }
693  }
694  }
695 }
696 
697 
698 
699 
700 // =====================================================================
701 // ---------------------------------------------------------------------
702 // ---------------------------------------------------------------------
703 // =====================================================================
704 /*
705  \fn DeallocateOutputImage()
706  \brief Free memory for the Image matrices dedicated to output writing on disk
707 */
709 {
710  if(m_verbose>=3) Cout("oImageSpace::DeallocateOutputImage ..."<< endl);
711 
712  if (m4p_outputImage)
713  {
714  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
715  {
716  if (m4p_outputImage[fr])
717  {
718  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
719  {
720  if (m4p_outputImage[fr][rg])
721  {
722  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
723  {
724  // Distinguish between intrinsic use of dynamic basis functions or not (in the latter,
725  // the m4p_outputImage points to the m4p_forwardImage, so we do not touch it)
726  if ( (!mp_ID->GetTimeStaticFlag() ||
727  !mp_ID->GetRespStaticFlag() ||
728  !mp_ID->GetCardStaticFlag()) &&
729  m4p_outputImage[fr][rg][cg] )
730  {
731  delete m4p_outputImage[fr][rg][cg];
732  }
733  }
734  delete[] m4p_outputImage[fr][rg];
735  }
736  }
737  delete[] m4p_outputImage[fr];
738  }
739  }
740  delete[] m4p_outputImage;
741  }
742  m4p_outputImage = NULL;
743 }
744 
745 
746 
747 // =====================================================================
748 // ---------------------------------------------------------------------
749 // ---------------------------------------------------------------------
750 // =====================================================================
751 /*
752  \fn InstantiateBwdImageForDeformation()
753  \brief Memory allocation for the buffer backward image required for image-based deformation
754 */
756 {
757  if(m_verbose>=3) Cout("oImageSpace::InstantiateBwdImageForDeformation ..."<< endl);
758 
760 
761  for (int img=0; img<m_nbBackwardImages; img++)
762  {
764  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
765  {
766  // The fourth pointer is for the number of respiratory basis functions
768  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
769  {
770  // The fifth pointer is for the number of cardiac basis functions
772  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
773  {
774  // The sixth pointer is for the 3D space
775  m5p_defTmpBackwardImage[img][tbf][rbf][cbf] = new FLTNB[mp_ID->GetNbVoxXYZ()];
776  }
777  }
778  }
779  }
780 
781 }
782 
783 
784 
785 
786 // =====================================================================
787 // ---------------------------------------------------------------------
788 // ---------------------------------------------------------------------
789 // =====================================================================
790 /*
791  \fn DeallocateBwdImageForDeformation()
792  \brief Free memory for the buffer backward image required for image-based deformation
793 */
795 {
796  if(m_verbose>=3) Cout("oImageSpace::DeallocateBwdImageForDeformation ..."<< endl);
797 
798  for (int img=0; img<m_nbBackwardImages; img++)
799  {
800  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
801  {
802  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
803  {
804  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
805  {
806  if (m5p_defTmpBackwardImage[img][tbf][rbf][cbf]) delete m5p_defTmpBackwardImage[img][tbf][rbf][cbf];
807  }
808  if (m5p_defTmpBackwardImage[img][tbf][rbf]) delete m5p_defTmpBackwardImage[img][tbf][rbf];
809  }
810  if (m5p_defTmpBackwardImage[img][tbf]) delete[] m5p_defTmpBackwardImage[img][tbf];
811  }
812  if (m5p_defTmpBackwardImage[img]) delete[] m5p_defTmpBackwardImage[img];
813  }
816 }
817 
818 
819 
820 
821 
822 // =====================================================================
823 // ---------------------------------------------------------------------
824 // ---------------------------------------------------------------------
825 // =====================================================================
826 /*
827  \fn InstantiateSensImageForDeformation()
828  \brief Memory allocation for the buffer sensitivity image required for image-based deformation (only for PET HISTOGRAM mode)
829 */
831 {
832  if(m_verbose>=3) Cout("oImageSpace::InstantiateSensImageForDeformation ..."<< endl);
833 
835 
836  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
837  {
839  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
840  {
842  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
843  m4p_defTmpSensitivityImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
844  }
845  }
846 }
847 
848 
849 
850 
851 // =====================================================================
852 // ---------------------------------------------------------------------
853 // ---------------------------------------------------------------------
854 // =====================================================================
855 /*
856  \fn DeallocateBwdImageForDeformation()
857  \brief Free memory for the buffer sensitivity image required for image-based deformation
858 */
860 {
861  if(m_verbose>=3) Cout("oImageSpace::DeallocateSensImageForDeformation ..."<< endl);
862 
863  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
864  {
865  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
866  {
867  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
868  if (m4p_defTmpSensitivityImage[fr][rg][cg]) delete[] m4p_defTmpSensitivityImage[fr][rg][cg];
869 
870  if (m4p_defTmpSensitivityImage[fr][rg]) delete[] m4p_defTmpSensitivityImage[fr][rg];
871  }
873  }
876 }
877 
878 
879 
880 
881 // =====================================================================
882 // ---------------------------------------------------------------------
883 // ---------------------------------------------------------------------
884 // =====================================================================
885 /*
886  \fn InstantiateVisitedVoxelsImage()
887  \brief Memory allocation and initialization for the image matrix containing binary information regarding which 3D voxels have been visited during the projection steps.
888 */
890 {
891  if(m_verbose>=3) Cout("oImageSpace::InstantiateVisitedVoxelsImage ..."<< endl);
892 
894  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
895  mp_visitedVoxelsImage[v] = 0.;
896 }
897 
898 
899 
900 
901 // =====================================================================
902 // ---------------------------------------------------------------------
903 // ---------------------------------------------------------------------
904 // =====================================================================
905 /*
906  \fn DeallocateVisitedVoxelsImage()
907  \brief Free memory for the image matrix containing binary information regarding which 3D voxels have been visited during the projection steps.
908 */
910 {
911  if(m_verbose>=3) Cout("oImageSpace::InstantiateVisitedVoxelsImage ..."<< endl);
912 
914  mp_visitedVoxelsImage = NULL;
915 }
916 
917 
918 
919 
920 // =====================================================================
921 // ---------------------------------------------------------------------
922 // ---------------------------------------------------------------------
923 // =====================================================================
924 /*
925  \fn LMS_DeallocateAttenuationImage()
926  \brief Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity generation)
927 */
929 {
930  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateAttenuationImage ..."<< endl);
931 
932  if(m4p_attenuation != NULL)
933  {
934  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
935  {
936  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates(); rg++)
937  {
938  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates(); cg++)
939  if(m4p_attenuation[fr][rg][cg] != NULL) delete m4p_attenuation[fr][rg][cg];
940 
941  if(m4p_attenuation[fr][rg] != NULL) delete m4p_attenuation[fr][rg];
942  }
943  if(m4p_attenuation[fr] != NULL) delete m4p_attenuation[fr];
944  }
945  if(m4p_attenuation != NULL) delete[] m4p_attenuation;
946  m4p_attenuation = NULL;
947  }
948 }
949 
950 
951 
952 
953 // =====================================================================
954 // ---------------------------------------------------------------------
955 // ---------------------------------------------------------------------
956 // =====================================================================
957 /*
958  \fn InitAttenuationImage()
959  \param a_pathToAtnImage : path to an existing image
960  \param a_value : value to initialize each voxel with, if an input image is not provided
961  \brief Memory allocation and initialisation for the attenuation image using either : - an image provided by the user (a_pathToAtnImage or a_pathToCTImage)
962  - the default value (=a_value)
963  \todo If fr/rg/cg > 1 and only one attenuation image is provided, initialize each redundant matrice with pointers to the first one (deal with that here)
964  \return 0 if success, positive value otherwise
965 */
966 int oImageSpace::InitAttenuationImage(const string& a_pathToAtnImage)
967 {
968  if(m_verbose>=3) Cout("oImageSpace::InitAttenuationImage ..."<< endl);
969 
970  // todo : read the number of dynamic images from Interfile header
971  // Allocate memory and initialize matrices according to the number of fr/rg/cg
972  // Provide this info to LMS_LoadAtnImage
973 
974 
975  // Initialization only if an image has been provided
976  if (!a_pathToAtnImage.empty()) // Image initiale
977  {
978  if(!m4p_attenuation)
979  {
981 
982  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
983  {
985 
986  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates(); rg++)
987  {
988  m4p_attenuation[fr][rg] = new FLTNB*[mp_ID->GetSensNbCardGates()];
989 
990  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates(); cg++)
991  m4p_attenuation[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
992  }
993  }
994  }
995 
996  if(LoadAttenuationImage(a_pathToAtnImage) )
997  {
998  Cerr("***** oImageSpace::LMS_InitAttenuationImage()-> Error while trying to read file at " << a_pathToAtnImage << endl);
999  return 1;
1000  }
1001  }
1002 
1003  return 0;
1004 }
1005 
1006 
1007 
1008 
1009 // =====================================================================
1010 // ---------------------------------------------------------------------
1011 // ---------------------------------------------------------------------
1012 // =====================================================================
1013 /*
1014  \fn LoadAttenuationImage()
1015  \param a_pathToAtnImage : path to an existing image
1016  \brief Load the attenuation image provided by the user in the m2p_attenuation matrix
1017  \todo If fr/rg/cg > 1 and only one attenuation image is provided, initialize each redundant matrice with pointers to the first one (deal with that here)
1018  \return 0 if success, positive value otherwise
1019 */
1020 int oImageSpace::LoadAttenuationImage(const string& a_pathToImage)
1021 {
1022  if(m_verbose>=3) Cout("oImageSpace::LMS_LoadAttenuationImage ..."<< endl);
1023 
1024  ifstream image_file(a_pathToImage.c_str(), ios::in|ios::binary); // Read the corresponding attenuation image
1025 
1026  if(!image_file.is_open())
1027  {
1028  Cerr("***** oImageSpace::LMS_LoadAttenuationImage()-> Error reading file !" << endl);
1029  return 1;
1030  }
1031  else
1032  {
1033  // Interfile image reading (INTF_LERP_ENABLED = interpolation allowed)
1035  {
1036  Cerr("***** oImageSpace::PROJ_LoadAttenuationImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
1037  return 1;
1038  }
1039 
1040  }
1041 
1042  image_file.close();
1043  return 0;
1044 }
1045 
1046 
1047 
1048 
1049 
1050 
1051 // =====================================================================
1052 // ---------------------------------------------------------------------
1053 // ---------------------------------------------------------------------
1054 // =====================================================================
1055 /*
1056  \fn InitImage
1057  \param a_pathToInitialImage : path to an existing image
1058  \param a_value : value to initialize each voxel with, if an input image is not provided
1059  \brief Initialize the main image, either using:
1060  - an existing image at path 'a_pathToInitialImage'
1061  - initialize each voxel with 'a_value'.
1062  \return 0 if success, positive value otherwise
1063 */
1064 int oImageSpace::InitImage(const string& a_pathToInitialImage, FLTNB a_value)
1065 {
1066  if(m_verbose>=3) Cout("oImageSpace::InitImage ..."<< endl);
1067 
1068  if (!a_pathToInitialImage.empty()) // Image initiale has been provided
1069  {
1070  if (LoadInitialImage(a_pathToInitialImage) )
1071  {
1072  Cerr("***** oImageSpace::InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
1073  return 1;
1074  }
1075  }
1076 
1077  else // Uniform initialization
1078  {
1079  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1080  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1081  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1082  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1083  m4p_image[tbf][rbf][cbf][v] = a_value;
1084  }
1085 
1086  return 0;
1087 }
1088 
1089 
1090 
1091 
1092 // =====================================================================
1093 // ---------------------------------------------------------------------
1094 // ---------------------------------------------------------------------
1095 // =====================================================================
1096 /*
1097  \fn LoadInitialImage()
1098  \param a_pathToImage : path to an existing image
1099  \brief Load the initial image provided by the user in the corresponding matrix
1100  \return 0 if success, positive value otherwise
1101 */
1102 int oImageSpace::LoadInitialImage(const string& a_pathToImage)
1103 {
1104  if(m_verbose>=3) Cout("oImageSpace::LoadInitialImage ..."<< endl);
1105 
1106  ifstream image_file;
1107  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
1108 
1109  if(!image_file.is_open())
1110  {
1111  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading file!" << endl);
1112  return 1;
1113  }
1114  else
1115  {
1116  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1117  // if(IntfReadImage(a_pathToImage, m4p_image, mp_ID, m_verbose, INTF_LERP_DISABLED))
1119  {
1120  Cerr("***** oImageSpace::LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
1121  return 1;
1122  }
1123  }
1124  image_file.close();
1125 
1126  return 0;
1127 }
1128 
1129 
1130 
1131 
1132 // =====================================================================
1133 // ---------------------------------------------------------------------
1134 // ---------------------------------------------------------------------
1135 // =====================================================================
1136 /*
1137  \fn InitializationBackwardImage()
1138  \brief Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on the fly).
1139 */
1141 {
1142  if(m_verbose>=3) Cout("oImageSpace::InitBackwardImage ..." << endl);
1143 
1144  // Reset backward images to 0.
1145  for (int img=0; img<m_nbBackwardImages; img++)
1146  {
1147  int th;
1148  #pragma omp parallel for private(th) schedule(static, 1)
1149  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1150  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1151  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1152  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1153  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1154  {
1155  m6p_backwardImage[img][th][tbf][rbf][cbf][v] = 0.;
1156  }
1157  }
1158 
1159  // If on-the-fly sensitivity, then reset to 0.
1160  if (!m_loadedSensitivity)
1161  {
1162  int th;
1163  #pragma omp parallel for private(th) schedule(static, 1)
1164  for (th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1165  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1166  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1167  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1168  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1169  m5p_sensitivity[th][fr][rg][cg][v] = 0.;
1170  }
1171 }
1172 
1173 
1174 
1175 
1176 
1177 
1178 // =====================================================================
1179 // ---------------------------------------------------------------------
1180 // ---------------------------------------------------------------------
1181 // =====================================================================
1182 /*
1183  \fn InitSensitivityImage()
1184  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
1185  \brief Initialization for the sensitivity image matrices
1186  \details Sensitivity image initialization depends on the reconstruction mode :
1187  - list-mode: Sensitivity image has been computed before reconstruction, it is loaded from the path provided in parameter
1188  - histogram: Sensitivity image is calculated on the fly during reconstruction.
1189  First dimension (thread) is only used in histogram mode, as the on-the-fly sensitivity image computation must be thread safe
1190  \return 0 if success, positive value otherwise
1191 */
1192 int oImageSpace::InitSensitivityImage(const string& a_pathToSensitivityImage)
1193 {
1194  if(m_verbose>=3) Cout("oImageSpace::InitSensitivityImage ..."<< endl);
1195 
1196  // =====================================================================================================
1197  // Case 1: a path to a sensitivity image is provided, meaning that we are in listmode and do not compute
1198  // the sensitivity on-the-fly. In that case, we do not take the multi-threading into account.
1199  // =====================================================================================================
1200 
1201  if (!a_pathToSensitivityImage.empty())
1202  {
1203  // Open file
1204  ifstream input_file;
1205  input_file.open(a_pathToSensitivityImage.c_str(), ios::binary | ios::in);
1206  // Read file
1207  if (input_file.is_open())
1208  {
1209  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
1210  if(IntfReadImage(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose, INTF_LERP_DISABLED) )
1211  {
1212  Cerr("***** oImageSpace::InitSensitivityImage()-> Error reading Interfile : " << a_pathToSensitivityImage << " !" << endl);
1213  return 1;
1214  }
1215  input_file.close();
1216  m_loadedSensitivity = true;
1217  }
1218  else
1219  {
1220  Cerr("***** oImageSpace::InitSensitivityImage() -> Input sensitivity file '" << a_pathToSensitivityImage << "' is missing or corrupted !" << endl);
1221  return 1;
1222  }
1223  }
1224 
1225  // =====================================================================================================
1226  // Case 2: no sensitivity provided, thus we are in histogram mode and will commpute it on-the-fly. We
1227  // thus need it to be thread-safe, so we allocated for all threads.
1228  // =====================================================================================================
1229 
1230  else
1231  {
1232 /*
1233  // Standard initialization
1234  for (int th=0; th<mp_ID->GetNbThreadsForProjection(); th++)
1235  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1236  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1237  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1238  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1239  m5p_sensitivity[th][fr][rg][cg][v] = -1.;
1240 * */
1241  // No need to initialize to any value, as it will be done by the InitBackwardImage function which is called at the beginning
1242  // of each subset, and deals with the sensitivity initialization to 0 when m_loadedSensitivity is false
1243  m_loadedSensitivity = false;
1244  }
1245 
1246  // End
1247  return 0;
1248 }
1249 
1250 
1251 
1252 
1253 
1254 // =====================================================================
1255 // ---------------------------------------------------------------------
1256 // ---------------------------------------------------------------------
1257 // =====================================================================
1258 /*
1259  \fn InitializationBwdImageForDeformation()
1260  \brief Initialize the buffer backward image dedicated to image-based deformation
1261 */
1263 {
1264  if(m_verbose>=3) Cout("oDeformationManager::InitBwdImageForDeformation ..." << endl);
1265 
1266  for (int img=0; img<m_nbBackwardImages; img++)
1267  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1268  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1269  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1270  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1271  m5p_defTmpBackwardImage[img][tbf][rbf][cbf][v] = 0.;
1272 }
1273 
1274 
1275 
1276 
1277 
1278 // =====================================================================
1279 // ---------------------------------------------------------------------
1280 // ---------------------------------------------------------------------
1281 // =====================================================================
1282 /*
1283  \fn InitializationSensImageForDeformation()
1284  \brief Initialize the buffer sensitivity image dedicated to image-based deformation, if required (histogram mode, as sensitivity is not loaded)
1285 */
1287 {
1288  if(m_verbose>=3) Cout("oImageSpace::InitSensImageForDeformation ..."<< endl);
1289 
1290  if(m_loadedSensitivity == false)
1291  {
1292  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1293  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1294  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1295  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
1296  m4p_defTmpSensitivityImage[fr][rg][cg][v] = 0.;
1297  }
1298 }
1299 
1300 
1301 
1302 
1303 
1304 // =====================================================================
1305 // ---------------------------------------------------------------------
1306 // ---------------------------------------------------------------------
1307 // =====================================================================
1308 /*
1309  \fn InstantiateOutputImage()
1310  \brief Compute output image using the m4p_image matrix and the time/respiratory/cardiac basis functions. Store the result in the m4p_outputImage matrix
1311  \details If time/respiratory/cardiac basis functions have been initialized, this function has no effect.
1312 */
1314 {
1315  if(m_verbose>=3) Cout("oImageSpace::ComputeOutputImage ..."<< endl);
1316 
1317  // First loop on frames
1318  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1319  {
1320  // Second loop on respiratory gates
1321  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1322  {
1323  // Third loop on cardiac gates
1324  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1325  {
1326  // Reset m4p_outputImage to 0
1327  int v;
1328  #pragma omp parallel for private(v) schedule(guided)
1329  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1330  m4p_outputImage[fr][rg][cg][v] = 0.;
1331  // First loop on time basis functions
1332  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1333  {
1334  // Get frame/basis coefficient and continue if null
1335  FLTNB time_basis_coef = mp_ID->GetTimeBasisCoefficient(tbf,fr);
1336  if (time_basis_coef==0.) continue;
1337  // Second loop on respiratory basis functions
1338  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1339  {
1340  // Get resp_gate/basis coefficient and continue if null
1341  FLTNB resp_basis_coef = mp_ID->GetRespBasisCoefficient(rbf,rg);
1342  if (resp_basis_coef==0.) continue;
1343  // Third loop on cardiac basis functions
1344  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1345  {
1346  // Get card_gate_basis coefficient and continue if null
1347  FLTNB card_basis_coef = mp_ID->GetCardBasisCoefficient(cbf,cg);
1348  if (card_basis_coef==0.) continue;
1349  // Compute global coefficient
1350  FLTNB global_basis_coeff = time_basis_coef * resp_basis_coef * card_basis_coef;
1351  // Loop on voxel with OpenMP (let the chunk size by default as all values are aligned in memory)
1352  #pragma omp parallel for private(v) schedule(guided)
1353  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1354  {
1355  // Add contribution from these basis functions
1356  m4p_outputImage[fr][rg][cg][v] += m4p_image[tbf][rbf][cbf][v] * global_basis_coeff;
1357  }
1358  }
1359  }
1360  }
1361  }
1362  }
1363  }
1364 
1365 }
1366 
1367 // =====================================================================
1368 // ---------------------------------------------------------------------
1369 // ---------------------------------------------------------------------
1370 // =====================================================================
1371 
1373 {
1374  // If no flip, then return
1375  if ( !mp_ID->GetFlipOutX() &&
1376  !mp_ID->GetFlipOutY() &&
1377  !mp_ID->GetFlipOutZ() ) return 0;
1378  // Verbose
1379  if (m_verbose>=2)
1380  {
1381  Cout("oImageSpace::ApplyOutputFlip() -> Flip image" << endl);
1382  if (mp_ID->GetFlipOutX()) Cout(" --> Over X" << endl);
1383  if (mp_ID->GetFlipOutY()) Cout(" --> Over Y" << endl);
1384  if (mp_ID->GetFlipOutZ()) Cout(" --> Over Z" << endl);
1385  }
1386  // First loop on frames
1387  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1388  {
1389  // Second loop on respiratory gates
1390  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1391  {
1392  // Third loop on cardiac gates
1393  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1394  {
1395  // Flip over Z
1396  if (mp_ID->GetFlipOutZ())
1397  {
1398  // Half loop over Z
1399  for (INTNB z_1=0; z_1<mp_ID->GetNbVoxZ()/2; z_1++)
1400  {
1401  // Compute opposite Z
1402  INTNB z_2 = mp_ID->GetNbVoxZ() - 1 - z_1;
1403  // For efficiency
1404  INTNB base_z1 = z_1 * mp_ID->GetNbVoxXY();
1405  INTNB base_z2 = z_2 * mp_ID->GetNbVoxXY();
1406  // Loop over Y
1407  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1408  {
1409  // For efficiency
1410  INTNB base_y = y * mp_ID->GetNbVoxX();
1411  // Loop over X
1412  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1413  {
1414  // Compute both indices
1415  INTNB indice_1 = base_z1 + base_y + x;
1416  INTNB indice_2 = base_z2 + base_y + x;
1417  // Switch voxels
1418  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1419  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1420  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1421  }
1422  }
1423  }
1424  }
1425  // Flip over Y
1426  if (mp_ID->GetFlipOutY())
1427  {
1428  // Loop over Z
1429  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1430  {
1431  // For efficiency
1432  INTNB base_z = z * mp_ID->GetNbVoxXY();
1433  // Half loop over Y
1434  for (INTNB y_1=0; y_1<mp_ID->GetNbVoxY()/2; y_1++)
1435  {
1436  // Compute opposite Y
1437  INTNB y_2 = mp_ID->GetNbVoxY() - 1 - y_1;
1438  // For efficiency
1439  INTNB base_y1 = y_1 * mp_ID->GetNbVoxX();
1440  INTNB base_y2 = y_2 * mp_ID->GetNbVoxX();
1441  // Loop over X
1442  for (INTNB x=0; x<mp_ID->GetNbVoxX(); x++)
1443  {
1444  // Compute both indices
1445  INTNB indice_1 = base_z + base_y1 + x;
1446  INTNB indice_2 = base_z + base_y2 + x;
1447  // Switch voxels
1448  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1449  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1450  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1451  }
1452  }
1453  }
1454  }
1455  // Flip over X
1456  if (mp_ID->GetFlipOutX())
1457  {
1458  // Loop over Z
1459  for (INTNB z=0; z<mp_ID->GetNbVoxZ(); z++)
1460  {
1461  // For efficiency
1462  INTNB base_z = z * mp_ID->GetNbVoxXY();
1463  // Loop over Y
1464  for (INTNB y=0; y<mp_ID->GetNbVoxY(); y++)
1465  {
1466  // For efficiency
1467  INTNB base_y = y * mp_ID->GetNbVoxX();
1468  // Half loop over X
1469  for (INTNB x_1=0; x_1<mp_ID->GetNbVoxX()/2; x_1++)
1470  {
1471  // Compute opposite X
1472  INTNB x_2 = mp_ID->GetNbVoxX() - 1 - x_1;
1473  // Compute both indices
1474  INTNB indice_1 = base_z + base_y + x_1;
1475  INTNB indice_2 = base_z + base_y + x_2;
1476  // Switch voxels
1477  FLTNB buffer = m4p_outputImage[fr][rg][cg][indice_1];
1478  m4p_outputImage[fr][rg][cg][indice_1] = m4p_outputImage[fr][rg][cg][indice_2];
1479  m4p_outputImage[fr][rg][cg][indice_2] = buffer;
1480  }
1481  }
1482  }
1483  }
1484  }
1485  }
1486  }
1487 
1488  // Normal end
1489  return 0;
1490 }
1491 
1492 // =====================================================================
1493 // ---------------------------------------------------------------------
1494 // ---------------------------------------------------------------------
1495 // =====================================================================
1496 
1498 {
1499  // If the output FOV percent is under 0 (the default value) and number of axial slices to be removed is 0, we do not mask anything
1500  if ( mp_ID->GetFOVOutPercent()<=0. &&
1501  mp_ID->GetNbSliceOutMask()==0 ) return 0;
1502  // Verbose
1503  if (m_verbose>=2)
1504  {
1505  Cout("oImageSpace::ApplyOutputFOVMasking() -> Mask output image" << endl);
1506  if (mp_ID->GetFOVOutPercent()>0.) Cout(" --> Mask transaxial FOV outside " << mp_ID->GetFOVOutPercent() << " %" << endl);
1507  if (mp_ID->GetNbSliceOutMask()>0) Cout(" --> Mask " << mp_ID->GetNbSliceOutMask() << " slices from both axial FOV limits" << endl);
1508  }
1509  // -----------------------------------------------
1510  // Transaxial FOV masking
1511  // -----------------------------------------------
1512  if (mp_ID->GetFOVOutPercent()>0.)
1513  {
1514  // Precast half the number of voxels over X and Y minus 1 (for efficiency)
1515  FLTNB flt_base_x = 0.5*((FLTNB)(mp_ID->GetNbVoxX()-1));
1516  FLTNB flt_base_y = 0.5*((FLTNB)(mp_ID->GetNbVoxY()-1));
1517  // Compute FOV elipse radius over X and Y, then squared
1518  FLTNB squared_radius_x = 0.5 * ((FLTNB)(mp_ID->GetNbVoxX())) * mp_ID->GetVoxSizeX()
1519  * mp_ID->GetFOVOutPercent() / 100.;
1520  squared_radius_x *= squared_radius_x;
1521  FLTNB squared_radius_y = 0.5 * ((FLTNB)(mp_ID->GetNbVoxY())) * mp_ID->GetVoxSizeY()
1522  * mp_ID->GetFOVOutPercent() / 100.;
1523  squared_radius_y *= squared_radius_y;
1524  // We assume that the computation of the distance from the center for a given
1525  // voxel and comparing it with the output FOV percentage costs more than performing
1526  // the loops in an inverse order compared to how the image is stored in memory.
1527  // Thus we begin the loops over X, then Y, then we test and if test passes, we
1528  // do the remaining loop over Z and over all dynamic dimensions.
1529  int x;
1530  #pragma omp parallel for private(x) schedule(guided)
1531  for (x=0; x<mp_ID->GetNbVoxX(); x++)
1532  {
1533  // Compute X distance from image center, then squared
1534  FLTNB squared_distance_x = (((FLTNB)x)-flt_base_x) * mp_ID->GetVoxSizeX();
1535  squared_distance_x *= squared_distance_x;
1536  // Loop over Y
1537  for (int y=0; y<mp_ID->GetNbVoxY(); y++)
1538  {
1539  // Compute Y distance from image center, then squared
1540  FLTNB squared_distance_y = (((FLTNB)y)-flt_base_y) * mp_ID->GetVoxSizeY();
1541  squared_distance_y *= squared_distance_y;
1542  // Test if the voxel is inside the FOV elipse, then we skip this voxel
1543  if ( squared_distance_x/squared_radius_x + squared_distance_y/squared_radius_y <= 1. ) continue;
1544  // Loop over Z
1545  for (int z=0; z<mp_ID->GetNbVoxZ(); z++)
1546  {
1547  // Compute global voxel index
1548  INTNB index = z*mp_ID->GetNbVoxXY() + y*mp_ID->GetNbVoxX() + x;
1549  // First loop on frames
1550  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1551  // Second loop on respiratory gates
1552  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1553  // Third loop on cardiac gates
1554  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1555  // Put 0 into this voxel
1556  m4p_outputImage[fr][rg][cg][index] = 0.;
1557  }
1558  }
1559  }
1560  }
1561  // -----------------------------------------------
1562  // Axial FOV masking
1563  // -----------------------------------------------
1564  if (mp_ID->GetNbSliceOutMask()>0)
1565  {
1566  INTNB removed_slices = mp_ID->GetNbSliceOutMask();
1567  // First loop on frames
1568  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1569  {
1570  // Second loop on respiratory gates
1571  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1572  {
1573  // Third loop on cardiac gates
1574  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1575  {
1576  // Mask slices
1577  for (int z=0; z<removed_slices; z++)
1578  {
1579  // First slices
1580  INTNB base_z_first = z*mp_ID->GetNbVoxXY();
1581  // Loop over Y and X
1582  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1583  {
1584  INTNB index = base_z_first + i;
1585  m4p_outputImage[fr][rg][cg][index] = 0.;
1586  }
1587  // Last slices
1588  INTNB base_z_last = (mp_ID->GetNbVoxZ()-1-z)*mp_ID->GetNbVoxXY();
1589  // Loop over Y and X
1590  for (int i=0; i<mp_ID->GetNbVoxXY(); i++)
1591  {
1592  INTNB index = base_z_last + i;
1593  m4p_outputImage[fr][rg][cg][index] = 0.;
1594  }
1595  }
1596  }
1597  }
1598  }
1599  }
1600  // End
1601  return 0;
1602 }
1603 
1604 
1605 // =====================================================================
1606 // ---------------------------------------------------------------------
1607 // ---------------------------------------------------------------------
1608 // =====================================================================
1609 /*
1610  \fn SaveOutputImage
1611  \param a_iteration : current iteration index
1612  \param a_subset : current number of subsets (or -1 by default)
1613  \brief Call the interfile function to write output image on disk
1614  \return 0 if success, positive value otherwise
1615 */
1616 int oImageSpace::SaveOutputImage(int a_iteration, int a_subset)
1617 {
1618  if(m_verbose>=3) Cout("oImageSpace::SaveOutputImage ..."<< endl);
1619 
1620  // Get the output manager
1621  sOutputManager* p_output_manager = sOutputManager::GetInstance();
1622 
1623  // Interfile writing
1624  // m4p_forwardImage is used here for the special case temporal regularization using intrinsic basis functions is used.
1625  // In this case, m4p_image contains basis function coefficients and not image values. These are computed in ComputeOutputImage()
1626  // In any other case, both matrices contain the image values.
1627 
1628  string path_to_img = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
1629 
1630  // Add a suffix for iteration
1631  if (a_iteration >= 0)
1632  {
1633  stringstream ss; ss << a_iteration + 1;
1634  path_to_img.append("_it").append(ss.str());
1635  }
1636 
1637  // Add a suffix for subset (if not negative by default), this means that we save a 'subset' image
1638  if (a_subset >= 0)
1639  {
1640  stringstream ss; ss << a_subset + 1;
1641  path_to_img.append("_ss").append(ss.str());
1642  }
1643 
1644  // We need one "mode" parameter which indicates if we would like to write 1 image file or several
1645  // Adjust functions according to that
1646  if(IntfWriteImgFile(path_to_img, m4p_outputImage, mp_ID, m_verbose) )
1647  {
1648  Cerr("***** oImageSpace::SaveOutputImage()-> Error writing Interfile of output image !" << endl);
1649  return 1;
1650  }
1651 
1652 
1653  // Normal end
1654  return 0;
1655 }
1656 
1657 
1658 
1659 
1660 
1661 // =====================================================================
1662 // ---------------------------------------------------------------------
1663 // ---------------------------------------------------------------------
1664 // =====================================================================
1665 /*
1666  \fn SaveDebugImage
1667  \param a_name : output name of the image
1668  \brief Just a debug function dedicated to write any kind of image on disk in raw format, for debugging purposes
1669 */
1670 void oImageSpace::SaveDebugImage(const string& a_name)
1671 {
1672  #ifdef CASTOR_DEBUG
1673  if(m_verbose>=3) Cout("oImageSpace::SaveDebugImage ..."<< endl);
1674  #endif
1675 
1676  ofstream output_file;
1677  output_file.open(a_name.c_str(), ios::binary | ios::out);
1678 
1679  //for (int fr=0; fr<ap_ImageDimensionsAndQuantification->GetNbTimeFrames(); fr++)
1680  // for (int rg=0; rg<ap_ImageDimensionsAndQuantification->GetNbRespGates(); rg++)
1681  // for (int cg=0; cg<ap_ImageDimensionsAndQuantification->GetNbCardGates(); cg++)
1682  // {
1683  // Write file
1684  output_file.write(reinterpret_cast<char*>(m6p_backwardImage[0][0][0][0][0]), mp_ID->GetNbVoxXYZ()*sizeof(FLTNB));
1685  // Close file
1686  output_file.close();
1687  // }
1688  // }
1689  //}
1690 }
1691 
1692 
1693 
1694 
1695 // =====================================================================
1696 // ---------------------------------------------------------------------
1697 // ---------------------------------------------------------------------
1698 // =====================================================================
1699 /*
1700  \fn PrepareForwardImage
1701  \brief Copy current image matrix in the forward-image buffer matrix
1702 */
1704 {
1705  if(m_verbose>=3) Cout("oImageSpace::PrepareForwardImage ..."<< endl);
1706 
1707  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1708  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1709  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1710  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1711  m4p_forwardImage[tbf][rbf][cbf][v] = m4p_image[tbf][rbf][cbf][v];
1712 }
1713 
1714 
1715 
1716 
1717 // =====================================================================
1718 // ---------------------------------------------------------------------
1719 // ---------------------------------------------------------------------
1720 // =====================================================================
1721 /*
1722  \fn Reduce
1723  \brief Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
1724  \todo Choose computation alternative (do not know yet which one is faster...)
1725 */
1727 {
1728  if(m_verbose>=3) Cout("oImageSpace::Reduce ..."<< endl);
1729 
1730  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
1731  if (m_verbose>=3) Cout("oImageSpace::Reduce() -> Merge parallel results" << endl);
1732  #endif
1733 
1734  // --------------------------------------------------------------------------------
1735  // Step 1: merge multi-threads results
1736  // --------------------------------------------------------------------------------
1737 
1738  #ifdef CASTOR_OMP
1739  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
1740 
1741  // Special case here where it appears to be always beneficial to use as many threads
1742  // as the number of threads used for projections. So we set it here, and set it back
1743  // at the end.
1744  omp_set_num_threads(mp_ID->GetNbThreadsForProjection());
1745 
1746  // todo : Choose computation alternative (do not know yet which one is faster...)
1747  int alternative = 2;
1748 
1749  // Alternative 1: standard loops based on the pointers hierarchy but mono-thread
1750  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
1751  {
1752  for (int img=0; img<m_nbBackwardImages; img++)
1753  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
1754  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1755  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1756  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1757  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1758  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
1759  }
1760  // Alternative 2: multi-thread merging with the voxels loop at first to be thread safe
1761  // Maybe faster than alternative 1, but the first loop on voxels breaks the memory order... have to try
1762  else if (alternative==2)
1763  {
1764  int v;
1765  #pragma omp parallel for private(v) schedule(guided)
1766  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1767  {
1768  for (int img=0; img<m_nbBackwardImages; img++)
1769  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
1770  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1771  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1772  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1773  m6p_backwardImage[img][0][tbf][rbf][cbf][v] += m6p_backwardImage[img][th][tbf][rbf][cbf][v];
1774  }
1775  }
1776  // If on-the-fly sensitivity, then do it also for it.
1777  if (!m_loadedSensitivity)
1778  {
1779  if (alternative==1 || mp_ID->GetNbThreadsForImageComputation()==1)
1780  {
1781  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
1782  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1783  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1784  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1785  for (int v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1786  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
1787  }
1788  else if (alternative==2)
1789  {
1790  int v;
1791  #pragma omp parallel for private(v) schedule(guided)
1792  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1793  {
1794  for (int th=1; th<mp_ID->GetNbThreadsForProjection(); th++)
1795  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1796  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1797  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1798  m5p_sensitivity[0][fr][rg][cg][v] += m5p_sensitivity[th][fr][rg][cg][v];
1799  }
1800  }
1801  }
1802  // Set back the number of threads to the one for image computation
1803  omp_set_num_threads(mp_ID->GetNbThreadsForImageComputation());
1804  #endif
1805 
1806  // --------------------------------------------------------------------------------
1807  // Step 2: merge multi-instance results (MPI)
1808  // --------------------------------------------------------------------------------
1809 
1810  #ifdef CASTOR_MPI
1811  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
1812  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
1813 
1814  // Merge backward images
1815  for (int img=0; img<m_nbBackwardImages; img++)
1816  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1817  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1818  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1819  MPI_Allreduce( MPI_IN_PLACE, &m6p_backwardImage[img][0][tbf][rbf][cbf][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
1820 
1821  // If on-the-fly sensitivity, then do it also for it
1822  if (!m_loadedSensitivity)
1823  {
1824  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
1825  for (int rg=0; rg<mp_ID->GetNbRespGates(); rg++)
1826  for (int cg=0; cg<mp_ID->GetNbCardGates(); cg++)
1827  MPI_Allreduce( MPI_IN_PLACE, &m5p_sensitivity[0][fr][rg][cg][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
1828  }
1829  #endif
1830 }
1831 
1832 
1833 
1834 
1835 
1836 // =====================================================================
1837 // ---------------------------------------------------------------------
1838 // ---------------------------------------------------------------------
1839 // =====================================================================
1840 /*
1841  \fn CleanNeverVisitedVoxels
1842  \brief Based on the visitedVoxelsImage, clean the never visited voxels in the image.
1843  This function must be called at the end of each iteration
1844 */
1846 {
1847  if(m_verbose>=3) Cout("oImageSpace::CleanNeverVisitedVoxels ..."<< endl);
1848 
1849  int v;
1850  #pragma omp parallel for private(v) schedule(guided)
1851  for (v=0; v<mp_ID->GetNbVoxXYZ(); v++)
1852  {
1853  // Clean this image voxel if the voxel was never visited
1854  if (mp_visitedVoxelsImage[v]==0.)
1855  {
1856  for (int tbf=0; tbf<mp_ID->GetNbTimeBasisFunctions(); tbf++)
1857  for (int rbf=0; rbf<mp_ID->GetNbRespBasisFunctions(); rbf++)
1858  for (int cbf=0; cbf<mp_ID->GetNbCardBasisFunctions(); cbf++)
1859  {
1860  m4p_image[tbf][rbf][cbf][v] = 0.;
1861  }
1862  }
1863  // Reset the image
1864  mp_visitedVoxelsImage[v] = 0.;
1865  }
1866 }
1867 
1868 
1869 
1870 
1871 
1872 
1873 
1874 
1875 
1876 // LIST-MODE SENSITIVITY GENERATION FUNCTIONS
1877 // =====================================================================
1878 // ---------------------------------------------------------------------
1879 // ---------------------------------------------------------------------
1880 // =====================================================================
1881 /*
1882  \fn LMS_InstantiateImage()
1883  \brief Allocate memory for the main image matrices (for list-mode sensitivity generation)
1884  \details This function is dedicated to list-mode sensitivity (LMS) generation.
1885 */
1887 {
1888  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateImage ..."<< endl);
1889 
1890  m4p_image = new FLTNB***[mp_ID->GetNbTimeFrames()];
1891 
1892  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1893  {
1894  m4p_image[fr] = new FLTNB**[mp_ID->GetSensNbRespGates()];
1895 
1896  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
1897  {
1898  m4p_image[fr][rg] = new FLTNB*[mp_ID->GetSensNbCardGates()];
1899 
1900  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
1901  m4p_image[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1902  }
1903  }
1904 }
1905 
1906 
1907 
1908 
1909 // =====================================================================
1910 // ---------------------------------------------------------------------
1911 // ---------------------------------------------------------------------
1912 // =====================================================================
1913 /*
1914  \fn LMS_DeallocateImage()
1915  \brief Free memory for the main image matrices (for list-mode sensitivity generation)
1916  \details This function is dedicated to list-mode sensitivity (LMS) generation.
1917 */
1919 {
1920  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateImage ..."<< endl);
1921 
1922  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1923  {
1924  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
1925  {
1926  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
1927  if(m4p_image[fr][rg][cg] != NULL) delete m4p_image[fr][rg][cg];
1928 
1929  if(m4p_image[fr][rg] != NULL) delete[] m4p_image[fr][rg];
1930  }
1931  if(m4p_image[fr] != NULL) delete[] m4p_image[fr] ;
1932  }
1933 
1934  if(m4p_image != NULL) delete[] m4p_image;
1935  m4p_image = NULL;
1936 }
1937 
1938 
1939 
1940 // =====================================================================
1941 // ---------------------------------------------------------------------
1942 // ---------------------------------------------------------------------
1943 // =====================================================================
1944 /*
1945  \fn LMS_InstantiateForwardImage()
1946  \brief Allocate memory for the forward image matrices (for list-mode sensitivity generation)
1947  \details This function is dedicated to list-mode sensitivity (LMS) generation.
1948 */
1950 {
1951  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateForwardImage ..."<< endl);
1952 
1954 
1955  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1956  {
1958 
1959  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
1960  {
1961  m4p_forwardImage[fr][rg] = new FLTNB*[mp_ID->GetSensNbCardGates()];
1962 
1963  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
1964  {
1965  m4p_forwardImage[fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
1966  }
1967  }
1968 
1969  }
1970 }
1971 
1972 
1973 
1974 // =====================================================================
1975 // ---------------------------------------------------------------------
1976 // ---------------------------------------------------------------------
1977 // =====================================================================
1978 /*
1979  \fn LMS_DeallocateForwardImage()
1980  \brief Free memory for the forward image matrices (for list-mode sensitivity generation)
1981  \details This function is dedicated to list-mode sensitivity (LMS) generation.
1982 */
1984 {
1985  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateForwardImage ..."<< endl);
1986 
1987  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
1988  {
1989  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
1990  {
1991  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
1992  {
1993  if(m4p_forwardImage[fr][rg][cg] != NULL) delete m4p_forwardImage[fr][rg][cg];
1994  }
1995 
1996  if(m4p_forwardImage[fr][rg] != NULL) delete[] m4p_forwardImage[fr][rg];
1997  }
1998 
1999  if(m4p_forwardImage[fr] != NULL) delete[] m4p_forwardImage[fr];
2000  }
2001 
2002  if(m4p_forwardImage != NULL) delete[] m4p_forwardImage;
2003  m4p_forwardImage = NULL;
2004 }
2005 
2006 // =====================================================================
2007 // ---------------------------------------------------------------------
2008 // ---------------------------------------------------------------------
2009 // =====================================================================
2010 /*
2011  \fn LMS_InstantiateSensitivityImage()
2012  \brief Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation)
2013  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2014 */
2016 {
2017  if(m_verbose>=3) Cout("oImageSpace::LMS_InstantiateSensitivityImage ..."<< endl);
2018 
2019  m5p_sensitivity = new FLTNB****[1];
2020  m5p_sensitivity[0] = new FLTNB***[mp_ID->GetNbTimeFrames()];
2021  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2022  {
2023  m5p_sensitivity[0][fr] = new FLTNB**[mp_ID->GetSensNbRespGates()];
2024  for (int rg=0; rg<mp_ID->GetSensNbRespGates(); rg++)
2025  {
2026  m5p_sensitivity[0][fr][rg] = new FLTNB*[mp_ID->GetSensNbCardGates()];
2027  for (int cg=0; cg<mp_ID->GetSensNbCardGates(); cg++)
2028  m5p_sensitivity[0][fr][rg][cg] = new FLTNB[mp_ID->GetNbVoxXYZ()];
2029  }
2030  }
2031 }
2032 
2033 
2034 
2035 
2036 // =====================================================================
2037 // ---------------------------------------------------------------------
2038 // ---------------------------------------------------------------------
2039 // =====================================================================
2040 /*
2041  \fn LMS_DeallocateSensitivityImage()
2042  \brief Free memory for the sensitivity image matrices (for list-mode sensitivity generation)
2043  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2044 */
2046 {
2047  if(m_verbose>=3) Cout("oImageSpace::LMS_DeallocateSensitivityImage ..."<< endl);
2048 
2049  for (int fr=0; fr<mp_ID->GetNbTimeFrames(); fr++)
2050  {
2051  for (int rg=0; rg<mp_ID->GetSensNbRespGates(); rg++)
2052  {
2053  for (int cg=0; cg<mp_ID->GetSensNbCardGates(); cg++)
2054  {
2055  if(m5p_sensitivity[0][fr][rg][cg] != NULL) delete m5p_sensitivity[0][fr][rg][cg];
2056  }
2057  if(m5p_sensitivity[0][fr][rg] != NULL)delete m5p_sensitivity[0][fr][rg];
2058  }
2059  if(m5p_sensitivity[0][fr] != NULL)delete m5p_sensitivity[0][fr];
2060  }
2061  if(m5p_sensitivity[0] != NULL) delete[] m5p_sensitivity[0];
2062 
2063  if(m5p_sensitivity != NULL) delete[] m5p_sensitivity;
2064  m5p_sensitivity = NULL;
2065 }
2066 
2067 
2068 
2069 
2070 // =====================================================================
2071 // ---------------------------------------------------------------------
2072 // ---------------------------------------------------------------------
2073 // =====================================================================
2074 /*
2075  \fn LMS_CopyAtnToImage()
2076  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_image matrix (for list-mode sensitivity generation)
2077  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2078 */
2080 {
2081  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToImage ..."<< endl);
2082 
2083  if(m4p_attenuation != NULL)
2084  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2085  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2086  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2087  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2088  {
2089  m4p_image[fr][rg][cg][v] = m4p_attenuation[fr][rg][cg][v];
2090  }
2091  else
2092  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2093  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2094  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2095  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2096  {
2097  m4p_image[fr][rg][cg][v] = 0;
2098  }
2099 
2100 }
2101 
2102 
2103 // =====================================================================
2104 // ---------------------------------------------------------------------
2105 // ---------------------------------------------------------------------
2106 // =====================================================================
2107 /*
2108  \fn LMS_CopyAtnToForwardImage()
2109  \brief Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_forwardImage matrix (for list-mode sensitivity generation)
2110  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2111 */
2113 {
2114  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyAtnToForwardImage ..."<< endl);
2115 
2116  if(m4p_attenuation != NULL)
2117  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2118  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2119  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2120  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2121  {
2122  m4p_forwardImage[fr][rg][cg][v] = m4p_attenuation[fr][rg][cg][v];
2123  }
2124  else
2125  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2126  for(int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2127  for(int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2128  for(int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2129  {
2130  m4p_forwardImage[fr][rg][cg][v] = 0;
2131  }
2132 }
2133 
2134 // =====================================================================
2135 // ---------------------------------------------------------------------
2136 // ---------------------------------------------------------------------
2137 // =====================================================================
2138 /*
2139  \fn LMS_CopyBackwardToSensitivityImage()
2140  \brief Copy the backward image containing the result of the sensitivity back-projection, into the sensitivity matrix.
2141  \details This function is dedicated to list-mode sensitivity (LMS) generation, it is useful to apply convolution and to save the image.
2142 */
2144 {
2145  if(m_verbose>=3) Cout("oImageSpace::LMS_CopyBackwardToSensitivity ..."<< endl);
2146  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2147  for (int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2148  for (int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2149  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2150  m5p_sensitivity[0][fr][rg][cg][v] = m6p_backwardImage[0][0][fr][rg][cg][v];
2151 }
2152 
2153 
2154 // =====================================================================
2155 // ---------------------------------------------------------------------
2156 // ---------------------------------------------------------------------
2157 // =====================================================================
2158 /*
2159  \fn LMS_PrepareForwardImage()
2160  \brief Copy current image in forward-image buffer (for list-mode sensitivity generation)
2161  \details This function is dedicated to list-mode sensitivity (LMS) generation.
2162 */
2164 {
2165  if(m_verbose>=3) Cout("oImageSpace::LMS_PrepareForwardImage ..."<< endl);
2166 
2167  for (int fr=0 ; fr<mp_ID->GetNbTimeFrames() ; fr++)
2168  for (int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2169  for (int cg=0 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2170  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2171  m4p_forwardImage[fr][rg][cg][v] = m4p_image[fr][rg][cg][v];
2172 }
2173 
2174 
2175 
2176 
2177 // =====================================================================
2178 // ---------------------------------------------------------------------
2179 // ---------------------------------------------------------------------
2180 // =====================================================================
2181 
2182 void oImageSpace::ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
2183 {
2184  #if defined(CASTOR_OMP) || defined(CASTOR_MPI)
2185  if (m_verbose>=3) Cout("oImageSpace::ReduceBackwardImage() -> Merge parallel results" << endl);
2186  #endif
2187 
2188  // --------------------------------------------------------------------------------
2189  // Step 1: merge multi-threads results
2190  // --------------------------------------------------------------------------------
2191 
2192  #ifdef CASTOR_OMP
2193  // Verbose
2194  if (m_verbose>=3) Cout(" --> Over threads ..." << endl);
2195  // Do it
2196  for (int th=1 ; th<mp_ID->GetNbThreadsForProjection() ; th++)
2197  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2198  m6p_backwardImage[a_imageIndex][0][a_timeIndex][a_respIndex][a_cardIndex][v] += m6p_backwardImage[a_imageIndex][th][a_timeIndex][a_respIndex][a_cardIndex][v];
2199  #endif
2200 
2201  // --------------------------------------------------------------------------------
2202  // Step 2: merge multi-instance results (MPI)
2203  // --------------------------------------------------------------------------------
2204 
2205  #ifdef CASTOR_MPI
2206  // We use the MPI_IN_PLACE as the send buffer so that it is the same as the result buffer
2207  if (m_verbose>=3) Cout(" --> Over instances ..." << endl);
2208  // Merge backward images
2209  MPI_Allreduce( MPI_IN_PLACE, &m6p_backwardImage[a_imageIndex][0][a_timeIndex][a_respIndex][a_cardIndex][0], mp_ID->GetNbVoxXYZ(), FLTNBMPI, MPI_SUM, MPI_COMM_WORLD );
2210  #endif
2211 }
2212 
2213 
2214 
2215 
2216 // =====================================================================
2217 // ---------------------------------------------------------------------
2218 // ---------------------------------------------------------------------
2219 // =====================================================================
2220 /*
2221  \fn LMS_SaveSensitivityImage
2222  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2223  \param ap_DeformationManager : Pointer to the deformation manager objet (required to retrieve the number of gates in the sensitivity image)
2224  \brief Call the interfile function to write the sensitivity image on disk
2225  \details If image deformation is enabled for respiratory/cardiac gated data, the gated images are summed up into one image and normalize
2226  \return 0 if success, positive value otherwise
2227  \todo Interfile management : if the number of sensitivity image to write is different to the actual number of image
2228  as it could be the case for dynamic imaging, if one sensitivity image is required for the whole dynamic serie,
2229  we will have to give this information to IntfWriteImgFile(), otherwise the Interfile functions will try to write several times the image (leading to segfault or error)
2230 */
2231 int oImageSpace::LMS_SaveSensitivityImage(const string& a_pathToSensitivityImage, oDeformationManager* ap_DeformationManager)
2232 {
2233  if(m_verbose>=3) Cout("oImageSpace::LMS_SaveSensitivityImage ..."<< endl);
2234 
2235  //sOutputManager* p_output_manager = sOutputManager::GetInstance();
2236 
2237  // According to the respiratory motion algorithm, images relatives to resp gates, cardiac gates, or both, have to be summed up and normalize into one image.
2238  // The resulting image will be summed up into the sensitivity matrix
2239 
2240  int nb_reco_card_images = ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetSensNbCardGates());
2241 
2242  if (ap_DeformationManager->GetNbSensImagesCardDeformation(mp_ID->GetSensNbCardGates()) == 0) // Summation of the card images
2243  {
2244  nb_reco_card_images = 1;
2245  FLTNB card_normalization = mp_ID->GetSensNbCardGates();
2246 
2247  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2248  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2249  for (int rg=0 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2250  {
2251  for (int cg=1 ; cg<mp_ID->GetSensNbCardGates() ; cg++)
2252  m5p_sensitivity[0][fr][rg][0][v] += m5p_sensitivity[0][fr][rg][cg][v];
2253  m5p_sensitivity[0][fr][rg][0][v] /= card_normalization;
2254  }
2255  }
2256 
2257  if (ap_DeformationManager->GetNbSensImagesRespDeformation(mp_ID->GetSensNbRespGates()) == 0) // Summation of the resp images
2258  {
2259  FLTNB resp_normalization = mp_ID->GetSensNbRespGates();
2260 
2261  for (int v=0 ; v<mp_ID->GetNbVoxXYZ() ; v++)
2262  for(int fr=0 ; fr<mp_ID->GetNbTimeFrames(); fr++)
2263  for (int cg=0 ; cg<nb_reco_card_images ; cg++) // Dual gating : Be sure to get the right number of card images
2264  {
2265  for (int rg=1 ; rg<mp_ID->GetSensNbRespGates() ; rg++)
2266  m5p_sensitivity[0][fr][0][cg][v] += m5p_sensitivity[0][fr][rg][cg][v];
2267  m5p_sensitivity[0][fr][0][cg][v] /= resp_normalization;
2268  }
2269  }
2270 
2271  if (SaveSensitivityImage(a_pathToSensitivityImage) )
2272  {
2273  Cerr("***** oImageSpace::LMS_SaveSensitivityImage()-> Error writing Sensitivity image !" << endl);
2274  return 1;
2275  }
2276 
2277  return 0;
2278 }
2279 
2280 
2281 
2282 
2283 // =====================================================================
2284 // ---------------------------------------------------------------------
2285 // ---------------------------------------------------------------------
2286 // =====================================================================
2287 /*
2288  \fn SaveSensitivityImage
2289  \param a_pathToSensitivityImage : path to the sensitivity image (should be provided only in list-mode reconstruction)
2290  \brief Call the interfile function to write the sensitivity image on disk
2291  \return 0 if success, positive value otherwise
2292 */
2293 int oImageSpace::SaveSensitivityImage(const string& a_pathToSensitivityImage)
2294 {
2295  if(m_verbose>=3) Cout("oImageSpace::SaveSensitivityImage ..."<< endl);
2296 
2297  if (IntfWriteImgFile(a_pathToSensitivityImage, m5p_sensitivity[0], mp_ID, m_verbose) )
2298  {
2299  Cerr("***** oImageSpace::SaveSensitivityImage()-> Error writing Interfile of output image !" << endl);
2300  return 1;
2301  }
2302 
2303  return 0;
2304 }
2305 
2306 
2307 
2308 
2309 // PROJECTION SCRIPT FUNCTIONS
2310 // =====================================================================
2311 // ---------------------------------------------------------------------
2312 // ---------------------------------------------------------------------
2313 // =====================================================================
2314 /*
2315  \fn PROJ_InstantiateProjectionImage()
2316  \param a_nbProjs : a number of projection slices in the projection
2317  \param a_nbPixels : a total number of pixels in the projection slices
2318  \brief Instanciate and initialize projection image for analytical projection
2319  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2320  \todo support for PET (requires projection format)
2321 */
2322 void oImageSpace::PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
2323 {
2324  if(m_verbose>=3) Cout("oImageSpace::PROJ_InstantiateProjectionImage ..."<< endl);
2325 
2326  m2p_projectionImage = new FLTNB*[a_nbProjs];
2327 
2328  for(int p=0 ; p<a_nbProjs ; p++)
2329  {
2330  m2p_projectionImage[p] = new FLTNB[a_nbPixels];
2331  for(int px=0 ; px<a_nbPixels ; px++)
2332  m2p_projectionImage[p][px] = 0;
2333  }
2334 }
2335 
2336 
2337 
2338 
2339 // =====================================================================
2340 // ---------------------------------------------------------------------
2341 // ---------------------------------------------------------------------
2342 // =====================================================================
2343 /*
2344  \fn PROJ_DeallocateProjectionImage()
2345  \param a_nbProjs : a number of projection slices in the projection
2346  \brief Free memory for the projection image for analytical projection
2347  \details This function is currently only dedicated to SPECT projection, and used by the analytical projection script
2348  \todo support for PET (requires projection format)
2349 */
2351 {
2352  if(m_verbose>=3) Cout("oImageSpace::PROJ_DeallocateProjectionImage ..."<< endl);
2353 
2354  for(int p=0 ; p<a_nbProjs ; p++)
2355  {
2356  delete[] m2p_projectionImage[p];
2357  }
2358 
2359  delete[] m2p_projectionImage;
2360  m2p_projectionImage = NULL;
2361 }
2362 
2363 
2364 
2365 
2366 // =====================================================================
2367 // ---------------------------------------------------------------------
2368 // ---------------------------------------------------------------------
2369 // =====================================================================
2370 /*
2371  \fn PROJ_InitImage()
2372  \param a_pathToInitialImage : path to the image to project
2373  \brief Load the initial image for the analytical projection
2374  \return 0 if success, positive value otherwise
2375 */
2376 int oImageSpace::PROJ_InitImage(const string& a_pathToInitialImage)
2377 {
2378  if(m_verbose>=3) Cout("oImageSpace::PROJ_InitImage ..."<< endl);
2379 
2380  if (!a_pathToInitialImage.empty()) // Image initiale
2381  {
2382  if(PROJ_LoadInitialImage(a_pathToInitialImage) )
2383  {
2384  Cerr("***** oImageSpace::PROJ_InitImage()-> Error while trying to read file at " << a_pathToInitialImage << endl);
2385  return 1;
2386  }
2387  }
2388  else
2389  {
2390  {
2391  Cerr("***** oImageSpace::PROJ_InitImage()-> No projected image provided ! " << endl);
2392  return 1;
2393  }
2394  }
2395  return 0;
2396 }
2397 
2398 
2399 
2400 
2401 // =====================================================================
2402 // ---------------------------------------------------------------------
2403 // ---------------------------------------------------------------------
2404 // =====================================================================
2405 /*
2406  \fn PROJ_LoadInitialImage()
2407  \param a_pathToImage : path to the image to project
2408  \brief Load the initial image for the analytical projection
2409  \return 0 if success, positive value otherwise
2410 */
2411 int oImageSpace::PROJ_LoadInitialImage(const string& a_pathToImage)
2412 {
2413  if(m_verbose>=3) Cout("oImageSpace::PROJ_LoadInitialImage ..."<< endl);
2414 
2415  ifstream image_file;
2416  image_file.open(a_pathToImage.c_str(), ios::in | ios::binary);
2417 
2418  if(!image_file.is_open())
2419  {
2420  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading file !" << endl);
2421  return 1;
2422  }
2423  else
2424  {
2425  // Interfile image reading (INTF_LERP_DISABLED = no interpolation allowed)
2427  {
2428  Cerr("***** oImageSpace::PROJ_LoadInitialImage()-> Error reading Interfile : " << a_pathToImage << " !" << endl);
2429  return 1;
2430  }
2431 
2432  }
2433  image_file.close();
2434  return 0;
2435 }
2436 
2437 
2438 
2439 
2440 
2441 // =====================================================================
2442 // ---------------------------------------------------------------------
2443 // ---------------------------------------------------------------------
2444 // =====================================================================
2445 /*
2446  \fn PROJ_SaveProjectionImage
2447  \brief Save an image of the projected data
2448  \return 0 if success, positive value otherwise
2449  \todo Support for PET projeted data
2450  implement interfile file recovery in a scanner class function
2451 */
2453 {
2454  if(m_verbose>=3) Cout("oImageSpace::PROJ_SaveProjectionImage ..."<< endl);
2455 
2456  string img_name = "_ProjectionImage";
2457 
2458  // Initialize interfile fields structure
2459  Intf_fields Img_fields;
2460  IntfKeySetFieldsOutput(&Img_fields , mp_ID);
2461 
2462  Img_fields.process_status = "acquired";
2463 
2464  // Get specific info about projection
2466 
2467  // Recover the values in the interfile structure.
2468  // todo : This step for PET projection
2469  // todo : Implement this step in scanner class functions
2470  // todo : For spect, deal with systems with several detector heads (currently assume one head)
2471  if(p_ScanMgr->GetScannerType() == 2)
2472  {
2473  uint16_t nb_projs, nb_heads;
2474  uint16_t nb_bins[2];
2475  FLTNB pix_size[2];
2476  FLTNB* p_angles;
2477  FLTNB* p_radius;
2478  int dir_rot;
2479  p_ScanMgr->GetSPECTSpecificParameters(&nb_projs,
2480  &nb_heads,
2481  nb_bins,
2482  pix_size,
2483  p_angles,
2484  p_radius,
2485  &dir_rot);
2486 
2487  Img_fields.nb_projections = nb_projs;
2488  Img_fields.nb_detector_heads = nb_heads;
2489  Img_fields.mtx_size[0] = nb_bins[0];
2490  Img_fields.mtx_size[1] = nb_bins[1];
2491  Img_fields.vox_size[0] = pix_size[0];
2492  Img_fields.vox_size[1] = pix_size[1];
2493  Img_fields.extent_rotation = 360; // 360 by default
2494  // Check rotation direction from the two first angles
2495  Img_fields.direction_rotation = (dir_rot == GEO_ROT_CCW)?
2496  "CCW" :
2497  "CW";
2498  Img_fields.first_angle = p_angles[0];
2499 
2500  // Note: In Interfile v3.3, doesn't seem to have a field to provide
2501  // individually each angle and Center of rotation.
2502  // Looks like it was planned for ulterior version, at least for the
2503  // center of rotation
2504  // For now, we just write each projection angle and radius as an array key
2505 
2506  string angles_str = "{";
2507  string radius_str = "{";
2508 
2509  // Flag to check if the radius is similar for each projection angle
2510  bool has_single_radius = true;
2511 
2512  for(uint16_t p=0 ; p<nb_projs ; p++)
2513  {
2514  stringstream ss_a, ss_r;
2515  // projection angles
2516  ss_a << p_angles[p];
2517  angles_str.append(ss_a.str());
2518  (p<nb_projs-1) ? angles_str.append(",") : angles_str.append("}");
2519 
2520  // radius
2521  ss_r << p_radius[p];
2522  radius_str.append(ss_r.str());
2523  (p<nb_projs-1) ? radius_str.append(",") : radius_str.append("}");
2524  if(p_radius[p] != p_radius[0])
2525  has_single_radius = false;
2526 
2527  // Some interfile editors struggle with long line, perhaps because
2528  // they are stored in char[256] or something like that
2529  // Hence, break line after a certain number of elements
2530  if((p+1)%30 == 0)
2531  {
2532  angles_str.append("\n");
2533  radius_str.append("\n");
2534  }
2535  }
2536 
2537  // If there is one common radius, write its value in the related string
2538  if(has_single_radius)
2539  {
2540  stringstream ss;
2541  ss << p_radius[0];
2542  radius_str = ss.str();
2543  }
2544 
2545  Img_fields.projection_angles = angles_str;
2546  Img_fields.radius = radius_str;
2547 
2548 
2549  // Common code to all modalities
2550  sOutputManager* p_output_manager = sOutputManager::GetInstance();
2551  string img_file_name = p_output_manager->GetPathName() + p_output_manager->GetBaseName();
2552  img_file_name.append("_ProjImage");
2553 
2554  if(IntfWriteProjFile(img_file_name, m2p_projectionImage, mp_ID, Img_fields, m_verbose) )
2555  {
2556  Cerr("***** oImageSpace::PROJ_SaveProjectionImage()-> Error writing Interfile of output image !" << endl);
2557  return 1;
2558  }
2559  }
2560  return 0;
2561 }
2562 
2563 
int GetNbSensImagesCardDeformation(int a_value)
return the required number of cardiac images in the sensitivity image depending on the cardiac deform...
#define INTF_LERP_DISABLED
Definition: oInterfileIO.hh:75
FLTNB **** m4p_forwardImage
Definition: oImageSpace.hh:68
void LMS_InstantiateImage()
Allocate memory for the main image matrices (for list-mode sensitivity generation) ...
void DeallocateBwdImageForDeformation()
Free memory for the buffer backward image required for image-based deformation.
Definition: oImageSpace.cc:794
int GetNbSensImagesRespDeformation(int a_value)
return the required number of respiratory images in the sensitivity image depending on the respirator...
static sScannerManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
FLTNB first_angle
FLTNB ** m2p_projectionImage
Definition: oImageSpace.hh:144
string direction_rotation
void DeallocateBackwardImageFromDynamicBasis()
Free memory for the backward image matrices.
Definition: oImageSpace.cc:247
int m_nbBackwardImages
Definition: oImageSpace.hh:596
FLTNB vox_size[3]
int GetNbCardBasisFunctions()
Get the number of cardiac basis functions.
#define FLTNB
Definition: gVariables.hh:55
#define INTF_LERP_ENABLED
Definition: oInterfileIO.hh:77
int InitMaskImage(const string &a_pathToImage)
Memory allocation and initialization for the mask image.
Definition: oImageSpace.cc:593
int PROJ_InitImage(const string &a_pathToInitialImage)
Load the initial image for the analytical projection.
FLTNB GetVoxSizeX()
Get the voxel's size along the X axis, in mm.
void DeallocateImage()
Free memory for the main image matrices.
Definition: oImageSpace.cc:100
#define FLTNBMPI
Definition: gVariables.hh:57
void InitSensImageForDeformation()
Initialize the buffer sensitivity image dedicated to image-based deformation, if required (histogram ...
int LMS_SaveSensitivityImage(const string &a_pathToSensitivityImage, oDeformationManager *ap_DeformationManager)
Call the interfile function to write the sensitivity image on disk.
void LMS_PrepareForwardImage()
Copy current image in forward-image buffer (for list-mode sensitivity generation) ...
FLTNB GetFOVOutPercent()
Get the percentage of transaxial unmasked FOV.
void Reduce()
Merge parallel results into the matrix of the backward image matrix of the first thread. Also for MPI.
FLTNB extent_rotation
FLTNB GetRespBasisCoefficient(int a_respBasisFunction, int a_respGate)
Get the respiratory basis coefficients for the provided respiratory gate and basis function...
int ApplyOutputFOVMasking()
Mask the outside of the transaxial FOV based on the m_fovOutPercent.
#define GEO_ROT_CCW
Definition: vScanner.hh:30
~oImageSpace()
oImageSpace destructor.
Definition: oImageSpace.cc:57
void DeallocateSensitivityImage()
Free memory for the sensitivity image matrices.
Definition: oImageSpace.cc:444
void InstantiateImage()
Allocate memory for the main image matrices.
Definition: oImageSpace.cc:64
int PROJ_SaveProjectionImage()
Save an image of the projected data (for analytic projection)
void InstantiateBackwardImageFromDynamicBasis(int a_nbBackwardImages)
Allocate memory for the backward image matrices and set the number of backward images for the whole c...
Definition: oImageSpace.cc:202
FLTNB **** m4p_anatomicalImage
Definition: oImageSpace.hh:93
int GetNbTimeBasisFunctions()
Get the number of time basis functions.
bool m_loadedSensitivity
Definition: oImageSpace.hh:593
void LMS_InstantiateForwardImage()
Allocate memory for the forward image matrices (for list-mode sensitivity generation) ...
void InstantiateSensitivityImage(const string &a_pathToSensitivityImage)
Allocate the sensitivity image matrices.
Definition: oImageSpace.cc:369
static sOutputManager * GetInstance()
Instanciate the singleton object and Initialize member variables if not already done, return a pointer to this object otherwise.
oImageDimensionsAndQuantification * mp_ID
Definition: oImageSpace.hh:591
void LMS_CopyAtnToImage()
Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m2p_image matrix...
int InitAttenuationImage(const string &a_pathToAtnImage)
Memory allocation and initialisation for the attenuation image using either :
Definition: oImageSpace.cc:966
int InitImage(const string &a_pathToInitialImage, FLTNB a_value)
Initialize the main image, either using:
void DeallocateVisitedVoxelsImage()
Free memory for the image matrix containing binary information regarding which 3D voxels have been vi...
Definition: oImageSpace.cc:909
void SaveDebugImage(const string &a_name)
Just a debug function dedicated to write any kind of image on disk in raw format, for debugging purpo...
void DeallocateAnatomicalImage()
Free memory for the anatomical image.
Definition: oImageSpace.cc:563
int GetSPECTSpecificParameters(uint16_t *ap_nbOfProjections, uint16_t *ap_nbHeads, uint16_t *ap_nbOfBins, FLTNB *ap_pixSizeXY, FLTNB *&ap_angles, FLTNB *&ap_CORtoDetectorDistance, int *ap_headRotDirection)
Transfer geometric information recovered from the datafile to the scanner object. ...
void DeallocateMaskImage()
Free memory for the mask image.
Definition: oImageSpace.cc:634
bool GetRespStaticFlag()
Get the respiratory static flag that says if the reconstruction has only one respiratory gate or not...
int PROJ_LoadInitialImage(const string &a_pathToImage)
Load the initial image for the analytical projection.
FLTNB * mp_visitedVoxelsImage
Definition: oImageSpace.hh:104
void PROJ_DeallocateProjectionImage(int a_nbProjs)
Free memory for the projection image for analytical projection.
void PrepareForwardImage()
Copy current image matrix in the forward-image buffer matrix.
FLTNB GetVoxSizeY()
Get the voxel's size along the Y axis, in mm.
#define Cerr(MESSAGE)
const string & GetPathName()
FLTNB ****** m6p_backwardImage
Definition: oImageSpace.hh:75
void InitBwdImageForDeformation()
Initialize the buffer backward image dedicated to image-based deformation.
bool GetCardStaticFlag()
Get the cardiac static flag that says if the reconstruction has only one cardiac gate or not...
FLTNB **** m4p_image
Definition: oImageSpace.hh:61
void InstantiateBwdImageForDeformation()
Memory allocation for the buffer backward image required for image-based deformation.
Definition: oImageSpace.cc:755
void InstantiateVisitedVoxelsImage()
Memory allocation and initialization for the image matrix containing binary information regarding whi...
Definition: oImageSpace.cc:889
void CleanNeverVisitedVoxels()
Based on the visitedVoxelsImage, clean the never visited voxels in the image. This function must be c...
Singleton class that manages output writing on disk (images, sinograms, etc). It also manages loggi...
bool m_loadedMask
Definition: oImageSpace.hh:595
int LoadAttenuationImage(const string &a_pathToImage)
Load the attenuation image provided by the user in the m2p_attenuation matrix.
void LMS_DeallocateSensitivityImage()
Free memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
int InitAnatomicalImage(const string &a_pathToAnatomicalImage)
Memory allocation and initialization for the anatomical image matrices.
Definition: oImageSpace.cc:510
Singleton class that Instantiate and initialize the scanner object.
int SaveOutputImage(int a_iteration, int a_subset=-1)
Call the interfile function to write output image on disk.
void LMS_DeallocateImage()
Free memory for the main image matrices (for list-mode sensitivity generation)
Declaration of class vOptimizer.
int GetSensNbRespGates()
call the eponym function from the oDynamicDataManager object
void ComputeOutputImage()
int LoadInitialImage(const string &a_pathToImage)
Load the initial image provided by the user in the corresponding matrix.
FLTNB **** m4p_outputImage
Definition: oImageSpace.hh:115
INTNB GetNbVoxXY()
Get the number of voxels in a slice.
FLTNB * mp_maskImage
Definition: oImageSpace.hh:100
uint32_t mtx_size[7]
FLTNB GetTimeBasisCoefficient(int a_timeBasisFunction, int a_timeFrame)
Get the time basis coefficients for the provided frame and basis function.
void InstantiateForwardImage()
Allocate memory for the forward image matrices.
Definition: oImageSpace.cc:133
bool GetTimeStaticFlag()
Get the time static flag that says if the reconstruction has only one frame or not.
uint16_t nb_projections
void PROJ_InstantiateProjectionImage(int a_nbProjs, int a_nbPixels)
Instanciate and initialize projection image for analytical projection.
bool GetFlipOutX()
Get the boolean saying if the output image must be flipped along the X axis.
#define INTNB
Definition: gVariables.hh:64
void ReduceBackwardImage(int a_imageIndex, int a_timeIndex, int a_respIndex, int a_cardIndex)
Merge parallel results into the backward image matrix of the first thread for the specific image / ti...
bool GetFlipOutY()
Get the boolean saying if the output image must be flipped along the Y axis.
const string & GetBaseName()
int ApplyOutputFlip()
Just flip the output image.
This class is designed to manage the image-based deformation part of the reconstruction.
bool m_loadedAnatomical
Definition: oImageSpace.hh:594
void InstantiateBackwardImageFromDynamicBins()
Allocate memory for the backward image matrices and initialize them.
Definition: oImageSpace.cc:290
Declaration of class oImageSpace.
int IntfWriteProjFile(const string &a_pathToImg, FLTNB **a2p_ImgMatrix, oImageDimensionsAndQuantification *ap_ID, Intf_fields a_Imgfields, int vb)
Function dedicated to Interfile image writing for projected data.
void LMS_DeallocateAttenuationImage()
Free memory for the Attenuation image matrices (for analytical projection or list-mode sensitivity ge...
Definition: oImageSpace.cc:928
int GetNbCardGates()
Get the number of cardiac gates.
void DeallocateOutputImage()
Free memory for the Image matrices dedicated to output writing on disk.
Definition: oImageSpace.cc:708
Interfile fields. This structure contains all the Interfile keys currently managed by CASToR Decl...
FLTNB ***** m5p_defTmpBackwardImage
Definition: oImageSpace.hh:124
uint16_t nb_detector_heads
void DeallocateSensImageForDeformation()
Definition: oImageSpace.cc:859
int GetNbThreadsForImageComputation()
Get the number of threads used for image operations.
int GetNbTimeFrames()
Get the number of time frames.
INTNB GetNbVoxXYZ()
Get the total number of voxels.
FLTNB **** m4p_defTmpSensitivityImage
Definition: oImageSpace.hh:135
int GetNbThreadsForProjection()
Get the number of threads used for projections.
oImageSpace()
oImageSpace constructor. Initialize the member variables to their default values. ...
Definition: oImageSpace.cc:28
void IntfKeySetFieldsOutput(Intf_fields *ap_IF, oImageDimensionsAndQuantification *ap_ID)
Init the keys of the Interfile header of an image to be written on disk.
INTNB GetNbVoxX()
Get the number of voxels along the X axis.
string projection_angles
INTNB GetNbVoxZ()
Get the number of voxels along the Z axis.
int GetNbRespGates()
Get the number of respiratory gates.
#define Cout(MESSAGE)
string process_status
FLTNB **** m4p_attenuation
Definition: oImageSpace.hh:108
void LMS_CopyAtnToForwardImage()
Copy the attenuation image contained in the 'm2p_attenuation' matrix inside the m4p_forwardImage matr...
void InstantiateSensImageForDeformation()
Memory allocation for the buffer sensitivity image required for image-based deformation (only for PET...
Definition: oImageSpace.cc:830
int SaveSensitivityImage(const string &a_pathToSensitivityImage)
Call the interfile function to write the sensitivity image on disk.
bool GetFlipOutZ()
Get the boolean saying if the output image must be flipped along the Z axis.
void DeallocateBackwardImageFromDynamicBins()
Free memory of the backward image matrices.
Definition: oImageSpace.cc:326
void InitBackwardImage()
Initialize each voxel of the backward images to 0, also for sensitivity if not loaded (estimated on t...
int InitSensitivityImage(const string &a_pathToSensitivityImage)
Initialization for the sensitivity image matrices.
void LMS_InstantiateSensitivityImage()
Allocate memory for the sensitivity image matrices (for list-mode sensitivity generation) ...
void DeallocateForwardImage()
Free memory for the forward image matrices.
Definition: oImageSpace.cc:169
int IntfReadImage(const string &a_pathToHeaderFile, FLTNB *ap_ImgMatrix, oImageDimensionsAndQuantification *ap_ID, int vb, bool a_lerpFlag)
Main function dedicated to Interfile 3D image loading.
int GetNbRespBasisFunctions()
Get the number of respiratory basis functions.
int IntfWriteImgFile(const string &a_pathToImg, FLTNB *ap_ImgMatrix, const Intf_fields &ap_IntfF, int vb)
Main function dedicated to Interfile 3D image writing. Recover image information from a provided In...
FLTNB ***** m5p_sensitivity
Definition: oImageSpace.hh:85
void LMS_CopyBackwardToSensitivity()
INTNB GetNbSliceOutMask()
Get the number of extrem slices that will be masked at each side.
INTNB GetNbVoxY()
Get the number of voxels along the Y axis.
void InstantiateOutputImage()
Instanciate Image matrices dedicated to output writing on disk.
Definition: oImageSpace.cc:661
void LMS_DeallocateForwardImage()
Free memory for the forward image matrices (for list-mode sensitivity generation) ...
FLTNB GetCardBasisCoefficient(int a_cardBasisFunction, int a_cardGate)
Get the cardiac basis coefficients for the provided cardiac gate and basis function.